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