1e884886fSBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> 3c6db04a5SJed Brown #include <../src/mat/impls/mffd/mffdimpl.h> /*I "petscmat.h" I*/ 4e884886fSBarry Smith 5f4259b30SLisandro Dalcin PetscFunctionList MatMFFDList = NULL; 6ace3abfcSBarry Smith PetscBool MatMFFDRegisterAllCalled = PETSC_FALSE; 7e884886fSBarry Smith 87087cfbeSBarry Smith PetscClassId MATMFFD_CLASSID; 9166c7f25SBarry Smith PetscLogEvent MATMFFD_Mult; 10e884886fSBarry Smith 11ace3abfcSBarry Smith static PetscBool MatMFFDPackageInitialized = PETSC_FALSE; 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 1811a5261eSBarry Smith .seealso: `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; 26*3ba16761SJacob 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 3511a5261eSBarry Smith .seealso: `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; 43*3ba16761SJacob 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)); 66*3ba16761SJacob 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); 77789d8953SBarry Smith PetscValidCharPointer(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)); 82*3ba16761SJacob 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)); 885f80ce2aSJacob Faibussowitsch PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown MatMFFD type %s given", ftype); 899566063dSJacob Faibussowitsch PetscCall((*r)(ctx)); 909566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)ctx, ftype)); 91*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 92789d8953SBarry Smith } 93789d8953SBarry Smith 94e884886fSBarry Smith /*@C 95e884886fSBarry Smith MatMFFDSetType - Sets the method that is used to compute the 96e884886fSBarry Smith differencing parameter for finite differene 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: 106e884886fSBarry Smith For example, such routines can compute h for use in 107e884886fSBarry Smith Jacobian-vector products of the form 108e884886fSBarry Smith 109e884886fSBarry Smith F(x+ha) - F(x) 110e884886fSBarry Smith F'(u)a ~= ---------------- 111e884886fSBarry Smith h 112e884886fSBarry Smith 11311a5261eSBarry Smith .seealso: `MATMFFD`, `MATMFFD_WP`, `MATMFFD_DS`, `MatCreateSNESMF()`, `MatMFFDRegister()`, `MatMFFDSetFunction()`, `MatCreateMFFD()` 114e884886fSBarry Smith @*/ 115d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetType(Mat mat, MatMFFDType ftype) 116d71ae5a4SJacob Faibussowitsch { 117e884886fSBarry Smith PetscFunctionBegin; 1180700a824SBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 119e884886fSBarry Smith PetscValidCharPointer(ftype, 2); 120cac4c232SBarry Smith PetscTryMethod(mat, "MatMFFDSetType_C", (Mat, MatMFFDType), (mat, ftype)); 121*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 122e884886fSBarry Smith } 123e884886fSBarry Smith 1245d21e1f8SStefano Zampini static PetscErrorCode MatGetDiagonal_MFFD(Mat, Vec); 1255d21e1f8SStefano Zampini 126e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void *, Vec); /* force argument to next function to not be extern C*/ 127d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunctioniBase_MFFD(Mat mat, FCN1 func) 128d71ae5a4SJacob Faibussowitsch { 129789d8953SBarry Smith MatMFFD ctx; 130e884886fSBarry Smith 131e884886fSBarry Smith PetscFunctionBegin; 1329566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 133e884886fSBarry Smith ctx->funcisetbase = func; 134*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 135e884886fSBarry Smith } 136e884886fSBarry Smith 137e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void *, PetscInt, Vec, PetscScalar *); /* force argument to next function to not be extern C*/ 138d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunctioni_MFFD(Mat mat, FCN2 funci) 139d71ae5a4SJacob Faibussowitsch { 140789d8953SBarry Smith MatMFFD ctx; 141e884886fSBarry Smith 142e884886fSBarry Smith PetscFunctionBegin; 1439566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 144e884886fSBarry Smith ctx->funci = funci; 1459566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(mat, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_MFFD)); 146*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1475d21e1f8SStefano Zampini } 148789d8953SBarry Smith 149d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDGetH_MFFD(Mat mat, PetscScalar *h) 150d71ae5a4SJacob Faibussowitsch { 151789d8953SBarry Smith MatMFFD ctx; 152789d8953SBarry Smith 153789d8953SBarry Smith PetscFunctionBegin; 1549566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 155789d8953SBarry Smith *h = ctx->currenth; 156*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 157e884886fSBarry Smith } 158e884886fSBarry Smith 159d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDResetHHistory_MFFD(Mat J) 160d71ae5a4SJacob Faibussowitsch { 161789d8953SBarry Smith MatMFFD ctx; 162bc0f08ceSBarry Smith 163bc0f08ceSBarry Smith PetscFunctionBegin; 1649566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J, &ctx)); 165bc0f08ceSBarry Smith ctx->ncurrenth = 0; 166*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 167bc0f08ceSBarry Smith } 168e884886fSBarry Smith 1691c84c290SBarry Smith /*@C 17011a5261eSBarry Smith MatMFFDRegister - Adds a method to the MATMFFD` registry. 1711c84c290SBarry Smith 1721c84c290SBarry Smith Not Collective 1731c84c290SBarry Smith 1741c84c290SBarry Smith Input Parameters: 1751c84c290SBarry Smith + name_solver - name of a new user-defined compute-h module 1761c84c290SBarry Smith - routine_create - routine to create method context 1771c84c290SBarry Smith 1781c84c290SBarry Smith Level: developer 1791c84c290SBarry Smith 18011a5261eSBarry Smith Note: 18111a5261eSBarry Smith `MatMFFDRegister()` may be called multiple times to add several user-defined solvers. 1821c84c290SBarry Smith 1831c84c290SBarry Smith Sample usage: 1841c84c290SBarry Smith .vb 185bdf89e91SBarry Smith MatMFFDRegister("my_h",MyHCreate); 1861c84c290SBarry Smith .ve 1871c84c290SBarry Smith 1881c84c290SBarry Smith Then, your solver can be chosen with the procedural interface via 18911a5261eSBarry Smith $ `MatMFFDSetType`(mfctx,"my_h") 1901c84c290SBarry Smith or at runtime via the option 191be7a6d03SBarry Smith $ -mat_mffd_type my_h 1921c84c290SBarry Smith 19311a5261eSBarry Smith .seealso: `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)); 200*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 201e884886fSBarry Smith } 202e884886fSBarry Smith 203e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/ 204d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_MFFD(Mat mat) 205d71ae5a4SJacob Faibussowitsch { 206789d8953SBarry Smith MatMFFD ctx; 207e884886fSBarry Smith 208e884886fSBarry Smith PetscFunctionBegin; 2099566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 2109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->w)); 2119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->current_u)); 21248a46eb9SPierre Jolivet if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f)); 213dbbe0bcdSBarry Smith PetscTryTypeMethod(ctx, destroy); 2149566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(&ctx)); 215e884886fSBarry Smith 2169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetBase_C", NULL)); 2179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctioniBase_C", NULL)); 2189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctioni_C", NULL)); 2199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunction_C", NULL)); 2209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctionError_C", NULL)); 2219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetCheckh_C", NULL)); 2229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetPeriod_C", NULL)); 2239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDResetHHistory_C", NULL)); 2249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetHHistory_C", NULL)); 2259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetType_C", NULL)); 2269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDGetH_C", NULL)); 2272e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatSNESMFSetReuseBase_C", NULL)); 2282e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatSNESMFGetReuseBase_C", NULL)); 229*3ba16761SJacob 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 } 272*3ba16761SJacob 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 278e884886fSBarry 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)); 292*3ba16761SJacob 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)); 31228b400f6SJacob 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"); 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)); 333*3ba16761SJacob 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)); 374*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 375e884886fSBarry Smith } 376e884886fSBarry Smith 377e884886fSBarry Smith /* 378e884886fSBarry 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 */ 385d71ae5a4SJacob Faibussowitsch 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)); 423*3ba16761SJacob 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; 450*3ba16761SJacob 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; 463*3ba16761SJacob 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: 47311a5261eSBarry Smith + A - 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 48211a5261eSBarry Smith .seealso: `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)); 493*3ba16761SJacob 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(); 520*3ba16761SJacob 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; 530*3ba16761SJacob 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; 541*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 542bc0f08ceSBarry Smith } 543bc0f08ceSBarry Smith 544d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat, PetscReal error) 545d71ae5a4SJacob Faibussowitsch { 546789d8953SBarry Smith MatMFFD ctx; 547bc0f08ceSBarry Smith 548bc0f08ceSBarry Smith PetscFunctionBegin; 5499566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 550bc0f08ceSBarry Smith if (error != PETSC_DEFAULT) ctx->error_rel = error; 551*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 552bc0f08ceSBarry Smith } 553bc0f08ceSBarry Smith 554d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetHHistory_MFFD(Mat J, PetscScalar history[], PetscInt nhistory) 555d71ae5a4SJacob Faibussowitsch { 556789d8953SBarry Smith MatMFFD ctx; 557789d8953SBarry Smith 5583b49f96aSBarry Smith PetscFunctionBegin; 5599566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J, &ctx)); 560789d8953SBarry Smith ctx->historyh = history; 561789d8953SBarry Smith ctx->maxcurrenth = nhistory; 562789d8953SBarry Smith ctx->currenth = 0.; 563*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5643b49f96aSBarry Smith } 5653b49f96aSBarry Smith 566e884886fSBarry Smith /*MC 567e884886fSBarry Smith MATMFFD - MATMFFD = "mffd" - A matrix free matrix type. 568e884886fSBarry Smith 569e884886fSBarry Smith Level: advanced 570e884886fSBarry Smith 57111a5261eSBarry Smith Developers Note: This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code 572789d8953SBarry Smith 573db781477SPatrick Sanan .seealso: `MatCreateMFFD()`, `MatCreateSNESMF()`, `MatMFFDSetFunction()`, `MatMFFDSetType()`, 574db781477SPatrick Sanan `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()` 575db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`, 576db781477SPatrick Sanan `MatMFFDGetH()`, 577e884886fSBarry Smith M*/ 578d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A) 579d71ae5a4SJacob Faibussowitsch { 580e884886fSBarry Smith MatMFFD mfctx; 581e884886fSBarry Smith 582e884886fSBarry Smith PetscFunctionBegin; 5839566063dSJacob Faibussowitsch PetscCall(MatMFFDInitializePackage()); 584e884886fSBarry Smith 5859566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(mfctx, MATMFFD_CLASSID, "MatMFFD", "Matrix-free Finite Differencing", "Mat", PetscObjectComm((PetscObject)A), NULL, NULL)); 5862205254eSKarl Rupp 587e884886fSBarry Smith mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 588e884886fSBarry Smith mfctx->recomputeperiod = 1; 589e884886fSBarry Smith mfctx->count = 0; 590e884886fSBarry Smith mfctx->currenth = 0.0; 5910298fd71SBarry Smith mfctx->historyh = NULL; 592e884886fSBarry Smith mfctx->ncurrenth = 0; 593e884886fSBarry Smith mfctx->maxcurrenth = 0; 594f4259b30SLisandro Dalcin ((PetscObject)mfctx)->type_name = NULL; 595e884886fSBarry Smith 596e884886fSBarry Smith /* 597e884886fSBarry Smith Create the empty data structure to contain compute-h routines. 598e884886fSBarry Smith These will be filled in below from the command line options or 599e884886fSBarry Smith a later call with MatMFFDSetType() or if that is not called 600e884886fSBarry Smith then it will default in the first use of MatMult_MFFD() 601e884886fSBarry Smith */ 602f4259b30SLisandro Dalcin mfctx->ops->compute = NULL; 603f4259b30SLisandro Dalcin mfctx->ops->destroy = NULL; 604f4259b30SLisandro Dalcin mfctx->ops->view = NULL; 605f4259b30SLisandro Dalcin mfctx->ops->setfromoptions = NULL; 606f4259b30SLisandro Dalcin mfctx->hctx = NULL; 607e884886fSBarry Smith 608f4259b30SLisandro Dalcin mfctx->func = NULL; 609f4259b30SLisandro Dalcin mfctx->funcctx = NULL; 6100298fd71SBarry Smith mfctx->w = NULL; 611789d8953SBarry Smith mfctx->mat = A; 612e884886fSBarry Smith 6139566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSHELL)); 6149566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(A, mfctx)); 6159566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))MatMult_MFFD)); 6169566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_DESTROY, (void (*)(void))MatDestroy_MFFD)); 6179566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_VIEW, (void (*)(void))MatView_MFFD)); 6189566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_ASSEMBLY_END, (void (*)(void))MatAssemblyEnd_MFFD)); 6199566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_SET_FROM_OPTIONS, (void (*)(void))MatSetFromOptions_MFFD)); 620e884886fSBarry Smith 6219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetBase_C", MatMFFDSetBase_MFFD)); 6229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctioniBase_C", MatMFFDSetFunctioniBase_MFFD)); 6239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctioni_C", MatMFFDSetFunctioni_MFFD)); 6249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunction_C", MatMFFDSetFunction_MFFD)); 6259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetCheckh_C", MatMFFDSetCheckh_MFFD)); 6269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetPeriod_C", MatMFFDSetPeriod_MFFD)); 6279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctionError_C", MatMFFDSetFunctionError_MFFD)); 6289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDResetHHistory_C", MatMFFDResetHHistory_MFFD)); 6299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetHHistory_C", MatMFFDSetHHistory_MFFD)); 6309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetType_C", MatMFFDSetType_MFFD)); 6319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDGetH_C", MatMFFDGetH_MFFD)); 6329566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATMFFD)); 633*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 634e884886fSBarry Smith } 635e884886fSBarry Smith 636e884886fSBarry Smith /*@ 63711a5261eSBarry Smith MatCreateMFFD - Creates a matrix-free matrix of type `MATMFFD`. See also `MatCreateSNESMF()` 638e884886fSBarry Smith 63911a5261eSBarry Smith Collective 640e884886fSBarry Smith 641e884886fSBarry Smith Input Parameters: 642fef1beadSBarry Smith + comm - MPI communicator 64311a5261eSBarry Smith . m - number of local rows (or `PETSC_DECIDE` to have calculated if M is given) 644fef1beadSBarry Smith This value should be the same as the local size used in creating the 645fef1beadSBarry Smith y vector for the matrix-vector product y = Ax. 646fef1beadSBarry Smith . n - This value should be the same as the local size used in creating the 64711a5261eSBarry Smith x vector for the matrix-vector product y = Ax. (or `PETSC_DECIDE` to have 648fef1beadSBarry Smith calculated if N is given) For square matrices n is almost always m. 64911a5261eSBarry Smith . M - number of global rows (or `PETSC_DETERMINE` to have calculated if m is given) 65011a5261eSBarry Smith - N - number of global columns (or `PETSC_DETERMINE` to have calculated if n is given) 651fef1beadSBarry Smith 652e884886fSBarry Smith Output Parameter: 653e884886fSBarry Smith . J - the matrix-free matrix 654e884886fSBarry Smith 65511a5261eSBarry Smith Options Database Keys: 65611a5261eSBarry Smith + -mat_mffd_type - wp or ds (see `MATMFFD_WP` or `MATMFFD_DS`) 657c92e3469SBarry Smith . -mat_mffd_err - square root of estimated relative error in function evaluation 658c92e3469SBarry Smith . -mat_mffd_period - how often h is recomputed, defaults to 1, every time 659c92e3469SBarry Smith . -mat_mffd_check_positivity - possibly decrease h until U + h*a has only positive values 66011a5261eSBarry Smith . -mat_mffd_umin <umin> - Sets umin (for default PETSc routine that computes h only) 661c92e3469SBarry Smith - -mat_mffd_complex - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing 662c92e3469SBarry Smith (requires real valued functions but that PETSc be configured for complex numbers) 6639c6ac3b3SBarry Smith 664e884886fSBarry Smith Level: advanced 665e884886fSBarry Smith 666e884886fSBarry Smith Notes: 667e884886fSBarry Smith The matrix-free matrix context merely contains the function pointers 668e884886fSBarry Smith and work space for performing finite difference approximations of 669e884886fSBarry Smith Jacobian-vector products, F'(u)*a, 670e884886fSBarry Smith 671e884886fSBarry Smith The default code uses the following approach to compute h 672e884886fSBarry Smith 673e884886fSBarry Smith .vb 674e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 675e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 676e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 677e884886fSBarry Smith where 678e884886fSBarry Smith error_rel = square root of relative error in function evaluation 679e884886fSBarry Smith umin = minimum iterate parameter 680e884886fSBarry Smith .ve 681e884886fSBarry Smith 68211a5261eSBarry Smith You can call `SNESSetJacobian()` with `MatMFFDComputeJacobian()` if you are using matrix and not a different 683ca93e954SBarry Smith preconditioner matrix 684ca93e954SBarry Smith 68511a5261eSBarry Smith The user can set the error_rel via `MatMFFDSetFunctionError()` and 68611a5261eSBarry Smith umin via `MatMFFDDSSetUmin()`; see Users-Manual: ch_snes for details. 687e884886fSBarry Smith 68811a5261eSBarry Smith The user should call `MatDestroy()` when finished with the matrix-free 689e884886fSBarry Smith matrix context. 690e884886fSBarry Smith 69111a5261eSBarry Smith .seealso: `MATMFFD`, `MatDestroy()`, `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()` 692db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`, 693db781477SPatrick Sanan `MatMFFDGetH()`, `MatMFFDRegister()`, `MatMFFDComputeJacobian()` 694e884886fSBarry Smith @*/ 695d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateMFFD(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, Mat *J) 696d71ae5a4SJacob Faibussowitsch { 697e884886fSBarry Smith PetscFunctionBegin; 6989566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, J)); 6999566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J, m, n, M, N)); 7009566063dSJacob Faibussowitsch PetscCall(MatSetType(*J, MATMFFD)); 7019566063dSJacob Faibussowitsch PetscCall(MatSetUp(*J)); 702*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 703e884886fSBarry Smith } 704e884886fSBarry Smith 705e884886fSBarry Smith /*@ 70611a5261eSBarry Smith MatMFFDGetH - Gets the last value that was used as the differencing for a `MATMFFD` matrix 707e884886fSBarry Smith parameter. 708e884886fSBarry Smith 709e884886fSBarry Smith Not Collective 710e884886fSBarry Smith 711e884886fSBarry Smith Input Parameters: 71211a5261eSBarry Smith . mat - the `MATMFFD` matrix 713e884886fSBarry Smith 714fd292e60Sprj- Output Parameter: 715e884886fSBarry Smith . h - the differencing step size 716e884886fSBarry Smith 717e884886fSBarry Smith Level: advanced 718e884886fSBarry Smith 71911a5261eSBarry Smith .seealso: `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDSetHHistory()`, `MatCreateMFFD()`, `MATMFFD`, `MatMFFDResetHHistory()` 720e884886fSBarry Smith @*/ 721d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDGetH(Mat mat, PetscScalar *h) 722d71ae5a4SJacob Faibussowitsch { 723e884886fSBarry Smith PetscFunctionBegin; 72488b4c220SStefano Zampini PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 725dadcf809SJacob Faibussowitsch PetscValidScalarPointer(h, 2); 726cac4c232SBarry Smith PetscUseMethod(mat, "MatMFFDGetH_C", (Mat, PetscScalar *), (mat, h)); 727*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 728e884886fSBarry Smith } 729e884886fSBarry Smith 730e884886fSBarry Smith /*@C 73111a5261eSBarry Smith MatMFFDSetFunction - Sets the function used in applying the matrix free `MATMFFD` matrix. 732e884886fSBarry Smith 733c3339decSBarry Smith Logically Collective 734e884886fSBarry Smith 735e884886fSBarry Smith Input Parameters: 73611a5261eSBarry Smith + mat - the matrix free matrix `MATMFFD` created via `MatCreateSNESMF()` or `MatCreateMFFD()` 737e884886fSBarry Smith . func - the function to use 738e884886fSBarry Smith - funcctx - optional function context passed to function 739e884886fSBarry Smith 74014f633e2SBarry Smith Calling Sequence of func: 74114f633e2SBarry Smith $ func (void *funcctx, Vec x, Vec f) 74214f633e2SBarry Smith 74314f633e2SBarry Smith + funcctx - user provided context 74414f633e2SBarry Smith . x - input vector 74514f633e2SBarry Smith - f - computed output function 74614f633e2SBarry Smith 747e884886fSBarry Smith Level: advanced 748e884886fSBarry Smith 749e884886fSBarry Smith Notes: 75011a5261eSBarry Smith If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix free 751e884886fSBarry Smith matrix inside your compute Jacobian routine 752e884886fSBarry Smith 75311a5261eSBarry Smith If this is not set then it will use the function set with `SNESSetFunction()` if `MatCreateSNESMF()` was used. 754e884886fSBarry Smith 75511a5261eSBarry Smith .seealso: `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, `MATMFFD`, 756db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()` 757e884886fSBarry Smith @*/ 758d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunction(Mat mat, PetscErrorCode (*func)(void *, Vec, Vec), void *funcctx) 759d71ae5a4SJacob Faibussowitsch { 760e884886fSBarry Smith PetscFunctionBegin; 76188b4c220SStefano Zampini PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 762cac4c232SBarry Smith PetscTryMethod(mat, "MatMFFDSetFunction_C", (Mat, PetscErrorCode(*)(void *, Vec, Vec), void *), (mat, func, funcctx)); 763*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 764e884886fSBarry Smith } 765e884886fSBarry Smith 766e884886fSBarry Smith /*@C 76711a5261eSBarry Smith MatMFFDSetFunctioni - Sets the function for a single component for a `MATMFFD` matrix 768e884886fSBarry Smith 769c3339decSBarry Smith Logically Collective 770e884886fSBarry Smith 771e884886fSBarry Smith Input Parameters: 77211a5261eSBarry Smith + mat - the matrix free matrix `MATMFFD` 773e884886fSBarry Smith - funci - the function to use 774e884886fSBarry Smith 775e884886fSBarry Smith Level: advanced 776e884886fSBarry Smith 777e884886fSBarry Smith Notes: 77811a5261eSBarry Smith If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix free 77994eb4328SStefano Zampini matrix inside your compute Jacobian routine. 78011a5261eSBarry Smith 78194eb4328SStefano Zampini This function is necessary to compute the diagonal of the matrix. 782486afcceSStefano Zampini funci must not contain any MPI call as it is called inside a loop on the local portion of the vector. 783e884886fSBarry Smith 78411a5261eSBarry Smith .seealso: `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()` 785e884886fSBarry Smith @*/ 786d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctioni(Mat mat, PetscErrorCode (*funci)(void *, PetscInt, Vec, PetscScalar *)) 787d71ae5a4SJacob Faibussowitsch { 788e884886fSBarry Smith PetscFunctionBegin; 7890700a824SBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 790cac4c232SBarry Smith PetscTryMethod(mat, "MatMFFDSetFunctioni_C", (Mat, PetscErrorCode(*)(void *, PetscInt, Vec, PetscScalar *)), (mat, funci)); 791*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 792e884886fSBarry Smith } 793e884886fSBarry Smith 794e884886fSBarry Smith /*@C 79511a5261eSBarry Smith MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation for a `MATMFFD` matrix 796e884886fSBarry Smith 797c3339decSBarry Smith Logically Collective 798e884886fSBarry Smith 799e884886fSBarry Smith Input Parameters: 80011a5261eSBarry Smith + mat - the `MATMFFD` matrix free matrix 801e884886fSBarry Smith - func - the function to use 802e884886fSBarry Smith 803e884886fSBarry Smith Level: advanced 804e884886fSBarry Smith 805e884886fSBarry Smith Notes: 80611a5261eSBarry Smith If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix free 80794eb4328SStefano Zampini matrix inside your compute Jacobian routine. 808e884886fSBarry Smith 80911a5261eSBarry Smith This function is necessary to compute the diagonal of the matrix, used for example with `PCJACOBI` 81011a5261eSBarry Smith 81111a5261eSBarry Smith .seealso: `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, `MATMFFD` 812db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()` 813e884886fSBarry Smith @*/ 814d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctioniBase(Mat mat, PetscErrorCode (*func)(void *, Vec)) 815d71ae5a4SJacob Faibussowitsch { 816e884886fSBarry Smith PetscFunctionBegin; 8170700a824SBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 818cac4c232SBarry Smith PetscTryMethod(mat, "MatMFFDSetFunctioniBase_C", (Mat, PetscErrorCode(*)(void *, Vec)), (mat, func)); 819*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 820e884886fSBarry Smith } 821e884886fSBarry Smith 822e884886fSBarry Smith /*@ 82311a5261eSBarry Smith MatMFFDSetPeriod - Sets how often h is recomputed for a `MATMFFD` matrix, by default it is every time 824e884886fSBarry Smith 825c3339decSBarry Smith Logically Collective 826e884886fSBarry Smith 827e884886fSBarry Smith Input Parameters: 82811a5261eSBarry Smith + mat - the `MATMFFD` matrix free matrix 829e884886fSBarry Smith - period - 1 for every time, 2 for every second etc 830e884886fSBarry Smith 831e884886fSBarry Smith Options Database Keys: 83267b8a455SSatish Balay . -mat_mffd_period <period> - Sets how often h is recomputed 833e884886fSBarry Smith 834e884886fSBarry Smith Level: advanced 835e884886fSBarry Smith 83611a5261eSBarry Smith .seealso: `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, 837db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()` 838e884886fSBarry Smith @*/ 839d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetPeriod(Mat mat, PetscInt period) 840d71ae5a4SJacob Faibussowitsch { 841e884886fSBarry Smith PetscFunctionBegin; 84288b4c220SStefano Zampini PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 84388b4c220SStefano Zampini PetscValidLogicalCollectiveInt(mat, period, 2); 844cac4c232SBarry Smith PetscTryMethod(mat, "MatMFFDSetPeriod_C", (Mat, PetscInt), (mat, period)); 845*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 846e884886fSBarry Smith } 847e884886fSBarry Smith 848e884886fSBarry Smith /*@ 84911a5261eSBarry Smith MatMFFDSetFunctionError - Sets the error_rel for the approximation of matrix-vector products using finite differences with the `MATMFFD` matrix 850e884886fSBarry Smith 851c3339decSBarry Smith Logically Collective 852e884886fSBarry Smith 853e884886fSBarry Smith Input Parameters: 85411a5261eSBarry Smith + mat - the `MATMFFD` matrix free matrix 85511a5261eSBarry Smith - error_rel - relative error (should be set to the square root of the relative error in the function evaluations) 856e884886fSBarry Smith 857e884886fSBarry Smith Options Database Keys: 858a2b725a8SWilliam Gropp . -mat_mffd_err <error_rel> - Sets error_rel 859e884886fSBarry Smith 860e884886fSBarry Smith Level: advanced 861e884886fSBarry Smith 86211a5261eSBarry Smith Note: 863e884886fSBarry Smith The default matrix-free matrix-vector product routine computes 864e884886fSBarry Smith .vb 865e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 866e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 867e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 868e884886fSBarry Smith .ve 869e884886fSBarry Smith 87011a5261eSBarry Smith .seealso: `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, `MATMFFD` 871db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()` 872e884886fSBarry Smith @*/ 873d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetFunctionError(Mat mat, PetscReal error) 874d71ae5a4SJacob Faibussowitsch { 875e884886fSBarry Smith PetscFunctionBegin; 87688b4c220SStefano Zampini PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 87788b4c220SStefano Zampini PetscValidLogicalCollectiveReal(mat, error, 2); 878cac4c232SBarry Smith PetscTryMethod(mat, "MatMFFDSetFunctionError_C", (Mat, PetscReal), (mat, error)); 879*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 880e884886fSBarry Smith } 881e884886fSBarry Smith 882e884886fSBarry Smith /*@ 883e884886fSBarry Smith MatMFFDSetHHistory - Sets an array to collect a history of the 88411a5261eSBarry Smith differencing values (h) computed for the matrix-free product `MATMFFD` matrix 885e884886fSBarry Smith 886c3339decSBarry Smith Logically Collective 887e884886fSBarry Smith 888e884886fSBarry Smith Input Parameters: 88911a5261eSBarry Smith + J - the `MATMFFD` matrix-free matrix 890a5b23f4aSJose E. Roman . history - space to hold the history 891e884886fSBarry Smith - nhistory - number of entries in history, if more entries are generated than 892e884886fSBarry Smith nhistory, then the later ones are discarded 893e884886fSBarry Smith 894e884886fSBarry Smith Level: advanced 895e884886fSBarry Smith 89611a5261eSBarry Smith Note: 89711a5261eSBarry Smith Use `MatMFFDResetHHistory()` to reset the history counter and collect 898e884886fSBarry Smith a new batch of differencing parameters, h. 899e884886fSBarry Smith 90011a5261eSBarry Smith .seealso: `MATMFFD`, `MatMFFDGetH()`, `MatCreateSNESMF()`, 901db781477SPatrick Sanan `MatMFFDResetHHistory()`, `MatMFFDSetFunctionError()` 902e884886fSBarry Smith @*/ 903d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetHHistory(Mat J, PetscScalar history[], PetscInt nhistory) 904d71ae5a4SJacob Faibussowitsch { 905e884886fSBarry Smith PetscFunctionBegin; 90688b4c220SStefano Zampini PetscValidHeaderSpecific(J, MAT_CLASSID, 1); 907dadcf809SJacob Faibussowitsch if (history) PetscValidScalarPointer(history, 2); 90888b4c220SStefano Zampini PetscValidLogicalCollectiveInt(J, nhistory, 3); 909cac4c232SBarry Smith PetscUseMethod(J, "MatMFFDSetHHistory_C", (Mat, PetscScalar[], PetscInt), (J, history, nhistory)); 910*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 911e884886fSBarry Smith } 912e884886fSBarry Smith 913e884886fSBarry Smith /*@ 914e884886fSBarry Smith MatMFFDResetHHistory - Resets the counter to zero to begin 91511a5261eSBarry Smith collecting a new set of differencing histories for the `MATMFFD` matrix 916e884886fSBarry Smith 917c3339decSBarry Smith Logically Collective 918e884886fSBarry Smith 919e884886fSBarry Smith Input Parameters: 920e884886fSBarry Smith . J - the matrix-free matrix context 921e884886fSBarry Smith 922e884886fSBarry Smith Level: advanced 923e884886fSBarry Smith 92411a5261eSBarry Smith Note: 92511a5261eSBarry Smith Use `MatMFFDSetHHistory()` to create the original history counter. 926e884886fSBarry Smith 92711a5261eSBarry Smith .seealso: `MATMFFD`, `MatMFFDGetH()`, `MatCreateSNESMF()`, 928db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDSetFunctionError()` 929e884886fSBarry Smith @*/ 930d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDResetHHistory(Mat J) 931d71ae5a4SJacob Faibussowitsch { 932e884886fSBarry Smith PetscFunctionBegin; 93388b4c220SStefano Zampini PetscValidHeaderSpecific(J, MAT_CLASSID, 1); 934cac4c232SBarry Smith PetscTryMethod(J, "MatMFFDResetHHistory_C", (Mat), (J)); 935*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 936e884886fSBarry Smith } 937e884886fSBarry Smith 938487a658cSBarry Smith /*@ 939e884886fSBarry Smith MatMFFDSetBase - Sets the vector U at which matrix vector products of the 94011a5261eSBarry Smith Jacobian are computed for the `MATMFFD` matrix 941e884886fSBarry Smith 942c3339decSBarry Smith Logically Collective 943e884886fSBarry Smith 944e884886fSBarry Smith Input Parameters: 94511a5261eSBarry Smith + J - the `MATMFFD` matrix 9463ec795f1SBarry Smith . U - the vector 947bcddec3dSBarry Smith - F - (optional) vector that contains F(u) if it has been already computed 948e884886fSBarry Smith 94995452b02SPatrick Sanan Notes: 95095452b02SPatrick Sanan This is rarely used directly 951e884886fSBarry Smith 9528af5ae88SBarry Smith If F is provided then it is not recomputed. Otherwise the function is evaluated at the base 95311a5261eSBarry Smith point during the first `MatMult()` after each call to `MatMFFDSetBase()`. 954dff2f722SBarry Smith 955e884886fSBarry Smith Level: advanced 956e884886fSBarry Smith 95711a5261eSBarry Smith .seealso: `MATMFFD`, `MatMult()`, `MatMFFDSetBase()` 958e884886fSBarry Smith @*/ 959d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetBase(Mat J, Vec U, Vec F) 960d71ae5a4SJacob Faibussowitsch { 961e884886fSBarry Smith PetscFunctionBegin; 9620700a824SBarry Smith PetscValidHeaderSpecific(J, MAT_CLASSID, 1); 9630700a824SBarry Smith PetscValidHeaderSpecific(U, VEC_CLASSID, 2); 9640700a824SBarry Smith if (F) PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 965cac4c232SBarry Smith PetscTryMethod(J, "MatMFFDSetBase_C", (Mat, Vec, Vec), (J, U, F)); 966*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 967e884886fSBarry Smith } 968e884886fSBarry Smith 969e884886fSBarry Smith /*@C 970e884886fSBarry Smith MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts 97111a5261eSBarry Smith it to satisfy some criteria for the `MATMFFD` matrix 972e884886fSBarry Smith 973c3339decSBarry Smith Logically Collective 974e884886fSBarry Smith 975e884886fSBarry Smith Input Parameters: 97611a5261eSBarry Smith + J - the `MATMFFD` matrix 977e884886fSBarry Smith . fun - the function that checks h 978e884886fSBarry Smith - ctx - any context needed by the function 979e884886fSBarry Smith 980e884886fSBarry Smith Options Database Keys: 98167b8a455SSatish Balay . -mat_mffd_check_positivity <bool> - Insure that U + h*a is non-negative 982e884886fSBarry Smith 983e884886fSBarry Smith Level: advanced 984e884886fSBarry Smith 98595452b02SPatrick Sanan Notes: 98611a5261eSBarry Smith For example, `MatMFFDCheckPositivity()` insures that all entries 987e884886fSBarry Smith of U + h*a are non-negative 988e884886fSBarry Smith 989755b7f64SBarry Smith The function you provide is called after the default h has been computed and allows you to 990755b7f64SBarry Smith modify it. 991755b7f64SBarry Smith 99211a5261eSBarry Smith .seealso: `MATMFFD`, `MatMFFDCheckPositivity()` 993e884886fSBarry Smith @*/ 994d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDSetCheckh(Mat J, PetscErrorCode (*fun)(void *, Vec, Vec, PetscScalar *), void *ctx) 995d71ae5a4SJacob Faibussowitsch { 996e884886fSBarry Smith PetscFunctionBegin; 9970700a824SBarry Smith PetscValidHeaderSpecific(J, MAT_CLASSID, 1); 998cac4c232SBarry Smith PetscTryMethod(J, "MatMFFDSetCheckh_C", (Mat, PetscErrorCode(*)(void *, Vec, Vec, PetscScalar *), void *), (J, fun, ctx)); 999*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1000e884886fSBarry Smith } 1001e884886fSBarry Smith 1002e884886fSBarry Smith /*@ 1003e884886fSBarry Smith MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or 100411a5261eSBarry Smith zero, decreases h until this is satisfied for a `MATMFFD` matrix 1005e884886fSBarry Smith 1006c3339decSBarry Smith Logically Collective 1007e884886fSBarry Smith 1008e884886fSBarry Smith Input Parameters: 1009e884886fSBarry Smith + U - base vector that is added to 1010e884886fSBarry Smith . a - vector that is added 1011e884886fSBarry Smith . h - scaling factor on a 1012e884886fSBarry Smith - dummy - context variable (unused) 1013e884886fSBarry Smith 1014e884886fSBarry Smith Options Database Keys: 101567b8a455SSatish Balay . -mat_mffd_check_positivity <bool> - Insure that U + h*a is nonnegative 1016e884886fSBarry Smith 1017e884886fSBarry Smith Level: advanced 1018e884886fSBarry Smith 101911a5261eSBarry Smith Note: 102095452b02SPatrick Sanan This is rarely used directly, rather it is passed as an argument to 102111a5261eSBarry Smith `MatMFFDSetCheckh()` 1022e884886fSBarry Smith 102311a5261eSBarry Smith .seealso: `MATMFFD`, `MatMFFDSetCheckh()` 1024e884886fSBarry Smith @*/ 1025d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDCheckPositivity(void *dummy, Vec U, Vec a, PetscScalar *h) 1026d71ae5a4SJacob Faibussowitsch { 1027e884886fSBarry Smith PetscReal val, minval; 1028e884886fSBarry Smith PetscScalar *u_vec, *a_vec; 1029e884886fSBarry Smith PetscInt i, n; 1030e884886fSBarry Smith MPI_Comm comm; 1031e884886fSBarry Smith 1032e884886fSBarry Smith PetscFunctionBegin; 103388b4c220SStefano Zampini PetscValidHeaderSpecific(U, VEC_CLASSID, 2); 103488b4c220SStefano Zampini PetscValidHeaderSpecific(a, VEC_CLASSID, 3); 1035dadcf809SJacob Faibussowitsch PetscValidScalarPointer(h, 4); 10369566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)U, &comm)); 10379566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u_vec)); 10389566063dSJacob Faibussowitsch PetscCall(VecGetArray(a, &a_vec)); 10399566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(U, &n)); 1040c068d9bbSLisandro Dalcin minval = PetscAbsScalar(*h) * PetscRealConstant(1.01); 1041e884886fSBarry Smith for (i = 0; i < n; i++) { 1042e884886fSBarry Smith if (PetscRealPart(u_vec[i] + *h * a_vec[i]) <= 0.0) { 1043e884886fSBarry Smith val = PetscAbsScalar(u_vec[i] / a_vec[i]); 1044e884886fSBarry Smith if (val < minval) minval = val; 1045e884886fSBarry Smith } 1046e884886fSBarry Smith } 10479566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u_vec)); 10489566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(a, &a_vec)); 10491c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&minval, &val, 1, MPIU_REAL, MPIU_MIN, comm)); 1050e884886fSBarry Smith if (val <= PetscAbsScalar(*h)) { 10519566063dSJacob Faibussowitsch PetscCall(PetscInfo(U, "Scaling back h from %g to %g\n", (double)PetscRealPart(*h), (double)(.99 * val))); 1052e884886fSBarry Smith if (PetscRealPart(*h) > 0.0) *h = 0.99 * val; 1053e884886fSBarry Smith else *h = -0.99 * val; 1054e884886fSBarry Smith } 1055*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1056e884886fSBarry Smith } 1057