xref: /petsc/src/snes/interface/noise/snesmfj2.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
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 
289371c9d4SSatish Balay PetscErrorCode SNESMatrixFreeDestroy2_Private(Mat mat) {
29d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
30d2dddef7SLois Curfman McInnes 
3107a3eed9SLois Curfman McInnes   PetscFunctionBegin;
329566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->w));
349566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceDestroy(&ctx->sp));
359566063dSJacob Faibussowitsch   if (ctx->jorge || ctx->compute_err) PetscCall(SNESDiffParameterDestroy_More(ctx->data));
369566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
3707a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
38d2dddef7SLois Curfman McInnes }
39d2dddef7SLois Curfman McInnes 
40d2dddef7SLois Curfman McInnes /*
41d2dddef7SLois Curfman McInnes    SNESMatrixFreeView2_Private - Views matrix-free parameters.
42d2dddef7SLois Curfman McInnes  */
439371c9d4SSatish Balay PetscErrorCode SNESMatrixFreeView2_Private(Mat J, PetscViewer viewer) {
44d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
45ace3abfcSBarry Smith   PetscBool      iascii;
46d2dddef7SLois Curfman McInnes 
4707a3eed9SLois Curfman McInnes   PetscFunctionBegin;
489566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J, &ctx));
499566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
5032077d6dSBarry Smith   if (iascii) {
519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  SNES matrix-free approximation:\n"));
5248a46eb9SPierre Jolivet     if (ctx->jorge) PetscCall(PetscViewerASCIIPrintf(viewer, "    using Jorge's method of determining differencing parameter\n"));
539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    err=%g (relative error in function evaluation)\n", (double)ctx->error_rel));
549566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    umin=%g (minimum iterate parameter)\n", (double)ctx->umin));
5548a46eb9SPierre Jolivet     if (ctx->compute_err) PetscCall(PetscViewerASCIIPrintf(viewer, "    freq_err=%" PetscInt_FMT " (frequency for computing err)\n", ctx->compute_err_freq));
56758f395bSBarry Smith   }
5707a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
58d2dddef7SLois Curfman McInnes }
59d2dddef7SLois Curfman McInnes 
60d2dddef7SLois Curfman McInnes /*
61d2dddef7SLois Curfman McInnes   SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector
62d2dddef7SLois Curfman McInnes   product, y = F'(u)*a:
63d2dddef7SLois Curfman McInnes         y = (F(u + ha) - F(u)) /h,
64d2dddef7SLois Curfman McInnes   where F = nonlinear function, as set by SNESSetFunction()
65d2dddef7SLois Curfman McInnes         u = current iterate
66d2dddef7SLois Curfman McInnes         h = difference interval
67d2dddef7SLois Curfman McInnes */
689371c9d4SSatish Balay PetscErrorCode SNESMatrixFreeMult2_Private(Mat mat, Vec a, Vec y) {
69d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
70d2dddef7SLois Curfman McInnes   SNES           snes;
71145595abSSatish Balay   PetscReal      h, norm, sum, umin, noise;
7279f36460SBarry Smith   PetscScalar    hs, dot;
73d2dddef7SLois Curfman McInnes   Vec            w, U, F;
74e6ed2b1fSLois Curfman McInnes   MPI_Comm       comm;
7577431f27SBarry Smith   PetscInt       iter;
765f80ce2aSJacob Faibussowitsch   PetscErrorCode (*eval_fct)(SNES, Vec, Vec);
77d2dddef7SLois Curfman McInnes 
7807a3eed9SLois Curfman McInnes   PetscFunctionBegin;
79d2dddef7SLois Curfman McInnes   /* We log matrix-free matrix-vector products separately, so that we can
80d2dddef7SLois Curfman McInnes      separate the performance monitoring from the cases that use conventional
81d2dddef7SLois Curfman McInnes      storage.  We may eventually modify event logging to associate events
82d2dddef7SLois Curfman McInnes      with particular objects, hence alleviating the more general problem. */
839566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MATMFFD_Mult, a, y, 0, 0));
84d2dddef7SLois Curfman McInnes 
859566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
869566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
87e6ed2b1fSLois Curfman McInnes   snes = ctx->snes;
88e6ed2b1fSLois Curfman McInnes   w    = ctx->w;
89e6ed2b1fSLois Curfman McInnes   umin = ctx->umin;
90e6ed2b1fSLois Curfman McInnes 
919566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &U));
92d2dddef7SLois Curfman McInnes   eval_fct = SNESComputeFunction;
939566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &F, NULL, NULL));
94295f7fbbSLois Curfman McInnes 
95d2dddef7SLois Curfman McInnes   /* Determine a "good" step size, h */
96d2dddef7SLois Curfman McInnes   if (ctx->need_h) {
97d2dddef7SLois Curfman McInnes     /* Use Jorge's method to compute h */
98d2dddef7SLois Curfman McInnes     if (ctx->jorge) {
999566063dSJacob Faibussowitsch       PetscCall(SNESDiffParameterCompute_More(snes, ctx->data, U, a, &noise, &h));
100d2dddef7SLois Curfman McInnes 
101d2dddef7SLois Curfman McInnes       /* Use the Brown/Saad method to compute h */
102d2dddef7SLois Curfman McInnes     } else {
103295f7fbbSLois Curfman McInnes       /* Compute error if desired */
1049566063dSJacob Faibussowitsch       PetscCall(SNESGetIterationNumber(snes, &iter));
105145595abSSatish Balay       if ((ctx->need_err) || ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter - 1) % ctx->compute_err_freq)))) {
106d2dddef7SLois Curfman McInnes         /* Use Jorge's method to compute noise */
1079566063dSJacob Faibussowitsch         PetscCall(SNESDiffParameterCompute_More(snes, ctx->data, U, a, &noise, &h));
1081aa26658SKarl Rupp 
109369cc0aeSBarry Smith         ctx->error_rel = PetscSqrtReal(noise);
1101aa26658SKarl Rupp 
1119566063dSJacob 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));
1121aa26658SKarl Rupp 
113295f7fbbSLois Curfman McInnes         ctx->compute_err_iter = iter;
1143b2a6d2fSBarry Smith         ctx->need_err         = PETSC_FALSE;
115d2dddef7SLois Curfman McInnes       }
116e6ed2b1fSLois Curfman McInnes 
1179566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(U, a, &dot));
1189566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(a, NORM_1, &sum));
1199566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(a, NORM_2, &norm));
1209566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(U, a, &dot));
1219566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(a, NORM_1, &sum));
1229566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(a, NORM_2, &norm));
123e6ed2b1fSLois Curfman McInnes 
124d2dddef7SLois Curfman McInnes       /* Safeguard for step sizes too small */
1251aa26658SKarl Rupp       if (sum == 0.0) {
1261aa26658SKarl Rupp         dot  = 1.0;
1271aa26658SKarl Rupp         norm = 1.0;
1281aa26658SKarl Rupp       } else if (PetscAbsScalar(dot) < umin * sum && PetscRealPart(dot) >= 0.0) dot = umin * sum;
129145595abSSatish Balay       else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin * sum) dot = -umin * sum;
130145595abSSatish Balay       h = PetscRealPart(ctx->error_rel * dot / (norm * norm));
131d2dddef7SLois Curfman McInnes     }
1321aa26658SKarl Rupp   } else h = ctx->h;
1331aa26658SKarl Rupp 
1349566063dSJacob Faibussowitsch   if (!ctx->jorge || !ctx->need_h) PetscCall(PetscInfo(snes, "h = %g\n", (double)h));
135d2dddef7SLois Curfman McInnes 
136d2dddef7SLois Curfman McInnes   /* Evaluate function at F(u + ha) */
137db0b4b35SLois Curfman McInnes   hs = h;
1389566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(w, hs, a, U));
1399566063dSJacob Faibussowitsch   PetscCall(eval_fct(snes, w, y));
1409566063dSJacob Faibussowitsch   PetscCall(VecAXPY(y, -1.0, F));
1419566063dSJacob Faibussowitsch   PetscCall(VecScale(y, 1.0 / hs));
1429566063dSJacob Faibussowitsch   if (mat->nullsp) PetscCall(MatNullSpaceRemove(mat->nullsp, y));
143d2dddef7SLois Curfman McInnes 
1449566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MATMFFD_Mult, a, y, 0, 0));
14507a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
146d2dddef7SLois Curfman McInnes }
147d2dddef7SLois Curfman McInnes 
148d2dddef7SLois Curfman McInnes /*@C
149f6dfbefdSBarry Smith    MatCreateSNESMFMore - Creates a matrix-free matrix
150f6dfbefdSBarry 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
151f6dfbefdSBarry Smith    the Jacobian argument for the routine `SNESSetJacobian()`.
152d2dddef7SLois Curfman McInnes 
153d2dddef7SLois Curfman McInnes    Input Parameters:
154f6dfbefdSBarry Smith +  snes - the `SNES` context
155f6dfbefdSBarry Smith -  x - vector where `SNES` solution is to be stored.
156d2dddef7SLois Curfman McInnes 
157d2dddef7SLois Curfman McInnes    Output Parameter:
158d2dddef7SLois Curfman McInnes .  J - the matrix-free matrix
159d2dddef7SLois Curfman McInnes 
160f6dfbefdSBarry Smith    Options Database Keys:
161f6dfbefdSBarry Smith +  -snes_mf_err <error_rel> - see `MatCreateSNESMF()`
162f6dfbefdSBarry Smith .  -snes_mf_umin <umin> - see `MatCreateSNESMF()`
163f6dfbefdSBarry Smith .  -snes_mf_compute_err - compute the square root or relative error in function
164f6dfbefdSBarry Smith .  -snes_mf_freq_err <freq> - set the frequency to recompute the parameters
165f6dfbefdSBarry Smith -  -snes_mf_jorge - use the method of Jorge More
166f6dfbefdSBarry Smith 
16790c26489SBarry Smith    Level: advanced
16890c26489SBarry Smith 
169d2dddef7SLois Curfman McInnes    Notes:
170f6dfbefdSBarry Smith    This is an experimental approach, use `MatCreateSNESMF()`.
171f6dfbefdSBarry Smith 
172d2dddef7SLois Curfman McInnes    The matrix-free matrix context merely contains the function pointers
173d2dddef7SLois Curfman McInnes    and work space for performing finite difference approximations of
174d2dddef7SLois Curfman McInnes    Jacobian-vector products, J(u)*a, via
175d2dddef7SLois Curfman McInnes 
176d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h,
177d2dddef7SLois Curfman McInnes $   where by default
178d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
179d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
180d2dddef7SLois Curfman McInnes $   where
181d2dddef7SLois Curfman McInnes $        error_rel = square root of relative error in
182d2dddef7SLois Curfman McInnes $                    function evaluation
183d2dddef7SLois Curfman McInnes $        umin = minimum iterate parameter
184d2dddef7SLois Curfman McInnes $   Alternatively, the differencing parameter, h, can be set using
185e6ed2b1fSLois Curfman McInnes $   Jorge's nifty new strategy if one specifies the option
186e6ed2b1fSLois Curfman McInnes $          -snes_mf_jorge
187d2dddef7SLois Curfman McInnes 
188f6dfbefdSBarry Smith    The user can set these parameters via `MatMFFDSetFunctionError()`.
189a7f22e61SSatish Balay    See Users-Manual: ch_snes for details
190d2dddef7SLois Curfman McInnes 
191f6dfbefdSBarry Smith    The user should call `MatDestroy()` when finished with the matrix-free
192d2dddef7SLois Curfman McInnes    matrix context.
193d2dddef7SLois Curfman McInnes 
194f6dfbefdSBarry Smith .seealso: `SNESCreateMF`(), `MatCreateMFFD()`, `MatDestroy()`, `MatMFFDSetFunctionError()`
195d2dddef7SLois Curfman McInnes @*/
196f6dfbefdSBarry Smith PetscErrorCode MatCreateSNESMFMore(SNES snes, Vec x, Mat *J) {
197d2dddef7SLois Curfman McInnes   MPI_Comm       comm;
198d2dddef7SLois Curfman McInnes   MFCtx_Private *mfctx;
199a7cc72afSBarry Smith   PetscInt       n, nloc;
200ace3abfcSBarry Smith   PetscBool      flg;
201d2dddef7SLois Curfman McInnes   char           p[64];
202d2dddef7SLois Curfman McInnes 
20307a3eed9SLois Curfman McInnes   PetscFunctionBegin;
204*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&mfctx));
2059e5d0892SLisandro Dalcin   mfctx->sp               = NULL;
206d2dddef7SLois Curfman McInnes   mfctx->snes             = snes;
20777d8c4bbSBarry Smith   mfctx->error_rel        = PETSC_SQRT_MACHINE_EPSILON;
208d2dddef7SLois Curfman McInnes   mfctx->umin             = 1.e-6;
209d2dddef7SLois Curfman McInnes   mfctx->h                = 0.0;
210a7cc72afSBarry Smith   mfctx->need_h           = PETSC_TRUE;
211a7cc72afSBarry Smith   mfctx->need_err         = PETSC_FALSE;
212a7cc72afSBarry Smith   mfctx->compute_err      = PETSC_FALSE;
213295f7fbbSLois Curfman McInnes   mfctx->compute_err_freq = 0;
214295f7fbbSLois Curfman McInnes   mfctx->compute_err_iter = -1;
21590d69ab7SBarry Smith   mfctx->compute_err      = PETSC_FALSE;
21690d69ab7SBarry Smith   mfctx->jorge            = PETSC_FALSE;
2171aa26658SKarl Rupp 
2189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_mf_err", &mfctx->error_rel, NULL));
2199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_mf_umin", &mfctx->umin, NULL));
2209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_mf_jorge", &mfctx->jorge, NULL));
2219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_mf_compute_err", &mfctx->compute_err, NULL));
2229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_mf_freq_err", &mfctx->compute_err_freq, &flg));
223295f7fbbSLois Curfman McInnes   if (flg) {
224295f7fbbSLois Curfman McInnes     if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0;
2253b2a6d2fSBarry Smith     mfctx->compute_err = PETSC_TRUE;
226295f7fbbSLois Curfman McInnes   }
2273b2a6d2fSBarry Smith   if (mfctx->compute_err) mfctx->need_err = PETSC_TRUE;
228295f7fbbSLois Curfman McInnes   if (mfctx->jorge || mfctx->compute_err) {
2299566063dSJacob Faibussowitsch     PetscCall(SNESDiffParameterCreate_More(snes, x, &mfctx->data));
2309e5d0892SLisandro Dalcin   } else mfctx->data = NULL;
231d2dddef7SLois Curfman McInnes 
2329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasHelp(((PetscObject)snes)->options, &flg));
2339566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(p, "-", sizeof(p)));
2349566063dSJacob Faibussowitsch   if (((PetscObject)snes)->prefix) PetscCall(PetscStrlcat(p, ((PetscObject)snes)->prefix, sizeof(p)));
235d2dddef7SLois Curfman McInnes   if (flg) {
2369566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), " Matrix-free Options (via SNES):\n"));
2379566063dSJacob 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));
2389566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), "   %ssnes_mf_umin <umin>: see users manual (default %g)\n", p, (double)mfctx->umin));
2399566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), "   %ssnes_mf_jorge: use Jorge More's method\n", p));
2409566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), "   %ssnes_mf_compute_err: compute sqrt or relative error in function\n", p));
2419566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), "   %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)\n", p));
2429566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), "   %ssnes_mf_noise_file <file>: set file for printing noise info\n", p));
243d2dddef7SLois Curfman McInnes   }
2449566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &mfctx->w));
2459566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)x, &comm));
2469566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &n));
2479566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x, &nloc));
2489566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, J));
2499566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J, nloc, n, n, n));
2509566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J, MATSHELL));
2519566063dSJacob Faibussowitsch   PetscCall(MatShellSetContext(*J, mfctx));
2529566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_MULT, (void (*)(void))SNESMatrixFreeMult2_Private));
2539566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_DESTROY, (void (*)(void))SNESMatrixFreeDestroy2_Private));
2549566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_VIEW, (void (*)(void))SNESMatrixFreeView2_Private));
2559566063dSJacob Faibussowitsch   PetscCall(MatSetUp(*J));
256d2dddef7SLois Curfman McInnes 
25707a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
258d2dddef7SLois Curfman McInnes }
259d2dddef7SLois Curfman McInnes 
26090c26489SBarry Smith /*@C
261f6dfbefdSBarry Smith    MatSNESMFMoreSetParameters - Sets the parameters for the approximation of
262f6dfbefdSBarry Smith    matrix-vector products using finite differences, see  `MatCreateSNESMFMore()`
263d2dddef7SLois Curfman McInnes 
264d2dddef7SLois Curfman McInnes    Input Parameters:
265758f395bSBarry Smith +  mat - the matrix
266f6dfbefdSBarry Smith .  error_rel - relative error (should be set to the square root of the relative error in the function evaluations)
267d2dddef7SLois Curfman McInnes .  umin - minimum allowable u-value
268758f395bSBarry Smith -  h - differencing parameter
269d2dddef7SLois Curfman McInnes 
270f6dfbefdSBarry Smith    Options Database Keys:
271f6dfbefdSBarry Smith +  -snes_mf_err <error_rel> - see `MatCreateSNESMF()`
272f6dfbefdSBarry Smith .  -snes_mf_umin <umin> - see `MatCreateSNESMF()`
273f6dfbefdSBarry Smith .  -snes_mf_compute_err - compute the square root or relative error in function
274f6dfbefdSBarry Smith .  -snes_mf_freq_err <freq> - set the frequency to recompute the parameters
275f6dfbefdSBarry Smith -  -snes_mf_jorge - use the method of Jorge More
276f6dfbefdSBarry Smith 
27790c26489SBarry Smith    Level: advanced
27890c26489SBarry Smith 
279f6dfbefdSBarry Smith    Note:
280d2dddef7SLois Curfman McInnes    If the user sets the parameter h directly, then this value will be used
281f6dfbefdSBarry Smith    instead of the default computation as discussed in `MatCreateSNESMFMore()`
282d2dddef7SLois Curfman McInnes 
283f6dfbefdSBarry Smith .seealso: `MatCreateSNESMF()`, `MatCreateSNESMFMore()`
284d2dddef7SLois Curfman McInnes @*/
285f6dfbefdSBarry Smith PetscErrorCode MatSNESMFMoreSetParameters(Mat mat, PetscReal error, PetscReal umin, PetscReal h) {
286d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
287d2dddef7SLois Curfman McInnes 
28807a3eed9SLois Curfman McInnes   PetscFunctionBegin;
2899566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
290d2dddef7SLois Curfman McInnes   if (ctx) {
291d2dddef7SLois Curfman McInnes     if (error != PETSC_DEFAULT) ctx->error_rel = error;
292d2dddef7SLois Curfman McInnes     if (umin != PETSC_DEFAULT) ctx->umin = umin;
293d2dddef7SLois Curfman McInnes     if (h != PETSC_DEFAULT) {
294d2dddef7SLois Curfman McInnes       ctx->h      = h;
295a7cc72afSBarry Smith       ctx->need_h = PETSC_FALSE;
296d2dddef7SLois Curfman McInnes     }
297d2dddef7SLois Curfman McInnes   }
29807a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
299d2dddef7SLois Curfman McInnes }
300d2dddef7SLois Curfman McInnes 
3019371c9d4SSatish Balay PetscErrorCode SNESUnSetMatrixFreeParameter(SNES snes) {
302d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
303d2dddef7SLois Curfman McInnes   Mat            mat;
304d2dddef7SLois Curfman McInnes 
30507a3eed9SLois Curfman McInnes   PetscFunctionBegin;
3069566063dSJacob Faibussowitsch   PetscCall(SNESGetJacobian(snes, &mat, NULL, NULL, NULL));
3079566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
308a7cc72afSBarry Smith   if (ctx) ctx->need_h = PETSC_TRUE;
30907a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
310d2dddef7SLois Curfman McInnes }
311