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