1 /*$Id: snesmfj.c,v 1.131 2001/09/05 18:45:40 bsmith Exp $*/ 2 3 #include "src/mat/matimpl.h" 4 #include "src/snes/mf/snesmfj.h" /*I "petscsnes.h" I*/ 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 (*)(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 231 PetscFunctionBegin; 232 ierr = MatSNESMFResetHHistory(J);CHKERRQ(ierr); 233 if (j->usesnes) { 234 ierr = SNESGetSolution(j->snes,&j->current_u);CHKERRQ(ierr); 235 ierr = SNESGetFunction(j->snes,&j->current_f,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 236 } 237 j->vshift = 0.0; 238 j->vscale = 1.0; 239 PetscFunctionReturn(0); 240 } 241 242 #undef __FUNCT__ 243 #define __FUNCT__ "MatMult_MFFD" 244 /* 245 MatSNESMFMult_Private - Default matrix-free form for Jacobian-vector 246 product, y = F'(u)*a: 247 248 y ~= (F(u + ha) - F(u))/h, 249 where F = nonlinear function, as set by SNESSetFunction() 250 u = current iterate 251 h = difference interval 252 */ 253 int MatMult_MFFD(Mat mat,Vec a,Vec y) 254 { 255 MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data; 256 SNES snes; 257 PetscScalar h,mone = -1.0; 258 Vec w,U,F; 259 int ierr,(*eval_fct)(SNES,Vec,Vec)=0; 260 261 PetscFunctionBegin; 262 /* We log matrix-free matrix-vector products separately, so that we can 263 separate the performance monitoring from the cases that use conventional 264 storage. We may eventually modify event logging to associate events 265 with particular objects, hence alleviating the more general problem. */ 266 ierr = PetscLogEventBegin(MAT_MultMatrixFree,a,y,0,0);CHKERRQ(ierr); 267 268 snes = ctx->snes; 269 w = ctx->w; 270 U = ctx->current_u; 271 272 /* 273 Compute differencing parameter 274 */ 275 if (!ctx->ops->compute) { 276 ierr = MatSNESMFSetType(mat,MATSNESMF_DEFAULT);CHKERRQ(ierr); 277 ierr = MatSNESMFSetFromOptions(mat);CHKERRQ(ierr); 278 } 279 ierr = (*ctx->ops->compute)(ctx,U,a,&h);CHKERRQ(ierr); 280 281 /* keep a record of the current differencing parameter h */ 282 ctx->currenth = h; 283 #if defined(PETSC_USE_COMPLEX) 284 PetscLogInfo(mat,"MatMult_MFFD:Current differencing parameter: %g + %g i\n",PetscRealPart(h),PetscImaginaryPart(h)); 285 #else 286 PetscLogInfo(mat,"MatMult_MFFD:Current differencing parameter: %15.12e\n",h); 287 #endif 288 if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 289 ctx->historyh[ctx->ncurrenth] = h; 290 } 291 ctx->ncurrenth++; 292 293 /* w = u + ha */ 294 ierr = VecWAXPY(&h,a,U,w);CHKERRQ(ierr); 295 296 if (ctx->usesnes) { 297 eval_fct = SNESComputeFunction; 298 F = ctx->current_f; 299 if (!F) SETERRQ(1,"You must call MatAssembly() even on matrix-free matrices"); 300 ierr = (*eval_fct)(snes,w,y);CHKERRQ(ierr); 301 } else { 302 F = ctx->funcvec; 303 /* compute func(U) as base for differencing */ 304 if (ctx->ncurrenth == 1) { 305 ierr = (*ctx->func)(snes,U,F,ctx->funcctx);CHKERRQ(ierr); 306 } 307 ierr = (*ctx->func)(snes,w,y,ctx->funcctx);CHKERRQ(ierr); 308 } 309 310 ierr = VecAXPY(&mone,F,y);CHKERRQ(ierr); 311 h = 1.0/h; 312 ierr = VecScale(&h,y);CHKERRQ(ierr); 313 314 315 if (ctx->vshift != 0.0 && ctx->vscale != 1.0) { 316 ierr = VecAXPBY(&ctx->vshift,&ctx->vscale,a,y);CHKERRQ(ierr); 317 } else if (ctx->vscale != 1.0) { 318 ierr = VecScale(&ctx->vscale,y);CHKERRQ(ierr); 319 } else if (ctx->vshift != 0.0) { 320 ierr = VecAXPY(&ctx->vshift,a,y);CHKERRQ(ierr); 321 } 322 323 if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);CHKERRQ(ierr);} 324 325 ierr = PetscLogEventEnd(MAT_MultMatrixFree,a,y,0,0);CHKERRQ(ierr); 326 PetscFunctionReturn(0); 327 } 328 329 #undef __FUNCT__ 330 #define __FUNCT__ "MatGetDiagonal_MFFD" 331 /* 332 MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix 333 334 y ~= (F(u + ha) - F(u))/h, 335 where F = nonlinear function, as set by SNESSetFunction() 336 u = current iterate 337 h = difference interval 338 */ 339 int MatGetDiagonal_MFFD(Mat mat,Vec a) 340 { 341 MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data; 342 PetscScalar h,*aa,*ww,v; 343 PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 344 Vec w,U; 345 int i,ierr,rstart,rend; 346 347 PetscFunctionBegin; 348 if (!ctx->funci) { 349 SETERRQ(1,"Requirers calling MatSNESMFSetFunctioni() first"); 350 } 351 352 w = ctx->w; 353 U = ctx->current_u; 354 ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr); 355 ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr); 356 ierr = VecCopy(U,w);CHKERRQ(ierr); 357 358 ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr); 359 ierr = VecGetArray(a,&aa);CHKERRQ(ierr); 360 for (i=rstart; i<rend; i++) { 361 ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 362 h = ww[i-rstart]; 363 if (h == 0.0) h = 1.0; 364 #if !defined(PETSC_USE_COMPLEX) 365 if (h < umin && h >= 0.0) h = umin; 366 else if (h < 0.0 && h > -umin) h = -umin; 367 #else 368 if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 369 else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 370 #endif 371 h *= epsilon; 372 373 ww[i-rstart] += h; 374 ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 375 ierr = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr); 376 aa[i-rstart] = (v - aa[i-rstart])/h; 377 378 /* possibly shift and scale result */ 379 aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart]; 380 381 ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 382 ww[i-rstart] -= h; 383 ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 384 } 385 ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr); 386 PetscFunctionReturn(0); 387 } 388 389 #undef __FUNCT__ 390 #define __FUNCT__ "MatShift_MFFD" 391 int MatShift_MFFD(PetscScalar *a,Mat Y) 392 { 393 MatSNESMFCtx shell = (MatSNESMFCtx)Y->data; 394 PetscFunctionBegin; 395 shell->vshift += *a; 396 PetscFunctionReturn(0); 397 } 398 399 #undef __FUNCT__ 400 #define __FUNCT__ "MatScale_MFFD" 401 int MatScale_MFFD(PetscScalar *a,Mat Y) 402 { 403 MatSNESMFCtx shell = (MatSNESMFCtx)Y->data; 404 PetscFunctionBegin; 405 shell->vscale *= *a; 406 PetscFunctionReturn(0); 407 } 408 409 410 #undef __FUNCT__ 411 #define __FUNCT__ "MatCreateSNESMF" 412 /*@C 413 MatCreateSNESMF - Creates a matrix-free matrix context for use with 414 a SNES solver. This matrix can be used as the Jacobian argument for 415 the routine SNESSetJacobian(). 416 417 Collective on SNES and Vec 418 419 Input Parameters: 420 + snes - the SNES context 421 - x - vector where SNES solution is to be stored. 422 423 Output Parameter: 424 . J - the matrix-free matrix 425 426 Level: advanced 427 428 Notes: 429 The matrix-free matrix context merely contains the function pointers 430 and work space for performing finite difference approximations of 431 Jacobian-vector products, F'(u)*a, 432 433 The default code uses the following approach to compute h 434 435 .vb 436 F'(u)*a = [F(u+h*a) - F(u)]/h where 437 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 438 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 439 where 440 error_rel = square root of relative error in function evaluation 441 umin = minimum iterate parameter 442 .ve 443 444 The user can set the error_rel via MatSNESMFSetFunctionError() and 445 umin via MatSNESMFDefaultSetUmin(); see the nonlinear solvers chapter 446 of the users manual for details. 447 448 The user should call MatDestroy() when finished with the matrix-free 449 matrix context. 450 451 Options Database Keys: 452 + -snes_mf_err <error_rel> - Sets error_rel 453 . -snes_mf_unim <umin> - Sets umin (for default PETSc routine that computes h only) 454 - -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h 455 456 .keywords: SNES, default, matrix-free, create, matrix 457 458 .seealso: MatDestroy(), MatSNESMFSetFunctionError(), MatSNESMFDefaultSetUmin() 459 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), MatCreateMF(), 460 MatSNESMFGetH(),MatSNESMFKSPMonitor(), MatSNESMFRegisterDynamic), MatSNESMFComputeJacobian() 461 462 @*/ 463 int MatCreateSNESMF(SNES snes,Vec x,Mat *J) 464 { 465 MatSNESMFCtx mfctx; 466 int ierr; 467 468 PetscFunctionBegin; 469 ierr = MatCreateMF(x,J);CHKERRQ(ierr); 470 471 mfctx = (MatSNESMFCtx)(*J)->data; 472 mfctx->snes = snes; 473 mfctx->usesnes = PETSC_TRUE; 474 PetscLogObjectParent(snes,*J); 475 PetscFunctionReturn(0); 476 } 477 478 EXTERN_C_BEGIN 479 #undef __FUNCT__ 480 #define __FUNCT__ "MatSNESMFSetBase_FD" 481 int MatSNESMFSetBase_FD(Mat J,Vec U) 482 { 483 int ierr; 484 MatSNESMFCtx ctx = (MatSNESMFCtx)J->data; 485 486 PetscFunctionBegin; 487 ierr = MatSNESMFResetHHistory(J);CHKERRQ(ierr); 488 ctx->current_u = U; 489 ctx->usesnes = PETSC_FALSE; 490 PetscFunctionReturn(0); 491 } 492 EXTERN_C_END 493 494 #undef __FUNCT__ 495 #define __FUNCT__ "MatSNESMFSetFromOptions" 496 /*@ 497 MatSNESMFSetFromOptions - Sets the MatSNESMF options from the command line 498 parameter. 499 500 Collective on Mat 501 502 Input Parameters: 503 . mat - the matrix obtained with MatCreateSNESMF() 504 505 Options Database Keys: 506 + -snes_mf_type - <default,wp> 507 - -snes_mf_err - square root of estimated relative error in function evaluation 508 - -snes_mf_period - how often h is recomputed, defaults to 1, everytime 509 510 Level: advanced 511 512 .keywords: SNES, matrix-free, parameters 513 514 .seealso: MatCreateSNESMF(),MatSNESMFSetHHistory(), 515 MatSNESMFResetHHistory(), MatSNESMFKSPMonitor() 516 @*/ 517 int MatSNESMFSetFromOptions(Mat mat) 518 { 519 MatSNESMFCtx mfctx = (MatSNESMFCtx)mat->data; 520 int ierr; 521 PetscTruth flg; 522 char ftype[256]; 523 524 PetscFunctionBegin; 525 if (!MatSNESMFRegisterAllCalled) {ierr = MatSNESMFRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 526 527 ierr = PetscOptionsBegin(mfctx->comm,mfctx->prefix,"Set matrix free computation parameters","MatSNESMF");CHKERRQ(ierr); 528 ierr = PetscOptionsList("-snes_mf_type","Matrix free type","MatSNESMFSetType",MatSNESMPetscFList,mfctx->type_name,ftype,256,&flg);CHKERRQ(ierr); 529 if (flg) { 530 ierr = MatSNESMFSetType(mat,ftype);CHKERRQ(ierr); 531 } 532 533 ierr = PetscOptionsReal("-snes_mf_err","set sqrt relative error in function","MatSNESMFSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr); 534 ierr = PetscOptionsInt("-snes_mf_period","how often h is recomputed","MatSNESMFSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr); 535 if (mfctx->snes) { 536 ierr = PetscOptionsName("-snes_mf_ksp_monitor","Monitor matrix-free parameters","MatSNESMFKSPMonitor",&flg);CHKERRQ(ierr); 537 if (flg) { 538 SLES sles; 539 KSP ksp; 540 ierr = SNESGetSLES(mfctx->snes,&sles);CHKERRQ(ierr); 541 ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr); 542 ierr = KSPSetMonitor(ksp,MatSNESMFKSPMonitor,PETSC_NULL,0);CHKERRQ(ierr); 543 } 544 } 545 if (mfctx->ops->setfromoptions) { 546 ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr); 547 } 548 ierr = PetscOptionsEnd();CHKERRQ(ierr); 549 PetscFunctionReturn(0); 550 } 551 552 #undef __FUNCT__ 553 #define __FUNCT__ "MatCreate_MFFD" 554 EXTERN_C_BEGIN 555 int MatCreate_MFFD(Mat A) 556 { 557 MatSNESMFCtx mfctx; 558 int ierr; 559 560 PetscFunctionBegin; 561 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 562 ierr = SNESInitializePackage(PETSC_NULL); CHKERRQ(ierr); 563 #endif 564 565 PetscHeaderCreate(mfctx,_p_MatSNESMFCtx,struct _MFOps,MATSNESMFCTX_COOKIE,0,"SNESMF",A->comm,MatDestroy_MFFD,MatView_MFFD); 566 PetscLogObjectCreate(mfctx); 567 mfctx->sp = 0; 568 mfctx->snes = 0; 569 mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 570 mfctx->recomputeperiod = 1; 571 mfctx->count = 0; 572 mfctx->currenth = 0.0; 573 mfctx->historyh = PETSC_NULL; 574 mfctx->ncurrenth = 0; 575 mfctx->maxcurrenth = 0; 576 mfctx->type_name = 0; 577 mfctx->usesnes = PETSC_FALSE; 578 579 mfctx->vshift = 0.0; 580 mfctx->vscale = 1.0; 581 582 /* 583 Create the empty data structure to contain compute-h routines. 584 These will be filled in below from the command line options or 585 a later call with MatSNESMFSetType() or if that is not called 586 then it will default in the first use of MatMult_MFFD() 587 */ 588 mfctx->ops->compute = 0; 589 mfctx->ops->destroy = 0; 590 mfctx->ops->view = 0; 591 mfctx->ops->setfromoptions = 0; 592 mfctx->hctx = 0; 593 594 mfctx->func = 0; 595 mfctx->funcctx = 0; 596 mfctx->funcvec = 0; 597 598 A->data = mfctx; 599 600 A->ops->mult = MatMult_MFFD; 601 A->ops->destroy = MatDestroy_MFFD; 602 A->ops->view = MatView_MFFD; 603 A->ops->assemblyend = MatAssemblyEnd_MFFD; 604 A->ops->getdiagonal = MatGetDiagonal_MFFD; 605 A->ops->scale = MatScale_MFFD; 606 A->ops->shift = MatShift_MFFD; 607 A->ops->setfromoptions = MatSNESMFSetFromOptions; 608 609 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetBase_C","MatSNESMFSetBase_FD",MatSNESMFSetBase_FD);CHKERRQ(ierr); 610 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetFunctioniBase_C","MatSNESMFSetFunctioniBase_FD",MatSNESMFSetFunctioniBase_FD);CHKERRQ(ierr); 611 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetFunctioni_C","MatSNESMFSetFunctioni_FD",MatSNESMFSetFunctioni_FD);CHKERRQ(ierr); 612 mfctx->mat = A; 613 ierr = VecCreateMPI(A->comm,A->n,A->N,&mfctx->w);CHKERRQ(ierr); 614 615 PetscFunctionReturn(0); 616 } 617 618 EXTERN_C_END 619 620 #undef __FUNCT__ 621 #define __FUNCT__ "MatCreateMF" 622 /*@C 623 MatCreateMF - Creates a matrix-free matrix. See also MatCreateSNESMF() 624 625 Collective on Vec 626 627 Input Parameters: 628 . x - vector that defines layout of the vectors and matrices 629 630 Output Parameter: 631 . J - the matrix-free matrix 632 633 Level: advanced 634 635 Notes: 636 The matrix-free matrix context merely contains the function pointers 637 and work space for performing finite difference approximations of 638 Jacobian-vector products, F'(u)*a, 639 640 The default code uses the following approach to compute h 641 642 .vb 643 F'(u)*a = [F(u+h*a) - F(u)]/h where 644 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 645 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 646 where 647 error_rel = square root of relative error in function evaluation 648 umin = minimum iterate parameter 649 .ve 650 651 The user can set the error_rel via MatSNESMFSetFunctionError() and 652 umin via MatSNESMFDefaultSetUmin(); see the nonlinear solvers chapter 653 of the users manual for details. 654 655 The user should call MatDestroy() when finished with the matrix-free 656 matrix context. 657 658 Options Database Keys: 659 + -snes_mf_err <error_rel> - Sets error_rel 660 . -snes_mf_unim <umin> - Sets umin (for default PETSc routine that computes h only) 661 - -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h 662 663 .keywords: default, matrix-free, create, matrix 664 665 .seealso: MatDestroy(), MatSNESMFSetFunctionError(), MatSNESMFDefaultSetUmin() 666 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), MatCreateSNESMF(), 667 MatSNESMFGetH(),MatSNESMFKSPMonitor(), MatSNESMFRegisterDynamic),, MatSNESMFComputeJacobian() 668 669 @*/ 670 int MatCreateMF(Vec x,Mat *J) 671 { 672 MPI_Comm comm; 673 int n,nloc,ierr; 674 675 PetscFunctionBegin; 676 ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr); 677 ierr = VecGetSize(x,&n);CHKERRQ(ierr); 678 ierr = VecGetLocalSize(x,&nloc);CHKERRQ(ierr); 679 ierr = MatCreate(comm,nloc,nloc,n,n,J);CHKERRQ(ierr); 680 ierr = MatRegister(MATMFFD,0,"MatCreate_MFFD",MatCreate_MFFD);CHKERRQ(ierr); 681 ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr); 682 PetscFunctionReturn(0); 683 } 684 685 686 #undef __FUNCT__ 687 #define __FUNCT__ "MatSNESMFGetH" 688 /*@ 689 MatSNESMFGetH - Gets the last value that was used as the differencing 690 parameter. 691 692 Not Collective 693 694 Input Parameters: 695 . mat - the matrix obtained with MatCreateSNESMF() 696 697 Output Paramter: 698 . h - the differencing step size 699 700 Level: advanced 701 702 .keywords: SNES, matrix-free, parameters 703 704 .seealso: MatCreateSNESMF(),MatSNESMFSetHHistory(), 705 MatSNESMFResetHHistory(),MatSNESMFKSPMonitor() 706 @*/ 707 int MatSNESMFGetH(Mat mat,PetscScalar *h) 708 { 709 MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data; 710 711 PetscFunctionBegin; 712 *h = ctx->currenth; 713 PetscFunctionReturn(0); 714 } 715 716 #undef __FUNCT__ 717 #define __FUNCT__ "MatSNESMFKSPMonitor" 718 /* 719 MatSNESMFKSPMonitor - A KSP monitor for use with the default PETSc 720 SNES matrix free routines. Prints the differencing parameter used at 721 each step. 722 */ 723 int MatSNESMFKSPMonitor(KSP ksp,int n,PetscReal rnorm,void *dummy) 724 { 725 PC pc; 726 MatSNESMFCtx ctx; 727 int ierr; 728 Mat mat; 729 MPI_Comm comm; 730 PetscTruth nonzeroinitialguess; 731 732 PetscFunctionBegin; 733 ierr = PetscObjectGetComm((PetscObject)ksp,&comm);CHKERRQ(ierr); 734 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 735 ierr = KSPGetInitialGuessNonzero(ksp,&nonzeroinitialguess);CHKERRQ(ierr); 736 ierr = PCGetOperators(pc,&mat,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 737 ctx = (MatSNESMFCtx)mat->data; 738 739 if (n > 0 || nonzeroinitialguess) { 740 #if defined(PETSC_USE_COMPLEX) 741 ierr = PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g + %g i\n",n,rnorm, 742 PetscRealPart(ctx->currenth),PetscImaginaryPart(ctx->currenth));CHKERRQ(ierr); 743 #else 744 ierr = PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g \n",n,rnorm,ctx->currenth);CHKERRQ(ierr); 745 #endif 746 } else { 747 ierr = PetscPrintf(comm,"%d KSP Residual norm %14.12e\n",n,rnorm);CHKERRQ(ierr); 748 } 749 PetscFunctionReturn(0); 750 } 751 752 #undef __FUNCT__ 753 #define __FUNCT__ "MatSNESMFSetFunction" 754 /*@C 755 MatSNESMFSetFunction - Sets the function used in applying the matrix free. 756 757 Collective on Mat 758 759 Input Parameters: 760 + mat - the matrix free matrix created via MatCreateSNESMF() 761 . v - workspace vector 762 . func - the function to use 763 - funcctx - optional function context passed to function 764 765 Level: advanced 766 767 Notes: 768 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 769 matrix inside your compute Jacobian routine 770 771 If this is not set then it will use the function set with SNESSetFunction() 772 773 .keywords: SNES, matrix-free, function 774 775 .seealso: MatCreateSNESMF(),MatSNESMFGetH(), 776 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 777 MatSNESMFKSPMonitor(), SNESetFunction() 778 @*/ 779 int MatSNESMFSetFunction(Mat mat,Vec v,int (*func)(SNES,Vec,Vec,void *),void *funcctx) 780 { 781 MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data; 782 783 PetscFunctionBegin; 784 ctx->func = func; 785 ctx->funcctx = funcctx; 786 ctx->funcvec = v; 787 PetscFunctionReturn(0); 788 } 789 790 #undef __FUNCT__ 791 #define __FUNCT__ "MatSNESMFSetFunctioni" 792 /*@C 793 MatSNESMFSetFunctioni - Sets the function for a single component 794 795 Collective on Mat 796 797 Input Parameters: 798 + mat - the matrix free matrix created via MatCreateSNESMF() 799 - funci - the function to use 800 801 Level: advanced 802 803 Notes: 804 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 805 matrix inside your compute Jacobian routine 806 807 808 .keywords: SNES, matrix-free, function 809 810 .seealso: MatCreateSNESMF(),MatSNESMFGetH(), 811 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 812 MatSNESMFKSPMonitor(), SNESetFunction() 813 @*/ 814 int MatSNESMFSetFunctioni(Mat mat,int (*funci)(int,Vec,PetscScalar*,void *)) 815 { 816 int ierr,(*f)(Mat,int (*)(int,Vec,PetscScalar*,void *)); 817 818 PetscFunctionBegin; 819 PetscValidHeaderSpecific(mat,MAT_COOKIE); 820 ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSNESMFSetFunctioni_C",(void (**)(void))&f);CHKERRQ(ierr); 821 if (f) { 822 ierr = (*f)(mat,funci);CHKERRQ(ierr); 823 } 824 PetscFunctionReturn(0); 825 } 826 827 828 #undef __FUNCT__ 829 #define __FUNCT__ "MatSNESMFSetFunctioniBase" 830 /*@C 831 MatSNESMFSetFunctioniBase - Sets the base vector for a single component function evaluation 832 833 Collective on Mat 834 835 Input Parameters: 836 + mat - the matrix free matrix created via MatCreateSNESMF() 837 - func - the function to use 838 839 Level: advanced 840 841 Notes: 842 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 843 matrix inside your compute Jacobian routine 844 845 846 .keywords: SNES, matrix-free, function 847 848 .seealso: MatCreateSNESMF(),MatSNESMFGetH(), 849 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 850 MatSNESMFKSPMonitor(), SNESetFunction() 851 @*/ 852 int MatSNESMFSetFunctioniBase(Mat mat,int (*func)(Vec,void *)) 853 { 854 int ierr,(*f)(Mat,int (*)(Vec,void *)); 855 856 PetscFunctionBegin; 857 PetscValidHeaderSpecific(mat,MAT_COOKIE); 858 ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSNESMFSetFunctioniBase_C",(void (**)(void))&f);CHKERRQ(ierr); 859 if (f) { 860 ierr = (*f)(mat,func);CHKERRQ(ierr); 861 } 862 PetscFunctionReturn(0); 863 } 864 865 866 #undef __FUNCT__ 867 #define __FUNCT__ "MatSNESMFSetPeriod" 868 /*@ 869 MatSNESMFSetPeriod - Sets how often h is recomputed, by default it is everytime 870 871 Collective on Mat 872 873 Input Parameters: 874 + mat - the matrix free matrix created via MatCreateSNESMF() 875 - period - 1 for everytime, 2 for every second etc 876 877 Options Database Keys: 878 + -snes_mf_period <period> 879 880 Level: advanced 881 882 883 .keywords: SNES, matrix-free, parameters 884 885 .seealso: MatCreateSNESMF(),MatSNESMFGetH(), 886 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 887 MatSNESMFKSPMonitor() 888 @*/ 889 int MatSNESMFSetPeriod(Mat mat,int period) 890 { 891 MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data; 892 893 PetscFunctionBegin; 894 ctx->recomputeperiod = period; 895 PetscFunctionReturn(0); 896 } 897 898 #undef __FUNCT__ 899 #define __FUNCT__ "MatSNESMFSetFunctionError" 900 /*@ 901 MatSNESMFSetFunctionError - Sets the error_rel for the approximation of 902 matrix-vector products using finite differences. 903 904 Collective on Mat 905 906 Input Parameters: 907 + mat - the matrix free matrix created via MatCreateSNESMF() 908 - error_rel - relative error (should be set to the square root of 909 the relative error in the function evaluations) 910 911 Options Database Keys: 912 + -snes_mf_err <error_rel> - Sets error_rel 913 914 Level: advanced 915 916 Notes: 917 The default matrix-free matrix-vector product routine computes 918 .vb 919 F'(u)*a = [F(u+h*a) - F(u)]/h where 920 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 921 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 922 .ve 923 924 .keywords: SNES, matrix-free, parameters 925 926 .seealso: MatCreateSNESMF(),MatSNESMFGetH(), 927 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 928 MatSNESMFKSPMonitor() 929 @*/ 930 int MatSNESMFSetFunctionError(Mat mat,PetscReal error) 931 { 932 MatSNESMFCtx ctx = (MatSNESMFCtx)mat->data; 933 934 PetscFunctionBegin; 935 if (error != PETSC_DEFAULT) ctx->error_rel = error; 936 PetscFunctionReturn(0); 937 } 938 939 #undef __FUNCT__ 940 #define __FUNCT__ "MatSNESMFAddNullSpace" 941 /*@ 942 MatSNESMFAddNullSpace - Provides a null space that an operator is 943 supposed to have. Since roundoff will create a small component in 944 the null space, if you know the null space you may have it 945 automatically removed. 946 947 Collective on Mat 948 949 Input Parameters: 950 + J - the matrix-free matrix context 951 - nullsp - object created with MatNullSpaceCreate() 952 953 Level: advanced 954 955 .keywords: SNES, matrix-free, null space 956 957 .seealso: MatNullSpaceCreate(), MatSNESMFGetH(), MatCreateSNESMF(), 958 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 959 MatSNESMFKSPMonitor(), MatSNESMFErrorRel() 960 @*/ 961 int MatSNESMFAddNullSpace(Mat J,MatNullSpace nullsp) 962 { 963 int ierr; 964 MatSNESMFCtx ctx = (MatSNESMFCtx)J->data; 965 MPI_Comm comm; 966 967 PetscFunctionBegin; 968 ierr = PetscObjectGetComm((PetscObject)J,&comm);CHKERRQ(ierr); 969 970 ctx->sp = nullsp; 971 ierr = PetscObjectReference((PetscObject)nullsp);CHKERRQ(ierr); 972 PetscFunctionReturn(0); 973 } 974 975 #undef __FUNCT__ 976 #define __FUNCT__ "MatSNESMFSetHHistory" 977 /*@ 978 MatSNESMFSetHHistory - Sets an array to collect a history of the 979 differencing values (h) computed for the matrix-free product. 980 981 Collective on Mat 982 983 Input Parameters: 984 + J - the matrix-free matrix context 985 . histroy - space to hold the history 986 - nhistory - number of entries in history, if more entries are generated than 987 nhistory, then the later ones are discarded 988 989 Level: advanced 990 991 Notes: 992 Use MatSNESMFResetHHistory() to reset the history counter and collect 993 a new batch of differencing parameters, h. 994 995 .keywords: SNES, matrix-free, h history, differencing history 996 997 .seealso: MatSNESMFGetH(), MatCreateSNESMF(), 998 MatSNESMFResetHHistory(), 999 MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError() 1000 1001 @*/ 1002 int MatSNESMFSetHHistory(Mat J,PetscScalar *history,int nhistory) 1003 { 1004 MatSNESMFCtx ctx = (MatSNESMFCtx)J->data; 1005 1006 PetscFunctionBegin; 1007 ctx->historyh = history; 1008 ctx->maxcurrenth = nhistory; 1009 ctx->currenth = 0; 1010 PetscFunctionReturn(0); 1011 } 1012 1013 #undef __FUNCT__ 1014 #define __FUNCT__ "MatSNESMFResetHHistory" 1015 /*@ 1016 MatSNESMFResetHHistory - Resets the counter to zero to begin 1017 collecting a new set of differencing histories. 1018 1019 Collective on Mat 1020 1021 Input Parameters: 1022 . J - the matrix-free matrix context 1023 1024 Level: advanced 1025 1026 Notes: 1027 Use MatSNESMFSetHHistory() to create the original history counter. 1028 1029 .keywords: SNES, matrix-free, h history, differencing history 1030 1031 .seealso: MatSNESMFGetH(), MatCreateSNESMF(), 1032 MatSNESMFSetHHistory(), 1033 MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError() 1034 1035 @*/ 1036 int MatSNESMFResetHHistory(Mat J) 1037 { 1038 MatSNESMFCtx ctx = (MatSNESMFCtx)J->data; 1039 1040 PetscFunctionBegin; 1041 ctx->ncurrenth = 0; 1042 PetscFunctionReturn(0); 1043 } 1044 1045 #undef __FUNCT__ 1046 #define __FUNCT__ "MatSNESMFComputeJacobian" 1047 int MatSNESMFComputeJacobian(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure *flag,void *dummy) 1048 { 1049 int ierr; 1050 PetscFunctionBegin; 1051 ierr = MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1052 ierr = MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1053 PetscFunctionReturn(0); 1054 } 1055 1056 #undef __FUNCT__ 1057 #define __FUNCT__ "MatSNESMFSetBase" 1058 int MatSNESMFSetBase(Mat J,Vec U) 1059 { 1060 int ierr,(*f)(Mat,Vec); 1061 1062 PetscFunctionBegin; 1063 PetscValidHeaderSpecific(J,MAT_COOKIE); 1064 PetscValidHeaderSpecific(U,VEC_COOKIE); 1065 ierr = PetscObjectQueryFunction((PetscObject)J,"MatSNESMFSetBase_C",(void (**)(void))&f);CHKERRQ(ierr); 1066 if (f) { 1067 ierr = (*f)(J,U);CHKERRQ(ierr); 1068 } 1069 PetscFunctionReturn(0); 1070 } 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080