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