1d2dddef7SLois Curfman McInnes 2af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 36370053bSBarry Smith /* matimpl.h is needed only for logging of matrix operation */ 4af0996ceSBarry Smith #include <petsc/private/matimpl.h> 5d2dddef7SLois Curfman McInnes 625acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESUnSetMatrixFreeParameter(SNES); 725acbd8eSLisandro Dalcin 825acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterCreate_More(SNES, Vec, void **); 925acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterCompute_More(SNES, void *, Vec, Vec, PetscReal *, PetscReal *); 1025acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterDestroy_More(void *); 11d2dddef7SLois Curfman McInnes 12d2dddef7SLois Curfman McInnes typedef struct { /* default context for matrix-free SNES */ 13d2dddef7SLois Curfman McInnes SNES snes; /* SNES context */ 14d2dddef7SLois Curfman McInnes Vec w; /* work vector */ 1574637425SBarry Smith MatNullSpace sp; /* null space context */ 16145595abSSatish Balay PetscReal error_rel; /* square root of relative error in computing function */ 17145595abSSatish Balay PetscReal umin; /* minimum allowable u'a value relative to |u|_1 */ 18f5af7f23SKarl Rupp PetscBool jorge; /* flag indicating use of Jorge's method for determining the differencing parameter */ 19145595abSSatish Balay PetscReal h; /* differencing parameter */ 20ace3abfcSBarry Smith PetscBool need_h; /* flag indicating whether we must compute h */ 21ace3abfcSBarry Smith PetscBool need_err; /* flag indicating whether we must currently compute error_rel */ 22ace3abfcSBarry Smith PetscBool compute_err; /* flag indicating whether we must ever compute error_rel */ 23a7cc72afSBarry Smith PetscInt compute_err_iter; /* last iter where we've computer error_rel */ 24a7cc72afSBarry Smith PetscInt compute_err_freq; /* frequency of computing error_rel */ 25d2dddef7SLois Curfman McInnes void *data; /* implementation-specific data */ 26d2dddef7SLois Curfman McInnes } MFCtx_Private; 27d2dddef7SLois Curfman McInnes 28d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMatrixFreeDestroy2_Private(Mat mat) 29d71ae5a4SJacob Faibussowitsch { 30d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 31d2dddef7SLois Curfman McInnes 3207a3eed9SLois Curfman McInnes PetscFunctionBegin; 339566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->w)); 359566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&ctx->sp)); 369566063dSJacob Faibussowitsch if (ctx->jorge || ctx->compute_err) PetscCall(SNESDiffParameterDestroy_More(ctx->data)); 379566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 39d2dddef7SLois Curfman McInnes } 40d2dddef7SLois Curfman McInnes 41d2dddef7SLois Curfman McInnes /* 42d2dddef7SLois Curfman McInnes SNESMatrixFreeView2_Private - Views matrix-free parameters. 43d2dddef7SLois Curfman McInnes */ 44d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMatrixFreeView2_Private(Mat J, PetscViewer viewer) 45d71ae5a4SJacob Faibussowitsch { 46d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 47ace3abfcSBarry Smith PetscBool iascii; 48d2dddef7SLois Curfman McInnes 4907a3eed9SLois Curfman McInnes PetscFunctionBegin; 509566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J, &ctx)); 519566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 5232077d6dSBarry Smith if (iascii) { 539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " SNES matrix-free approximation:\n")); 5448a46eb9SPierre Jolivet if (ctx->jorge) PetscCall(PetscViewerASCIIPrintf(viewer, " using Jorge's method of determining differencing parameter\n")); 559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " err=%g (relative error in function evaluation)\n", (double)ctx->error_rel)); 569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " umin=%g (minimum iterate parameter)\n", (double)ctx->umin)); 5748a46eb9SPierre Jolivet if (ctx->compute_err) PetscCall(PetscViewerASCIIPrintf(viewer, " freq_err=%" PetscInt_FMT " (frequency for computing err)\n", ctx->compute_err_freq)); 58758f395bSBarry Smith } 593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 60d2dddef7SLois Curfman McInnes } 61d2dddef7SLois Curfman McInnes 62d2dddef7SLois Curfman McInnes /* 63d2dddef7SLois Curfman McInnes SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector 64d2dddef7SLois Curfman McInnes product, y = F'(u)*a: 65d2dddef7SLois Curfman McInnes y = (F(u + ha) - F(u)) /h, 66d2dddef7SLois Curfman McInnes where F = nonlinear function, as set by SNESSetFunction() 67d2dddef7SLois Curfman McInnes u = current iterate 68d2dddef7SLois Curfman McInnes h = difference interval 69d2dddef7SLois Curfman McInnes */ 70d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMatrixFreeMult2_Private(Mat mat, Vec a, Vec y) 71d71ae5a4SJacob Faibussowitsch { 72d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 73d2dddef7SLois Curfman McInnes SNES snes; 74145595abSSatish Balay PetscReal h, norm, sum, umin, noise; 7579f36460SBarry Smith PetscScalar hs, dot; 76d2dddef7SLois Curfman McInnes Vec w, U, F; 77e6ed2b1fSLois Curfman McInnes MPI_Comm comm; 7877431f27SBarry Smith PetscInt iter; 795f80ce2aSJacob Faibussowitsch PetscErrorCode (*eval_fct)(SNES, Vec, Vec); 80d2dddef7SLois Curfman McInnes 8107a3eed9SLois Curfman McInnes PetscFunctionBegin; 82d2dddef7SLois Curfman McInnes /* We log matrix-free matrix-vector products separately, so that we can 83d2dddef7SLois Curfman McInnes separate the performance monitoring from the cases that use conventional 84d2dddef7SLois Curfman McInnes storage. We may eventually modify event logging to associate events 85d2dddef7SLois Curfman McInnes with particular objects, hence alleviating the more general problem. */ 869566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MATMFFD_Mult, a, y, 0, 0)); 87d2dddef7SLois Curfman McInnes 889566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)mat, &comm)); 899566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 90e6ed2b1fSLois Curfman McInnes snes = ctx->snes; 91e6ed2b1fSLois Curfman McInnes w = ctx->w; 92e6ed2b1fSLois Curfman McInnes umin = ctx->umin; 93e6ed2b1fSLois Curfman McInnes 949566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &U)); 95d2dddef7SLois Curfman McInnes eval_fct = SNESComputeFunction; 969566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &F, NULL, NULL)); 97295f7fbbSLois Curfman McInnes 98d2dddef7SLois Curfman McInnes /* Determine a "good" step size, h */ 99d2dddef7SLois Curfman McInnes if (ctx->need_h) { 100d2dddef7SLois Curfman McInnes /* Use Jorge's method to compute h */ 101d2dddef7SLois Curfman McInnes if (ctx->jorge) { 1029566063dSJacob Faibussowitsch PetscCall(SNESDiffParameterCompute_More(snes, ctx->data, U, a, &noise, &h)); 103d2dddef7SLois Curfman McInnes 104d2dddef7SLois Curfman McInnes /* Use the Brown/Saad method to compute h */ 105d2dddef7SLois Curfman McInnes } else { 106295f7fbbSLois Curfman McInnes /* Compute error if desired */ 1079566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &iter)); 108145595abSSatish Balay if ((ctx->need_err) || ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter - 1) % ctx->compute_err_freq)))) { 109d2dddef7SLois Curfman McInnes /* Use Jorge's method to compute noise */ 1109566063dSJacob Faibussowitsch PetscCall(SNESDiffParameterCompute_More(snes, ctx->data, U, a, &noise, &h)); 1111aa26658SKarl Rupp 112369cc0aeSBarry Smith ctx->error_rel = PetscSqrtReal(noise); 1131aa26658SKarl Rupp 1149566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Using Jorge's noise: noise=%g, sqrt(noise)=%g, h_more=%g\n", (double)noise, (double)ctx->error_rel, (double)h)); 1151aa26658SKarl Rupp 116295f7fbbSLois Curfman McInnes ctx->compute_err_iter = iter; 1173b2a6d2fSBarry Smith ctx->need_err = PETSC_FALSE; 118d2dddef7SLois Curfman McInnes } 119e6ed2b1fSLois Curfman McInnes 1209566063dSJacob Faibussowitsch PetscCall(VecDotBegin(U, a, &dot)); 1219566063dSJacob Faibussowitsch PetscCall(VecNormBegin(a, NORM_1, &sum)); 1229566063dSJacob Faibussowitsch PetscCall(VecNormBegin(a, NORM_2, &norm)); 1239566063dSJacob Faibussowitsch PetscCall(VecDotEnd(U, a, &dot)); 1249566063dSJacob Faibussowitsch PetscCall(VecNormEnd(a, NORM_1, &sum)); 1259566063dSJacob Faibussowitsch PetscCall(VecNormEnd(a, NORM_2, &norm)); 126e6ed2b1fSLois Curfman McInnes 127d2dddef7SLois Curfman McInnes /* Safeguard for step sizes too small */ 1281aa26658SKarl Rupp if (sum == 0.0) { 1291aa26658SKarl Rupp dot = 1.0; 1301aa26658SKarl Rupp norm = 1.0; 1311aa26658SKarl Rupp } else if (PetscAbsScalar(dot) < umin * sum && PetscRealPart(dot) >= 0.0) dot = umin * sum; 132145595abSSatish Balay else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin * sum) dot = -umin * sum; 133145595abSSatish Balay h = PetscRealPart(ctx->error_rel * dot / (norm * norm)); 134d2dddef7SLois Curfman McInnes } 1351aa26658SKarl Rupp } else h = ctx->h; 1361aa26658SKarl Rupp 1379566063dSJacob Faibussowitsch if (!ctx->jorge || !ctx->need_h) PetscCall(PetscInfo(snes, "h = %g\n", (double)h)); 138d2dddef7SLois Curfman McInnes 139d2dddef7SLois Curfman McInnes /* Evaluate function at F(u + ha) */ 140db0b4b35SLois Curfman McInnes hs = h; 1419566063dSJacob Faibussowitsch PetscCall(VecWAXPY(w, hs, a, U)); 1429566063dSJacob Faibussowitsch PetscCall(eval_fct(snes, w, y)); 1439566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, F)); 1449566063dSJacob Faibussowitsch PetscCall(VecScale(y, 1.0 / hs)); 1459566063dSJacob Faibussowitsch if (mat->nullsp) PetscCall(MatNullSpaceRemove(mat->nullsp, y)); 146d2dddef7SLois Curfman McInnes 1479566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MATMFFD_Mult, a, y, 0, 0)); 1483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 149d2dddef7SLois Curfman McInnes } 150d2dddef7SLois Curfman McInnes 151d2dddef7SLois Curfman McInnes /*@C 152f6dfbefdSBarry Smith MatCreateSNESMFMore - Creates a matrix-free matrix 153f6dfbefdSBarry Smith context for use with a `SNES` solver that uses the More method to compute an optimal h based on the noise of the function. This matrix can be used as 154f6dfbefdSBarry Smith the Jacobian argument for the routine `SNESSetJacobian()`. 155d2dddef7SLois Curfman McInnes 156d2dddef7SLois Curfman McInnes Input Parameters: 157f6dfbefdSBarry Smith + snes - the `SNES` context 158f6dfbefdSBarry Smith - x - vector where `SNES` solution is to be stored. 159d2dddef7SLois Curfman McInnes 160d2dddef7SLois Curfman McInnes Output Parameter: 161d2dddef7SLois Curfman McInnes . J - the matrix-free matrix 162d2dddef7SLois Curfman McInnes 163f6dfbefdSBarry Smith Options Database Keys: 164f6dfbefdSBarry Smith + -snes_mf_err <error_rel> - see `MatCreateSNESMF()` 165f6dfbefdSBarry Smith . -snes_mf_umin <umin> - see `MatCreateSNESMF()` 166f6dfbefdSBarry Smith . -snes_mf_compute_err - compute the square root or relative error in function 167f6dfbefdSBarry Smith . -snes_mf_freq_err <freq> - set the frequency to recompute the parameters 168f6dfbefdSBarry Smith - -snes_mf_jorge - use the method of Jorge More 169f6dfbefdSBarry Smith 17090c26489SBarry Smith Level: advanced 17190c26489SBarry Smith 172d2dddef7SLois Curfman McInnes Notes: 173f6dfbefdSBarry Smith This is an experimental approach, use `MatCreateSNESMF()`. 174f6dfbefdSBarry Smith 175d2dddef7SLois Curfman McInnes The matrix-free matrix context merely contains the function pointers 176d2dddef7SLois Curfman McInnes and work space for performing finite difference approximations of 177d2dddef7SLois Curfman McInnes Jacobian-vector products, J(u)*a, via 178d2dddef7SLois Curfman McInnes 17937fdd005SBarry Smith .vb 18037fdd005SBarry Smith J(u)*a = [J(u+h*a) - J(u)]/h, 18137fdd005SBarry Smith where by default 18237fdd005SBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 18337fdd005SBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 18437fdd005SBarry Smith where 18537fdd005SBarry Smith error_rel = square root of relative error in function evaluation 18637fdd005SBarry Smith umin = minimum iterate parameter 18737fdd005SBarry Smith Alternatively, the differencing parameter, h, can be set using 18837fdd005SBarry Smith Jorge's nifty new strategy if one specifies the option 18937fdd005SBarry Smith -snes_mf_jorge 19037fdd005SBarry Smith .ve 191d2dddef7SLois Curfman McInnes 192f6dfbefdSBarry Smith The user can set these parameters via `MatMFFDSetFunctionError()`. 193d2dddef7SLois Curfman McInnes 194f6dfbefdSBarry Smith The user should call `MatDestroy()` when finished with the matrix-free 195d2dddef7SLois Curfman McInnes matrix context. 196d2dddef7SLois Curfman McInnes 19776fbde31SPierre Jolivet .seealso: `SNESCreateMF()`, `MatCreateMFFD()`, `MatDestroy()`, `MatMFFDSetFunctionError()` 198d2dddef7SLois Curfman McInnes @*/ 199d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSNESMFMore(SNES snes, Vec x, Mat *J) 200d71ae5a4SJacob Faibussowitsch { 201d2dddef7SLois Curfman McInnes MPI_Comm comm; 202d2dddef7SLois Curfman McInnes MFCtx_Private *mfctx; 203a7cc72afSBarry Smith PetscInt n, nloc; 204ace3abfcSBarry Smith PetscBool flg; 205d2dddef7SLois Curfman McInnes char p[64]; 206d2dddef7SLois Curfman McInnes 20707a3eed9SLois Curfman McInnes PetscFunctionBegin; 2084dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&mfctx)); 2099e5d0892SLisandro Dalcin mfctx->sp = NULL; 210d2dddef7SLois Curfman McInnes mfctx->snes = snes; 21177d8c4bbSBarry Smith mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 212d2dddef7SLois Curfman McInnes mfctx->umin = 1.e-6; 213d2dddef7SLois Curfman McInnes mfctx->h = 0.0; 214a7cc72afSBarry Smith mfctx->need_h = PETSC_TRUE; 215a7cc72afSBarry Smith mfctx->need_err = PETSC_FALSE; 216a7cc72afSBarry Smith mfctx->compute_err = PETSC_FALSE; 217295f7fbbSLois Curfman McInnes mfctx->compute_err_freq = 0; 218295f7fbbSLois Curfman McInnes mfctx->compute_err_iter = -1; 21990d69ab7SBarry Smith mfctx->compute_err = PETSC_FALSE; 22090d69ab7SBarry Smith mfctx->jorge = PETSC_FALSE; 2211aa26658SKarl Rupp 2229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_mf_err", &mfctx->error_rel, NULL)); 2239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_mf_umin", &mfctx->umin, NULL)); 2249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_mf_jorge", &mfctx->jorge, NULL)); 2259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_mf_compute_err", &mfctx->compute_err, NULL)); 2269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_mf_freq_err", &mfctx->compute_err_freq, &flg)); 227295f7fbbSLois Curfman McInnes if (flg) { 228295f7fbbSLois Curfman McInnes if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0; 2293b2a6d2fSBarry Smith mfctx->compute_err = PETSC_TRUE; 230295f7fbbSLois Curfman McInnes } 2313b2a6d2fSBarry Smith if (mfctx->compute_err) mfctx->need_err = PETSC_TRUE; 232295f7fbbSLois Curfman McInnes if (mfctx->jorge || mfctx->compute_err) { 2339566063dSJacob Faibussowitsch PetscCall(SNESDiffParameterCreate_More(snes, x, &mfctx->data)); 2349e5d0892SLisandro Dalcin } else mfctx->data = NULL; 235d2dddef7SLois Curfman McInnes 2369566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasHelp(((PetscObject)snes)->options, &flg)); 2379566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(p, "-", sizeof(p))); 2389566063dSJacob Faibussowitsch if (((PetscObject)snes)->prefix) PetscCall(PetscStrlcat(p, ((PetscObject)snes)->prefix, sizeof(p))); 239d2dddef7SLois Curfman McInnes if (flg) { 2409566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), " Matrix-free Options (via SNES):\n")); 2419566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), " %ssnes_mf_err <err>: set sqrt of relative error in function (default %g)\n", p, (double)mfctx->error_rel)); 2429566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), " %ssnes_mf_umin <umin>: see users manual (default %g)\n", p, (double)mfctx->umin)); 2439566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), " %ssnes_mf_jorge: use Jorge More's method\n", p)); 2449566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), " %ssnes_mf_compute_err: compute sqrt or relative error in function\n", p)); 2459566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), " %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)\n", p)); 2469566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), " %ssnes_mf_noise_file <file>: set file for printing noise info\n", p)); 247d2dddef7SLois Curfman McInnes } 2489566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &mfctx->w)); 2499566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)x, &comm)); 2509566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &n)); 2519566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x, &nloc)); 2529566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, J)); 2539566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J, nloc, n, n, n)); 2549566063dSJacob Faibussowitsch PetscCall(MatSetType(*J, MATSHELL)); 2559566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(*J, mfctx)); 2569566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_MULT, (void (*)(void))SNESMatrixFreeMult2_Private)); 2579566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_DESTROY, (void (*)(void))SNESMatrixFreeDestroy2_Private)); 2589566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_VIEW, (void (*)(void))SNESMatrixFreeView2_Private)); 2599566063dSJacob Faibussowitsch PetscCall(MatSetUp(*J)); 260d2dddef7SLois Curfman McInnes 2613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 262d2dddef7SLois Curfman McInnes } 263d2dddef7SLois Curfman McInnes 26490c26489SBarry Smith /*@C 265f6dfbefdSBarry Smith MatSNESMFMoreSetParameters - Sets the parameters for the approximation of 266f6dfbefdSBarry Smith matrix-vector products using finite differences, see `MatCreateSNESMFMore()` 267d2dddef7SLois Curfman McInnes 268d2dddef7SLois Curfman McInnes Input Parameters: 269758f395bSBarry Smith + mat - the matrix 270*e4094ef1SJacob Faibussowitsch . error - relative error (should be set to the square root of the relative error in the function evaluations) 271d2dddef7SLois Curfman McInnes . umin - minimum allowable u-value 272758f395bSBarry Smith - h - differencing parameter 273d2dddef7SLois Curfman McInnes 274f6dfbefdSBarry Smith Options Database Keys: 275f6dfbefdSBarry Smith + -snes_mf_err <error_rel> - see `MatCreateSNESMF()` 276f6dfbefdSBarry Smith . -snes_mf_umin <umin> - see `MatCreateSNESMF()` 277f6dfbefdSBarry Smith . -snes_mf_compute_err - compute the square root or relative error in function 278f6dfbefdSBarry Smith . -snes_mf_freq_err <freq> - set the frequency to recompute the parameters 279f6dfbefdSBarry Smith - -snes_mf_jorge - use the method of Jorge More 280f6dfbefdSBarry Smith 28190c26489SBarry Smith Level: advanced 28290c26489SBarry Smith 283f6dfbefdSBarry Smith Note: 284d2dddef7SLois Curfman McInnes If the user sets the parameter h directly, then this value will be used 285f6dfbefdSBarry Smith instead of the default computation as discussed in `MatCreateSNESMFMore()` 286d2dddef7SLois Curfman McInnes 287f6dfbefdSBarry Smith .seealso: `MatCreateSNESMF()`, `MatCreateSNESMFMore()` 288d2dddef7SLois Curfman McInnes @*/ 289d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSNESMFMoreSetParameters(Mat mat, PetscReal error, PetscReal umin, PetscReal h) 290d71ae5a4SJacob Faibussowitsch { 291d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 292d2dddef7SLois Curfman McInnes 29307a3eed9SLois Curfman McInnes PetscFunctionBegin; 2949566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 295d2dddef7SLois Curfman McInnes if (ctx) { 29613bcc0bdSJacob Faibussowitsch if (error != (PetscReal)PETSC_DEFAULT) ctx->error_rel = error; 29713bcc0bdSJacob Faibussowitsch if (umin != (PetscReal)PETSC_DEFAULT) ctx->umin = umin; 29813bcc0bdSJacob Faibussowitsch if (h != (PetscReal)PETSC_DEFAULT) { 299d2dddef7SLois Curfman McInnes ctx->h = h; 300a7cc72afSBarry Smith ctx->need_h = PETSC_FALSE; 301d2dddef7SLois Curfman McInnes } 302d2dddef7SLois Curfman McInnes } 3033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 304d2dddef7SLois Curfman McInnes } 305d2dddef7SLois Curfman McInnes 306d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESUnSetMatrixFreeParameter(SNES snes) 307d71ae5a4SJacob Faibussowitsch { 308d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 309d2dddef7SLois Curfman McInnes Mat mat; 310d2dddef7SLois Curfman McInnes 31107a3eed9SLois Curfman McInnes PetscFunctionBegin; 3129566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &mat, NULL, NULL, NULL)); 3139566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 314a7cc72afSBarry Smith if (ctx) ctx->need_h = PETSC_TRUE; 3153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 316d2dddef7SLois Curfman McInnes } 317