1*a1f56445SPierre Jolivet #include <../src/mat/impls/shell/shell.h> /*I "petscmat.h" I*/ 2*a1f56445SPierre Jolivet #include <../src/mat/impls/mffd/mffdimpl.h> 3e884886fSBarry Smith 4f4259b30SLisandro Dalcin PetscFunctionList MatMFFDList = NULL; 5ace3abfcSBarry Smith PetscBool MatMFFDRegisterAllCalled = PETSC_FALSE; 6e884886fSBarry Smith 77087cfbeSBarry Smith PetscClassId MATMFFD_CLASSID; 8166c7f25SBarry Smith PetscLogEvent MATMFFD_Mult; 9e884886fSBarry Smith 10ace3abfcSBarry Smith static PetscBool MatMFFDPackageInitialized = PETSC_FALSE; 1166976f2fSJacob Faibussowitsch 12b022a5c1SBarry Smith /*@C 1311a5261eSBarry Smith MatMFFDFinalizePackage - This function destroys everything in the MATMFFD` package. It is 1411a5261eSBarry Smith called from `PetscFinalize()`. 15b022a5c1SBarry Smith 16b022a5c1SBarry Smith Level: developer 17b022a5c1SBarry Smith 181cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `PetscFinalize()`, `MatCreateMFFD()`, `MatCreateSNESMF()` 19b022a5c1SBarry Smith @*/ 20d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDFinalizePackage(void) 21d71ae5a4SJacob Faibussowitsch { 22b022a5c1SBarry Smith PetscFunctionBegin; 239566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&MatMFFDList)); 24b022a5c1SBarry Smith MatMFFDPackageInitialized = PETSC_FALSE; 25b022a5c1SBarry Smith MatMFFDRegisterAllCalled = PETSC_FALSE; 263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 27b022a5c1SBarry Smith } 28b022a5c1SBarry Smith 293ec795f1SBarry Smith /*@C 3011a5261eSBarry Smith MatMFFDInitializePackage - This function initializes everything in the MATMFFD` package. It is called 3111a5261eSBarry Smith from `MatInitializePackage()`. 323ec795f1SBarry Smith 333ec795f1SBarry Smith Level: developer 343ec795f1SBarry Smith 351cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `PetscInitialize()` 363ec795f1SBarry Smith @*/ 37d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDInitializePackage(void) 38d71ae5a4SJacob Faibussowitsch { 393ec795f1SBarry Smith char logList[256]; 408e81d068SLisandro Dalcin PetscBool opt, pkg; 413ec795f1SBarry Smith 423ec795f1SBarry Smith PetscFunctionBegin; 433ba16761SJacob Faibussowitsch if (MatMFFDPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 44b022a5c1SBarry Smith MatMFFDPackageInitialized = PETSC_TRUE; 453ec795f1SBarry Smith /* Register Classes */ 469566063dSJacob Faibussowitsch PetscCall(PetscClassIdRegister("MatMFFD", &MATMFFD_CLASSID)); 473ec795f1SBarry Smith /* Register Constructors */ 489566063dSJacob Faibussowitsch PetscCall(MatMFFDRegisterAll()); 493ec795f1SBarry Smith /* Register Events */ 509566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("MatMult MF", MATMFFD_CLASSID, &MATMFFD_Mult)); 51e94e781bSJacob Faibussowitsch /* Process Info */ 52e94e781bSJacob Faibussowitsch { 53e94e781bSJacob Faibussowitsch PetscClassId classids[1]; 54e94e781bSJacob Faibussowitsch 55e94e781bSJacob Faibussowitsch classids[0] = MATMFFD_CLASSID; 569566063dSJacob Faibussowitsch PetscCall(PetscInfoProcessClass("matmffd", 1, classids)); 573ec795f1SBarry Smith } 583ec795f1SBarry Smith /* Process summary exclusions */ 599566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log_exclude", logList, sizeof(logList), &opt)); 603ec795f1SBarry Smith if (opt) { 619566063dSJacob Faibussowitsch PetscCall(PetscStrInList("matmffd", logList, ',', &pkg)); 629566063dSJacob Faibussowitsch if (pkg) PetscCall(PetscLogEventExcludeClass(MATMFFD_CLASSID)); 633ec795f1SBarry Smith } 648e81d068SLisandro Dalcin /* Register package finalizer */ 659566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(MatMFFDFinalizePackage)); 663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 673ec795f1SBarry Smith } 683ec795f1SBarry Smith 69d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetType_MFFD(Mat mat, MatMFFDType ftype) 70d71ae5a4SJacob Faibussowitsch { 71789d8953SBarry Smith MatMFFD ctx; 72789d8953SBarry Smith PetscBool match; 735f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(MatMFFD); 74789d8953SBarry Smith 75789d8953SBarry Smith PetscFunctionBegin; 76789d8953SBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 774f572ea9SToby Isaac PetscAssertPointer(ftype, 2); 789566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 79789d8953SBarry Smith 80789d8953SBarry Smith /* already set, so just return */ 819566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ctx, ftype, &match)); 823ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 83789d8953SBarry Smith 84789d8953SBarry Smith /* destroy the old one if it exists */ 85dbbe0bcdSBarry Smith PetscTryTypeMethod(ctx, destroy); 86789d8953SBarry Smith 879566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(MatMFFDList, ftype, &r)); 886adde796SStefano Zampini PetscCheck(r, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown MatMFFD type %s given", ftype); 899566063dSJacob Faibussowitsch PetscCall((*r)(ctx)); 909566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)ctx, ftype)); 913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 92789d8953SBarry Smith } 93789d8953SBarry Smith 94e884886fSBarry Smith /*@C 95e884886fSBarry Smith MatMFFDSetType - Sets the method that is used to compute the 96da81f932SPierre Jolivet differencing parameter for finite difference matrix-free formulations. 97e884886fSBarry Smith 98e884886fSBarry Smith Input Parameters: 9911a5261eSBarry Smith + mat - the "matrix-free" matrix created via `MatCreateSNESMF()`, or `MatCreateMFFD()` 10011a5261eSBarry Smith or `MatSetType`(mat,`MATMFFD`); 10111a5261eSBarry Smith - ftype - the type requested, either `MATMFFD_WP` or `MATMFFD_DS` 102e884886fSBarry Smith 103e884886fSBarry Smith Level: advanced 104e884886fSBarry Smith 10511a5261eSBarry Smith Note: 1062ef1f0ffSBarry Smith For example, such routines can compute `h` for use in 107e884886fSBarry Smith Jacobian-vector products of the form 1082ef1f0ffSBarry Smith .vb 109e884886fSBarry Smith 110e884886fSBarry Smith F(x+ha) - F(x) 111e884886fSBarry Smith F'(u)a ~= ---------------- 112e884886fSBarry Smith h 1132ef1f0ffSBarry Smith .ve 114e884886fSBarry Smith 1151cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MATMFFD_WP`, `MATMFFD_DS`, `MatCreateSNESMF()`, `MatMFFDRegister()`, `MatMFFDSetFunction()`, `MatCreateMFFD()` 116e884886fSBarry Smith @*/ 117d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetType(Mat mat, MatMFFDType ftype) 118d71ae5a4SJacob Faibussowitsch { 119e884886fSBarry Smith PetscFunctionBegin; 1200700a824SBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 1214f572ea9SToby Isaac PetscAssertPointer(ftype, 2); 122cac4c232SBarry Smith PetscTryMethod(mat, "MatMFFDSetType_C", (Mat, MatMFFDType), (mat, ftype)); 1233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 124e884886fSBarry Smith } 125e884886fSBarry Smith 1265d21e1f8SStefano Zampini static PetscErrorCode MatGetDiagonal_MFFD(Mat, Vec); 1275d21e1f8SStefano Zampini 128e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void *, Vec); /* force argument to next function to not be extern C*/ 129d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunctioniBase_MFFD(Mat mat, FCN1 func) 130d71ae5a4SJacob Faibussowitsch { 131789d8953SBarry Smith MatMFFD ctx; 132e884886fSBarry Smith 133e884886fSBarry Smith PetscFunctionBegin; 1349566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 135e884886fSBarry Smith ctx->funcisetbase = func; 1363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 137e884886fSBarry Smith } 138e884886fSBarry Smith 139e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void *, PetscInt, Vec, PetscScalar *); /* force argument to next function to not be extern C*/ 140d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunctioni_MFFD(Mat mat, FCN2 funci) 141d71ae5a4SJacob Faibussowitsch { 142789d8953SBarry Smith MatMFFD ctx; 143e884886fSBarry Smith 144e884886fSBarry Smith PetscFunctionBegin; 1459566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 146e884886fSBarry Smith ctx->funci = funci; 1479566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(mat, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_MFFD)); 1483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1495d21e1f8SStefano Zampini } 150789d8953SBarry Smith 151d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDGetH_MFFD(Mat mat, PetscScalar *h) 152d71ae5a4SJacob Faibussowitsch { 153789d8953SBarry Smith MatMFFD ctx; 154789d8953SBarry Smith 155789d8953SBarry Smith PetscFunctionBegin; 1569566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 157789d8953SBarry Smith *h = ctx->currenth; 1583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 159e884886fSBarry Smith } 160e884886fSBarry Smith 161d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDResetHHistory_MFFD(Mat J) 162d71ae5a4SJacob Faibussowitsch { 163789d8953SBarry Smith MatMFFD ctx; 164bc0f08ceSBarry Smith 165bc0f08ceSBarry Smith PetscFunctionBegin; 1669566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J, &ctx)); 167bc0f08ceSBarry Smith ctx->ncurrenth = 0; 1683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 169bc0f08ceSBarry Smith } 170e884886fSBarry Smith 1711c84c290SBarry Smith /*@C 17220f4b53cSBarry Smith MatMFFDRegister - Adds a method to the `MATMFFD` registry. 1731c84c290SBarry Smith 1741c84c290SBarry Smith Not Collective 1751c84c290SBarry Smith 1761c84c290SBarry Smith Input Parameters: 17720f4b53cSBarry Smith + sname - name of a new user-defined compute-h module 17820f4b53cSBarry Smith - function - routine to create method context 1791c84c290SBarry Smith 1801c84c290SBarry Smith Level: developer 1811c84c290SBarry Smith 18211a5261eSBarry Smith Note: 18311a5261eSBarry Smith `MatMFFDRegister()` may be called multiple times to add several user-defined solvers. 1841c84c290SBarry Smith 185fe59aa6dSJacob Faibussowitsch Example Usage: 1861c84c290SBarry Smith .vb 187bdf89e91SBarry Smith MatMFFDRegister("my_h", MyHCreate); 1881c84c290SBarry Smith .ve 1891c84c290SBarry Smith 1901c84c290SBarry Smith Then, your solver can be chosen with the procedural interface via 19111a5261eSBarry Smith $ `MatMFFDSetType`(mfctx, "my_h") 1921c84c290SBarry Smith or at runtime via the option 193be7a6d03SBarry Smith $ -mat_mffd_type my_h 1941c84c290SBarry Smith 1951cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDRegisterAll()`, `MatMFFDRegisterDestroy()` 1961c84c290SBarry Smith @*/ 197d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDRegister(const char sname[], PetscErrorCode (*function)(MatMFFD)) 198d71ae5a4SJacob Faibussowitsch { 199e884886fSBarry Smith PetscFunctionBegin; 2009566063dSJacob Faibussowitsch PetscCall(MatInitializePackage()); 2019566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&MatMFFDList, sname, function)); 2023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 203e884886fSBarry Smith } 204e884886fSBarry Smith 205d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_MFFD(Mat mat) 206d71ae5a4SJacob Faibussowitsch { 207789d8953SBarry Smith MatMFFD ctx; 208e884886fSBarry Smith 209e884886fSBarry Smith PetscFunctionBegin; 2109566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 2119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->w)); 2129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->current_u)); 21348a46eb9SPierre Jolivet if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f)); 214dbbe0bcdSBarry Smith PetscTryTypeMethod(ctx, destroy); 2159566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(&ctx)); 216e884886fSBarry Smith 2179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetBase_C", NULL)); 2189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctioniBase_C", NULL)); 2199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctioni_C", NULL)); 2209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunction_C", NULL)); 2219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctionError_C", NULL)); 2229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetCheckh_C", NULL)); 2239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetPeriod_C", NULL)); 2249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDResetHHistory_C", NULL)); 2259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetHHistory_C", NULL)); 2269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetType_C", NULL)); 2279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDGetH_C", NULL)); 2282e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatSNESMFSetReuseBase_C", NULL)); 2292e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatSNESMFGetReuseBase_C", NULL)); 230ee912ba9SBarry Smith PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatShellSetContext_C", NULL)); 2313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 232e884886fSBarry Smith } 233e884886fSBarry Smith 234e884886fSBarry Smith /* 235e884886fSBarry Smith MatMFFDView_MFFD - Views matrix-free parameters. 236e884886fSBarry Smith 237e884886fSBarry Smith */ 238d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_MFFD(Mat J, PetscViewer viewer) 239d71ae5a4SJacob Faibussowitsch { 240789d8953SBarry Smith MatMFFD ctx; 241a283635eSDmitry Karpeev PetscBool iascii, viewbase, viewfunction; 242a283635eSDmitry Karpeev const char *prefix; 243e884886fSBarry Smith 244e884886fSBarry Smith PetscFunctionBegin; 2459566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J, &ctx)); 2469566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 247e884886fSBarry Smith if (iascii) { 2489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Matrix-free approximation:\n")); 2499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "err=%g (relative error in function evaluation)\n", (double)ctx->error_rel)); 2517adad957SLisandro Dalcin if (!((PetscObject)ctx)->type_name) { 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "The compute h routine has not yet been set\n")); 253e884886fSBarry Smith } else { 2549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Using %s compute h routine\n", ((PetscObject)ctx)->type_name)); 255e884886fSBarry Smith } 256c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 25748a46eb9SPierre Jolivet if (ctx->usecomplex) PetscCall(PetscViewerASCIIPrintf(viewer, "Using Lyness complex number trick to compute the matrix-vector product\n")); 258c92e3469SBarry Smith #endif 259dbbe0bcdSBarry Smith PetscTryTypeMethod(ctx, view, viewer); 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)J, &prefix)); 261a283635eSDmitry Karpeev 2629566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(((PetscObject)J)->options, prefix, "-mat_mffd_view_base", &viewbase)); 263a283635eSDmitry Karpeev if (viewbase) { 2649566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Base:\n")); 2659566063dSJacob Faibussowitsch PetscCall(VecView(ctx->current_u, viewer)); 266a283635eSDmitry Karpeev } 2679566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(((PetscObject)J)->options, prefix, "-mat_mffd_view_function", &viewfunction)); 268a283635eSDmitry Karpeev if (viewfunction) { 2699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Function:\n")); 2709566063dSJacob Faibussowitsch PetscCall(VecView(ctx->current_f, viewer)); 271a283635eSDmitry Karpeev } 2729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 27311aeaf0aSBarry Smith } 2743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 275e884886fSBarry Smith } 276e884886fSBarry Smith 277e884886fSBarry Smith /* 278e884886fSBarry Smith MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This 279e884886fSBarry Smith allows the user to indicate the beginning of a new linear solve by calling 28001c1178eSBarry Smith MatAssemblyXXX() on the matrix-free matrix. This then allows the 2811d0fab5eSBarry Smith MatCreateMFFD_WP() to properly compute ||U|| only the first time 282e884886fSBarry Smith in the linear solver rather than every time. 2835a576424SJed Brown 284cc2e6a90SBarry Smith This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence 285cc2e6a90SBarry Smith it must be labeled as PETSC_EXTERN 286e884886fSBarry Smith */ 287d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J, MatAssemblyType mt) 288d71ae5a4SJacob Faibussowitsch { 289789d8953SBarry Smith MatMFFD j; 290e884886fSBarry Smith 291e884886fSBarry Smith PetscFunctionBegin; 2929566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J, &j)); 2939566063dSJacob Faibussowitsch PetscCall(MatMFFDResetHHistory(J)); 2943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 295e884886fSBarry Smith } 296e884886fSBarry Smith 297e884886fSBarry Smith /* 298e884886fSBarry Smith MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a: 299e884886fSBarry Smith 300e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 301e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 302e884886fSBarry Smith u = current iterate 303e884886fSBarry Smith h = difference interval 304e884886fSBarry Smith */ 305d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_MFFD(Mat mat, Vec a, Vec y) 306d71ae5a4SJacob Faibussowitsch { 307789d8953SBarry Smith MatMFFD ctx; 308e884886fSBarry Smith PetscScalar h; 309e884886fSBarry Smith Vec w, U, F; 310ace3abfcSBarry Smith PetscBool zeroa; 311e884886fSBarry Smith 312e884886fSBarry Smith PetscFunctionBegin; 3139566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 31428b400f6SJacob Faibussowitsch 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"); 315e884886fSBarry Smith /* We log matrix-free matrix-vector products separately, so that we can 316e884886fSBarry Smith separate the performance monitoring from the cases that use conventional 317e884886fSBarry Smith storage. We may eventually modify event logging to associate events 318e884886fSBarry Smith with particular objects, hence alleviating the more general problem. */ 3199566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MATMFFD_Mult, a, y, 0, 0)); 320e884886fSBarry Smith 321e884886fSBarry Smith w = ctx->w; 322e884886fSBarry Smith U = ctx->current_u; 3233ec795f1SBarry Smith F = ctx->current_f; 324e884886fSBarry Smith /* 325e884886fSBarry Smith Compute differencing parameter 326e884886fSBarry Smith */ 3274a2f8832SBarry Smith if (!((PetscObject)ctx)->type_name) { 3289566063dSJacob Faibussowitsch PetscCall(MatMFFDSetType(mat, MATMFFD_WP)); 3299566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(mat)); 330e884886fSBarry Smith } 331dbbe0bcdSBarry Smith PetscUseTypeMethod(ctx, compute, U, a, &h, &zeroa); 332e884886fSBarry Smith if (zeroa) { 3339566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 33485d14658SLisandro Dalcin PetscCall(PetscLogEventEnd(MATMFFD_Mult, a, y, 0, 0)); 3353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 336e884886fSBarry Smith } 337e884886fSBarry Smith 33808401ef6SPierre Jolivet PetscCheck(!mat->erroriffailure || !PetscIsInfOrNanScalar(h), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Computed Nan differencing parameter h"); 33948a46eb9SPierre Jolivet if (ctx->checkh) PetscCall((*ctx->checkh)(ctx->checkhctx, U, a, &h)); 340e884886fSBarry Smith 341e884886fSBarry Smith /* keep a record of the current differencing parameter h */ 342e884886fSBarry Smith ctx->currenth = h; 343e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX) 3449566063dSJacob Faibussowitsch PetscCall(PetscInfo(mat, "Current differencing parameter: %g + %g i\n", (double)PetscRealPart(h), (double)PetscImaginaryPart(h))); 345e884886fSBarry Smith #else 3469566063dSJacob Faibussowitsch PetscCall(PetscInfo(mat, "Current differencing parameter: %15.12e\n", (double)PetscRealPart(h))); 347e884886fSBarry Smith #endif 348ad540459SPierre Jolivet if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) ctx->historyh[ctx->ncurrenth] = h; 349e884886fSBarry Smith ctx->ncurrenth++; 350e884886fSBarry Smith 351c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 352c92e3469SBarry Smith if (ctx->usecomplex) h = PETSC_i * h; 353c92e3469SBarry Smith #endif 354c92e3469SBarry Smith 355e884886fSBarry Smith /* w = u + ha */ 3569566063dSJacob Faibussowitsch PetscCall(VecWAXPY(w, h, a, U)); 357e884886fSBarry Smith 3581b797266SDmitry Karpeev /* compute func(U) as base for differencing; only needed first time in and not when provided by user */ 35948a46eb9SPierre Jolivet if (ctx->ncurrenth == 1 && ctx->current_f_allocated) PetscCall((*ctx->func)(ctx->funcctx, U, F)); 3609566063dSJacob Faibussowitsch PetscCall((*ctx->func)(ctx->funcctx, w, y)); 361e884886fSBarry Smith 362c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 363c92e3469SBarry Smith if (ctx->usecomplex) { 3649566063dSJacob Faibussowitsch PetscCall(VecImaginaryPart(y)); 365c92e3469SBarry Smith h = PetscImaginaryPart(h); 366c92e3469SBarry Smith } else { 3679566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, F)); 368c92e3469SBarry Smith } 369c92e3469SBarry Smith #else 3709566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, F)); 371c92e3469SBarry Smith #endif 3729566063dSJacob Faibussowitsch PetscCall(VecScale(y, 1.0 / h)); 3739566063dSJacob Faibussowitsch if (mat->nullsp) PetscCall(MatNullSpaceRemove(mat->nullsp, y)); 374e884886fSBarry Smith 3759566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MATMFFD_Mult, a, y, 0, 0)); 3763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 377e884886fSBarry Smith } 378e884886fSBarry Smith 379e884886fSBarry Smith /* 38001c1178eSBarry Smith MatGetDiagonal_MFFD - Gets the diagonal for a matrix-free matrix 381e884886fSBarry Smith 382e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 383e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 384e884886fSBarry Smith u = current iterate 385e884886fSBarry Smith h = difference interval 386e884886fSBarry Smith */ 38766976f2fSJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_MFFD(Mat mat, Vec a) 388d71ae5a4SJacob Faibussowitsch { 389789d8953SBarry Smith MatMFFD ctx; 390e884886fSBarry Smith PetscScalar h, *aa, *ww, v; 391e884886fSBarry Smith PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON, umin = 100.0 * PETSC_SQRT_MACHINE_EPSILON; 392e884886fSBarry Smith Vec w, U; 393e884886fSBarry Smith PetscInt i, rstart, rend; 394e884886fSBarry Smith 395e884886fSBarry Smith PetscFunctionBegin; 3969566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 39728b400f6SJacob Faibussowitsch PetscCheck(ctx->func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Requires calling MatMFFDSetFunction() first"); 39828b400f6SJacob Faibussowitsch PetscCheck(ctx->funci, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Requires calling MatMFFDSetFunctioni() first"); 399e884886fSBarry Smith w = ctx->w; 400e884886fSBarry Smith U = ctx->current_u; 4019566063dSJacob Faibussowitsch PetscCall((*ctx->func)(ctx->funcctx, U, a)); 4021baa6e33SBarry Smith if (ctx->funcisetbase) PetscCall((*ctx->funcisetbase)(ctx->funcctx, U)); 4039566063dSJacob Faibussowitsch PetscCall(VecCopy(U, w)); 404e884886fSBarry Smith 4059566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(a, &rstart, &rend)); 4069566063dSJacob Faibussowitsch PetscCall(VecGetArray(a, &aa)); 407e884886fSBarry Smith for (i = rstart; i < rend; i++) { 4089566063dSJacob Faibussowitsch PetscCall(VecGetArray(w, &ww)); 409e884886fSBarry Smith h = ww[i - rstart]; 410e884886fSBarry Smith if (h == 0.0) h = 1.0; 411e884886fSBarry Smith if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 412e884886fSBarry Smith else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 413e884886fSBarry Smith h *= epsilon; 414e884886fSBarry Smith 415e884886fSBarry Smith ww[i - rstart] += h; 4169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w, &ww)); 4179566063dSJacob Faibussowitsch PetscCall((*ctx->funci)(ctx->funcctx, i, w, &v)); 418e884886fSBarry Smith aa[i - rstart] = (v - aa[i - rstart]) / h; 419e884886fSBarry Smith 4209566063dSJacob Faibussowitsch PetscCall(VecGetArray(w, &ww)); 421e884886fSBarry Smith ww[i - rstart] -= h; 4229566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w, &ww)); 423e884886fSBarry Smith } 4249566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(a, &aa)); 4253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 426e884886fSBarry Smith } 427e884886fSBarry Smith 428d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J, Vec U, Vec F) 429d71ae5a4SJacob Faibussowitsch { 430789d8953SBarry Smith MatMFFD ctx; 431e884886fSBarry Smith 432e884886fSBarry Smith PetscFunctionBegin; 4339566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J, &ctx)); 4349566063dSJacob Faibussowitsch PetscCall(MatMFFDResetHHistory(J)); 435c51ad4d4SStefano Zampini if (!ctx->current_u) { 4369566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U, &ctx->current_u)); 4379566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(ctx->current_u)); 438c51ad4d4SStefano Zampini } 4399566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(ctx->current_u)); 4409566063dSJacob Faibussowitsch PetscCall(VecCopy(U, ctx->current_u)); 4419566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(ctx->current_u)); 44252121784SBarry Smith if (F) { 4439566063dSJacob Faibussowitsch if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f)); 4443ec795f1SBarry Smith ctx->current_f = F; 445cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_FALSE; 446cfe22f5eSBarry Smith } else if (!ctx->current_f_allocated) { 4479566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(J, NULL, &ctx->current_f)); 448cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_TRUE; 44952121784SBarry Smith } 45048a46eb9SPierre Jolivet if (!ctx->w) PetscCall(VecDuplicate(ctx->current_u, &ctx->w)); 451e884886fSBarry Smith J->assembled = PETSC_TRUE; 4523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 453e884886fSBarry Smith } 4545a576424SJed Brown 455e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void *, Vec, Vec, PetscScalar *); /* force argument to next function to not be extern C*/ 456b2573a8aSBarry Smith 457d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J, FCN3 fun, void *ectx) 458d71ae5a4SJacob Faibussowitsch { 459789d8953SBarry Smith MatMFFD ctx; 460e884886fSBarry Smith 461e884886fSBarry Smith PetscFunctionBegin; 4629566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J, &ctx)); 463e884886fSBarry Smith ctx->checkh = fun; 464e884886fSBarry Smith ctx->checkhctx = ectx; 4653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 466e884886fSBarry Smith } 467e884886fSBarry Smith 4686aa9148fSLisandro Dalcin /*@C 4696aa9148fSLisandro Dalcin MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all 47011a5261eSBarry Smith MATMFFD` options in the database. 4716aa9148fSLisandro Dalcin 472c3339decSBarry Smith Collective 4736aa9148fSLisandro Dalcin 474d8d19677SJose E. Roman Input Parameters: 475fe59aa6dSJacob Faibussowitsch + mat - the `MATMFFD` context 4766aa9148fSLisandro Dalcin - prefix - the prefix to prepend to all option names 4776aa9148fSLisandro Dalcin 47811a5261eSBarry Smith Note: 4796aa9148fSLisandro Dalcin A hyphen (-) must NOT be given at the beginning of the prefix name. 4806aa9148fSLisandro Dalcin The first character of all runtime options is AUTOMATICALLY the hyphen. 4816aa9148fSLisandro Dalcin 4826aa9148fSLisandro Dalcin Level: advanced 4836aa9148fSLisandro Dalcin 4841cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatSetFromOptions()`, `MatCreateSNESMF()`, `MatCreateMFFD()` 4856aa9148fSLisandro Dalcin @*/ 486d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat, const char prefix[]) 487d71ae5a4SJacob Faibussowitsch { 488789d8953SBarry Smith MatMFFD mfctx; 4895fd66863SKarl Rupp 4906aa9148fSLisandro Dalcin PetscFunctionBegin; 4910700a824SBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 4929566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &mfctx)); 4930700a824SBarry Smith PetscValidHeaderSpecific(mfctx, MATMFFD_CLASSID, 1); 4949566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mfctx, prefix)); 4953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4966aa9148fSLisandro Dalcin } 4976aa9148fSLisandro Dalcin 498d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetFromOptions_MFFD(Mat mat, PetscOptionItems *PetscOptionsObject) 499d71ae5a4SJacob Faibussowitsch { 500789d8953SBarry Smith MatMFFD mfctx; 501ace3abfcSBarry Smith PetscBool flg; 502e884886fSBarry Smith char ftype[256]; 503e884886fSBarry Smith 504e884886fSBarry Smith PetscFunctionBegin; 5059566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &mfctx)); 506dbbe0bcdSBarry Smith PetscValidHeaderSpecific(mfctx, MATMFFD_CLASSID, 1); 507d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)mfctx); 5089566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-mat_mffd_type", "Matrix free type", "MatMFFDSetType", MatMFFDList, ((PetscObject)mfctx)->type_name, ftype, 256, &flg)); 5091baa6e33SBarry Smith if (flg) PetscCall(MatMFFDSetType(mat, ftype)); 510e884886fSBarry Smith 5119566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_mffd_err", "set sqrt relative error in function", "MatMFFDSetFunctionError", mfctx->error_rel, &mfctx->error_rel, NULL)); 5129566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_mffd_period", "how often h is recomputed", "MatMFFDSetPeriod", mfctx->recomputeperiod, &mfctx->recomputeperiod, NULL)); 513e884886fSBarry Smith 51490d69ab7SBarry Smith flg = PETSC_FALSE; 5159566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_mffd_check_positivity", "Insure that U + h*a is nonnegative", "MatMFFDSetCheckh", flg, &flg, NULL)); 5161baa6e33SBarry Smith if (flg) PetscCall(MatMFFDSetCheckh(mat, MatMFFDCheckPositivity, NULL)); 517c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 5189566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_mffd_complex", "Use Lyness complex number trick to compute the matrix-vector product", "None", mfctx->usecomplex, &mfctx->usecomplex, NULL)); 519c92e3469SBarry Smith #endif 520dbbe0bcdSBarry Smith PetscTryTypeMethod(mfctx, setfromoptions, PetscOptionsObject); 521d0609cedSBarry Smith PetscOptionsEnd(); 5223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 523e884886fSBarry Smith } 524e884886fSBarry Smith 525d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat, PetscInt period) 526d71ae5a4SJacob Faibussowitsch { 527789d8953SBarry Smith MatMFFD ctx; 528bc0f08ceSBarry Smith 529bc0f08ceSBarry Smith PetscFunctionBegin; 5309566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 531bc0f08ceSBarry Smith ctx->recomputeperiod = period; 5323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 533bc0f08ceSBarry Smith } 534bc0f08ceSBarry Smith 535d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat, PetscErrorCode (*func)(void *, Vec, Vec), void *funcctx) 536d71ae5a4SJacob Faibussowitsch { 537789d8953SBarry Smith MatMFFD ctx; 538bc0f08ceSBarry Smith 539bc0f08ceSBarry Smith PetscFunctionBegin; 5409566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 541bc0f08ceSBarry Smith ctx->func = func; 542bc0f08ceSBarry Smith ctx->funcctx = funcctx; 5433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 544bc0f08ceSBarry Smith } 545bc0f08ceSBarry Smith 546d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat, PetscReal error) 547d71ae5a4SJacob Faibussowitsch { 54813bcc0bdSJacob Faibussowitsch PetscFunctionBegin; 54913bcc0bdSJacob Faibussowitsch if (error != (PetscReal)PETSC_DEFAULT) { 550789d8953SBarry Smith MatMFFD ctx; 551bc0f08ceSBarry Smith 5529566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 55313bcc0bdSJacob Faibussowitsch ctx->error_rel = error; 55413bcc0bdSJacob Faibussowitsch } 5553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 556bc0f08ceSBarry Smith } 557bc0f08ceSBarry Smith 55866976f2fSJacob Faibussowitsch static PetscErrorCode MatMFFDSetHHistory_MFFD(Mat J, PetscScalar history[], PetscInt nhistory) 559d71ae5a4SJacob Faibussowitsch { 560789d8953SBarry Smith MatMFFD ctx; 561789d8953SBarry Smith 5623b49f96aSBarry Smith PetscFunctionBegin; 5639566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J, &ctx)); 564789d8953SBarry Smith ctx->historyh = history; 565789d8953SBarry Smith ctx->maxcurrenth = nhistory; 566789d8953SBarry Smith ctx->currenth = 0.; 5673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5683b49f96aSBarry Smith } 5693b49f96aSBarry Smith 570e884886fSBarry Smith /*MC 57101c1178eSBarry Smith MATMFFD - "mffd" - A matrix-free matrix type. 572e884886fSBarry Smith 573e884886fSBarry Smith Level: advanced 574e884886fSBarry Smith 575*a1f56445SPierre Jolivet Developer Notes: 5762ef1f0ffSBarry Smith This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code 577789d8953SBarry Smith 578ee912ba9SBarry Smith Users should not MatShell... operations on this class, there is some error checking for that incorrect usage 579ee912ba9SBarry Smith 5801cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateMFFD()`, `MatCreateSNESMF()`, `MatMFFDSetFunction()`, `MatMFFDSetType()`, 581db781477SPatrick Sanan `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()` 582db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`, 583db781477SPatrick Sanan `MatMFFDGetH()`, 584e884886fSBarry Smith M*/ 585d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A) 586d71ae5a4SJacob Faibussowitsch { 587e884886fSBarry Smith MatMFFD mfctx; 588e884886fSBarry Smith 589e884886fSBarry Smith PetscFunctionBegin; 5909566063dSJacob Faibussowitsch PetscCall(MatMFFDInitializePackage()); 591e884886fSBarry Smith 5929566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(mfctx, MATMFFD_CLASSID, "MatMFFD", "Matrix-free Finite Differencing", "Mat", PetscObjectComm((PetscObject)A), NULL, NULL)); 5932205254eSKarl Rupp 594e884886fSBarry Smith mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 595e884886fSBarry Smith mfctx->recomputeperiod = 1; 596e884886fSBarry Smith mfctx->count = 0; 597e884886fSBarry Smith mfctx->currenth = 0.0; 5980298fd71SBarry Smith mfctx->historyh = NULL; 599e884886fSBarry Smith mfctx->ncurrenth = 0; 600e884886fSBarry Smith mfctx->maxcurrenth = 0; 601f4259b30SLisandro Dalcin ((PetscObject)mfctx)->type_name = NULL; 602e884886fSBarry Smith 603e884886fSBarry Smith /* 604e884886fSBarry Smith Create the empty data structure to contain compute-h routines. 605e884886fSBarry Smith These will be filled in below from the command line options or 606e884886fSBarry Smith a later call with MatMFFDSetType() or if that is not called 607e884886fSBarry Smith then it will default in the first use of MatMult_MFFD() 608e884886fSBarry Smith */ 609f4259b30SLisandro Dalcin mfctx->ops->compute = NULL; 610f4259b30SLisandro Dalcin mfctx->ops->destroy = NULL; 611f4259b30SLisandro Dalcin mfctx->ops->view = NULL; 612f4259b30SLisandro Dalcin mfctx->ops->setfromoptions = NULL; 613f4259b30SLisandro Dalcin mfctx->hctx = NULL; 614e884886fSBarry Smith 615f4259b30SLisandro Dalcin mfctx->func = NULL; 616f4259b30SLisandro Dalcin mfctx->funcctx = NULL; 6170298fd71SBarry Smith mfctx->w = NULL; 618789d8953SBarry Smith mfctx->mat = A; 619e884886fSBarry Smith 6209566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSHELL)); 6219566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(A, mfctx)); 6229566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))MatMult_MFFD)); 6239566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_DESTROY, (void (*)(void))MatDestroy_MFFD)); 6249566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_VIEW, (void (*)(void))MatView_MFFD)); 6259566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_ASSEMBLY_END, (void (*)(void))MatAssemblyEnd_MFFD)); 6269566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_SET_FROM_OPTIONS, (void (*)(void))MatSetFromOptions_MFFD)); 627*a1f56445SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", MatShellSetContext_Immutable)); 628*a1f56445SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable)); 629*a1f56445SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable)); 630e884886fSBarry Smith 6319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetBase_C", MatMFFDSetBase_MFFD)); 6329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctioniBase_C", MatMFFDSetFunctioniBase_MFFD)); 6339566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctioni_C", MatMFFDSetFunctioni_MFFD)); 6349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunction_C", MatMFFDSetFunction_MFFD)); 6359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetCheckh_C", MatMFFDSetCheckh_MFFD)); 6369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetPeriod_C", MatMFFDSetPeriod_MFFD)); 6379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctionError_C", MatMFFDSetFunctionError_MFFD)); 6389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDResetHHistory_C", MatMFFDResetHHistory_MFFD)); 6399566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetHHistory_C", MatMFFDSetHHistory_MFFD)); 6409566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetType_C", MatMFFDSetType_MFFD)); 6419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDGetH_C", MatMFFDGetH_MFFD)); 6429566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATMFFD)); 6433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 644e884886fSBarry Smith } 645e884886fSBarry Smith 646e884886fSBarry Smith /*@ 64711a5261eSBarry Smith MatCreateMFFD - Creates a matrix-free matrix of type `MATMFFD`. See also `MatCreateSNESMF()` 648e884886fSBarry Smith 64911a5261eSBarry Smith Collective 650e884886fSBarry Smith 651e884886fSBarry Smith Input Parameters: 652fef1beadSBarry Smith + comm - MPI communicator 6532ef1f0ffSBarry Smith . m - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given) 654fef1beadSBarry Smith This value should be the same as the local size used in creating the 655fef1beadSBarry Smith y vector for the matrix-vector product y = Ax. 656fef1beadSBarry Smith . n - This value should be the same as the local size used in creating the 65711a5261eSBarry Smith x vector for the matrix-vector product y = Ax. (or `PETSC_DECIDE` to have 6582ef1f0ffSBarry Smith calculated if `N` is given) For square matrices `n` is almost always `m`. 6592ef1f0ffSBarry Smith . M - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given) 6602ef1f0ffSBarry Smith - N - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given) 661fef1beadSBarry Smith 662e884886fSBarry Smith Output Parameter: 663e884886fSBarry Smith . J - the matrix-free matrix 664e884886fSBarry Smith 66511a5261eSBarry Smith Options Database Keys: 66611a5261eSBarry Smith + -mat_mffd_type - wp or ds (see `MATMFFD_WP` or `MATMFFD_DS`) 667c92e3469SBarry Smith . -mat_mffd_err - square root of estimated relative error in function evaluation 668c92e3469SBarry Smith . -mat_mffd_period - how often h is recomputed, defaults to 1, every time 6692ef1f0ffSBarry Smith . -mat_mffd_check_positivity - possibly decrease `h` until U + h*a has only positive values 67011a5261eSBarry Smith . -mat_mffd_umin <umin> - Sets umin (for default PETSc routine that computes h only) 671c92e3469SBarry Smith - -mat_mffd_complex - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing 672c92e3469SBarry Smith (requires real valued functions but that PETSc be configured for complex numbers) 6739c6ac3b3SBarry Smith 674e884886fSBarry Smith Level: advanced 675e884886fSBarry Smith 676e884886fSBarry Smith Notes: 677e884886fSBarry Smith The matrix-free matrix context merely contains the function pointers 678e884886fSBarry Smith and work space for performing finite difference approximations of 679e884886fSBarry Smith Jacobian-vector products, F'(u)*a, 680e884886fSBarry Smith 681e884886fSBarry Smith The default code uses the following approach to compute h 682e884886fSBarry Smith 683e884886fSBarry Smith .vb 684e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 685e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 686e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 687e884886fSBarry Smith where 688e884886fSBarry Smith error_rel = square root of relative error in function evaluation 689e884886fSBarry Smith umin = minimum iterate parameter 690e884886fSBarry Smith .ve 691e884886fSBarry Smith 69211a5261eSBarry Smith You can call `SNESSetJacobian()` with `MatMFFDComputeJacobian()` if you are using matrix and not a different 693ca93e954SBarry Smith preconditioner matrix 694ca93e954SBarry Smith 69511a5261eSBarry Smith The user can set the error_rel via `MatMFFDSetFunctionError()` and 6962ef1f0ffSBarry Smith umin via `MatMFFDDSSetUmin()`. 697e884886fSBarry Smith 69811a5261eSBarry Smith The user should call `MatDestroy()` when finished with the matrix-free 699e884886fSBarry Smith matrix context. 700e884886fSBarry Smith 7011cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatDestroy()`, `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()` 702db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`, 703db781477SPatrick Sanan `MatMFFDGetH()`, `MatMFFDRegister()`, `MatMFFDComputeJacobian()` 704e884886fSBarry Smith @*/ 705d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateMFFD(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, Mat *J) 706d71ae5a4SJacob Faibussowitsch { 707e884886fSBarry Smith PetscFunctionBegin; 7089566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, J)); 7099566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J, m, n, M, N)); 7109566063dSJacob Faibussowitsch PetscCall(MatSetType(*J, MATMFFD)); 7119566063dSJacob Faibussowitsch PetscCall(MatSetUp(*J)); 7123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 713e884886fSBarry Smith } 714e884886fSBarry Smith 715e884886fSBarry Smith /*@ 71611a5261eSBarry Smith MatMFFDGetH - Gets the last value that was used as the differencing for a `MATMFFD` matrix 717e884886fSBarry Smith parameter. 718e884886fSBarry Smith 719e884886fSBarry Smith Not Collective 720e884886fSBarry Smith 721e884886fSBarry Smith Input Parameters: 72211a5261eSBarry Smith . mat - the `MATMFFD` matrix 723e884886fSBarry Smith 724fd292e60Sprj- Output Parameter: 725e884886fSBarry Smith . h - the differencing step size 726e884886fSBarry Smith 727e884886fSBarry Smith Level: advanced 728e884886fSBarry Smith 729fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDSetHHistory()`, `MatCreateMFFD()`, `MatMFFDResetHHistory()` 730e884886fSBarry Smith @*/ 731d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDGetH(Mat mat, PetscScalar *h) 732d71ae5a4SJacob Faibussowitsch { 733e884886fSBarry Smith PetscFunctionBegin; 73488b4c220SStefano Zampini PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7354f572ea9SToby Isaac PetscAssertPointer(h, 2); 736cac4c232SBarry Smith PetscUseMethod(mat, "MatMFFDGetH_C", (Mat, PetscScalar *), (mat, h)); 7373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 738e884886fSBarry Smith } 739e884886fSBarry Smith 740e884886fSBarry Smith /*@C 74101c1178eSBarry Smith MatMFFDSetFunction - Sets the function used in applying the matrix-free `MATMFFD` matrix. 742e884886fSBarry Smith 743c3339decSBarry Smith Logically Collective 744e884886fSBarry Smith 745e884886fSBarry Smith Input Parameters: 74601c1178eSBarry Smith + mat - the matrix-free matrix `MATMFFD` created via `MatCreateSNESMF()` or `MatCreateMFFD()` 747e884886fSBarry Smith . func - the function to use 748e884886fSBarry Smith - funcctx - optional function context passed to function 749e884886fSBarry Smith 7502920cce0SJacob Faibussowitsch Calling sequence of `func`: 75114f633e2SBarry Smith + funcctx - user provided context 75214f633e2SBarry Smith . x - input vector 75314f633e2SBarry Smith - f - computed output function 75414f633e2SBarry Smith 755e884886fSBarry Smith Level: advanced 756e884886fSBarry Smith 757e884886fSBarry Smith Notes: 75801c1178eSBarry Smith If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix-free 759e884886fSBarry Smith matrix inside your compute Jacobian routine 760e884886fSBarry Smith 76111a5261eSBarry Smith If this is not set then it will use the function set with `SNESSetFunction()` if `MatCreateSNESMF()` was used. 762e884886fSBarry Smith 763fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, 7642920cce0SJacob Faibussowitsch `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESSetFunction()` 765e884886fSBarry Smith @*/ 7662920cce0SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunction(Mat mat, PetscErrorCode (*func)(void *funcctx, Vec x, Vec f), void *funcctx) 767d71ae5a4SJacob Faibussowitsch { 768e884886fSBarry Smith PetscFunctionBegin; 76988b4c220SStefano Zampini PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 770cac4c232SBarry Smith PetscTryMethod(mat, "MatMFFDSetFunction_C", (Mat, PetscErrorCode(*)(void *, Vec, Vec), void *), (mat, func, funcctx)); 7713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 772e884886fSBarry Smith } 773e884886fSBarry Smith 774e884886fSBarry Smith /*@C 77511a5261eSBarry Smith MatMFFDSetFunctioni - Sets the function for a single component for a `MATMFFD` matrix 776e884886fSBarry Smith 777c3339decSBarry Smith Logically Collective 778e884886fSBarry Smith 779e884886fSBarry Smith Input Parameters: 78001c1178eSBarry Smith + mat - the matrix-free matrix `MATMFFD` 781e884886fSBarry Smith - funci - the function to use 782e884886fSBarry Smith 783e884886fSBarry Smith Level: advanced 784e884886fSBarry Smith 785e884886fSBarry Smith Notes: 78601c1178eSBarry Smith If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix-free 78794eb4328SStefano Zampini matrix inside your compute Jacobian routine. 78811a5261eSBarry Smith 78994eb4328SStefano Zampini This function is necessary to compute the diagonal of the matrix. 790486afcceSStefano Zampini funci must not contain any MPI call as it is called inside a loop on the local portion of the vector. 791e884886fSBarry Smith 7921cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()` 793e884886fSBarry Smith @*/ 794d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctioni(Mat mat, PetscErrorCode (*funci)(void *, PetscInt, Vec, PetscScalar *)) 795d71ae5a4SJacob Faibussowitsch { 796e884886fSBarry Smith PetscFunctionBegin; 7970700a824SBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 798cac4c232SBarry Smith PetscTryMethod(mat, "MatMFFDSetFunctioni_C", (Mat, PetscErrorCode(*)(void *, PetscInt, Vec, PetscScalar *)), (mat, funci)); 7993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 800e884886fSBarry Smith } 801e884886fSBarry Smith 802e884886fSBarry Smith /*@C 80311a5261eSBarry Smith MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation for a `MATMFFD` matrix 804e884886fSBarry Smith 805c3339decSBarry Smith Logically Collective 806e884886fSBarry Smith 807e884886fSBarry Smith Input Parameters: 80801c1178eSBarry Smith + mat - the `MATMFFD` matrix-free matrix 809e884886fSBarry Smith - func - the function to use 810e884886fSBarry Smith 811e884886fSBarry Smith Level: advanced 812e884886fSBarry Smith 813e884886fSBarry Smith Notes: 81401c1178eSBarry Smith If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix-free 81594eb4328SStefano Zampini matrix inside your compute Jacobian routine. 816e884886fSBarry Smith 81711a5261eSBarry Smith This function is necessary to compute the diagonal of the matrix, used for example with `PCJACOBI` 81811a5261eSBarry Smith 819fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, 820db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()` 821e884886fSBarry Smith @*/ 822d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctioniBase(Mat mat, PetscErrorCode (*func)(void *, Vec)) 823d71ae5a4SJacob Faibussowitsch { 824e884886fSBarry Smith PetscFunctionBegin; 8250700a824SBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 826cac4c232SBarry Smith PetscTryMethod(mat, "MatMFFDSetFunctioniBase_C", (Mat, PetscErrorCode(*)(void *, Vec)), (mat, func)); 8273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 828e884886fSBarry Smith } 829e884886fSBarry Smith 830e884886fSBarry Smith /*@ 83111a5261eSBarry Smith MatMFFDSetPeriod - Sets how often h is recomputed for a `MATMFFD` matrix, by default it is every time 832e884886fSBarry Smith 833c3339decSBarry Smith Logically Collective 834e884886fSBarry Smith 835e884886fSBarry Smith Input Parameters: 83601c1178eSBarry Smith + mat - the `MATMFFD` matrix-free matrix 837e884886fSBarry Smith - period - 1 for every time, 2 for every second etc 838e884886fSBarry Smith 8392ef1f0ffSBarry Smith Options Database Key: 8402ef1f0ffSBarry Smith . -mat_mffd_period <period> - Sets how often `h` is recomputed 841e884886fSBarry Smith 842e884886fSBarry Smith Level: advanced 843e884886fSBarry Smith 8441cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, 845db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()` 846e884886fSBarry Smith @*/ 847d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetPeriod(Mat mat, PetscInt period) 848d71ae5a4SJacob Faibussowitsch { 849e884886fSBarry Smith PetscFunctionBegin; 85088b4c220SStefano Zampini PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 85188b4c220SStefano Zampini PetscValidLogicalCollectiveInt(mat, period, 2); 852cac4c232SBarry Smith PetscTryMethod(mat, "MatMFFDSetPeriod_C", (Mat, PetscInt), (mat, period)); 8533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 854e884886fSBarry Smith } 855e884886fSBarry Smith 856e884886fSBarry Smith /*@ 85711a5261eSBarry Smith MatMFFDSetFunctionError - Sets the error_rel for the approximation of matrix-vector products using finite differences with the `MATMFFD` matrix 858e884886fSBarry Smith 859c3339decSBarry Smith Logically Collective 860e884886fSBarry Smith 861e884886fSBarry Smith Input Parameters: 86201c1178eSBarry Smith + mat - the `MATMFFD` matrix-free matrix 863fe59aa6dSJacob Faibussowitsch - error - relative error (should be set to the square root of the relative error in the function evaluations) 864e884886fSBarry Smith 8652ef1f0ffSBarry Smith Options Database Key: 866a2b725a8SWilliam Gropp . -mat_mffd_err <error_rel> - Sets error_rel 867e884886fSBarry Smith 868e884886fSBarry Smith Level: advanced 869e884886fSBarry Smith 87011a5261eSBarry Smith Note: 871e884886fSBarry Smith The default matrix-free matrix-vector product routine computes 872e884886fSBarry Smith .vb 873e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 874e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 875e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 876e884886fSBarry Smith .ve 877e884886fSBarry Smith 878fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, 879db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()` 880e884886fSBarry Smith @*/ 881d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctionError(Mat mat, PetscReal error) 882d71ae5a4SJacob Faibussowitsch { 883e884886fSBarry Smith PetscFunctionBegin; 88488b4c220SStefano Zampini PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 88588b4c220SStefano Zampini PetscValidLogicalCollectiveReal(mat, error, 2); 886cac4c232SBarry Smith PetscTryMethod(mat, "MatMFFDSetFunctionError_C", (Mat, PetscReal), (mat, error)); 8873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 888e884886fSBarry Smith } 889e884886fSBarry Smith 890e884886fSBarry Smith /*@ 891e884886fSBarry Smith MatMFFDSetHHistory - Sets an array to collect a history of the 89211a5261eSBarry Smith differencing values (h) computed for the matrix-free product `MATMFFD` matrix 893e884886fSBarry Smith 894c3339decSBarry Smith Logically Collective 895e884886fSBarry Smith 896e884886fSBarry Smith Input Parameters: 89711a5261eSBarry Smith + J - the `MATMFFD` matrix-free matrix 898a5b23f4aSJose E. Roman . history - space to hold the history 899e884886fSBarry Smith - nhistory - number of entries in history, if more entries are generated than 900e884886fSBarry Smith nhistory, then the later ones are discarded 901e884886fSBarry Smith 902e884886fSBarry Smith Level: advanced 903e884886fSBarry Smith 90411a5261eSBarry Smith Note: 90511a5261eSBarry Smith Use `MatMFFDResetHHistory()` to reset the history counter and collect 906e884886fSBarry Smith a new batch of differencing parameters, h. 907e884886fSBarry Smith 9081cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDGetH()`, `MatCreateSNESMF()`, 909db781477SPatrick Sanan `MatMFFDResetHHistory()`, `MatMFFDSetFunctionError()` 910e884886fSBarry Smith @*/ 911d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetHHistory(Mat J, PetscScalar history[], PetscInt nhistory) 912d71ae5a4SJacob Faibussowitsch { 913e884886fSBarry Smith PetscFunctionBegin; 91488b4c220SStefano Zampini PetscValidHeaderSpecific(J, MAT_CLASSID, 1); 9154f572ea9SToby Isaac if (history) PetscAssertPointer(history, 2); 91688b4c220SStefano Zampini PetscValidLogicalCollectiveInt(J, nhistory, 3); 917cac4c232SBarry Smith PetscUseMethod(J, "MatMFFDSetHHistory_C", (Mat, PetscScalar[], PetscInt), (J, history, nhistory)); 9183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 919e884886fSBarry Smith } 920e884886fSBarry Smith 921e884886fSBarry Smith /*@ 922e884886fSBarry Smith MatMFFDResetHHistory - Resets the counter to zero to begin 92311a5261eSBarry Smith collecting a new set of differencing histories for the `MATMFFD` matrix 924e884886fSBarry Smith 925c3339decSBarry Smith Logically Collective 926e884886fSBarry Smith 9272fe279fdSBarry Smith Input Parameter: 928e884886fSBarry Smith . J - the matrix-free matrix context 929e884886fSBarry Smith 930e884886fSBarry Smith Level: advanced 931e884886fSBarry Smith 93211a5261eSBarry Smith Note: 93311a5261eSBarry Smith Use `MatMFFDSetHHistory()` to create the original history counter. 934e884886fSBarry Smith 9351cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDGetH()`, `MatCreateSNESMF()`, 936db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDSetFunctionError()` 937e884886fSBarry Smith @*/ 938d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDResetHHistory(Mat J) 939d71ae5a4SJacob Faibussowitsch { 940e884886fSBarry Smith PetscFunctionBegin; 94188b4c220SStefano Zampini PetscValidHeaderSpecific(J, MAT_CLASSID, 1); 942cac4c232SBarry Smith PetscTryMethod(J, "MatMFFDResetHHistory_C", (Mat), (J)); 9433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 944e884886fSBarry Smith } 945e884886fSBarry Smith 946487a658cSBarry Smith /*@ 9472ef1f0ffSBarry Smith MatMFFDSetBase - Sets the vector `U` at which matrix vector products of the 94811a5261eSBarry Smith Jacobian are computed for the `MATMFFD` matrix 949e884886fSBarry Smith 950c3339decSBarry Smith Logically Collective 951e884886fSBarry Smith 952e884886fSBarry Smith Input Parameters: 95311a5261eSBarry Smith + J - the `MATMFFD` matrix 9543ec795f1SBarry Smith . U - the vector 955bcddec3dSBarry Smith - F - (optional) vector that contains F(u) if it has been already computed 956e884886fSBarry Smith 9572ef1f0ffSBarry Smith Level: advanced 9582ef1f0ffSBarry Smith 95995452b02SPatrick Sanan Notes: 96095452b02SPatrick Sanan This is rarely used directly 961e884886fSBarry Smith 9622ef1f0ffSBarry Smith If `F` is provided then it is not recomputed. Otherwise the function is evaluated at the base 96311a5261eSBarry Smith point during the first `MatMult()` after each call to `MatMFFDSetBase()`. 964dff2f722SBarry Smith 96542747ad1SJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMult()` 966e884886fSBarry Smith @*/ 967d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetBase(Mat J, Vec U, Vec F) 968d71ae5a4SJacob Faibussowitsch { 969e884886fSBarry Smith PetscFunctionBegin; 9700700a824SBarry Smith PetscValidHeaderSpecific(J, MAT_CLASSID, 1); 9710700a824SBarry Smith PetscValidHeaderSpecific(U, VEC_CLASSID, 2); 9720700a824SBarry Smith if (F) PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 973cac4c232SBarry Smith PetscTryMethod(J, "MatMFFDSetBase_C", (Mat, Vec, Vec), (J, U, F)); 9743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 975e884886fSBarry Smith } 976e884886fSBarry Smith 977e884886fSBarry Smith /*@C 978e884886fSBarry Smith MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts 97911a5261eSBarry Smith it to satisfy some criteria for the `MATMFFD` matrix 980e884886fSBarry Smith 981c3339decSBarry Smith Logically Collective 982e884886fSBarry Smith 983e884886fSBarry Smith Input Parameters: 98411a5261eSBarry Smith + J - the `MATMFFD` matrix 9852ef1f0ffSBarry Smith . fun - the function that checks `h` 986e884886fSBarry Smith - ctx - any context needed by the function 987e884886fSBarry Smith 988e884886fSBarry Smith Options Database Keys: 98967b8a455SSatish Balay . -mat_mffd_check_positivity <bool> - Insure that U + h*a is non-negative 990e884886fSBarry Smith 991e884886fSBarry Smith Level: advanced 992e884886fSBarry Smith 99395452b02SPatrick Sanan Notes: 9942ef1f0ffSBarry Smith For example, `MatMFFDCheckPositivity()` insures that all entries of U + h*a are non-negative 995e884886fSBarry Smith 9962ef1f0ffSBarry Smith The function you provide is called after the default `h` has been computed and allows you to 997755b7f64SBarry Smith modify it. 998755b7f64SBarry Smith 9991cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDCheckPositivity()` 1000e884886fSBarry Smith @*/ 1001d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetCheckh(Mat J, PetscErrorCode (*fun)(void *, Vec, Vec, PetscScalar *), void *ctx) 1002d71ae5a4SJacob Faibussowitsch { 1003e884886fSBarry Smith PetscFunctionBegin; 10040700a824SBarry Smith PetscValidHeaderSpecific(J, MAT_CLASSID, 1); 1005cac4c232SBarry Smith PetscTryMethod(J, "MatMFFDSetCheckh_C", (Mat, PetscErrorCode(*)(void *, Vec, Vec, PetscScalar *), void *), (J, fun, ctx)); 10063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1007e884886fSBarry Smith } 1008e884886fSBarry Smith 1009e884886fSBarry Smith /*@ 1010e884886fSBarry Smith MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or 101111a5261eSBarry Smith zero, decreases h until this is satisfied for a `MATMFFD` matrix 1012e884886fSBarry Smith 1013c3339decSBarry Smith Logically Collective 1014e884886fSBarry Smith 1015e884886fSBarry Smith Input Parameters: 1016e884886fSBarry Smith + U - base vector that is added to 1017e884886fSBarry Smith . a - vector that is added 1018e884886fSBarry Smith . h - scaling factor on a 1019e884886fSBarry Smith - dummy - context variable (unused) 1020e884886fSBarry Smith 1021e884886fSBarry Smith Options Database Keys: 102267b8a455SSatish Balay . -mat_mffd_check_positivity <bool> - Insure that U + h*a is nonnegative 1023e884886fSBarry Smith 1024e884886fSBarry Smith Level: advanced 1025e884886fSBarry Smith 102611a5261eSBarry Smith Note: 10272ef1f0ffSBarry Smith This is rarely used directly, rather it is passed as an argument to `MatMFFDSetCheckh()` 1028e884886fSBarry Smith 10291cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDSetCheckh()` 1030e884886fSBarry Smith @*/ 1031d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDCheckPositivity(void *dummy, Vec U, Vec a, PetscScalar *h) 1032d71ae5a4SJacob Faibussowitsch { 1033e884886fSBarry Smith PetscReal val, minval; 1034e884886fSBarry Smith PetscScalar *u_vec, *a_vec; 1035e884886fSBarry Smith PetscInt i, n; 1036e884886fSBarry Smith MPI_Comm comm; 1037e884886fSBarry Smith 1038e884886fSBarry Smith PetscFunctionBegin; 103988b4c220SStefano Zampini PetscValidHeaderSpecific(U, VEC_CLASSID, 2); 104088b4c220SStefano Zampini PetscValidHeaderSpecific(a, VEC_CLASSID, 3); 10414f572ea9SToby Isaac PetscAssertPointer(h, 4); 10429566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)U, &comm)); 10439566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u_vec)); 10449566063dSJacob Faibussowitsch PetscCall(VecGetArray(a, &a_vec)); 10459566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(U, &n)); 1046c068d9bbSLisandro Dalcin minval = PetscAbsScalar(*h) * PetscRealConstant(1.01); 1047e884886fSBarry Smith for (i = 0; i < n; i++) { 1048e884886fSBarry Smith if (PetscRealPart(u_vec[i] + *h * a_vec[i]) <= 0.0) { 1049e884886fSBarry Smith val = PetscAbsScalar(u_vec[i] / a_vec[i]); 1050e884886fSBarry Smith if (val < minval) minval = val; 1051e884886fSBarry Smith } 1052e884886fSBarry Smith } 10539566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u_vec)); 10549566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(a, &a_vec)); 10551c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&minval, &val, 1, MPIU_REAL, MPIU_MIN, comm)); 1056e884886fSBarry Smith if (val <= PetscAbsScalar(*h)) { 10579566063dSJacob Faibussowitsch PetscCall(PetscInfo(U, "Scaling back h from %g to %g\n", (double)PetscRealPart(*h), (double)(.99 * val))); 1058e884886fSBarry Smith if (PetscRealPart(*h) > 0.0) *h = 0.99 * val; 1059e884886fSBarry Smith else *h = -0.99 * val; 1060e884886fSBarry Smith } 10613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1062e884886fSBarry Smith } 1063