1 2 #ifdef PETSC_RCS_HEADER 3 static char vcid[] = "$Id: snesmfj.c,v 1.79 1999/03/07 17:29:26 bsmith Exp bsmith $"; 4 #endif 5 6 #include "src/snes/snesimpl.h" 7 #include "src/snes/mf/snesmfj.h" /*I "snes.h" I*/ 8 9 FList MatSNESFDMFList = 0; 10 int MatSNESFDMFRegisterAllCalled = 0; 11 12 13 #undef __FUNC__ 14 #define __FUNC__ "MatSNESFDMFSetType" 15 /*@ 16 MatSNESFDMFSetType - Sets the method that is used to compute the h in the 17 finite difference matrix free formulation. 18 19 Input Parameters: 20 + mat - the matrix free matrix created via MatCreateSNESFDMF() 21 - ftype - the type requested 22 23 Level: advanced 24 25 .seealso: MatCreateSNESFDMF(), MatSNESFDMFRegister() 26 @*/ 27 int MatSNESFDMFSetType(Mat mat,char *ftype) 28 { 29 int ierr, (*r)(MatSNESFDMFCtx); 30 MatSNESFDMFCtx ctx; 31 32 PetscFunctionBegin; 33 ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 34 35 /* already set, so just return */ 36 if (PetscTypeCompare(ctx->type_name,ftype)) PetscFunctionReturn(0); 37 38 /* destroy the old one if it exists */ 39 if (ctx->ops->destroy) { 40 ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr); 41 } 42 43 /* Get the function pointers for the iterative method requested */ 44 if (!MatSNESFDMFRegisterAllCalled) {ierr = MatSNESFDMFRegisterAll(PETSC_NULL); CHKERRQ(ierr);} 45 46 ierr = FListFind(ctx->comm, MatSNESFDMFList, ftype,(int (**)(void *)) &r );CHKERRQ(ierr); 47 48 if (!r) SETERRQ(1,1,"Unknown MatSNESFDMF type given"); 49 50 ierr = (*r)(ctx); CHKERRQ(ierr); 51 ierr = PetscStrncpy(ctx->type_name,ftype,256);CHKERRQ(ierr); 52 53 PetscFunctionReturn(0); 54 } 55 56 /*MC 57 MatSNESFDMFRegister - Adds a method to the MatSNESFDMF registry 58 59 Synopsis: 60 MatSNESFDMFRegister(char *name_solver,char *path,char *name_create,int (*routine_create)(MatSNESFDMF)) 61 62 Not Collective 63 64 Input Parameters: 65 + name_solver - name of a new user-defined compute-h module 66 . path - path (either absolute or relative) the library containing this solver 67 . name_create - name of routine to create method context 68 - routine_create - routine to create method context 69 70 Level: developer 71 72 Notes: 73 MatSNESFDMFRegister() may be called multiple times to add several user-defined solvers. 74 75 If dynamic libraries are used, then the fourth input argument (routine_create) 76 is ignored. 77 78 Sample usage: 79 .vb 80 MatSNESFDMFRegister("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a, 81 "MyHCreate",MyHCreate); 82 .ve 83 84 Then, your solver can be chosen with the procedural interface via 85 $ MatSNESFDMFSetType(mfctx,"my_h") 86 or at runtime via the option 87 $ -snes_mf_type my_h 88 89 .keywords: MatSNESFDMF, register 90 91 .seealso: MatSNESFDMFRegisterAll(), MatSNESFDMFRegisterDestroy() 92 M*/ 93 94 #undef __FUNC__ 95 #define __FUNC__ "MatSNESFDMFRegister_Private" 96 int MatSNESFDMFRegister_Private(char *sname,char *path,char *name,int (*function)(MatSNESFDMFCtx)) 97 { 98 int ierr; 99 char fullname[256]; 100 101 PetscFunctionBegin; 102 PetscStrcpy(fullname,path); PetscStrcat(fullname,":");PetscStrcat(fullname,name); 103 ierr = FListAdd_Private(&MatSNESFDMFList,sname,fullname,(int (*)(void*))function);CHKERRQ(ierr); 104 PetscFunctionReturn(0); 105 } 106 107 108 #undef __FUNC__ 109 #define __FUNC__ "MatSNESFDMFRegisterDestroy" 110 /*@C 111 MatSNESFDMFRegisterDestroy - Frees the list of MatSNESFDMF methods that were 112 registered by MatSNESFDMFRegister(). 113 114 Not Collective 115 116 Level: developer 117 118 .keywords: MatSNESFDMF, register, destroy 119 120 .seealso: MatSNESFDMFRegister(), MatSNESFDMFRegisterAll() 121 @*/ 122 int MatSNESFDMFRegisterDestroy(void) 123 { 124 int ierr; 125 126 PetscFunctionBegin; 127 if (MatSNESFDMFList) { 128 ierr = FListDestroy( MatSNESFDMFList );CHKERRQ(ierr); 129 MatSNESFDMFList = 0; 130 } 131 MatSNESFDMFRegisterAllCalled = 0; 132 PetscFunctionReturn(0); 133 } 134 135 /* ----------------------------------------------------------------------------------------*/ 136 #undef __FUNC__ 137 #define __FUNC__ "MatSNESFDMFDestroy_Private" 138 int MatSNESFDMFDestroy_Private(Mat mat) 139 { 140 int ierr; 141 MatSNESFDMFCtx ctx; 142 143 PetscFunctionBegin; 144 ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 145 ierr = VecDestroy(ctx->w); CHKERRQ(ierr); 146 if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);} 147 if (ctx->sp) {ierr = PCNullSpaceDestroy(ctx->sp);CHKERRQ(ierr);} 148 PetscFree(ctx->ops); 149 PetscFree(ctx); 150 PetscFunctionReturn(0); 151 } 152 153 #undef __FUNC__ 154 #define __FUNC__ "MatSNESFDMFView_Private" 155 /* 156 MatSNESFDMFView_Private - Views matrix-free parameters. 157 158 */ 159 int MatSNESFDMFView_Private(Mat J,Viewer viewer) 160 { 161 int ierr; 162 MatSNESFDMFCtx ctx; 163 MPI_Comm comm; 164 FILE *fd; 165 ViewerType vtype; 166 167 PetscFunctionBegin; 168 ierr = PetscObjectGetComm((PetscObject)J,&comm);CHKERRQ(ierr); 169 ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 170 ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 171 ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 172 if (PetscTypeCompare(vtype,ASCII_VIEWER)) { 173 PetscFPrintf(comm,fd," SNES matrix-free approximation:\n"); 174 PetscFPrintf(comm,fd," err=%g (relative error in function evaluation)\n",ctx->error_rel); 175 PetscFPrintf(ctx->comm,fd," Using %s compute h routine\n",ctx->type_name); 176 if (ctx->ops->view) { 177 ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr); 178 } 179 } else { 180 SETERRQ(1,1,"Viewer type not supported for this object"); 181 } 182 PetscFunctionReturn(0); 183 } 184 185 #undef __FUNC__ 186 #define __FUNC__ "MatSNESFDMFAssemblyEnd_Private" 187 /* 188 MatSNESFDMFAssemblyEnd_Private - Resets the ctx->ncurrenth to zero. This 189 allows the user to indicate the beginning of a new linear solve by call 190 MatAssemblyXXX() on the matrix free matrix. This then allows the 191 MatSNESFDMFCreate_WP() to properly compute the || U|| only the first 192 time in the linear solver rather than every time 193 194 */ 195 int MatSNESFDMFAssemblyEnd_Private(Mat J) 196 { 197 int ierr; 198 199 PetscFunctionBegin; 200 ierr = MatSNESFDMFResetHHistory(J);CHKERRQ(ierr); 201 PetscFunctionReturn(0); 202 } 203 204 205 #undef __FUNC__ 206 #define __FUNC__ "MatSNESFDMFMult_Private" 207 /* 208 MatSNESFDMFMult_Private - Default matrix-free form for Jacobian-vector 209 product, y = F'(u)*a: 210 211 y ~= ( F(u + ha) - F(u) )/h, 212 where F = nonlinear function, as set by SNESSetFunction() 213 u = current iterate 214 h = difference interval 215 */ 216 int MatSNESFDMFMult_Private(Mat mat,Vec a,Vec y) 217 { 218 MatSNESFDMFCtx ctx; 219 SNES snes; 220 Scalar h, mone = -1.0; 221 Vec w,U,F; 222 int ierr, (*eval_fct)(SNES,Vec,Vec)=0; 223 224 PetscFunctionBegin; 225 /* We log matrix-free matrix-vector products separately, so that we can 226 separate the performance monitoring from the cases that use conventional 227 storage. We may eventually modify event logging to associate events 228 with particular objects, hence alleviating the more general problem. */ 229 PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0); 230 231 ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 232 snes = ctx->snes; 233 w = ctx->w; 234 235 ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr); 236 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 237 eval_fct = SNESComputeFunction; 238 ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr); 239 } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 240 eval_fct = SNESComputeGradient; 241 ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr); 242 } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid method class"); 243 244 if (!ctx->ops->compute) { 245 ierr = MatSNESFDMFSetType(mat,"default");CHKERRQ(ierr); 246 ierr = MatSNESFDMFSetFromOptions(mat);CHKERRQ(ierr); 247 } 248 ierr = (*ctx->ops->compute)(ctx,U,a,&h); CHKERRQ(ierr); 249 250 /* keep a record of the current differencing parameter h */ 251 ctx->currenth = h; 252 #if defined(USE_PETSC_COMPLEX) 253 PLogInfo(mat,"Current differencing parameter: %g + %g i\n",PetscReal(h),PetscImaginary(h)); 254 #else 255 PLogInfo(mat,"Current differencing parameter: %g\n",h); 256 #endif 257 if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 258 ctx->historyh[ctx->ncurrenth++] = h; 259 } else { 260 ctx->ncurrenth++; 261 } 262 263 /* Evaluate function at F(u + ha) */ 264 ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr); 265 ierr = eval_fct(snes,w,y); CHKERRQ(ierr); 266 267 ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr); 268 h = 1.0/h; 269 ierr = VecScale(&h,y); CHKERRQ(ierr); 270 if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);} 271 272 PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0); 273 PetscFunctionReturn(0); 274 } 275 276 #undef __FUNC__ 277 #define __FUNC__ "MatCreateSNESFDMF" 278 /*@C 279 MatCreateSNESFDMF - Creates a matrix-free matrix 280 context for use with a SNES solver. This matrix can be used as 281 the Jacobian argument for the routine SNESSetJacobian(). 282 283 Collective on SNES and Vec 284 285 Input Parameters: 286 + snes - the SNES context 287 - x - vector where SNES solution is to be stored. 288 289 Output Parameter: 290 . J - the matrix-free matrix 291 292 Level: advanced 293 294 Notes: 295 The matrix-free matrix context merely contains the function pointers 296 and work space for performing finite difference approximations of 297 Jacobian-vector products, J(u)*a, 298 299 The default code uses the following approach to compute h 300 301 .vb 302 J(u)*a = [J(u+h*a) - J(u)]/h where 303 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 304 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 305 where 306 error_rel = square root of relative error in function evaluation 307 umin = minimum iterate parameter 308 .ve 309 310 The user can set the error_rel via MatSNESFDMFSetFunctionError() and 311 umin via MatSNESFDMFDefaultSetUmin() 312 See the nonlinear solvers chapter of the users manual for details. 313 314 The user should call MatDestroy() when finished with the matrix-free 315 matrix context. 316 317 Options Database Keys: 318 + -snes_mf_err <error_rel> - Sets error_rel 319 . -snes_mf_unim <umin> - Sets umin (for default PETSc routine that computes h only) 320 - -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h 321 322 .keywords: SNES, default, matrix-free, create, matrix 323 324 .seealso: MatDestroy(), MatSNESFDMFSetFunctionError(), MatSNESFDMFDefaultSetUmin() 325 MatSNESFDMFSetHHistory(), MatSNESFDMFResetHHistory(), 326 MatSNESFDMFGetH(),MatSNESFDMFKSPMonitor(), MatSNESFDMFRegister() 327 328 @*/ 329 int MatCreateSNESFDMF(SNES snes,Vec x, Mat *J) 330 { 331 MPI_Comm comm; 332 MatSNESFDMFCtx mfctx; 333 int n, nloc, ierr; 334 335 PetscFunctionBegin; 336 mfctx = (MatSNESFDMFCtx) PetscMalloc(sizeof(struct _p_MatSNESFDMFCtx)); CHKPTRQ(mfctx); 337 PLogObjectMemory(snes,sizeof(MatSNESFDMFCtx)); 338 mfctx->comm = snes->comm; 339 mfctx->sp = 0; 340 mfctx->snes = snes; 341 mfctx->error_rel = 1.e-8; /* assumes double precision */ 342 mfctx->currenth = 0.0; 343 mfctx->historyh = PETSC_NULL; 344 mfctx->ncurrenth = 0; 345 mfctx->maxcurrenth = 0; 346 347 /* 348 Create the empty data structure to contain compute-h routines. 349 These will be filled in below from the command line options or 350 a later call with MatSNESFDMFSetType() or if that is not called 351 then it will default in the first use of MatSNESFDMFMult_private() 352 */ 353 mfctx->ops = (MFOps *)PetscMalloc(sizeof(MFOps)); CHKPTRQ(mfctx->ops); 354 mfctx->ops->compute = 0; 355 mfctx->ops->destroy = 0; 356 mfctx->ops->view = 0; 357 mfctx->ops->printhelp = 0; 358 mfctx->ops->setfromoptions = 0; 359 mfctx->hctx = 0; 360 361 ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr); 362 ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr); 363 ierr = VecGetSize(x,&n); CHKERRQ(ierr); 364 ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr); 365 ierr = MatCreateShell(comm,nloc,nloc,n,n,mfctx,J); CHKERRQ(ierr); 366 ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)MatSNESFDMFMult_Private);CHKERRQ(ierr); 367 ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)MatSNESFDMFDestroy_Private);CHKERRQ(ierr); 368 ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)MatSNESFDMFView_Private); CHKERRQ(ierr); 369 ierr = MatShellSetOperation(*J,MATOP_ASSEMBLY_END,(void *)MatSNESFDMFAssemblyEnd_Private);CHKERRQ(ierr); 370 PLogObjectParent(*J,mfctx->w); 371 PLogObjectParent(snes,*J); 372 373 mfctx->mat = *J; 374 375 376 PetscFunctionReturn(0); 377 } 378 379 #undef __FUNC__ 380 #define __FUNC__ "MatSNESFDMFGetH" 381 /*@ 382 MatSNESFDMFSetFromOptions - Sets the MatSNESFDMF options from the command line 383 parameter. 384 385 Collective on Mat 386 387 Input Parameters: 388 . mat - the matrix obtained with MatCreateSNESFDMF() 389 390 Level: advanced 391 392 .keywords: SNES, matrix-free, parameters 393 394 .seealso: MatCreateSNESFDMF(),MatSNESFDMFSetHHistory(), 395 MatSNESFDMFResetHHistory(),MatSNESFDMFKSPMonitor() 396 @*/ 397 int MatSNESFDMFSetFromOptions(Mat mat) 398 { 399 MatSNESFDMFCtx mfctx; 400 int ierr,flg; 401 char ftype[256],p[64]; 402 403 PetscFunctionBegin; 404 ierr = MatShellGetContext(mat,(void **)&mfctx); CHKERRQ(ierr); 405 if (mfctx) { 406 /* allow user to set the type */ 407 ierr = OptionsGetString(mfctx->snes->prefix,"-snes_mf_type",ftype,256,&flg);CHKERRQ(ierr); 408 if (flg) { 409 ierr = MatSNESFDMFSetType(mat,ftype);CHKERRQ(ierr); 410 } 411 412 ierr = OptionsGetDouble(mfctx->snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg);CHKERRQ(ierr); 413 if (mfctx->ops->setfromoptions) { 414 ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr); 415 } 416 417 ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 418 PetscStrcpy(p,"-"); 419 if (mfctx->snes->prefix) PetscStrcat(p,mfctx->snes->prefix); 420 if (flg) { 421 (*PetscHelpPrintf)(mfctx->snes->comm," %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel); 422 if (mfctx->ops->printhelp) { 423 (*mfctx->ops->printhelp)(mfctx);CHKERRQ(ierr); 424 } 425 } 426 } 427 PetscFunctionReturn(0); 428 } 429 430 #undef __FUNC__ 431 #define __FUNC__ "MatSNESFDMFGetH" 432 /*@ 433 MatSNESFDMFGetH - Gets the last h that was used as the differencing 434 parameter. 435 436 Not Collective 437 438 Input Parameters: 439 . mat - the matrix obtained with MatCreateSNESFDMF() 440 441 Output Paramter: 442 . h - the differencing step size 443 444 Level: advanced 445 446 .keywords: SNES, matrix-free, parameters 447 448 .seealso: MatCreateSNESFDMF(),MatSNESFDMFSetHHistory(), 449 MatSNESFDMFResetHHistory(),MatSNESFDMFKSPMonitor() 450 @*/ 451 int MatSNESFDMFGetH(Mat mat,Scalar *h) 452 { 453 MatSNESFDMFCtx ctx; 454 int ierr; 455 456 PetscFunctionBegin; 457 ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 458 if (ctx) { 459 *h = ctx->currenth; 460 } 461 PetscFunctionReturn(0); 462 } 463 464 #undef __FUNC__ 465 #define __FUNC__ "MatSNESFDMFKSPMonitor" 466 /* 467 MatSNESFDMFKSPMonitor - A KSP monitor for use with the default PETSc 468 SNES matrix free routines. Prints the h differencing parameter used at each 469 timestep. 470 471 */ 472 int MatSNESFDMFKSPMonitor(KSP ksp,int n,double rnorm,void *dummy) 473 { 474 PC pc; 475 MatSNESFDMFCtx ctx; 476 int ierr; 477 Mat mat; 478 MPI_Comm comm; 479 PetscTruth nonzeroinitialguess; 480 481 PetscFunctionBegin; 482 ierr = PetscObjectGetComm((PetscObject)ksp,&comm);CHKERRQ(ierr); 483 ierr = KSPGetPC(ksp,&pc); CHKERRQ(ierr); 484 ierr = KSPGetInitialGuessNonzero(ksp,&nonzeroinitialguess);CHKERRQ(ierr); 485 ierr = PCGetOperators(pc,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 486 ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 487 if (!ctx) { 488 SETERRQ(1,1,"Matrix is not a matrix free shell matrix"); 489 } 490 if (n > 0 || nonzeroinitialguess) { 491 #if defined(USE_PETSC_COMPLEX) 492 PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g + %g i\n",n,rnorm, 493 PetscReal(ctx->currenth),PetscImaginary(ctx->currenth)); 494 #else 495 PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g \n",n,rnorm,ctx->currenth); 496 #endif 497 } else { 498 PetscPrintf(comm,"%d KSP Residual norm %14.12e\n",n,rnorm); 499 } 500 PetscFunctionReturn(0); 501 } 502 503 #undef __FUNC__ 504 #define __FUNC__ "MatSNESFDMFSetFunctionError" 505 /*@ 506 MatSNESFDMFSetFunctionError - Sets the error_rel for the approximation of 507 matrix-vector products using finite differences. 508 509 Collective on Mat 510 511 Input Parameters: 512 + mat - the matrix free matrix created via MatCreateSNESFDMF() 513 - error_rel - relative error (should be set to the square root of 514 the relative error in the function evaluations) 515 516 Options Database Keys: 517 + -snes_mf_err <error_rel> - Sets error_rel 518 519 Level: advanced 520 521 Notes: 522 The default matrix-free matrix-vector product routine computes 523 .vb 524 J(u)*a = [J(u+h*a) - J(u)]/h where 525 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 526 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 527 .ve 528 529 530 .keywords: SNES, matrix-free, parameters 531 532 .seealso: MatCreateSNESFDMF(),MatSNESFDMFGetH(), 533 MatSNESFDMFSetHHistory(), MatSNESFDMFResetHHistory(), 534 MatSNESFDMFKSPMonitor() 535 @*/ 536 int MatSNESFDMFSetFunctionError(Mat mat,double error) 537 { 538 MatSNESFDMFCtx ctx; 539 int ierr; 540 541 PetscFunctionBegin; 542 ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 543 if (ctx) { 544 if (error != PETSC_DEFAULT) ctx->error_rel = error; 545 } 546 PetscFunctionReturn(0); 547 } 548 549 #undef __FUNC__ 550 #define __FUNC__ "MatSNESFDMFAddNullSpace" 551 /*@ 552 MatSNESFDMFAddNullSpace - Provides a null space that 553 an operator is supposed to have. Since roundoff will create a 554 small component in the null space, if you know the null space 555 you may have it automatically removed. 556 557 Collective on Mat 558 559 Input Parameters: 560 + J - the matrix-free matrix context 561 . has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants 562 . n - number of vectors (excluding constant vector) in null space 563 - vecs - the vectors that span the null space (excluding the constant vector); 564 these vectors must be orthonormal 565 566 Level: advanced 567 568 .keywords: SNES, matrix-free, null space 569 570 .seealso: MatSNESFDMFGetH(), MatCreateSNESFDMF(), 571 MatSNESFDMFSetHHistory(), MatSNESFDMFResetHHistory(), 572 MatSNESFDMFKSPMonitor(), MatSNESFDMFErrorRel() 573 574 @*/ 575 int MatSNESFDMFAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs) 576 { 577 int ierr; 578 MatSNESFDMFCtx ctx; 579 MPI_Comm comm; 580 581 PetscFunctionBegin; 582 ierr = PetscObjectGetComm((PetscObject)J,&comm);CHKERRQ(ierr); 583 584 ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 585 /* no context indicates that it is not the "matrix free" matrix type */ 586 if (!ctx) PetscFunctionReturn(0); 587 ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr); 588 PetscFunctionReturn(0); 589 } 590 591 #undef __FUNC__ 592 #define __FUNC__ "MatSNESFDMFSetHHistory" 593 /*@ 594 MatSNESFDMFSetHHistory - Sets an array to collect a history 595 of the differencing values h computed for the matrix free product 596 597 Collective on Mat 598 599 Input Parameters: 600 + J - the matrix-free matrix context 601 . histroy - space to hold the h history 602 - nhistory - number of entries in history, if more h are generated than 603 nhistory the later ones are discarded 604 605 Level: advanced 606 607 Notes: 608 Use MatSNESFDMFResetHHistory() to reset the history counter 609 and collect a new batch of h. 610 611 .keywords: SNES, matrix-free, h history, differencing history 612 613 .seealso: MatSNESFDMFGetH(), MatCreateSNESFDMF(), 614 MatSNESFDMFResetHHistory(), 615 MatSNESFDMFKSPMonitor(), MatSNESFDMFSetFunctionError() 616 617 @*/ 618 int MatSNESFDMFSetHHistory(Mat J,Scalar *history,int nhistory) 619 { 620 int ierr; 621 MatSNESFDMFCtx ctx; 622 623 PetscFunctionBegin; 624 625 ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 626 /* no context indicates that it is not the "matrix free" matrix type */ 627 if (!ctx) PetscFunctionReturn(0); 628 ctx->historyh = history; 629 ctx->maxcurrenth = nhistory; 630 ctx->currenth = 0; 631 632 PetscFunctionReturn(0); 633 } 634 635 #undef __FUNC__ 636 #define __FUNC__ "MatSNESFDMFResetHHistory" 637 /*@ 638 MatSNESFDMFResetHHistory - Resets the counter to zero to begin 639 collecting a new set of differencing histories. 640 641 Collective on Mat 642 643 Input Parameters: 644 . J - the matrix-free matrix context 645 646 Level: advanced 647 648 Notes: 649 Use MatSNESFDMFSetHHistory() to create the original history counter 650 651 .keywords: SNES, matrix-free, h history, differencing history 652 653 .seealso: MatSNESFDMFGetH(), MatCreateSNESFDMF(), 654 MatSNESFDMFSetHHistory(), 655 MatSNESFDMFKSPMonitor(), MatSNESFDMFSetFunctionError() 656 657 @*/ 658 int MatSNESFDMFResetHHistory(Mat J) 659 { 660 int ierr; 661 MatSNESFDMFCtx ctx; 662 663 PetscFunctionBegin; 664 665 ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 666 /* no context indicates that it is not the "matrix free" matrix type */ 667 if (!ctx) PetscFunctionReturn(0); 668 ctx->ncurrenth = 0; 669 670 PetscFunctionReturn(0); 671 } 672 673