1 2 #include <petsc/private/matimpl.h> 3 #include <../src/mat/impls/mffd/mffdimpl.h> /*I "petscmat.h" I*/ 4 5 PetscFunctionList MatMFFDList = NULL; 6 PetscBool MatMFFDRegisterAllCalled = PETSC_FALSE; 7 8 PetscClassId MATMFFD_CLASSID; 9 PetscLogEvent MATMFFD_Mult; 10 11 static PetscBool MatMFFDPackageInitialized = PETSC_FALSE; 12 /*@C 13 MatMFFDFinalizePackage - This function destroys everything in the MatMFFD package. It is 14 called from PetscFinalize(). 15 16 Level: developer 17 18 .seealso: PetscFinalize(), MatCreateMFFD(), MatCreateSNESMF() 19 @*/ 20 PetscErrorCode MatMFFDFinalizePackage(void) 21 { 22 PetscFunctionBegin; 23 PetscCall(PetscFunctionListDestroy(&MatMFFDList)); 24 MatMFFDPackageInitialized = PETSC_FALSE; 25 MatMFFDRegisterAllCalled = PETSC_FALSE; 26 PetscFunctionReturn(0); 27 } 28 29 /*@C 30 MatMFFDInitializePackage - This function initializes everything in the MatMFFD package. It is called 31 from MatInitializePackage(). 32 33 Level: developer 34 35 .seealso: PetscInitialize() 36 @*/ 37 PetscErrorCode MatMFFDInitializePackage(void) 38 { 39 char logList[256]; 40 PetscBool opt,pkg; 41 42 PetscFunctionBegin; 43 if (MatMFFDPackageInitialized) PetscFunctionReturn(0); 44 MatMFFDPackageInitialized = PETSC_TRUE; 45 /* Register Classes */ 46 PetscCall(PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID)); 47 /* Register Constructors */ 48 PetscCall(MatMFFDRegisterAll()); 49 /* Register Events */ 50 PetscCall(PetscLogEventRegister("MatMult MF",MATMFFD_CLASSID,&MATMFFD_Mult)); 51 /* Process Info */ 52 { 53 PetscClassId classids[1]; 54 55 classids[0] = MATMFFD_CLASSID; 56 PetscCall(PetscInfoProcessClass("matmffd", 1, classids)); 57 } 58 /* Process summary exclusions */ 59 PetscCall(PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt)); 60 if (opt) { 61 PetscCall(PetscStrInList("matmffd",logList,',',&pkg)); 62 if (pkg) PetscCall(PetscLogEventExcludeClass(MATMFFD_CLASSID)); 63 } 64 /* Register package finalizer */ 65 PetscCall(PetscRegisterFinalize(MatMFFDFinalizePackage)); 66 PetscFunctionReturn(0); 67 } 68 69 static PetscErrorCode MatMFFDSetType_MFFD(Mat mat,MatMFFDType ftype) 70 { 71 MatMFFD ctx; 72 PetscBool match; 73 PetscErrorCode (*r)(MatMFFD); 74 75 PetscFunctionBegin; 76 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 77 PetscValidCharPointer(ftype,2); 78 PetscCall(MatShellGetContext(mat,&ctx)); 79 80 /* already set, so just return */ 81 PetscCall(PetscObjectTypeCompare((PetscObject)ctx,ftype,&match)); 82 if (match) PetscFunctionReturn(0); 83 84 /* destroy the old one if it exists */ 85 if (ctx->ops->destroy) PetscCall((*ctx->ops->destroy)(ctx)); 86 87 PetscCall(PetscFunctionListFind(MatMFFDList,ftype,&r)); 88 PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype); 89 PetscCall((*r)(ctx)); 90 PetscCall(PetscObjectChangeTypeName((PetscObject)ctx,ftype)); 91 PetscFunctionReturn(0); 92 } 93 94 /*@C 95 MatMFFDSetType - Sets the method that is used to compute the 96 differencing parameter for finite differene matrix-free formulations. 97 98 Input Parameters: 99 + mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD() 100 or MatSetType(mat,MATMFFD); 101 - ftype - the type requested, either MATMFFD_WP or MATMFFD_DS 102 103 Level: advanced 104 105 Notes: 106 For example, such routines can compute h for use in 107 Jacobian-vector products of the form 108 109 F(x+ha) - F(x) 110 F'(u)a ~= ---------------- 111 h 112 113 .seealso: MatCreateSNESMF(), MatMFFDRegister(), MatMFFDSetFunction(), MatCreateMFFD() 114 @*/ 115 PetscErrorCode MatMFFDSetType(Mat mat,MatMFFDType ftype) 116 { 117 PetscFunctionBegin; 118 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 119 PetscValidCharPointer(ftype,2); 120 PetscTryMethod(mat,"MatMFFDSetType_C",(Mat,MatMFFDType),(mat,ftype)); 121 PetscFunctionReturn(0); 122 } 123 124 static PetscErrorCode MatGetDiagonal_MFFD(Mat,Vec); 125 126 typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/ 127 static PetscErrorCode MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func) 128 { 129 MatMFFD ctx; 130 131 PetscFunctionBegin; 132 PetscCall(MatShellGetContext(mat,&ctx)); 133 ctx->funcisetbase = func; 134 PetscFunctionReturn(0); 135 } 136 137 typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 138 static PetscErrorCode MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci) 139 { 140 MatMFFD ctx; 141 142 PetscFunctionBegin; 143 PetscCall(MatShellGetContext(mat,&ctx)); 144 ctx->funci = funci; 145 PetscCall(MatShellSetOperation(mat,MATOP_GET_DIAGONAL,(void (*)(void))MatGetDiagonal_MFFD)); 146 PetscFunctionReturn(0); 147 } 148 149 static PetscErrorCode MatMFFDGetH_MFFD(Mat mat,PetscScalar *h) 150 { 151 MatMFFD ctx; 152 153 PetscFunctionBegin; 154 PetscCall(MatShellGetContext(mat,&ctx)); 155 *h = ctx->currenth; 156 PetscFunctionReturn(0); 157 } 158 159 static PetscErrorCode MatMFFDResetHHistory_MFFD(Mat J) 160 { 161 MatMFFD ctx; 162 163 PetscFunctionBegin; 164 PetscCall(MatShellGetContext(J,&ctx)); 165 ctx->ncurrenth = 0; 166 PetscFunctionReturn(0); 167 } 168 169 /*@C 170 MatMFFDRegister - Adds a method to the MatMFFD registry. 171 172 Not Collective 173 174 Input Parameters: 175 + name_solver - name of a new user-defined compute-h module 176 - routine_create - routine to create method context 177 178 Level: developer 179 180 Notes: 181 MatMFFDRegister() may be called multiple times to add several user-defined solvers. 182 183 Sample usage: 184 .vb 185 MatMFFDRegister("my_h",MyHCreate); 186 .ve 187 188 Then, your solver can be chosen with the procedural interface via 189 $ MatMFFDSetType(mfctx,"my_h") 190 or at runtime via the option 191 $ -mat_mffd_type my_h 192 193 .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 194 @*/ 195 PetscErrorCode MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD)) 196 { 197 PetscFunctionBegin; 198 PetscCall(MatInitializePackage()); 199 PetscCall(PetscFunctionListAdd(&MatMFFDList,sname,function)); 200 PetscFunctionReturn(0); 201 } 202 203 /* ----------------------------------------------------------------------------------------*/ 204 static PetscErrorCode MatDestroy_MFFD(Mat mat) 205 { 206 MatMFFD ctx; 207 208 PetscFunctionBegin; 209 PetscCall(MatShellGetContext(mat,&ctx)); 210 PetscCall(VecDestroy(&ctx->w)); 211 PetscCall(VecDestroy(&ctx->current_u)); 212 if (ctx->current_f_allocated) { 213 PetscCall(VecDestroy(&ctx->current_f)); 214 } 215 if (ctx->ops->destroy) PetscCall((*ctx->ops->destroy)(ctx)); 216 PetscCall(PetscHeaderDestroy(&ctx)); 217 218 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C",NULL)); 219 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",NULL)); 220 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",NULL)); 221 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C",NULL)); 222 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C",NULL)); 223 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C",NULL)); 224 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C",NULL)); 225 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C",NULL)); 226 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetHHistory_C",NULL)); 227 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetType_C",NULL)); 228 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDGetH_C",NULL)); 229 PetscFunctionReturn(0); 230 } 231 232 /* 233 MatMFFDView_MFFD - Views matrix-free parameters. 234 235 */ 236 static PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer) 237 { 238 MatMFFD ctx; 239 PetscBool iascii, viewbase, viewfunction; 240 const char *prefix; 241 242 PetscFunctionBegin; 243 PetscCall(MatShellGetContext(J,&ctx)); 244 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 245 if (iascii) { 246 PetscCall(PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n")); 247 PetscCall(PetscViewerASCIIPushTab(viewer)); 248 PetscCall(PetscViewerASCIIPrintf(viewer,"err=%g (relative error in function evaluation)\n",(double)ctx->error_rel)); 249 if (!((PetscObject)ctx)->type_name) { 250 PetscCall(PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n")); 251 } else { 252 PetscCall(PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name)); 253 } 254 #if defined(PETSC_USE_COMPLEX) 255 if (ctx->usecomplex) { 256 PetscCall(PetscViewerASCIIPrintf(viewer,"Using Lyness complex number trick to compute the matrix-vector product\n")); 257 } 258 #endif 259 if (ctx->ops->view) { 260 PetscCall((*ctx->ops->view)(ctx,viewer)); 261 } 262 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)J, &prefix)); 263 264 PetscCall(PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_base", &viewbase)); 265 if (viewbase) { 266 PetscCall(PetscViewerASCIIPrintf(viewer, "Base:\n")); 267 PetscCall(VecView(ctx->current_u, viewer)); 268 } 269 PetscCall(PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_function", &viewfunction)); 270 if (viewfunction) { 271 PetscCall(PetscViewerASCIIPrintf(viewer, "Function:\n")); 272 PetscCall(VecView(ctx->current_f, viewer)); 273 } 274 PetscCall(PetscViewerASCIIPopTab(viewer)); 275 } 276 PetscFunctionReturn(0); 277 } 278 279 /* 280 MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This 281 allows the user to indicate the beginning of a new linear solve by calling 282 MatAssemblyXXX() on the matrix free matrix. This then allows the 283 MatCreateMFFD_WP() to properly compute ||U|| only the first time 284 in the linear solver rather than every time. 285 286 This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence 287 it must be labeled as PETSC_EXTERN 288 */ 289 PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt) 290 { 291 MatMFFD j; 292 293 PetscFunctionBegin; 294 PetscCall(MatShellGetContext(J,&j)); 295 PetscCall(MatMFFDResetHHistory(J)); 296 PetscFunctionReturn(0); 297 } 298 299 /* 300 MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a: 301 302 y ~= (F(u + ha) - F(u))/h, 303 where F = nonlinear function, as set by SNESSetFunction() 304 u = current iterate 305 h = difference interval 306 */ 307 static PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y) 308 { 309 MatMFFD ctx; 310 PetscScalar h; 311 Vec w,U,F; 312 PetscBool zeroa; 313 314 PetscFunctionBegin; 315 PetscCall(MatShellGetContext(mat,&ctx)); 316 PetscCheck(ctx->current_u,PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MatMFFDSetBase() has not been called, this is often caused by forgetting to call \n\t\tMatAssemblyBegin/End on the first Mat in the SNES compute function"); 317 /* We log matrix-free matrix-vector products separately, so that we can 318 separate the performance monitoring from the cases that use conventional 319 storage. We may eventually modify event logging to associate events 320 with particular objects, hence alleviating the more general problem. */ 321 PetscCall(PetscLogEventBegin(MATMFFD_Mult,a,y,0,0)); 322 323 w = ctx->w; 324 U = ctx->current_u; 325 F = ctx->current_f; 326 /* 327 Compute differencing parameter 328 */ 329 if (!((PetscObject)ctx)->type_name) { 330 PetscCall(MatMFFDSetType(mat,MATMFFD_WP)); 331 PetscCall(MatSetFromOptions(mat)); 332 } 333 PetscCall((*ctx->ops->compute)(ctx,U,a,&h,&zeroa)); 334 if (zeroa) { 335 PetscCall(VecSet(y,0.0)); 336 PetscFunctionReturn(0); 337 } 338 339 PetscCheck(!mat->erroriffailure || !PetscIsInfOrNanScalar(h),PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h"); 340 if (ctx->checkh) { 341 PetscCall((*ctx->checkh)(ctx->checkhctx,U,a,&h)); 342 } 343 344 /* keep a record of the current differencing parameter h */ 345 ctx->currenth = h; 346 #if defined(PETSC_USE_COMPLEX) 347 PetscCall(PetscInfo(mat,"Current differencing parameter: %g + %g i\n",(double)PetscRealPart(h),(double)PetscImaginaryPart(h))); 348 #else 349 PetscCall(PetscInfo(mat,"Current differencing parameter: %15.12e\n",(double)PetscRealPart(h))); 350 #endif 351 if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 352 ctx->historyh[ctx->ncurrenth] = h; 353 } 354 ctx->ncurrenth++; 355 356 #if defined(PETSC_USE_COMPLEX) 357 if (ctx->usecomplex) h = PETSC_i*h; 358 #endif 359 360 /* w = u + ha */ 361 PetscCall(VecWAXPY(w,h,a,U)); 362 363 /* compute func(U) as base for differencing; only needed first time in and not when provided by user */ 364 if (ctx->ncurrenth == 1 && ctx->current_f_allocated) { 365 PetscCall((*ctx->func)(ctx->funcctx,U,F)); 366 } 367 PetscCall((*ctx->func)(ctx->funcctx,w,y)); 368 369 #if defined(PETSC_USE_COMPLEX) 370 if (ctx->usecomplex) { 371 PetscCall(VecImaginaryPart(y)); 372 h = PetscImaginaryPart(h); 373 } else { 374 PetscCall(VecAXPY(y,-1.0,F)); 375 } 376 #else 377 PetscCall(VecAXPY(y,-1.0,F)); 378 #endif 379 PetscCall(VecScale(y,1.0/h)); 380 if (mat->nullsp) PetscCall(MatNullSpaceRemove(mat->nullsp,y)); 381 382 PetscCall(PetscLogEventEnd(MATMFFD_Mult,a,y,0,0)); 383 PetscFunctionReturn(0); 384 } 385 386 /* 387 MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix 388 389 y ~= (F(u + ha) - F(u))/h, 390 where F = nonlinear function, as set by SNESSetFunction() 391 u = current iterate 392 h = difference interval 393 */ 394 PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a) 395 { 396 MatMFFD ctx; 397 PetscScalar h,*aa,*ww,v; 398 PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 399 Vec w,U; 400 PetscInt i,rstart,rend; 401 402 PetscFunctionBegin; 403 PetscCall(MatShellGetContext(mat,&ctx)); 404 PetscCheck(ctx->func,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunction() first"); 405 PetscCheck(ctx->funci,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunctioni() first"); 406 w = ctx->w; 407 U = ctx->current_u; 408 PetscCall((*ctx->func)(ctx->funcctx,U,a)); 409 if (ctx->funcisetbase) { 410 PetscCall((*ctx->funcisetbase)(ctx->funcctx,U)); 411 } 412 PetscCall(VecCopy(U,w)); 413 414 PetscCall(VecGetOwnershipRange(a,&rstart,&rend)); 415 PetscCall(VecGetArray(a,&aa)); 416 for (i=rstart; i<rend; i++) { 417 PetscCall(VecGetArray(w,&ww)); 418 h = ww[i-rstart]; 419 if (h == 0.0) h = 1.0; 420 if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 421 else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 422 h *= epsilon; 423 424 ww[i-rstart] += h; 425 PetscCall(VecRestoreArray(w,&ww)); 426 PetscCall((*ctx->funci)(ctx->funcctx,i,w,&v)); 427 aa[i-rstart] = (v - aa[i-rstart])/h; 428 429 PetscCall(VecGetArray(w,&ww)); 430 ww[i-rstart] -= h; 431 PetscCall(VecRestoreArray(w,&ww)); 432 } 433 PetscCall(VecRestoreArray(a,&aa)); 434 PetscFunctionReturn(0); 435 } 436 437 PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F) 438 { 439 MatMFFD ctx; 440 441 PetscFunctionBegin; 442 PetscCall(MatShellGetContext(J,&ctx)); 443 PetscCall(MatMFFDResetHHistory(J)); 444 if (!ctx->current_u) { 445 PetscCall(VecDuplicate(U,&ctx->current_u)); 446 PetscCall(VecLockReadPush(ctx->current_u)); 447 } 448 PetscCall(VecLockReadPop(ctx->current_u)); 449 PetscCall(VecCopy(U,ctx->current_u)); 450 PetscCall(VecLockReadPush(ctx->current_u)); 451 if (F) { 452 if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f)); 453 ctx->current_f = F; 454 ctx->current_f_allocated = PETSC_FALSE; 455 } else if (!ctx->current_f_allocated) { 456 PetscCall(MatCreateVecs(J,NULL,&ctx->current_f)); 457 ctx->current_f_allocated = PETSC_TRUE; 458 } 459 if (!ctx->w) { 460 PetscCall(VecDuplicate(ctx->current_u,&ctx->w)); 461 } 462 J->assembled = PETSC_TRUE; 463 PetscFunctionReturn(0); 464 } 465 466 typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 467 468 static PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx) 469 { 470 MatMFFD ctx; 471 472 PetscFunctionBegin; 473 PetscCall(MatShellGetContext(J,&ctx)); 474 ctx->checkh = fun; 475 ctx->checkhctx = ectx; 476 PetscFunctionReturn(0); 477 } 478 479 /*@C 480 MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all 481 MatMFFD options in the database. 482 483 Collective on Mat 484 485 Input Parameters: 486 + A - the Mat context 487 - prefix - the prefix to prepend to all option names 488 489 Notes: 490 A hyphen (-) must NOT be given at the beginning of the prefix name. 491 The first character of all runtime options is AUTOMATICALLY the hyphen. 492 493 Level: advanced 494 495 .seealso: MatSetFromOptions(), MatCreateSNESMF(), MatCreateMFFD() 496 @*/ 497 PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat,const char prefix[]) 498 { 499 MatMFFD mfctx; 500 501 PetscFunctionBegin; 502 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 503 PetscCall(MatShellGetContext(mat,&mfctx)); 504 PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 505 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix)); 506 PetscFunctionReturn(0); 507 } 508 509 static PetscErrorCode MatSetFromOptions_MFFD(PetscOptionItems *PetscOptionsObject,Mat mat) 510 { 511 MatMFFD mfctx; 512 PetscBool flg; 513 char ftype[256]; 514 515 PetscFunctionBegin; 516 PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 517 PetscCall(MatShellGetContext(mat,&mfctx)); 518 PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,2); 519 PetscObjectOptionsBegin((PetscObject)mfctx); 520 PetscCall(PetscOptionsFList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg)); 521 if (flg) { 522 PetscCall(MatMFFDSetType(mat,ftype)); 523 } 524 525 PetscCall(PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,NULL)); 526 PetscCall(PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,NULL)); 527 528 flg = PETSC_FALSE; 529 PetscCall(PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL)); 530 if (flg) { 531 PetscCall(MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,NULL)); 532 } 533 #if defined(PETSC_USE_COMPLEX) 534 PetscCall(PetscOptionsBool("-mat_mffd_complex","Use Lyness complex number trick to compute the matrix-vector product","None",mfctx->usecomplex,&mfctx->usecomplex,NULL)); 535 #endif 536 if (mfctx->ops->setfromoptions) { 537 PetscCall((*mfctx->ops->setfromoptions)(PetscOptionsObject,mfctx)); 538 } 539 PetscOptionsEnd(); 540 PetscFunctionReturn(0); 541 } 542 543 static PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period) 544 { 545 MatMFFD ctx; 546 547 PetscFunctionBegin; 548 PetscCall(MatShellGetContext(mat,&ctx)); 549 ctx->recomputeperiod = period; 550 PetscFunctionReturn(0); 551 } 552 553 static PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 554 { 555 MatMFFD ctx; 556 557 PetscFunctionBegin; 558 PetscCall(MatShellGetContext(mat,&ctx)); 559 ctx->func = func; 560 ctx->funcctx = funcctx; 561 PetscFunctionReturn(0); 562 } 563 564 static PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error) 565 { 566 MatMFFD ctx; 567 568 PetscFunctionBegin; 569 PetscCall(MatShellGetContext(mat,&ctx)); 570 if (error != PETSC_DEFAULT) ctx->error_rel = error; 571 PetscFunctionReturn(0); 572 } 573 574 PetscErrorCode MatMFFDSetHHistory_MFFD(Mat J,PetscScalar history[],PetscInt nhistory) 575 { 576 MatMFFD ctx; 577 578 PetscFunctionBegin; 579 PetscCall(MatShellGetContext(J,&ctx)); 580 ctx->historyh = history; 581 ctx->maxcurrenth = nhistory; 582 ctx->currenth = 0.; 583 PetscFunctionReturn(0); 584 } 585 586 /*MC 587 MATMFFD - MATMFFD = "mffd" - A matrix free matrix type. 588 589 Level: advanced 590 591 Developers Note: This is implemented on top of MATSHELL to get support for scaling and shifting without requiring duplicate code 592 593 .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction(), MatMFFDSetType(), 594 MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 595 MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 596 MatMFFDGetH(), 597 M*/ 598 PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A) 599 { 600 MatMFFD mfctx; 601 602 PetscFunctionBegin; 603 PetscCall(MatMFFDInitializePackage()); 604 605 PetscCall(PetscHeaderCreate(mfctx,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),NULL,NULL)); 606 607 mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 608 mfctx->recomputeperiod = 1; 609 mfctx->count = 0; 610 mfctx->currenth = 0.0; 611 mfctx->historyh = NULL; 612 mfctx->ncurrenth = 0; 613 mfctx->maxcurrenth = 0; 614 ((PetscObject)mfctx)->type_name = NULL; 615 616 /* 617 Create the empty data structure to contain compute-h routines. 618 These will be filled in below from the command line options or 619 a later call with MatMFFDSetType() or if that is not called 620 then it will default in the first use of MatMult_MFFD() 621 */ 622 mfctx->ops->compute = NULL; 623 mfctx->ops->destroy = NULL; 624 mfctx->ops->view = NULL; 625 mfctx->ops->setfromoptions = NULL; 626 mfctx->hctx = NULL; 627 628 mfctx->func = NULL; 629 mfctx->funcctx = NULL; 630 mfctx->w = NULL; 631 mfctx->mat = A; 632 633 PetscCall(MatSetType(A,MATSHELL)); 634 PetscCall(MatShellSetContext(A,mfctx)); 635 PetscCall(MatShellSetOperation(A,MATOP_MULT,(void (*)(void))MatMult_MFFD)); 636 PetscCall(MatShellSetOperation(A,MATOP_DESTROY,(void (*)(void))MatDestroy_MFFD)); 637 PetscCall(MatShellSetOperation(A,MATOP_VIEW,(void (*)(void))MatView_MFFD)); 638 PetscCall(MatShellSetOperation(A,MATOP_ASSEMBLY_END,(void (*)(void))MatAssemblyEnd_MFFD)); 639 PetscCall(MatShellSetOperation(A,MATOP_SET_FROM_OPTIONS,(void (*)(void))MatSetFromOptions_MFFD)); 640 641 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD)); 642 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD)); 643 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD)); 644 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD)); 645 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD)); 646 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD)); 647 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD)); 648 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD)); 649 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetHHistory_C",MatMFFDSetHHistory_MFFD)); 650 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetType_C",MatMFFDSetType_MFFD)); 651 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDGetH_C",MatMFFDGetH_MFFD)); 652 PetscCall(PetscObjectChangeTypeName((PetscObject)A,MATMFFD)); 653 PetscFunctionReturn(0); 654 } 655 656 /*@ 657 MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF() 658 659 Collective on Vec 660 661 Input Parameters: 662 + comm - MPI communicator 663 . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 664 This value should be the same as the local size used in creating the 665 y vector for the matrix-vector product y = Ax. 666 . n - This value should be the same as the local size used in creating the 667 x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 668 calculated if N is given) For square matrices n is almost always m. 669 . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 670 - N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 671 672 Output Parameter: 673 . J - the matrix-free matrix 674 675 Options Database Keys: call MatSetFromOptions() to trigger these 676 + -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS) 677 . -mat_mffd_err - square root of estimated relative error in function evaluation 678 . -mat_mffd_period - how often h is recomputed, defaults to 1, everytime 679 . -mat_mffd_check_positivity - possibly decrease h until U + h*a has only positive values 680 - -mat_mffd_complex - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing 681 (requires real valued functions but that PETSc be configured for complex numbers) 682 683 Level: advanced 684 685 Notes: 686 The matrix-free matrix context merely contains the function pointers 687 and work space for performing finite difference approximations of 688 Jacobian-vector products, F'(u)*a, 689 690 The default code uses the following approach to compute h 691 692 .vb 693 F'(u)*a = [F(u+h*a) - F(u)]/h where 694 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 695 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 696 where 697 error_rel = square root of relative error in function evaluation 698 umin = minimum iterate parameter 699 .ve 700 701 You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different 702 preconditioner matrix 703 704 The user can set the error_rel via MatMFFDSetFunctionError() and 705 umin via MatMFFDDSSetUmin(); see Users-Manual: ch_snes for details. 706 707 The user should call MatDestroy() when finished with the matrix-free 708 matrix context. 709 710 Options Database Keys: 711 + -mat_mffd_err <error_rel> - Sets error_rel 712 . -mat_mffd_umin <umin> - Sets umin (for default PETSc routine that computes h only) 713 - -mat_mffd_check_positivity - check positivity 714 715 .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 716 MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 717 MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian() 718 719 @*/ 720 PetscErrorCode MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J) 721 { 722 PetscFunctionBegin; 723 PetscCall(MatCreate(comm,J)); 724 PetscCall(MatSetSizes(*J,m,n,M,N)); 725 PetscCall(MatSetType(*J,MATMFFD)); 726 PetscCall(MatSetUp(*J)); 727 PetscFunctionReturn(0); 728 } 729 730 /*@ 731 MatMFFDGetH - Gets the last value that was used as the differencing 732 parameter. 733 734 Not Collective 735 736 Input Parameters: 737 . mat - the matrix obtained with MatCreateSNESMF() 738 739 Output Parameter: 740 . h - the differencing step size 741 742 Level: advanced 743 744 .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory() 745 @*/ 746 PetscErrorCode MatMFFDGetH(Mat mat,PetscScalar *h) 747 { 748 PetscFunctionBegin; 749 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 750 PetscValidScalarPointer(h,2); 751 PetscUseMethod(mat,"MatMFFDGetH_C",(Mat,PetscScalar*),(mat,h)); 752 PetscFunctionReturn(0); 753 } 754 755 /*@C 756 MatMFFDSetFunction - Sets the function used in applying the matrix free. 757 758 Logically Collective on Mat 759 760 Input Parameters: 761 + mat - the matrix free matrix created via MatCreateSNESMF() or MatCreateMFFD() 762 . func - the function to use 763 - funcctx - optional function context passed to function 764 765 Calling Sequence of func: 766 $ func (void *funcctx, Vec x, Vec f) 767 768 + funcctx - user provided context 769 . x - input vector 770 - f - computed output function 771 772 Level: advanced 773 774 Notes: 775 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 776 matrix inside your compute Jacobian routine 777 778 If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used. 779 780 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD, 781 MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 782 @*/ 783 PetscErrorCode MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 784 { 785 PetscFunctionBegin; 786 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 787 PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx)); 788 PetscFunctionReturn(0); 789 } 790 791 /*@C 792 MatMFFDSetFunctioni - Sets the function for a single component 793 794 Logically Collective on Mat 795 796 Input Parameters: 797 + mat - the matrix free matrix created via MatCreateSNESMF() 798 - funci - the function to use 799 800 Level: advanced 801 802 Notes: 803 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 804 matrix inside your compute Jacobian routine. 805 This function is necessary to compute the diagonal of the matrix. 806 funci must not contain any MPI call as it is called inside a loop on the local portion of the vector. 807 808 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal() 809 810 @*/ 811 PetscErrorCode MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*)) 812 { 813 PetscFunctionBegin; 814 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 815 PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci)); 816 PetscFunctionReturn(0); 817 } 818 819 /*@C 820 MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation 821 822 Logically Collective on Mat 823 824 Input Parameters: 825 + mat - the matrix free matrix created via MatCreateSNESMF() 826 - func - the function to use 827 828 Level: advanced 829 830 Notes: 831 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 832 matrix inside your compute Jacobian routine. 833 This function is necessary to compute the diagonal of the matrix. 834 835 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 836 MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal() 837 @*/ 838 PetscErrorCode MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec)) 839 { 840 PetscFunctionBegin; 841 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 842 PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func)); 843 PetscFunctionReturn(0); 844 } 845 846 /*@ 847 MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime 848 849 Logically Collective on Mat 850 851 Input Parameters: 852 + mat - the matrix free matrix created via MatCreateSNESMF() 853 - period - 1 for everytime, 2 for every second etc 854 855 Options Database Keys: 856 . -mat_mffd_period <period> - Sets how often h is recomputed 857 858 Level: advanced 859 860 .seealso: MatCreateSNESMF(),MatMFFDGetH(), 861 MatMFFDSetHHistory(), MatMFFDResetHHistory() 862 @*/ 863 PetscErrorCode MatMFFDSetPeriod(Mat mat,PetscInt period) 864 { 865 PetscFunctionBegin; 866 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 867 PetscValidLogicalCollectiveInt(mat,period,2); 868 PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period)); 869 PetscFunctionReturn(0); 870 } 871 872 /*@ 873 MatMFFDSetFunctionError - Sets the error_rel for the approximation of 874 matrix-vector products using finite differences. 875 876 Logically Collective on Mat 877 878 Input Parameters: 879 + mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF() 880 - error_rel - relative error (should be set to the square root of 881 the relative error in the function evaluations) 882 883 Options Database Keys: 884 . -mat_mffd_err <error_rel> - Sets error_rel 885 886 Level: advanced 887 888 Notes: 889 The default matrix-free matrix-vector product routine computes 890 .vb 891 F'(u)*a = [F(u+h*a) - F(u)]/h where 892 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 893 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 894 .ve 895 896 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 897 MatMFFDSetHHistory(), MatMFFDResetHHistory() 898 @*/ 899 PetscErrorCode MatMFFDSetFunctionError(Mat mat,PetscReal error) 900 { 901 PetscFunctionBegin; 902 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 903 PetscValidLogicalCollectiveReal(mat,error,2); 904 PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error)); 905 PetscFunctionReturn(0); 906 } 907 908 /*@ 909 MatMFFDSetHHistory - Sets an array to collect a history of the 910 differencing values (h) computed for the matrix-free product. 911 912 Logically Collective on Mat 913 914 Input Parameters: 915 + J - the matrix-free matrix context 916 . history - space to hold the history 917 - nhistory - number of entries in history, if more entries are generated than 918 nhistory, then the later ones are discarded 919 920 Level: advanced 921 922 Notes: 923 Use MatMFFDResetHHistory() to reset the history counter and collect 924 a new batch of differencing parameters, h. 925 926 .seealso: MatMFFDGetH(), MatCreateSNESMF(), 927 MatMFFDResetHHistory(), MatMFFDSetFunctionError() 928 929 @*/ 930 PetscErrorCode MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory) 931 { 932 PetscFunctionBegin; 933 PetscValidHeaderSpecific(J,MAT_CLASSID,1); 934 if (history) PetscValidScalarPointer(history,2); 935 PetscValidLogicalCollectiveInt(J,nhistory,3); 936 PetscUseMethod(J,"MatMFFDSetHHistory_C",(Mat,PetscScalar[],PetscInt),(J,history,nhistory)); 937 PetscFunctionReturn(0); 938 } 939 940 /*@ 941 MatMFFDResetHHistory - Resets the counter to zero to begin 942 collecting a new set of differencing histories. 943 944 Logically Collective on Mat 945 946 Input Parameters: 947 . J - the matrix-free matrix context 948 949 Level: advanced 950 951 Notes: 952 Use MatMFFDSetHHistory() to create the original history counter. 953 954 .seealso: MatMFFDGetH(), MatCreateSNESMF(), 955 MatMFFDSetHHistory(), MatMFFDSetFunctionError() 956 957 @*/ 958 PetscErrorCode MatMFFDResetHHistory(Mat J) 959 { 960 PetscFunctionBegin; 961 PetscValidHeaderSpecific(J,MAT_CLASSID,1); 962 PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J)); 963 PetscFunctionReturn(0); 964 } 965 966 /*@ 967 MatMFFDSetBase - Sets the vector U at which matrix vector products of the 968 Jacobian are computed 969 970 Logically Collective on Mat 971 972 Input Parameters: 973 + J - the MatMFFD matrix 974 . U - the vector 975 - F - (optional) vector that contains F(u) if it has been already computed 976 977 Notes: 978 This is rarely used directly 979 980 If F is provided then it is not recomputed. Otherwise the function is evaluated at the base 981 point during the first MatMult() after each call to MatMFFDSetBase(). 982 983 Level: advanced 984 985 @*/ 986 PetscErrorCode MatMFFDSetBase(Mat J,Vec U,Vec F) 987 { 988 PetscFunctionBegin; 989 PetscValidHeaderSpecific(J,MAT_CLASSID,1); 990 PetscValidHeaderSpecific(U,VEC_CLASSID,2); 991 if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3); 992 PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F)); 993 PetscFunctionReturn(0); 994 } 995 996 /*@C 997 MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts 998 it to satisfy some criteria 999 1000 Logically Collective on Mat 1001 1002 Input Parameters: 1003 + J - the MatMFFD matrix 1004 . fun - the function that checks h 1005 - ctx - any context needed by the function 1006 1007 Options Database Keys: 1008 . -mat_mffd_check_positivity <bool> - Insure that U + h*a is non-negative 1009 1010 Level: advanced 1011 1012 Notes: 1013 For example, MatMFFDCheckPositivity() insures that all entries 1014 of U + h*a are non-negative 1015 1016 The function you provide is called after the default h has been computed and allows you to 1017 modify it. 1018 1019 .seealso: MatMFFDCheckPositivity() 1020 @*/ 1021 PetscErrorCode MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx) 1022 { 1023 PetscFunctionBegin; 1024 PetscValidHeaderSpecific(J,MAT_CLASSID,1); 1025 PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx)); 1026 PetscFunctionReturn(0); 1027 } 1028 1029 /*@ 1030 MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or 1031 zero, decreases h until this is satisfied. 1032 1033 Logically Collective on Vec 1034 1035 Input Parameters: 1036 + U - base vector that is added to 1037 . a - vector that is added 1038 . h - scaling factor on a 1039 - dummy - context variable (unused) 1040 1041 Options Database Keys: 1042 . -mat_mffd_check_positivity <bool> - Insure that U + h*a is nonnegative 1043 1044 Level: advanced 1045 1046 Notes: 1047 This is rarely used directly, rather it is passed as an argument to 1048 MatMFFDSetCheckh() 1049 1050 .seealso: MatMFFDSetCheckh() 1051 @*/ 1052 PetscErrorCode MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h) 1053 { 1054 PetscReal val, minval; 1055 PetscScalar *u_vec, *a_vec; 1056 PetscInt i,n; 1057 MPI_Comm comm; 1058 1059 PetscFunctionBegin; 1060 PetscValidHeaderSpecific(U,VEC_CLASSID,2); 1061 PetscValidHeaderSpecific(a,VEC_CLASSID,3); 1062 PetscValidScalarPointer(h,4); 1063 PetscCall(PetscObjectGetComm((PetscObject)U,&comm)); 1064 PetscCall(VecGetArray(U,&u_vec)); 1065 PetscCall(VecGetArray(a,&a_vec)); 1066 PetscCall(VecGetLocalSize(U,&n)); 1067 minval = PetscAbsScalar(*h)*PetscRealConstant(1.01); 1068 for (i=0; i<n; i++) { 1069 if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) { 1070 val = PetscAbsScalar(u_vec[i]/a_vec[i]); 1071 if (val < minval) minval = val; 1072 } 1073 } 1074 PetscCall(VecRestoreArray(U,&u_vec)); 1075 PetscCall(VecRestoreArray(a,&a_vec)); 1076 PetscCall(MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm)); 1077 if (val <= PetscAbsScalar(*h)) { 1078 PetscCall(PetscInfo(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val))); 1079 if (PetscRealPart(*h) > 0.0) *h = 0.99*val; 1080 else *h = -0.99*val; 1081 } 1082 PetscFunctionReturn(0); 1083 } 1084