1a1f56445SPierre Jolivet #include <../src/mat/impls/shell/shell.h> /*I "petscmat.h" I*/ 2a1f56445SPierre 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 190*a3b724e8SBarry Smith Then, your solver can be chosen with the procedural interface via `MatMFFDSetType`(mfctx, "my_h")` or at runtime via the option 191*a3b724e8SBarry Smith `-mat_mffd_type my_h` 1921c84c290SBarry Smith 1931cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDRegisterAll()`, `MatMFFDRegisterDestroy()` 1941c84c290SBarry Smith @*/ 195d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDRegister(const char sname[], PetscErrorCode (*function)(MatMFFD)) 196d71ae5a4SJacob Faibussowitsch { 197e884886fSBarry Smith PetscFunctionBegin; 1989566063dSJacob Faibussowitsch PetscCall(MatInitializePackage()); 1999566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&MatMFFDList, sname, function)); 2003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 201e884886fSBarry Smith } 202e884886fSBarry Smith 203d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_MFFD(Mat mat) 204d71ae5a4SJacob Faibussowitsch { 205789d8953SBarry Smith MatMFFD ctx; 206e884886fSBarry Smith 207e884886fSBarry Smith PetscFunctionBegin; 2089566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 2099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->w)); 2109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->current_u)); 21148a46eb9SPierre Jolivet if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f)); 212dbbe0bcdSBarry Smith PetscTryTypeMethod(ctx, destroy); 2139566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(&ctx)); 214e884886fSBarry Smith 2159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetBase_C", NULL)); 2169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctioniBase_C", NULL)); 2179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctioni_C", NULL)); 2189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunction_C", NULL)); 2199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctionError_C", NULL)); 2209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetCheckh_C", NULL)); 2219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetPeriod_C", NULL)); 2229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDResetHHistory_C", NULL)); 2239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetHHistory_C", NULL)); 2249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetType_C", NULL)); 2259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDGetH_C", NULL)); 2262e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatSNESMFSetReuseBase_C", NULL)); 2272e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatSNESMFGetReuseBase_C", NULL)); 228ee912ba9SBarry Smith PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatShellSetContext_C", NULL)); 2293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 230e884886fSBarry Smith } 231e884886fSBarry Smith 232e884886fSBarry Smith /* 233e884886fSBarry Smith MatMFFDView_MFFD - Views matrix-free parameters. 234e884886fSBarry Smith 235e884886fSBarry Smith */ 236d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_MFFD(Mat J, PetscViewer viewer) 237d71ae5a4SJacob Faibussowitsch { 238789d8953SBarry Smith MatMFFD ctx; 239a283635eSDmitry Karpeev PetscBool iascii, viewbase, viewfunction; 240a283635eSDmitry Karpeev const char *prefix; 241e884886fSBarry Smith 242e884886fSBarry Smith PetscFunctionBegin; 2439566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J, &ctx)); 2449566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 245e884886fSBarry Smith if (iascii) { 2469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Matrix-free approximation:\n")); 2479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "err=%g (relative error in function evaluation)\n", (double)ctx->error_rel)); 2497adad957SLisandro Dalcin if (!((PetscObject)ctx)->type_name) { 2509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "The compute h routine has not yet been set\n")); 251e884886fSBarry Smith } else { 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Using %s compute h routine\n", ((PetscObject)ctx)->type_name)); 253e884886fSBarry Smith } 254c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 25548a46eb9SPierre Jolivet if (ctx->usecomplex) PetscCall(PetscViewerASCIIPrintf(viewer, "Using Lyness complex number trick to compute the matrix-vector product\n")); 256c92e3469SBarry Smith #endif 257dbbe0bcdSBarry Smith PetscTryTypeMethod(ctx, view, viewer); 2589566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)J, &prefix)); 259a283635eSDmitry Karpeev 2609566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(((PetscObject)J)->options, prefix, "-mat_mffd_view_base", &viewbase)); 261a283635eSDmitry Karpeev if (viewbase) { 2629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Base:\n")); 2639566063dSJacob Faibussowitsch PetscCall(VecView(ctx->current_u, viewer)); 264a283635eSDmitry Karpeev } 2659566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(((PetscObject)J)->options, prefix, "-mat_mffd_view_function", &viewfunction)); 266a283635eSDmitry Karpeev if (viewfunction) { 2679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Function:\n")); 2689566063dSJacob Faibussowitsch PetscCall(VecView(ctx->current_f, viewer)); 269a283635eSDmitry Karpeev } 2709566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 27111aeaf0aSBarry Smith } 2723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 273e884886fSBarry Smith } 274e884886fSBarry Smith 275e884886fSBarry Smith /* 276e884886fSBarry Smith MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This 277e884886fSBarry Smith allows the user to indicate the beginning of a new linear solve by calling 27801c1178eSBarry Smith MatAssemblyXXX() on the matrix-free matrix. This then allows the 2791d0fab5eSBarry Smith MatCreateMFFD_WP() to properly compute ||U|| only the first time 280e884886fSBarry Smith in the linear solver rather than every time. 2815a576424SJed Brown 282cc2e6a90SBarry Smith This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence 283cc2e6a90SBarry Smith it must be labeled as PETSC_EXTERN 284e884886fSBarry Smith */ 285d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J, MatAssemblyType mt) 286d71ae5a4SJacob Faibussowitsch { 287789d8953SBarry Smith MatMFFD j; 288e884886fSBarry Smith 289e884886fSBarry Smith PetscFunctionBegin; 2909566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J, &j)); 2919566063dSJacob Faibussowitsch PetscCall(MatMFFDResetHHistory(J)); 2923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 293e884886fSBarry Smith } 294e884886fSBarry Smith 295e884886fSBarry Smith /* 296e884886fSBarry Smith MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a: 297e884886fSBarry Smith 298e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 299e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 300e884886fSBarry Smith u = current iterate 301e884886fSBarry Smith h = difference interval 302e884886fSBarry Smith */ 303d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_MFFD(Mat mat, Vec a, Vec y) 304d71ae5a4SJacob Faibussowitsch { 305789d8953SBarry Smith MatMFFD ctx; 306e884886fSBarry Smith PetscScalar h; 307e884886fSBarry Smith Vec w, U, F; 308ace3abfcSBarry Smith PetscBool zeroa; 309e884886fSBarry Smith 310e884886fSBarry Smith PetscFunctionBegin; 3119566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 31200045ab3SPierre Jolivet PetscCheck(ctx->current_u, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "MatMFFDSetBase() has not been called, this is often caused by forgetting to call MatAssemblyBegin/End on the first Mat in the SNES compute function"); 313e884886fSBarry Smith /* We log matrix-free matrix-vector products separately, so that we can 314e884886fSBarry Smith separate the performance monitoring from the cases that use conventional 315e884886fSBarry Smith storage. We may eventually modify event logging to associate events 316e884886fSBarry Smith with particular objects, hence alleviating the more general problem. */ 3179566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MATMFFD_Mult, a, y, 0, 0)); 318e884886fSBarry Smith 319e884886fSBarry Smith w = ctx->w; 320e884886fSBarry Smith U = ctx->current_u; 3213ec795f1SBarry Smith F = ctx->current_f; 322e884886fSBarry Smith /* 323e884886fSBarry Smith Compute differencing parameter 324e884886fSBarry Smith */ 3254a2f8832SBarry Smith if (!((PetscObject)ctx)->type_name) { 3269566063dSJacob Faibussowitsch PetscCall(MatMFFDSetType(mat, MATMFFD_WP)); 3279566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(mat)); 328e884886fSBarry Smith } 329dbbe0bcdSBarry Smith PetscUseTypeMethod(ctx, compute, U, a, &h, &zeroa); 330e884886fSBarry Smith if (zeroa) { 3319566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 33285d14658SLisandro Dalcin PetscCall(PetscLogEventEnd(MATMFFD_Mult, a, y, 0, 0)); 3333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 334e884886fSBarry Smith } 335e884886fSBarry Smith 33608401ef6SPierre Jolivet PetscCheck(!mat->erroriffailure || !PetscIsInfOrNanScalar(h), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Computed Nan differencing parameter h"); 33748a46eb9SPierre Jolivet if (ctx->checkh) PetscCall((*ctx->checkh)(ctx->checkhctx, U, a, &h)); 338e884886fSBarry Smith 339e884886fSBarry Smith /* keep a record of the current differencing parameter h */ 340e884886fSBarry Smith ctx->currenth = h; 341e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX) 3429566063dSJacob Faibussowitsch PetscCall(PetscInfo(mat, "Current differencing parameter: %g + %g i\n", (double)PetscRealPart(h), (double)PetscImaginaryPart(h))); 343e884886fSBarry Smith #else 3449566063dSJacob Faibussowitsch PetscCall(PetscInfo(mat, "Current differencing parameter: %15.12e\n", (double)PetscRealPart(h))); 345e884886fSBarry Smith #endif 346ad540459SPierre Jolivet if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) ctx->historyh[ctx->ncurrenth] = h; 347e884886fSBarry Smith ctx->ncurrenth++; 348e884886fSBarry Smith 349c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 350c92e3469SBarry Smith if (ctx->usecomplex) h = PETSC_i * h; 351c92e3469SBarry Smith #endif 352c92e3469SBarry Smith 353e884886fSBarry Smith /* w = u + ha */ 3549566063dSJacob Faibussowitsch PetscCall(VecWAXPY(w, h, a, U)); 355e884886fSBarry Smith 3561b797266SDmitry Karpeev /* compute func(U) as base for differencing; only needed first time in and not when provided by user */ 35748a46eb9SPierre Jolivet if (ctx->ncurrenth == 1 && ctx->current_f_allocated) PetscCall((*ctx->func)(ctx->funcctx, U, F)); 3589566063dSJacob Faibussowitsch PetscCall((*ctx->func)(ctx->funcctx, w, y)); 359e884886fSBarry Smith 360c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 361c92e3469SBarry Smith if (ctx->usecomplex) { 3629566063dSJacob Faibussowitsch PetscCall(VecImaginaryPart(y)); 363c92e3469SBarry Smith h = PetscImaginaryPart(h); 364c92e3469SBarry Smith } else { 3659566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, F)); 366c92e3469SBarry Smith } 367c92e3469SBarry Smith #else 3689566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, F)); 369c92e3469SBarry Smith #endif 3709566063dSJacob Faibussowitsch PetscCall(VecScale(y, 1.0 / h)); 3719566063dSJacob Faibussowitsch if (mat->nullsp) PetscCall(MatNullSpaceRemove(mat->nullsp, y)); 372e884886fSBarry Smith 3739566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MATMFFD_Mult, a, y, 0, 0)); 3743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 375e884886fSBarry Smith } 376e884886fSBarry Smith 377e884886fSBarry Smith /* 37801c1178eSBarry Smith MatGetDiagonal_MFFD - Gets the diagonal for a matrix-free matrix 379e884886fSBarry Smith 380e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 381e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 382e884886fSBarry Smith u = current iterate 383e884886fSBarry Smith h = difference interval 384e884886fSBarry Smith */ 38566976f2fSJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_MFFD(Mat mat, Vec a) 386d71ae5a4SJacob Faibussowitsch { 387789d8953SBarry Smith MatMFFD ctx; 388e884886fSBarry Smith PetscScalar h, *aa, *ww, v; 389e884886fSBarry Smith PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON, umin = 100.0 * PETSC_SQRT_MACHINE_EPSILON; 390e884886fSBarry Smith Vec w, U; 391e884886fSBarry Smith PetscInt i, rstart, rend; 392e884886fSBarry Smith 393e884886fSBarry Smith PetscFunctionBegin; 3949566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 39528b400f6SJacob Faibussowitsch PetscCheck(ctx->func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Requires calling MatMFFDSetFunction() first"); 39628b400f6SJacob Faibussowitsch PetscCheck(ctx->funci, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Requires calling MatMFFDSetFunctioni() first"); 397e884886fSBarry Smith w = ctx->w; 398e884886fSBarry Smith U = ctx->current_u; 3999566063dSJacob Faibussowitsch PetscCall((*ctx->func)(ctx->funcctx, U, a)); 4001baa6e33SBarry Smith if (ctx->funcisetbase) PetscCall((*ctx->funcisetbase)(ctx->funcctx, U)); 4019566063dSJacob Faibussowitsch PetscCall(VecCopy(U, w)); 402e884886fSBarry Smith 4039566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(a, &rstart, &rend)); 4049566063dSJacob Faibussowitsch PetscCall(VecGetArray(a, &aa)); 405e884886fSBarry Smith for (i = rstart; i < rend; i++) { 4069566063dSJacob Faibussowitsch PetscCall(VecGetArray(w, &ww)); 407e884886fSBarry Smith h = ww[i - rstart]; 408e884886fSBarry Smith if (h == 0.0) h = 1.0; 409e884886fSBarry Smith if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 410e884886fSBarry Smith else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 411e884886fSBarry Smith h *= epsilon; 412e884886fSBarry Smith 413e884886fSBarry Smith ww[i - rstart] += h; 4149566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w, &ww)); 4159566063dSJacob Faibussowitsch PetscCall((*ctx->funci)(ctx->funcctx, i, w, &v)); 416e884886fSBarry Smith aa[i - rstart] = (v - aa[i - rstart]) / h; 417e884886fSBarry Smith 4189566063dSJacob Faibussowitsch PetscCall(VecGetArray(w, &ww)); 419e884886fSBarry Smith ww[i - rstart] -= h; 4209566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w, &ww)); 421e884886fSBarry Smith } 4229566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(a, &aa)); 4233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 424e884886fSBarry Smith } 425e884886fSBarry Smith 426d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J, Vec U, Vec F) 427d71ae5a4SJacob Faibussowitsch { 428789d8953SBarry Smith MatMFFD ctx; 429e884886fSBarry Smith 430e884886fSBarry Smith PetscFunctionBegin; 4319566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J, &ctx)); 4329566063dSJacob Faibussowitsch PetscCall(MatMFFDResetHHistory(J)); 433c51ad4d4SStefano Zampini if (!ctx->current_u) { 4349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U, &ctx->current_u)); 4359566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(ctx->current_u)); 436c51ad4d4SStefano Zampini } 4379566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(ctx->current_u)); 4389566063dSJacob Faibussowitsch PetscCall(VecCopy(U, ctx->current_u)); 4399566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(ctx->current_u)); 44052121784SBarry Smith if (F) { 4419566063dSJacob Faibussowitsch if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f)); 4423ec795f1SBarry Smith ctx->current_f = F; 443cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_FALSE; 444cfe22f5eSBarry Smith } else if (!ctx->current_f_allocated) { 4459566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(J, NULL, &ctx->current_f)); 446cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_TRUE; 44752121784SBarry Smith } 44848a46eb9SPierre Jolivet if (!ctx->w) PetscCall(VecDuplicate(ctx->current_u, &ctx->w)); 449e884886fSBarry Smith J->assembled = PETSC_TRUE; 4503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 451e884886fSBarry Smith } 4525a576424SJed Brown 453e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void *, Vec, Vec, PetscScalar *); /* force argument to next function to not be extern C*/ 454b2573a8aSBarry Smith 455d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J, FCN3 fun, void *ectx) 456d71ae5a4SJacob Faibussowitsch { 457789d8953SBarry Smith MatMFFD ctx; 458e884886fSBarry Smith 459e884886fSBarry Smith PetscFunctionBegin; 4609566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J, &ctx)); 461e884886fSBarry Smith ctx->checkh = fun; 462e884886fSBarry Smith ctx->checkhctx = ectx; 4633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 464e884886fSBarry Smith } 465e884886fSBarry Smith 4666aa9148fSLisandro Dalcin /*@C 4676aa9148fSLisandro Dalcin MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all 46811a5261eSBarry Smith MATMFFD` options in the database. 4696aa9148fSLisandro Dalcin 470c3339decSBarry Smith Collective 4716aa9148fSLisandro Dalcin 472d8d19677SJose E. Roman Input Parameters: 473fe59aa6dSJacob Faibussowitsch + mat - the `MATMFFD` context 4746aa9148fSLisandro Dalcin - prefix - the prefix to prepend to all option names 4756aa9148fSLisandro Dalcin 47611a5261eSBarry Smith Note: 4776aa9148fSLisandro Dalcin A hyphen (-) must NOT be given at the beginning of the prefix name. 4786aa9148fSLisandro Dalcin The first character of all runtime options is AUTOMATICALLY the hyphen. 4796aa9148fSLisandro Dalcin 4806aa9148fSLisandro Dalcin Level: advanced 4816aa9148fSLisandro Dalcin 4821cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatSetFromOptions()`, `MatCreateSNESMF()`, `MatCreateMFFD()` 4836aa9148fSLisandro Dalcin @*/ 484d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat, const char prefix[]) 485d71ae5a4SJacob Faibussowitsch { 486789d8953SBarry Smith MatMFFD mfctx; 4875fd66863SKarl Rupp 4886aa9148fSLisandro Dalcin PetscFunctionBegin; 4890700a824SBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 4909566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &mfctx)); 4910700a824SBarry Smith PetscValidHeaderSpecific(mfctx, MATMFFD_CLASSID, 1); 4929566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mfctx, prefix)); 4933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4946aa9148fSLisandro Dalcin } 4956aa9148fSLisandro Dalcin 496d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetFromOptions_MFFD(Mat mat, PetscOptionItems *PetscOptionsObject) 497d71ae5a4SJacob Faibussowitsch { 498789d8953SBarry Smith MatMFFD mfctx; 499ace3abfcSBarry Smith PetscBool flg; 500e884886fSBarry Smith char ftype[256]; 501e884886fSBarry Smith 502e884886fSBarry Smith PetscFunctionBegin; 5039566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &mfctx)); 504dbbe0bcdSBarry Smith PetscValidHeaderSpecific(mfctx, MATMFFD_CLASSID, 1); 505d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)mfctx); 5069566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-mat_mffd_type", "Matrix free type", "MatMFFDSetType", MatMFFDList, ((PetscObject)mfctx)->type_name, ftype, 256, &flg)); 5071baa6e33SBarry Smith if (flg) PetscCall(MatMFFDSetType(mat, ftype)); 508e884886fSBarry Smith 5099566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_mffd_err", "set sqrt relative error in function", "MatMFFDSetFunctionError", mfctx->error_rel, &mfctx->error_rel, NULL)); 5109566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_mffd_period", "how often h is recomputed", "MatMFFDSetPeriod", mfctx->recomputeperiod, &mfctx->recomputeperiod, NULL)); 511e884886fSBarry Smith 51290d69ab7SBarry Smith flg = PETSC_FALSE; 5139566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_mffd_check_positivity", "Insure that U + h*a is nonnegative", "MatMFFDSetCheckh", flg, &flg, NULL)); 5141baa6e33SBarry Smith if (flg) PetscCall(MatMFFDSetCheckh(mat, MatMFFDCheckPositivity, NULL)); 515c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 5169566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_mffd_complex", "Use Lyness complex number trick to compute the matrix-vector product", "None", mfctx->usecomplex, &mfctx->usecomplex, NULL)); 517c92e3469SBarry Smith #endif 518dbbe0bcdSBarry Smith PetscTryTypeMethod(mfctx, setfromoptions, PetscOptionsObject); 519d0609cedSBarry Smith PetscOptionsEnd(); 5203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 521e884886fSBarry Smith } 522e884886fSBarry Smith 523d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat, PetscInt period) 524d71ae5a4SJacob Faibussowitsch { 525789d8953SBarry Smith MatMFFD ctx; 526bc0f08ceSBarry Smith 527bc0f08ceSBarry Smith PetscFunctionBegin; 5289566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 529bc0f08ceSBarry Smith ctx->recomputeperiod = period; 5303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 531bc0f08ceSBarry Smith } 532bc0f08ceSBarry Smith 533d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat, PetscErrorCode (*func)(void *, Vec, Vec), void *funcctx) 534d71ae5a4SJacob Faibussowitsch { 535789d8953SBarry Smith MatMFFD ctx; 536bc0f08ceSBarry Smith 537bc0f08ceSBarry Smith PetscFunctionBegin; 5389566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 539bc0f08ceSBarry Smith ctx->func = func; 540bc0f08ceSBarry Smith ctx->funcctx = funcctx; 5413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 542bc0f08ceSBarry Smith } 543bc0f08ceSBarry Smith 544d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat, PetscReal error) 545d71ae5a4SJacob Faibussowitsch { 54613bcc0bdSJacob Faibussowitsch PetscFunctionBegin; 54713bcc0bdSJacob Faibussowitsch if (error != (PetscReal)PETSC_DEFAULT) { 548789d8953SBarry Smith MatMFFD ctx; 549bc0f08ceSBarry Smith 5509566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 55113bcc0bdSJacob Faibussowitsch ctx->error_rel = error; 55213bcc0bdSJacob Faibussowitsch } 5533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 554bc0f08ceSBarry Smith } 555bc0f08ceSBarry Smith 55666976f2fSJacob Faibussowitsch static PetscErrorCode MatMFFDSetHHistory_MFFD(Mat J, PetscScalar history[], PetscInt nhistory) 557d71ae5a4SJacob Faibussowitsch { 558789d8953SBarry Smith MatMFFD ctx; 559789d8953SBarry Smith 5603b49f96aSBarry Smith PetscFunctionBegin; 5619566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J, &ctx)); 562789d8953SBarry Smith ctx->historyh = history; 563789d8953SBarry Smith ctx->maxcurrenth = nhistory; 564789d8953SBarry Smith ctx->currenth = 0.; 5653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5663b49f96aSBarry Smith } 5673b49f96aSBarry Smith 568e884886fSBarry Smith /*MC 56901c1178eSBarry Smith MATMFFD - "mffd" - A matrix-free matrix type. 570e884886fSBarry Smith 571e884886fSBarry Smith Level: advanced 572e884886fSBarry Smith 573a1f56445SPierre Jolivet Developer Notes: 5742ef1f0ffSBarry Smith This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code 575789d8953SBarry Smith 576ee912ba9SBarry Smith Users should not MatShell... operations on this class, there is some error checking for that incorrect usage 577ee912ba9SBarry Smith 5781cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateMFFD()`, `MatCreateSNESMF()`, `MatMFFDSetFunction()`, `MatMFFDSetType()`, 579db781477SPatrick Sanan `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()` 580db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`, 581db781477SPatrick Sanan `MatMFFDGetH()`, 582e884886fSBarry Smith M*/ 583d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A) 584d71ae5a4SJacob Faibussowitsch { 585e884886fSBarry Smith MatMFFD mfctx; 586e884886fSBarry Smith 587e884886fSBarry Smith PetscFunctionBegin; 5889566063dSJacob Faibussowitsch PetscCall(MatMFFDInitializePackage()); 589e884886fSBarry Smith 5909566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(mfctx, MATMFFD_CLASSID, "MatMFFD", "Matrix-free Finite Differencing", "Mat", PetscObjectComm((PetscObject)A), NULL, NULL)); 5912205254eSKarl Rupp 592e884886fSBarry Smith mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 593e884886fSBarry Smith mfctx->recomputeperiod = 1; 594e884886fSBarry Smith mfctx->count = 0; 595e884886fSBarry Smith mfctx->currenth = 0.0; 5960298fd71SBarry Smith mfctx->historyh = NULL; 597e884886fSBarry Smith mfctx->ncurrenth = 0; 598e884886fSBarry Smith mfctx->maxcurrenth = 0; 599f4259b30SLisandro Dalcin ((PetscObject)mfctx)->type_name = NULL; 600e884886fSBarry Smith 601e884886fSBarry Smith /* 602e884886fSBarry Smith Create the empty data structure to contain compute-h routines. 603e884886fSBarry Smith These will be filled in below from the command line options or 604e884886fSBarry Smith a later call with MatMFFDSetType() or if that is not called 605e884886fSBarry Smith then it will default in the first use of MatMult_MFFD() 606e884886fSBarry Smith */ 607f4259b30SLisandro Dalcin mfctx->ops->compute = NULL; 608f4259b30SLisandro Dalcin mfctx->ops->destroy = NULL; 609f4259b30SLisandro Dalcin mfctx->ops->view = NULL; 610f4259b30SLisandro Dalcin mfctx->ops->setfromoptions = NULL; 611f4259b30SLisandro Dalcin mfctx->hctx = NULL; 612e884886fSBarry Smith 613f4259b30SLisandro Dalcin mfctx->func = NULL; 614f4259b30SLisandro Dalcin mfctx->funcctx = NULL; 6150298fd71SBarry Smith mfctx->w = NULL; 616789d8953SBarry Smith mfctx->mat = A; 617e884886fSBarry Smith 6189566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSHELL)); 6199566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(A, mfctx)); 6209566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))MatMult_MFFD)); 6219566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_DESTROY, (void (*)(void))MatDestroy_MFFD)); 6229566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_VIEW, (void (*)(void))MatView_MFFD)); 6239566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_ASSEMBLY_END, (void (*)(void))MatAssemblyEnd_MFFD)); 6249566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_SET_FROM_OPTIONS, (void (*)(void))MatSetFromOptions_MFFD)); 625a1f56445SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", MatShellSetContext_Immutable)); 626a1f56445SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable)); 627a1f56445SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable)); 628e884886fSBarry Smith 6299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetBase_C", MatMFFDSetBase_MFFD)); 6309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctioniBase_C", MatMFFDSetFunctioniBase_MFFD)); 6319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctioni_C", MatMFFDSetFunctioni_MFFD)); 6329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunction_C", MatMFFDSetFunction_MFFD)); 6339566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetCheckh_C", MatMFFDSetCheckh_MFFD)); 6349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetPeriod_C", MatMFFDSetPeriod_MFFD)); 6359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctionError_C", MatMFFDSetFunctionError_MFFD)); 6369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDResetHHistory_C", MatMFFDResetHHistory_MFFD)); 6379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetHHistory_C", MatMFFDSetHHistory_MFFD)); 6389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetType_C", MatMFFDSetType_MFFD)); 6399566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDGetH_C", MatMFFDGetH_MFFD)); 6409566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATMFFD)); 6413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 642e884886fSBarry Smith } 643e884886fSBarry Smith 644e884886fSBarry Smith /*@ 64511a5261eSBarry Smith MatCreateMFFD - Creates a matrix-free matrix of type `MATMFFD`. See also `MatCreateSNESMF()` 646e884886fSBarry Smith 64711a5261eSBarry Smith Collective 648e884886fSBarry Smith 649e884886fSBarry Smith Input Parameters: 650fef1beadSBarry Smith + comm - MPI communicator 6512ef1f0ffSBarry Smith . m - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given) 652fef1beadSBarry Smith This value should be the same as the local size used in creating the 653fef1beadSBarry Smith y vector for the matrix-vector product y = Ax. 654fef1beadSBarry Smith . n - This value should be the same as the local size used in creating the 65511a5261eSBarry Smith x vector for the matrix-vector product y = Ax. (or `PETSC_DECIDE` to have 6562ef1f0ffSBarry Smith calculated if `N` is given) For square matrices `n` is almost always `m`. 6572ef1f0ffSBarry Smith . M - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given) 6582ef1f0ffSBarry Smith - N - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given) 659fef1beadSBarry Smith 660e884886fSBarry Smith Output Parameter: 661e884886fSBarry Smith . J - the matrix-free matrix 662e884886fSBarry Smith 66311a5261eSBarry Smith Options Database Keys: 66411a5261eSBarry Smith + -mat_mffd_type - wp or ds (see `MATMFFD_WP` or `MATMFFD_DS`) 665c92e3469SBarry Smith . -mat_mffd_err - square root of estimated relative error in function evaluation 666c92e3469SBarry Smith . -mat_mffd_period - how often h is recomputed, defaults to 1, every time 6672ef1f0ffSBarry Smith . -mat_mffd_check_positivity - possibly decrease `h` until U + h*a has only positive values 66811a5261eSBarry Smith . -mat_mffd_umin <umin> - Sets umin (for default PETSc routine that computes h only) 669c92e3469SBarry Smith - -mat_mffd_complex - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing 670c92e3469SBarry Smith (requires real valued functions but that PETSc be configured for complex numbers) 6719c6ac3b3SBarry Smith 672e884886fSBarry Smith Level: advanced 673e884886fSBarry Smith 674e884886fSBarry Smith Notes: 675e884886fSBarry Smith The matrix-free matrix context merely contains the function pointers 676e884886fSBarry Smith and work space for performing finite difference approximations of 677e884886fSBarry Smith Jacobian-vector products, F'(u)*a, 678e884886fSBarry Smith 679e884886fSBarry Smith The default code uses the following approach to compute h 680e884886fSBarry Smith 681e884886fSBarry Smith .vb 682e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 683e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 684e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 685e884886fSBarry Smith where 686e884886fSBarry Smith error_rel = square root of relative error in function evaluation 687e884886fSBarry Smith umin = minimum iterate parameter 688e884886fSBarry Smith .ve 689e884886fSBarry Smith 69011a5261eSBarry Smith You can call `SNESSetJacobian()` with `MatMFFDComputeJacobian()` if you are using matrix and not a different 691ca93e954SBarry Smith preconditioner matrix 692ca93e954SBarry Smith 69311a5261eSBarry Smith The user can set the error_rel via `MatMFFDSetFunctionError()` and 6942ef1f0ffSBarry Smith umin via `MatMFFDDSSetUmin()`. 695e884886fSBarry Smith 69611a5261eSBarry Smith The user should call `MatDestroy()` when finished with the matrix-free 697e884886fSBarry Smith matrix context. 698e884886fSBarry Smith 6991cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatDestroy()`, `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()` 700db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`, 701db781477SPatrick Sanan `MatMFFDGetH()`, `MatMFFDRegister()`, `MatMFFDComputeJacobian()` 702e884886fSBarry Smith @*/ 703d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateMFFD(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, Mat *J) 704d71ae5a4SJacob Faibussowitsch { 705e884886fSBarry Smith PetscFunctionBegin; 7069566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, J)); 7079566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J, m, n, M, N)); 7089566063dSJacob Faibussowitsch PetscCall(MatSetType(*J, MATMFFD)); 7099566063dSJacob Faibussowitsch PetscCall(MatSetUp(*J)); 7103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 711e884886fSBarry Smith } 712e884886fSBarry Smith 713e884886fSBarry Smith /*@ 71411a5261eSBarry Smith MatMFFDGetH - Gets the last value that was used as the differencing for a `MATMFFD` matrix 715e884886fSBarry Smith parameter. 716e884886fSBarry Smith 717e884886fSBarry Smith Not Collective 718e884886fSBarry Smith 719e884886fSBarry Smith Input Parameters: 72011a5261eSBarry Smith . mat - the `MATMFFD` matrix 721e884886fSBarry Smith 722fd292e60Sprj- Output Parameter: 723e884886fSBarry Smith . h - the differencing step size 724e884886fSBarry Smith 725e884886fSBarry Smith Level: advanced 726e884886fSBarry Smith 727fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDSetHHistory()`, `MatCreateMFFD()`, `MatMFFDResetHHistory()` 728e884886fSBarry Smith @*/ 729d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDGetH(Mat mat, PetscScalar *h) 730d71ae5a4SJacob Faibussowitsch { 731e884886fSBarry Smith PetscFunctionBegin; 73288b4c220SStefano Zampini PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7334f572ea9SToby Isaac PetscAssertPointer(h, 2); 734cac4c232SBarry Smith PetscUseMethod(mat, "MatMFFDGetH_C", (Mat, PetscScalar *), (mat, h)); 7353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 736e884886fSBarry Smith } 737e884886fSBarry Smith 738e884886fSBarry Smith /*@C 73901c1178eSBarry Smith MatMFFDSetFunction - Sets the function used in applying the matrix-free `MATMFFD` matrix. 740e884886fSBarry Smith 741c3339decSBarry Smith Logically Collective 742e884886fSBarry Smith 743e884886fSBarry Smith Input Parameters: 74401c1178eSBarry Smith + mat - the matrix-free matrix `MATMFFD` created via `MatCreateSNESMF()` or `MatCreateMFFD()` 745e884886fSBarry Smith . func - the function to use 746e884886fSBarry Smith - funcctx - optional function context passed to function 747e884886fSBarry Smith 7482920cce0SJacob Faibussowitsch Calling sequence of `func`: 74914f633e2SBarry Smith + funcctx - user provided context 75014f633e2SBarry Smith . x - input vector 75114f633e2SBarry Smith - f - computed output function 75214f633e2SBarry Smith 753e884886fSBarry Smith Level: advanced 754e884886fSBarry Smith 755e884886fSBarry Smith Notes: 75601c1178eSBarry Smith If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix-free 757e884886fSBarry Smith matrix inside your compute Jacobian routine 758e884886fSBarry Smith 75911a5261eSBarry Smith If this is not set then it will use the function set with `SNESSetFunction()` if `MatCreateSNESMF()` was used. 760e884886fSBarry Smith 761fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, 7622920cce0SJacob Faibussowitsch `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESSetFunction()` 763e884886fSBarry Smith @*/ 7642920cce0SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunction(Mat mat, PetscErrorCode (*func)(void *funcctx, Vec x, Vec f), void *funcctx) 765d71ae5a4SJacob Faibussowitsch { 766e884886fSBarry Smith PetscFunctionBegin; 76788b4c220SStefano Zampini PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 768cac4c232SBarry Smith PetscTryMethod(mat, "MatMFFDSetFunction_C", (Mat, PetscErrorCode(*)(void *, Vec, Vec), void *), (mat, func, funcctx)); 7693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 770e884886fSBarry Smith } 771e884886fSBarry Smith 772e884886fSBarry Smith /*@C 77311a5261eSBarry Smith MatMFFDSetFunctioni - Sets the function for a single component for a `MATMFFD` matrix 774e884886fSBarry Smith 775c3339decSBarry Smith Logically Collective 776e884886fSBarry Smith 777e884886fSBarry Smith Input Parameters: 77801c1178eSBarry Smith + mat - the matrix-free matrix `MATMFFD` 779e884886fSBarry Smith - funci - the function to use 780e884886fSBarry Smith 781e884886fSBarry Smith Level: advanced 782e884886fSBarry Smith 783e884886fSBarry Smith Notes: 78401c1178eSBarry Smith If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix-free 78594eb4328SStefano Zampini matrix inside your compute Jacobian routine. 78611a5261eSBarry Smith 78794eb4328SStefano Zampini This function is necessary to compute the diagonal of the matrix. 788486afcceSStefano Zampini funci must not contain any MPI call as it is called inside a loop on the local portion of the vector. 789e884886fSBarry Smith 7901cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()` 791e884886fSBarry Smith @*/ 792d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctioni(Mat mat, PetscErrorCode (*funci)(void *, PetscInt, Vec, PetscScalar *)) 793d71ae5a4SJacob Faibussowitsch { 794e884886fSBarry Smith PetscFunctionBegin; 7950700a824SBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 796cac4c232SBarry Smith PetscTryMethod(mat, "MatMFFDSetFunctioni_C", (Mat, PetscErrorCode(*)(void *, PetscInt, Vec, PetscScalar *)), (mat, funci)); 7973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 798e884886fSBarry Smith } 799e884886fSBarry Smith 800e884886fSBarry Smith /*@C 80111a5261eSBarry Smith MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation for a `MATMFFD` matrix 802e884886fSBarry Smith 803c3339decSBarry Smith Logically Collective 804e884886fSBarry Smith 805e884886fSBarry Smith Input Parameters: 80601c1178eSBarry Smith + mat - the `MATMFFD` matrix-free matrix 807e884886fSBarry Smith - func - the function to use 808e884886fSBarry Smith 809e884886fSBarry Smith Level: advanced 810e884886fSBarry Smith 811e884886fSBarry Smith Notes: 81201c1178eSBarry Smith If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix-free 81394eb4328SStefano Zampini matrix inside your compute Jacobian routine. 814e884886fSBarry Smith 81511a5261eSBarry Smith This function is necessary to compute the diagonal of the matrix, used for example with `PCJACOBI` 81611a5261eSBarry Smith 817fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, 818db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()` 819e884886fSBarry Smith @*/ 820d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctioniBase(Mat mat, PetscErrorCode (*func)(void *, Vec)) 821d71ae5a4SJacob Faibussowitsch { 822e884886fSBarry Smith PetscFunctionBegin; 8230700a824SBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 824cac4c232SBarry Smith PetscTryMethod(mat, "MatMFFDSetFunctioniBase_C", (Mat, PetscErrorCode(*)(void *, Vec)), (mat, func)); 8253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 826e884886fSBarry Smith } 827e884886fSBarry Smith 828e884886fSBarry Smith /*@ 82911a5261eSBarry Smith MatMFFDSetPeriod - Sets how often h is recomputed for a `MATMFFD` matrix, by default it is every time 830e884886fSBarry Smith 831c3339decSBarry Smith Logically Collective 832e884886fSBarry Smith 833e884886fSBarry Smith Input Parameters: 83401c1178eSBarry Smith + mat - the `MATMFFD` matrix-free matrix 835e884886fSBarry Smith - period - 1 for every time, 2 for every second etc 836e884886fSBarry Smith 8372ef1f0ffSBarry Smith Options Database Key: 8382ef1f0ffSBarry Smith . -mat_mffd_period <period> - Sets how often `h` is recomputed 839e884886fSBarry Smith 840e884886fSBarry Smith Level: advanced 841e884886fSBarry Smith 8421cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, 843db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()` 844e884886fSBarry Smith @*/ 845d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetPeriod(Mat mat, PetscInt period) 846d71ae5a4SJacob Faibussowitsch { 847e884886fSBarry Smith PetscFunctionBegin; 84888b4c220SStefano Zampini PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 84988b4c220SStefano Zampini PetscValidLogicalCollectiveInt(mat, period, 2); 850cac4c232SBarry Smith PetscTryMethod(mat, "MatMFFDSetPeriod_C", (Mat, PetscInt), (mat, period)); 8513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 852e884886fSBarry Smith } 853e884886fSBarry Smith 854e884886fSBarry Smith /*@ 85511a5261eSBarry Smith MatMFFDSetFunctionError - Sets the error_rel for the approximation of matrix-vector products using finite differences with the `MATMFFD` matrix 856e884886fSBarry Smith 857c3339decSBarry Smith Logically Collective 858e884886fSBarry Smith 859e884886fSBarry Smith Input Parameters: 86001c1178eSBarry Smith + mat - the `MATMFFD` matrix-free matrix 861fe59aa6dSJacob Faibussowitsch - error - relative error (should be set to the square root of the relative error in the function evaluations) 862e884886fSBarry Smith 8632ef1f0ffSBarry Smith Options Database Key: 864a2b725a8SWilliam Gropp . -mat_mffd_err <error_rel> - Sets error_rel 865e884886fSBarry Smith 866e884886fSBarry Smith Level: advanced 867e884886fSBarry Smith 86811a5261eSBarry Smith Note: 869e884886fSBarry Smith The default matrix-free matrix-vector product routine computes 870e884886fSBarry Smith .vb 871e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 872e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 873e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 874e884886fSBarry Smith .ve 875e884886fSBarry Smith 876fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, 877db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()` 878e884886fSBarry Smith @*/ 879d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctionError(Mat mat, PetscReal error) 880d71ae5a4SJacob Faibussowitsch { 881e884886fSBarry Smith PetscFunctionBegin; 88288b4c220SStefano Zampini PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 88388b4c220SStefano Zampini PetscValidLogicalCollectiveReal(mat, error, 2); 884cac4c232SBarry Smith PetscTryMethod(mat, "MatMFFDSetFunctionError_C", (Mat, PetscReal), (mat, error)); 8853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 886e884886fSBarry Smith } 887e884886fSBarry Smith 888e884886fSBarry Smith /*@ 889e884886fSBarry Smith MatMFFDSetHHistory - Sets an array to collect a history of the 89011a5261eSBarry Smith differencing values (h) computed for the matrix-free product `MATMFFD` matrix 891e884886fSBarry Smith 892c3339decSBarry Smith Logically Collective 893e884886fSBarry Smith 894e884886fSBarry Smith Input Parameters: 89511a5261eSBarry Smith + J - the `MATMFFD` matrix-free matrix 896a5b23f4aSJose E. Roman . history - space to hold the history 897e884886fSBarry Smith - nhistory - number of entries in history, if more entries are generated than 898e884886fSBarry Smith nhistory, then the later ones are discarded 899e884886fSBarry Smith 900e884886fSBarry Smith Level: advanced 901e884886fSBarry Smith 90211a5261eSBarry Smith Note: 90311a5261eSBarry Smith Use `MatMFFDResetHHistory()` to reset the history counter and collect 904e884886fSBarry Smith a new batch of differencing parameters, h. 905e884886fSBarry Smith 9061cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDGetH()`, `MatCreateSNESMF()`, 907db781477SPatrick Sanan `MatMFFDResetHHistory()`, `MatMFFDSetFunctionError()` 908e884886fSBarry Smith @*/ 909d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetHHistory(Mat J, PetscScalar history[], PetscInt nhistory) 910d71ae5a4SJacob Faibussowitsch { 911e884886fSBarry Smith PetscFunctionBegin; 91288b4c220SStefano Zampini PetscValidHeaderSpecific(J, MAT_CLASSID, 1); 9134f572ea9SToby Isaac if (history) PetscAssertPointer(history, 2); 91488b4c220SStefano Zampini PetscValidLogicalCollectiveInt(J, nhistory, 3); 915cac4c232SBarry Smith PetscUseMethod(J, "MatMFFDSetHHistory_C", (Mat, PetscScalar[], PetscInt), (J, history, nhistory)); 9163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 917e884886fSBarry Smith } 918e884886fSBarry Smith 919e884886fSBarry Smith /*@ 920e884886fSBarry Smith MatMFFDResetHHistory - Resets the counter to zero to begin 92111a5261eSBarry Smith collecting a new set of differencing histories for the `MATMFFD` matrix 922e884886fSBarry Smith 923c3339decSBarry Smith Logically Collective 924e884886fSBarry Smith 9252fe279fdSBarry Smith Input Parameter: 926e884886fSBarry Smith . J - the matrix-free matrix context 927e884886fSBarry Smith 928e884886fSBarry Smith Level: advanced 929e884886fSBarry Smith 93011a5261eSBarry Smith Note: 93111a5261eSBarry Smith Use `MatMFFDSetHHistory()` to create the original history counter. 932e884886fSBarry Smith 9331cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDGetH()`, `MatCreateSNESMF()`, 934db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDSetFunctionError()` 935e884886fSBarry Smith @*/ 936d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDResetHHistory(Mat J) 937d71ae5a4SJacob Faibussowitsch { 938e884886fSBarry Smith PetscFunctionBegin; 93988b4c220SStefano Zampini PetscValidHeaderSpecific(J, MAT_CLASSID, 1); 940cac4c232SBarry Smith PetscTryMethod(J, "MatMFFDResetHHistory_C", (Mat), (J)); 9413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 942e884886fSBarry Smith } 943e884886fSBarry Smith 944487a658cSBarry Smith /*@ 9452ef1f0ffSBarry Smith MatMFFDSetBase - Sets the vector `U` at which matrix vector products of the 94611a5261eSBarry Smith Jacobian are computed for the `MATMFFD` matrix 947e884886fSBarry Smith 948c3339decSBarry Smith Logically Collective 949e884886fSBarry Smith 950e884886fSBarry Smith Input Parameters: 95111a5261eSBarry Smith + J - the `MATMFFD` matrix 9523ec795f1SBarry Smith . U - the vector 953bcddec3dSBarry Smith - F - (optional) vector that contains F(u) if it has been already computed 954e884886fSBarry Smith 9552ef1f0ffSBarry Smith Level: advanced 9562ef1f0ffSBarry Smith 95795452b02SPatrick Sanan Notes: 95895452b02SPatrick Sanan This is rarely used directly 959e884886fSBarry Smith 9602ef1f0ffSBarry Smith If `F` is provided then it is not recomputed. Otherwise the function is evaluated at the base 96111a5261eSBarry Smith point during the first `MatMult()` after each call to `MatMFFDSetBase()`. 962dff2f722SBarry Smith 96342747ad1SJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMult()` 964e884886fSBarry Smith @*/ 965d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetBase(Mat J, Vec U, Vec F) 966d71ae5a4SJacob Faibussowitsch { 967e884886fSBarry Smith PetscFunctionBegin; 9680700a824SBarry Smith PetscValidHeaderSpecific(J, MAT_CLASSID, 1); 9690700a824SBarry Smith PetscValidHeaderSpecific(U, VEC_CLASSID, 2); 9700700a824SBarry Smith if (F) PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 971cac4c232SBarry Smith PetscTryMethod(J, "MatMFFDSetBase_C", (Mat, Vec, Vec), (J, U, F)); 9723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 973e884886fSBarry Smith } 974e884886fSBarry Smith 975e884886fSBarry Smith /*@C 976e884886fSBarry Smith MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts 97711a5261eSBarry Smith it to satisfy some criteria for the `MATMFFD` matrix 978e884886fSBarry Smith 979c3339decSBarry Smith Logically Collective 980e884886fSBarry Smith 981e884886fSBarry Smith Input Parameters: 98211a5261eSBarry Smith + J - the `MATMFFD` matrix 9832ef1f0ffSBarry Smith . fun - the function that checks `h` 984e884886fSBarry Smith - ctx - any context needed by the function 985e884886fSBarry Smith 986e884886fSBarry Smith Options Database Keys: 98767b8a455SSatish Balay . -mat_mffd_check_positivity <bool> - Insure that U + h*a is non-negative 988e884886fSBarry Smith 989e884886fSBarry Smith Level: advanced 990e884886fSBarry Smith 99195452b02SPatrick Sanan Notes: 9922ef1f0ffSBarry Smith For example, `MatMFFDCheckPositivity()` insures that all entries of U + h*a are non-negative 993e884886fSBarry Smith 9942ef1f0ffSBarry Smith The function you provide is called after the default `h` has been computed and allows you to 995755b7f64SBarry Smith modify it. 996755b7f64SBarry Smith 9971cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDCheckPositivity()` 998e884886fSBarry Smith @*/ 999d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetCheckh(Mat J, PetscErrorCode (*fun)(void *, Vec, Vec, PetscScalar *), void *ctx) 1000d71ae5a4SJacob Faibussowitsch { 1001e884886fSBarry Smith PetscFunctionBegin; 10020700a824SBarry Smith PetscValidHeaderSpecific(J, MAT_CLASSID, 1); 1003cac4c232SBarry Smith PetscTryMethod(J, "MatMFFDSetCheckh_C", (Mat, PetscErrorCode(*)(void *, Vec, Vec, PetscScalar *), void *), (J, fun, ctx)); 10043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1005e884886fSBarry Smith } 1006e884886fSBarry Smith 1007e884886fSBarry Smith /*@ 1008e884886fSBarry Smith MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or 100911a5261eSBarry Smith zero, decreases h until this is satisfied for a `MATMFFD` matrix 1010e884886fSBarry Smith 1011c3339decSBarry Smith Logically Collective 1012e884886fSBarry Smith 1013e884886fSBarry Smith Input Parameters: 1014e884886fSBarry Smith + U - base vector that is added to 1015e884886fSBarry Smith . a - vector that is added 1016e884886fSBarry Smith . h - scaling factor on a 1017e884886fSBarry Smith - dummy - context variable (unused) 1018e884886fSBarry Smith 1019e884886fSBarry Smith Options Database Keys: 102067b8a455SSatish Balay . -mat_mffd_check_positivity <bool> - Insure that U + h*a is nonnegative 1021e884886fSBarry Smith 1022e884886fSBarry Smith Level: advanced 1023e884886fSBarry Smith 102411a5261eSBarry Smith Note: 10252ef1f0ffSBarry Smith This is rarely used directly, rather it is passed as an argument to `MatMFFDSetCheckh()` 1026e884886fSBarry Smith 10271cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDSetCheckh()` 1028e884886fSBarry Smith @*/ 1029d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDCheckPositivity(void *dummy, Vec U, Vec a, PetscScalar *h) 1030d71ae5a4SJacob Faibussowitsch { 1031e884886fSBarry Smith PetscReal val, minval; 1032e884886fSBarry Smith PetscScalar *u_vec, *a_vec; 1033e884886fSBarry Smith PetscInt i, n; 1034e884886fSBarry Smith MPI_Comm comm; 1035e884886fSBarry Smith 1036e884886fSBarry Smith PetscFunctionBegin; 103788b4c220SStefano Zampini PetscValidHeaderSpecific(U, VEC_CLASSID, 2); 103888b4c220SStefano Zampini PetscValidHeaderSpecific(a, VEC_CLASSID, 3); 10394f572ea9SToby Isaac PetscAssertPointer(h, 4); 10409566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)U, &comm)); 10419566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u_vec)); 10429566063dSJacob Faibussowitsch PetscCall(VecGetArray(a, &a_vec)); 10439566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(U, &n)); 1044c068d9bbSLisandro Dalcin minval = PetscAbsScalar(*h) * PetscRealConstant(1.01); 1045e884886fSBarry Smith for (i = 0; i < n; i++) { 1046e884886fSBarry Smith if (PetscRealPart(u_vec[i] + *h * a_vec[i]) <= 0.0) { 1047e884886fSBarry Smith val = PetscAbsScalar(u_vec[i] / a_vec[i]); 1048e884886fSBarry Smith if (val < minval) minval = val; 1049e884886fSBarry Smith } 1050e884886fSBarry Smith } 10519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u_vec)); 10529566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(a, &a_vec)); 10531c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&minval, &val, 1, MPIU_REAL, MPIU_MIN, comm)); 1054e884886fSBarry Smith if (val <= PetscAbsScalar(*h)) { 10559566063dSJacob Faibussowitsch PetscCall(PetscInfo(U, "Scaling back h from %g to %g\n", (double)PetscRealPart(*h), (double)(.99 * val))); 1056e884886fSBarry Smith if (PetscRealPart(*h) > 0.0) *h = 0.99 * val; 1057e884886fSBarry Smith else *h = -0.99 * val; 1058e884886fSBarry Smith } 10593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1060e884886fSBarry Smith } 1061