1 /*$Id: snesmfj.c,v 1.130 2001/08/21 21:03:56 bsmith Exp bsmith $*/ 2 3 #include "src/snes/mf/snesmfj.h" /*I "petscsnes.h" I*/ 4 #include "src/mat/matimpl.h" 5 6 PetscFList MatSNESMPetscFList = 0; 7 PetscTruth MatSNESMFRegisterAllCalled = PETSC_FALSE; 8 9 #undef __FUNCT__ 10 #define __FUNCT__ "MatSNESMFSetType" 11 /*@C 12 MatSNESMFSetType - Sets the method that is used to compute the 13 differencing parameter for finite differene matrix-free formulations. 14 15 Input Parameters: 16 + mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMF() 17 or MatSetType(mat,MATMFFD); 18 - ftype - the type requested 19 20 Level: advanced 21 22 Notes: 23 For example, such routines can compute h for use in 24 Jacobian-vector products of the form 25 26 F(x+ha) - F(x) 27 F'(u)a ~= ---------------- 28 h 29 30 .seealso: MatCreateSNESMF(), MatSNESMFRegisterDynamic) 31 @*/ 32 int MatSNESMFSetType(Mat mat,MatSNESMFType ftype) 33 { 34 int ierr,(*r)(MatSNESMFCtx); 35 MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data; 36 PetscTruth match; 37 38 PetscFunctionBegin; 39 PetscValidHeaderSpecific(mat,MAT_COOKIE); 40 PetscValidCharPointer(ftype); 41 42 /* already set, so just return */ 43 ierr = PetscTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr); 44 if (match) PetscFunctionReturn(0); 45 46 /* destroy the old one if it exists */ 47 if (ctx->ops->destroy) { 48 ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr); 49 } 50 51 /* Get the function pointers for the requrested method */ 52 if (!MatSNESMFRegisterAllCalled) {ierr = MatSNESMFRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 53 54 ierr = PetscFListFind(ctx->comm,MatSNESMPetscFList,ftype,(void (**)(void)) &r);CHKERRQ(ierr); 55 56 if (!r) SETERRQ(1,"Unknown MatSNESMF type given"); 57 58 ierr = (*r)(ctx);CHKERRQ(ierr); 59 60 ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr); 61 62 PetscFunctionReturn(0); 63 } 64 65 EXTERN_C_BEGIN 66 #undef __FUNCT__ 67 #define __FUNCT__ "MatSNESMFSetFunctioniBase_FD" 68 int MatSNESMFSetFunctioniBase_FD(Mat mat,int (*func)(Vec,void *)) 69 { 70 MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data; 71 72 PetscFunctionBegin; 73 ctx->funcisetbase = func; 74 PetscFunctionReturn(0); 75 } 76 EXTERN_C_END 77 78 EXTERN_C_BEGIN 79 #undef __FUNCT__ 80 #define __FUNCT__ "MatSNESMFSetFunctioni_FD" 81 int MatSNESMFSetFunctioni_FD(Mat mat,int (*funci)(int,Vec,PetscScalar*,void *)) 82 { 83 MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data; 84 85 PetscFunctionBegin; 86 ctx->funci = funci; 87 PetscFunctionReturn(0); 88 } 89 EXTERN_C_END 90 91 /*MC 92 MatSNESMFRegisterDynamic - Adds a method to the MatSNESMF registry. 93 94 Synopsis: 95 int MatSNESMFRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(MatSNESMF)) 96 97 Not Collective 98 99 Input Parameters: 100 + name_solver - name of a new user-defined compute-h module 101 . path - path (either absolute or relative) the library containing this solver 102 . name_create - name of routine to create method context 103 - routine_create - routine to create method context 104 105 Level: developer 106 107 Notes: 108 MatSNESMFRegisterDynamic) may be called multiple times to add several user-defined solvers. 109 110 If dynamic libraries are used, then the fourth input argument (routine_create) 111 is ignored. 112 113 Sample usage: 114 .vb 115 MatSNESMFRegisterDynamic"my_h",/home/username/my_lib/lib/libO/solaris/mylib.a, 116 "MyHCreate",MyHCreate); 117 .ve 118 119 Then, your solver can be chosen with the procedural interface via 120 $ MatSNESMFSetType(mfctx,"my_h") 121 or at runtime via the option 122 $ -snes_mf_type my_h 123 124 .keywords: MatSNESMF, register 125 126 .seealso: MatSNESMFRegisterAll(), MatSNESMFRegisterDestroy() 127 M*/ 128 129 #undef __FUNCT__ 130 #define __FUNCT__ "MatSNESMFRegister" 131 int MatSNESMFRegister(char *sname,char *path,char *name,int (*function)(MatSNESMFCtx)) 132 { 133 int ierr; 134 char fullname[256]; 135 136 PetscFunctionBegin; 137 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 138 ierr = PetscFListAdd(&MatSNESMPetscFList,sname,fullname,(void (*)())function);CHKERRQ(ierr); 139 PetscFunctionReturn(0); 140 } 141 142 143 #undef __FUNCT__ 144 #define __FUNCT__ "MatSNESMFRegisterDestroy" 145 /*@C 146 MatSNESMFRegisterDestroy - Frees the list of MatSNESMF methods that were 147 registered by MatSNESMFRegisterDynamic). 148 149 Not Collective 150 151 Level: developer 152 153 .keywords: MatSNESMF, register, destroy 154 155 .seealso: MatSNESMFRegisterDynamic), MatSNESMFRegisterAll() 156 @*/ 157 int MatSNESMFRegisterDestroy(void) 158 { 159 int ierr; 160 161 PetscFunctionBegin; 162 if (MatSNESMPetscFList) { 163 ierr = PetscFListDestroy(&MatSNESMPetscFList);CHKERRQ(ierr); 164 MatSNESMPetscFList = 0; 165 } 166 MatSNESMFRegisterAllCalled = PETSC_FALSE; 167 PetscFunctionReturn(0); 168 } 169 170 /* ----------------------------------------------------------------------------------------*/ 171 #undef __FUNCT__ 172 #define __FUNCT__ "MatDestroy_MFFD" 173 int MatDestroy_MFFD(Mat mat) 174 { 175 int ierr; 176 MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data; 177 178 PetscFunctionBegin; 179 ierr = VecDestroy(ctx->w);CHKERRQ(ierr); 180 if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);} 181 if (ctx->sp) {ierr = MatNullSpaceDestroy(ctx->sp);CHKERRQ(ierr);} 182 PetscHeaderDestroy(ctx); 183 PetscFunctionReturn(0); 184 } 185 186 #undef __FUNCT__ 187 #define __FUNCT__ "MatView_MFFD" 188 /* 189 MatSNESMFView_MFFD - Views matrix-free parameters. 190 191 */ 192 int MatView_MFFD(Mat J,PetscViewer viewer) 193 { 194 int ierr; 195 MatSNESMFCtx ctx = (MatSNESMFCtx)J->data; 196 PetscTruth isascii; 197 198 PetscFunctionBegin; 199 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 200 if (isascii) { 201 ierr = PetscViewerASCIIPrintf(viewer," SNES matrix-free approximation:\n");CHKERRQ(ierr); 202 ierr = PetscViewerASCIIPrintf(viewer," err=%g (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr); 203 if (!ctx->type_name) { 204 ierr = PetscViewerASCIIPrintf(viewer," The compute h routine has not yet been set\n");CHKERRQ(ierr); 205 } else { 206 ierr = PetscViewerASCIIPrintf(viewer," Using %s compute h routine\n",ctx->type_name);CHKERRQ(ierr); 207 } 208 if (ctx->ops->view) { 209 ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr); 210 } 211 } else { 212 SETERRQ1(1,"Viewer type %s not supported for SNES matrix free matrix",((PetscObject)viewer)->type_name); 213 } 214 PetscFunctionReturn(0); 215 } 216 217 #undef __FUNCT__ 218 #define __FUNCT__ "MatAssemblyEnd_MFFD" 219 /* 220 MatSNESMFAssemblyEnd_Private - Resets the ctx->ncurrenth to zero. This 221 allows the user to indicate the beginning of a new linear solve by calling 222 MatAssemblyXXX() on the matrix free matrix. This then allows the 223 MatSNESMFCreate_WP() to properly compute ||U|| only the first time 224 in the linear solver rather than every time. 225 */ 226 int MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt) 227 { 228 int ierr; 229 MatSNESMFCtx j = (MatSNESMFCtx)J->data; 230 SNESProblemType type; 231 232 PetscFunctionBegin; 233 ierr = MatSNESMFResetHHistory(J);CHKERRQ(ierr); 234 if (j->usesnes) { 235 ierr = SNESGetSolution(j->snes,&j->current_u);CHKERRQ(ierr); 236 ierr = SNESGetProblemType(j->snes,&type);CHKERRQ(ierr); 237 if (type == SNES_NONLINEAR_EQUATIONS) { 238 ierr = SNESGetFunction(j->snes,&j->current_f,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 239 } else if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 240 ierr = SNESGetGradient(j->snes,&j->current_f,PETSC_NULL);CHKERRQ(ierr); 241 } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid method class"); 242 } 243 j->vshift = 0.0; 244 j->vscale = 1.0; 245 PetscFunctionReturn(0); 246 } 247 248 #undef __FUNCT__ 249 #define __FUNCT__ "MatMult_MFFD" 250 /* 251 MatSNESMFMult_Private - Default matrix-free form for Jacobian-vector 252 product, y = F'(u)*a: 253 254 y ~= (F(u + ha) - F(u))/h, 255 where F = nonlinear function, as set by SNESSetFunction() 256 u = current iterate 257 h = difference interval 258 */ 259 int MatMult_MFFD(Mat mat,Vec a,Vec y) 260 { 261 MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data; 262 SNES snes; 263 PetscScalar h,mone = -1.0; 264 Vec w,U,F; 265 int ierr,(*eval_fct)(SNES,Vec,Vec)=0; 266 SNESProblemType type; 267 268 PetscFunctionBegin; 269 /* We log matrix-free matrix-vector products separately, so that we can 270 separate the performance monitoring from the cases that use conventional 271 storage. We may eventually modify event logging to associate events 272 with particular objects, hence alleviating the more general problem. */ 273 ierr = PetscLogEventBegin(MAT_MatrixFreeMult,a,y,0,0);CHKERRQ(ierr); 274 275 snes = ctx->snes; 276 w = ctx->w; 277 U = ctx->current_u; 278 279 /* 280 Compute differencing parameter 281 */ 282 if (!ctx->ops->compute) { 283 ierr = MatSNESMFSetType(mat,MATSNESMF_DEFAULT);CHKERRQ(ierr); 284 ierr = MatSNESMFSetFromOptions(mat);CHKERRQ(ierr); 285 } 286 ierr = (*ctx->ops->compute)(ctx,U,a,&h);CHKERRQ(ierr); 287 288 /* keep a record of the current differencing parameter h */ 289 ctx->currenth = h; 290 #if defined(PETSC_USE_COMPLEX) 291 PetscLogInfo(mat,"MatMult_MFFD:Current differencing parameter: %g + %g i\n",PetscRealPart(h),PetscImaginaryPart(h)); 292 #else 293 PetscLogInfo(mat,"MatMult_MFFD:Current differencing parameter: %15.12e\n",h); 294 #endif 295 if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 296 ctx->historyh[ctx->ncurrenth] = h; 297 } 298 ctx->ncurrenth++; 299 300 /* w = u + ha */ 301 ierr = VecWAXPY(&h,a,U,w);CHKERRQ(ierr); 302 303 if (ctx->usesnes) { 304 ierr = SNESGetProblemType(snes,&type);CHKERRQ(ierr); 305 if (type == SNES_NONLINEAR_EQUATIONS) { 306 eval_fct = SNESComputeFunction; 307 } else if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 308 eval_fct = SNESComputeGradient; 309 } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid method class"); 310 F = ctx->current_f; 311 if (!F) SETERRQ(1,"You must call MatAssembly() even on matrix-free matrices"); 312 ierr = (*eval_fct)(snes,w,y);CHKERRQ(ierr); 313 } else { 314 F = ctx->funcvec; 315 /* compute func(U) as base for differencing */ 316 if (ctx->ncurrenth == 1) { 317 ierr = (*ctx->func)(snes,U,F,ctx->funcctx);CHKERRQ(ierr); 318 } 319 ierr = (*ctx->func)(snes,w,y,ctx->funcctx);CHKERRQ(ierr); 320 } 321 322 ierr = VecAXPY(&mone,F,y);CHKERRQ(ierr); 323 h = 1.0/h; 324 ierr = VecScale(&h,y);CHKERRQ(ierr); 325 326 327 if (ctx->vshift != 0.0 && ctx->vscale != 1.0) { 328 ierr = VecAXPBY(&ctx->vshift,&ctx->vscale,a,y);CHKERRQ(ierr); 329 } else if (ctx->vscale != 1.0) { 330 ierr = VecScale(&ctx->vscale,y);CHKERRQ(ierr); 331 } else if (ctx->vshift != 0.0) { 332 ierr = VecAXPY(&ctx->vshift,a,y);CHKERRQ(ierr); 333 } 334 335 if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);CHKERRQ(ierr);} 336 337 ierr = PetscLogEventEnd(MAT_MatrixFreeMult,a,y,0,0);CHKERRQ(ierr); 338 PetscFunctionReturn(0); 339 } 340 341 #undef __FUNCT__ 342 #define __FUNCT__ "MatGetDiagonal_MFFD" 343 /* 344 MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix 345 346 y ~= (F(u + ha) - F(u))/h, 347 where F = nonlinear function, as set by SNESSetFunction() 348 u = current iterate 349 h = difference interval 350 */ 351 int MatGetDiagonal_MFFD(Mat mat,Vec a) 352 { 353 MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data; 354 PetscScalar h,*aa,*ww,v; 355 PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 356 Vec w,U; 357 int i,ierr,rstart,rend; 358 359 PetscFunctionBegin; 360 if (!ctx->funci) { 361 SETERRQ(1,"Requirers calling MatSNESMFSetFunctioni() first"); 362 } 363 364 w = ctx->w; 365 U = ctx->current_u; 366 ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr); 367 ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr); 368 ierr = VecCopy(U,w);CHKERRQ(ierr); 369 370 ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr); 371 ierr = VecGetArray(a,&aa);CHKERRQ(ierr); 372 for (i=rstart; i<rend; i++) { 373 ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 374 h = ww[i-rstart]; 375 if (h == 0.0) h = 1.0; 376 #if !defined(PETSC_USE_COMPLEX) 377 if (h < umin && h >= 0.0) h = umin; 378 else if (h < 0.0 && h > -umin) h = -umin; 379 #else 380 if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 381 else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 382 #endif 383 h *= epsilon; 384 385 ww[i-rstart] += h; 386 ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 387 ierr = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr); 388 aa[i-rstart] = (v - aa[i-rstart])/h; 389 390 /* possibly shift and scale result */ 391 aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart]; 392 393 ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 394 ww[i-rstart] -= h; 395 ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 396 } 397 ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr); 398 PetscFunctionReturn(0); 399 } 400 401 #undef __FUNCT__ 402 #define __FUNCT__ "MatShift_MFFD" 403 int MatShift_MFFD(PetscScalar *a,Mat Y) 404 { 405 MatSNESMFCtx shell = (MatSNESMFCtx)Y->data; 406 PetscFunctionBegin; 407 shell->vshift += *a; 408 PetscFunctionReturn(0); 409 } 410 411 #undef __FUNCT__ 412 #define __FUNCT__ "MatScale_MFFD" 413 int MatScale_MFFD(PetscScalar *a,Mat Y) 414 { 415 MatSNESMFCtx shell = (MatSNESMFCtx)Y->data; 416 PetscFunctionBegin; 417 shell->vscale *= *a; 418 PetscFunctionReturn(0); 419 } 420 421 422 #undef __FUNCT__ 423 #define __FUNCT__ "MatCreateSNESMF" 424 /*@C 425 MatCreateSNESMF - Creates a matrix-free matrix context for use with 426 a SNES solver. This matrix can be used as the Jacobian argument for 427 the routine SNESSetJacobian(). 428 429 Collective on SNES and Vec 430 431 Input Parameters: 432 + snes - the SNES context 433 - x - vector where SNES solution is to be stored. 434 435 Output Parameter: 436 . J - the matrix-free matrix 437 438 Level: advanced 439 440 Notes: 441 The matrix-free matrix context merely contains the function pointers 442 and work space for performing finite difference approximations of 443 Jacobian-vector products, F'(u)*a, 444 445 The default code uses the following approach to compute h 446 447 .vb 448 F'(u)*a = [F(u+h*a) - F(u)]/h where 449 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 450 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 451 where 452 error_rel = square root of relative error in function evaluation 453 umin = minimum iterate parameter 454 .ve 455 456 The user can set the error_rel via MatSNESMFSetFunctionError() and 457 umin via MatSNESMFDefaultSetUmin(); see the nonlinear solvers chapter 458 of the users manual for details. 459 460 The user should call MatDestroy() when finished with the matrix-free 461 matrix context. 462 463 Options Database Keys: 464 + -snes_mf_err <error_rel> - Sets error_rel 465 . -snes_mf_unim <umin> - Sets umin (for default PETSc routine that computes h only) 466 - -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h 467 468 .keywords: SNES, default, matrix-free, create, matrix 469 470 .seealso: MatDestroy(), MatSNESMFSetFunctionError(), MatSNESMFDefaultSetUmin() 471 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), MatCreateMF(), 472 MatSNESMFGetH(),MatSNESMFKSPMonitor(), MatSNESMFRegisterDynamic), MatSNESMFComputeJacobian() 473 474 @*/ 475 int MatCreateSNESMF(SNES snes,Vec x,Mat *J) 476 { 477 MatSNESMFCtx mfctx; 478 int ierr; 479 480 PetscFunctionBegin; 481 ierr = MatCreateMF(x,J);CHKERRQ(ierr); 482 483 mfctx = (MatSNESMFCtx)(*J)->data; 484 mfctx->snes = snes; 485 mfctx->usesnes = PETSC_TRUE; 486 PetscLogObjectParent(snes,*J); 487 PetscFunctionReturn(0); 488 } 489 490 EXTERN_C_BEGIN 491 #undef __FUNCT__ 492 #define __FUNCT__ "MatSNESMFSetBase_FD" 493 int MatSNESMFSetBase_FD(Mat J,Vec U) 494 { 495 int ierr; 496 MatSNESMFCtx ctx = (MatSNESMFCtx)J->data; 497 498 PetscFunctionBegin; 499 ierr = MatSNESMFResetHHistory(J);CHKERRQ(ierr); 500 ctx->current_u = U; 501 ctx->usesnes = PETSC_FALSE; 502 PetscFunctionReturn(0); 503 } 504 EXTERN_C_END 505 506 #undef __FUNCT__ 507 #define __FUNCT__ "MatSNESMFSetFromOptions" 508 /*@ 509 MatSNESMFSetFromOptions - Sets the MatSNESMF options from the command line 510 parameter. 511 512 Collective on Mat 513 514 Input Parameters: 515 . mat - the matrix obtained with MatCreateSNESMF() 516 517 Options Database Keys: 518 + -snes_mf_type - <default,wp> 519 - -snes_mf_err - square root of estimated relative error in function evaluation 520 - -snes_mf_period - how often h is recomputed, defaults to 1, everytime 521 522 Level: advanced 523 524 .keywords: SNES, matrix-free, parameters 525 526 .seealso: MatCreateSNESMF(),MatSNESMFSetHHistory(), 527 MatSNESMFResetHHistory(), MatSNESMFKSPMonitor() 528 @*/ 529 int MatSNESMFSetFromOptions(Mat mat) 530 { 531 MatSNESMFCtx mfctx = (MatSNESMFCtx)mat->data; 532 int ierr; 533 PetscTruth flg; 534 char ftype[256]; 535 536 PetscFunctionBegin; 537 if (!MatSNESMFRegisterAllCalled) {ierr = MatSNESMFRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 538 539 ierr = PetscOptionsBegin(mfctx->comm,mfctx->prefix,"Set matrix free computation parameters","MatSNESMF");CHKERRQ(ierr); 540 ierr = PetscOptionsList("-snes_mf_type","Matrix free type","MatSNESMFSetType",MatSNESMPetscFList,mfctx->type_name,ftype,256,&flg);CHKERRQ(ierr); 541 if (flg) { 542 ierr = MatSNESMFSetType(mat,ftype);CHKERRQ(ierr); 543 } 544 545 ierr = PetscOptionsReal("-snes_mf_err","set sqrt relative error in function","MatSNESMFSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr); 546 ierr = PetscOptionsInt("-snes_mf_period","how often h is recomputed","MatSNESMFSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr); 547 if (mfctx->snes) { 548 ierr = PetscOptionsName("-snes_mf_ksp_monitor","Monitor matrix-free parameters","MatSNESMFKSPMonitor",&flg);CHKERRQ(ierr); 549 if (flg) { 550 SLES sles; 551 KSP ksp; 552 ierr = SNESGetSLES(mfctx->snes,&sles);CHKERRQ(ierr); 553 ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr); 554 ierr = KSPSetMonitor(ksp,MatSNESMFKSPMonitor,PETSC_NULL,0);CHKERRQ(ierr); 555 } 556 } 557 if (mfctx->ops->setfromoptions) { 558 ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr); 559 } 560 ierr = PetscOptionsEnd();CHKERRQ(ierr); 561 PetscFunctionReturn(0); 562 } 563 564 #undef __FUNCT__ 565 #define __FUNCT__ "MatCreate_MFFD" 566 EXTERN_C_BEGIN 567 int MatCreate_MFFD(Mat A) 568 { 569 MatSNESMFCtx mfctx; 570 int ierr; 571 572 PetscFunctionBegin; 573 PetscHeaderCreate(mfctx,_p_MatSNESMFCtx,struct _MFOps,MATSNESMFCTX_COOKIE,0,"SNESMF",A->comm,MatDestroy_MFFD,MatView_MFFD); 574 PetscLogObjectCreate(mfctx); 575 mfctx->sp = 0; 576 mfctx->snes = 0; 577 mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 578 mfctx->recomputeperiod = 1; 579 mfctx->count = 0; 580 mfctx->currenth = 0.0; 581 mfctx->historyh = PETSC_NULL; 582 mfctx->ncurrenth = 0; 583 mfctx->maxcurrenth = 0; 584 mfctx->type_name = 0; 585 mfctx->usesnes = PETSC_FALSE; 586 587 mfctx->vshift = 0.0; 588 mfctx->vscale = 1.0; 589 590 /* 591 Create the empty data structure to contain compute-h routines. 592 These will be filled in below from the command line options or 593 a later call with MatSNESMFSetType() or if that is not called 594 then it will default in the first use of MatMult_MFFD() 595 */ 596 mfctx->ops->compute = 0; 597 mfctx->ops->destroy = 0; 598 mfctx->ops->view = 0; 599 mfctx->ops->setfromoptions = 0; 600 mfctx->hctx = 0; 601 602 mfctx->func = 0; 603 mfctx->funcctx = 0; 604 mfctx->funcvec = 0; 605 606 A->data = mfctx; 607 608 A->ops->mult = MatMult_MFFD; 609 A->ops->destroy = MatDestroy_MFFD; 610 A->ops->view = MatView_MFFD; 611 A->ops->assemblyend = MatAssemblyEnd_MFFD; 612 A->ops->getdiagonal = MatGetDiagonal_MFFD; 613 A->ops->scale = MatScale_MFFD; 614 A->ops->shift = MatShift_MFFD; 615 A->ops->setfromoptions = MatSNESMFSetFromOptions; 616 617 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetBase_C","MatSNESMFSetBase_FD",MatSNESMFSetBase_FD);CHKERRQ(ierr); 618 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetFunctioniBase_C","MatSNESMFSetFunctioniBase_FD",MatSNESMFSetFunctioniBase_FD);CHKERRQ(ierr); 619 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetFunctioni_C","MatSNESMFSetFunctioni_FD",MatSNESMFSetFunctioni_FD);CHKERRQ(ierr); 620 mfctx->mat = A; 621 ierr = VecCreateMPI(A->comm,A->n,A->N,&mfctx->w);CHKERRQ(ierr); 622 623 PetscFunctionReturn(0); 624 } 625 626 EXTERN_C_END 627 628 #undef __FUNCT__ 629 #define __FUNCT__ "MatCreateMF" 630 /*@C 631 MatCreateMF - Creates a matrix-free matrix. See also MatCreateSNESMF() 632 633 Collective on Vec 634 635 Input Parameters: 636 . x - vector that defines layout of the vectors and matrices 637 638 Output Parameter: 639 . J - the matrix-free matrix 640 641 Level: advanced 642 643 Notes: 644 The matrix-free matrix context merely contains the function pointers 645 and work space for performing finite difference approximations of 646 Jacobian-vector products, F'(u)*a, 647 648 The default code uses the following approach to compute h 649 650 .vb 651 F'(u)*a = [F(u+h*a) - F(u)]/h where 652 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 653 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 654 where 655 error_rel = square root of relative error in function evaluation 656 umin = minimum iterate parameter 657 .ve 658 659 The user can set the error_rel via MatSNESMFSetFunctionError() and 660 umin via MatSNESMFDefaultSetUmin(); see the nonlinear solvers chapter 661 of the users manual for details. 662 663 The user should call MatDestroy() when finished with the matrix-free 664 matrix context. 665 666 Options Database Keys: 667 + -snes_mf_err <error_rel> - Sets error_rel 668 . -snes_mf_unim <umin> - Sets umin (for default PETSc routine that computes h only) 669 - -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h 670 671 .keywords: default, matrix-free, create, matrix 672 673 .seealso: MatDestroy(), MatSNESMFSetFunctionError(), MatSNESMFDefaultSetUmin() 674 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), MatCreateSNESMF(), 675 MatSNESMFGetH(),MatSNESMFKSPMonitor(), MatSNESMFRegisterDynamic),, MatSNESMFComputeJacobian() 676 677 @*/ 678 int MatCreateMF(Vec x,Mat *J) 679 { 680 MPI_Comm comm; 681 int n,nloc,ierr; 682 683 PetscFunctionBegin; 684 ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr); 685 ierr = VecGetSize(x,&n);CHKERRQ(ierr); 686 ierr = VecGetLocalSize(x,&nloc);CHKERRQ(ierr); 687 ierr = MatCreate(comm,nloc,nloc,n,n,J);CHKERRQ(ierr); 688 ierr = MatRegister(MATMFFD,0,"MatCreate_MFFD",MatCreate_MFFD);CHKERRQ(ierr); 689 ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr); 690 PetscFunctionReturn(0); 691 } 692 693 694 #undef __FUNCT__ 695 #define __FUNCT__ "MatSNESMFGetH" 696 /*@ 697 MatSNESMFGetH - Gets the last value that was used as the differencing 698 parameter. 699 700 Not Collective 701 702 Input Parameters: 703 . mat - the matrix obtained with MatCreateSNESMF() 704 705 Output Paramter: 706 . h - the differencing step size 707 708 Level: advanced 709 710 .keywords: SNES, matrix-free, parameters 711 712 .seealso: MatCreateSNESMF(),MatSNESMFSetHHistory(), 713 MatSNESMFResetHHistory(),MatSNESMFKSPMonitor() 714 @*/ 715 int MatSNESMFGetH(Mat mat,PetscScalar *h) 716 { 717 MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data; 718 719 PetscFunctionBegin; 720 *h = ctx->currenth; 721 PetscFunctionReturn(0); 722 } 723 724 #undef __FUNCT__ 725 #define __FUNCT__ "MatSNESMFKSPMonitor" 726 /* 727 MatSNESMFKSPMonitor - A KSP monitor for use with the default PETSc 728 SNES matrix free routines. Prints the differencing parameter used at 729 each step. 730 */ 731 int MatSNESMFKSPMonitor(KSP ksp,int n,PetscReal rnorm,void *dummy) 732 { 733 PC pc; 734 MatSNESMFCtx ctx; 735 int ierr; 736 Mat mat; 737 MPI_Comm comm; 738 PetscTruth nonzeroinitialguess; 739 740 PetscFunctionBegin; 741 ierr = PetscObjectGetComm((PetscObject)ksp,&comm);CHKERRQ(ierr); 742 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 743 ierr = KSPGetInitialGuessNonzero(ksp,&nonzeroinitialguess);CHKERRQ(ierr); 744 ierr = PCGetOperators(pc,&mat,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 745 ctx = (MatSNESMFCtx)mat->data; 746 747 if (n > 0 || nonzeroinitialguess) { 748 #if defined(PETSC_USE_COMPLEX) 749 ierr = PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g + %g i\n",n,rnorm, 750 PetscRealPart(ctx->currenth),PetscImaginaryPart(ctx->currenth));CHKERRQ(ierr); 751 #else 752 ierr = PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g \n",n,rnorm,ctx->currenth);CHKERRQ(ierr); 753 #endif 754 } else { 755 ierr = PetscPrintf(comm,"%d KSP Residual norm %14.12e\n",n,rnorm);CHKERRQ(ierr); 756 } 757 PetscFunctionReturn(0); 758 } 759 760 #undef __FUNCT__ 761 #define __FUNCT__ "MatSNESMFSetFunction" 762 /*@C 763 MatSNESMFSetFunction - Sets the function used in applying the matrix free. 764 765 Collective on Mat 766 767 Input Parameters: 768 + mat - the matrix free matrix created via MatCreateSNESMF() 769 . v - workspace vector 770 . func - the function to use 771 - funcctx - optional function context passed to function 772 773 Level: advanced 774 775 Notes: 776 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 777 matrix inside your compute Jacobian routine 778 779 If this is not set then it will use the function set with SNESSetFunction() 780 781 .keywords: SNES, matrix-free, function 782 783 .seealso: MatCreateSNESMF(),MatSNESMFGetH(), 784 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 785 MatSNESMFKSPMonitor(), SNESetFunction() 786 @*/ 787 int MatSNESMFSetFunction(Mat mat,Vec v,int (*func)(SNES,Vec,Vec,void *),void *funcctx) 788 { 789 MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data; 790 791 PetscFunctionBegin; 792 ctx->func = func; 793 ctx->funcctx = funcctx; 794 ctx->funcvec = v; 795 PetscFunctionReturn(0); 796 } 797 798 #undef __FUNCT__ 799 #define __FUNCT__ "MatSNESMFSetFunctioni" 800 /*@C 801 MatSNESMFSetFunctioni - Sets the function for a single component 802 803 Collective on Mat 804 805 Input Parameters: 806 + mat - the matrix free matrix created via MatCreateSNESMF() 807 - funci - the function to use 808 809 Level: advanced 810 811 Notes: 812 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 813 matrix inside your compute Jacobian routine 814 815 816 .keywords: SNES, matrix-free, function 817 818 .seealso: MatCreateSNESMF(),MatSNESMFGetH(), 819 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 820 MatSNESMFKSPMonitor(), SNESetFunction() 821 @*/ 822 int MatSNESMFSetFunctioni(Mat mat,int (*funci)(int,Vec,PetscScalar*,void *)) 823 { 824 int ierr,(*f)(Mat,int (*)(int,Vec,PetscScalar*,void *)); 825 826 PetscFunctionBegin; 827 PetscValidHeaderSpecific(mat,MAT_COOKIE); 828 ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSNESMFSetFunctioni_C",(void (**)())&f);CHKERRQ(ierr); 829 if (f) { 830 ierr = (*f)(mat,funci);CHKERRQ(ierr); 831 } 832 PetscFunctionReturn(0); 833 } 834 835 836 #undef __FUNCT__ 837 #define __FUNCT__ "MatSNESMFSetFunctioniBase" 838 /*@C 839 MatSNESMFSetFunctioniBase - Sets the base vector for a single component function evaluation 840 841 Collective on Mat 842 843 Input Parameters: 844 + mat - the matrix free matrix created via MatCreateSNESMF() 845 - func - the function to use 846 847 Level: advanced 848 849 Notes: 850 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 851 matrix inside your compute Jacobian routine 852 853 854 .keywords: SNES, matrix-free, function 855 856 .seealso: MatCreateSNESMF(),MatSNESMFGetH(), 857 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 858 MatSNESMFKSPMonitor(), SNESetFunction() 859 @*/ 860 int MatSNESMFSetFunctioniBase(Mat mat,int (*func)(Vec,void *)) 861 { 862 int ierr,(*f)(Mat,int (*)(Vec,void *)); 863 864 PetscFunctionBegin; 865 PetscValidHeaderSpecific(mat,MAT_COOKIE); 866 ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSNESMFSetFunctioniBase_C",(void (**)())&f);CHKERRQ(ierr); 867 if (f) { 868 ierr = (*f)(mat,func);CHKERRQ(ierr); 869 } 870 PetscFunctionReturn(0); 871 } 872 873 874 #undef __FUNCT__ 875 #define __FUNCT__ "MatSNESMFSetPeriod" 876 /*@ 877 MatSNESMFSetPeriod - Sets how often h is recomputed, by default it is everytime 878 879 Collective on Mat 880 881 Input Parameters: 882 + mat - the matrix free matrix created via MatCreateSNESMF() 883 - period - 1 for everytime, 2 for every second etc 884 885 Options Database Keys: 886 + -snes_mf_period <period> 887 888 Level: advanced 889 890 891 .keywords: SNES, matrix-free, parameters 892 893 .seealso: MatCreateSNESMF(),MatSNESMFGetH(), 894 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 895 MatSNESMFKSPMonitor() 896 @*/ 897 int MatSNESMFSetPeriod(Mat mat,int period) 898 { 899 MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data; 900 901 PetscFunctionBegin; 902 ctx->recomputeperiod = period; 903 PetscFunctionReturn(0); 904 } 905 906 #undef __FUNCT__ 907 #define __FUNCT__ "MatSNESMFSetFunctionError" 908 /*@ 909 MatSNESMFSetFunctionError - Sets the error_rel for the approximation of 910 matrix-vector products using finite differences. 911 912 Collective on Mat 913 914 Input Parameters: 915 + mat - the matrix free matrix created via MatCreateSNESMF() 916 - error_rel - relative error (should be set to the square root of 917 the relative error in the function evaluations) 918 919 Options Database Keys: 920 + -snes_mf_err <error_rel> - Sets error_rel 921 922 Level: advanced 923 924 Notes: 925 The default matrix-free matrix-vector product routine computes 926 .vb 927 F'(u)*a = [F(u+h*a) - F(u)]/h where 928 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 929 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 930 .ve 931 932 .keywords: SNES, matrix-free, parameters 933 934 .seealso: MatCreateSNESMF(),MatSNESMFGetH(), 935 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 936 MatSNESMFKSPMonitor() 937 @*/ 938 int MatSNESMFSetFunctionError(Mat mat,PetscReal error) 939 { 940 MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data; 941 942 PetscFunctionBegin; 943 if (error != PETSC_DEFAULT) ctx->error_rel = error; 944 PetscFunctionReturn(0); 945 } 946 947 #undef __FUNCT__ 948 #define __FUNCT__ "MatSNESMFAddNullSpace" 949 /*@ 950 MatSNESMFAddNullSpace - Provides a null space that an operator is 951 supposed to have. Since roundoff will create a small component in 952 the null space, if you know the null space you may have it 953 automatically removed. 954 955 Collective on Mat 956 957 Input Parameters: 958 + J - the matrix-free matrix context 959 - nullsp - object created with MatNullSpaceCreate() 960 961 Level: advanced 962 963 .keywords: SNES, matrix-free, null space 964 965 .seealso: MatNullSpaceCreate(), MatSNESMFGetH(), MatCreateSNESMF(), 966 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 967 MatSNESMFKSPMonitor(), MatSNESMFErrorRel() 968 @*/ 969 int MatSNESMFAddNullSpace(Mat J,MatNullSpace nullsp) 970 { 971 int ierr; 972 MatSNESMFCtx ctx = (MatSNESMFCtx)J->data; 973 MPI_Comm comm; 974 975 PetscFunctionBegin; 976 ierr = PetscObjectGetComm((PetscObject)J,&comm);CHKERRQ(ierr); 977 978 ctx->sp = nullsp; 979 ierr = PetscObjectReference((PetscObject)nullsp);CHKERRQ(ierr); 980 PetscFunctionReturn(0); 981 } 982 983 #undef __FUNCT__ 984 #define __FUNCT__ "MatSNESMFSetHHistory" 985 /*@ 986 MatSNESMFSetHHistory - Sets an array to collect a history of the 987 differencing values (h) computed for the matrix-free product. 988 989 Collective on Mat 990 991 Input Parameters: 992 + J - the matrix-free matrix context 993 . histroy - space to hold the history 994 - nhistory - number of entries in history, if more entries are generated than 995 nhistory, then the later ones are discarded 996 997 Level: advanced 998 999 Notes: 1000 Use MatSNESMFResetHHistory() to reset the history counter and collect 1001 a new batch of differencing parameters, h. 1002 1003 .keywords: SNES, matrix-free, h history, differencing history 1004 1005 .seealso: MatSNESMFGetH(), MatCreateSNESMF(), 1006 MatSNESMFResetHHistory(), 1007 MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError() 1008 1009 @*/ 1010 int MatSNESMFSetHHistory(Mat J,PetscScalar *history,int nhistory) 1011 { 1012 MatSNESMFCtx ctx = (MatSNESMFCtx)J->data; 1013 1014 PetscFunctionBegin; 1015 ctx->historyh = history; 1016 ctx->maxcurrenth = nhistory; 1017 ctx->currenth = 0; 1018 PetscFunctionReturn(0); 1019 } 1020 1021 #undef __FUNCT__ 1022 #define __FUNCT__ "MatSNESMFResetHHistory" 1023 /*@ 1024 MatSNESMFResetHHistory - Resets the counter to zero to begin 1025 collecting a new set of differencing histories. 1026 1027 Collective on Mat 1028 1029 Input Parameters: 1030 . J - the matrix-free matrix context 1031 1032 Level: advanced 1033 1034 Notes: 1035 Use MatSNESMFSetHHistory() to create the original history counter. 1036 1037 .keywords: SNES, matrix-free, h history, differencing history 1038 1039 .seealso: MatSNESMFGetH(), MatCreateSNESMF(), 1040 MatSNESMFSetHHistory(), 1041 MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError() 1042 1043 @*/ 1044 int MatSNESMFResetHHistory(Mat J) 1045 { 1046 MatSNESMFCtx ctx = (MatSNESMFCtx)J->data; 1047 1048 PetscFunctionBegin; 1049 ctx->ncurrenth = 0; 1050 PetscFunctionReturn(0); 1051 } 1052 1053 #undef __FUNCT__ 1054 #define __FUNCT__ "MatSNESMFComputeJacobian" 1055 int MatSNESMFComputeJacobian(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure *flag,void *dummy) 1056 { 1057 int ierr; 1058 PetscFunctionBegin; 1059 ierr = MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1060 ierr = MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1061 PetscFunctionReturn(0); 1062 } 1063 1064 #undef __FUNCT__ 1065 #define __FUNCT__ "MatSNESMFSetBase" 1066 int MatSNESMFSetBase(Mat J,Vec U) 1067 { 1068 int ierr,(*f)(Mat,Vec); 1069 1070 PetscFunctionBegin; 1071 PetscValidHeaderSpecific(J,MAT_COOKIE); 1072 PetscValidHeaderSpecific(U,VEC_COOKIE); 1073 ierr = PetscObjectQueryFunction((PetscObject)J,"MatSNESMFSetBase_C",(void (**)())&f);CHKERRQ(ierr); 1074 if (f) { 1075 ierr = (*f)(J,U);CHKERRQ(ierr); 1076 } 1077 PetscFunctionReturn(0); 1078 } 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088