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