xref: /petsc/src/snes/interface/noise/snesmfj2.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 PETSC_INTERN PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES, Vec, Mat *);
825acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDefaultMatrixFreeSetParameters2(Mat, PetscReal, PetscReal, PetscReal);
925acbd8eSLisandro Dalcin 
1025acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterCreate_More(SNES, Vec, void **);
1125acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterCompute_More(SNES, void *, Vec, Vec, PetscReal *, PetscReal *);
1225acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterDestroy_More(void *);
13d2dddef7SLois Curfman McInnes 
14d2dddef7SLois Curfman McInnes typedef struct {                 /* default context for matrix-free SNES */
15d2dddef7SLois Curfman McInnes   SNES         snes;             /* SNES context */
16d2dddef7SLois Curfman McInnes   Vec          w;                /* work vector */
1774637425SBarry Smith   MatNullSpace sp;               /* null space context */
18145595abSSatish Balay   PetscReal    error_rel;        /* square root of relative error in computing function */
19145595abSSatish Balay   PetscReal    umin;             /* minimum allowable u'a value relative to |u|_1 */
20f5af7f23SKarl Rupp   PetscBool    jorge;            /* flag indicating use of Jorge's method for determining the differencing parameter */
21145595abSSatish Balay   PetscReal    h;                /* differencing parameter */
22ace3abfcSBarry Smith   PetscBool    need_h;           /* flag indicating whether we must compute h */
23ace3abfcSBarry Smith   PetscBool    need_err;         /* flag indicating whether we must currently compute error_rel */
24ace3abfcSBarry Smith   PetscBool    compute_err;      /* flag indicating whether we must ever compute error_rel */
25a7cc72afSBarry Smith   PetscInt     compute_err_iter; /* last iter where we've computer error_rel */
26a7cc72afSBarry Smith   PetscInt     compute_err_freq; /* frequency of computing error_rel */
27d2dddef7SLois Curfman McInnes   void        *data;             /* implementation-specific data */
28d2dddef7SLois Curfman McInnes } MFCtx_Private;
29d2dddef7SLois Curfman McInnes 
309371c9d4SSatish Balay PetscErrorCode SNESMatrixFreeDestroy2_Private(Mat mat) {
31d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
32d2dddef7SLois Curfman McInnes 
3307a3eed9SLois Curfman McInnes   PetscFunctionBegin;
349566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->w));
369566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceDestroy(&ctx->sp));
379566063dSJacob Faibussowitsch   if (ctx->jorge || ctx->compute_err) PetscCall(SNESDiffParameterDestroy_More(ctx->data));
389566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
3907a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
40d2dddef7SLois Curfman McInnes }
41d2dddef7SLois Curfman McInnes 
42d2dddef7SLois Curfman McInnes /*
43d2dddef7SLois Curfman McInnes    SNESMatrixFreeView2_Private - Views matrix-free parameters.
44d2dddef7SLois Curfman McInnes  */
459371c9d4SSatish Balay PetscErrorCode SNESMatrixFreeView2_Private(Mat J, PetscViewer viewer) {
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"));
54*48a46eb9SPierre 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));
57*48a46eb9SPierre Jolivet     if (ctx->compute_err) PetscCall(PetscViewerASCIIPrintf(viewer, "    freq_err=%" PetscInt_FMT " (frequency for computing err)\n", ctx->compute_err_freq));
58758f395bSBarry Smith   }
5907a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
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 */
709371c9d4SSatish Balay PetscErrorCode SNESMatrixFreeMult2_Private(Mat mat, Vec a, Vec y) {
71d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
72d2dddef7SLois Curfman McInnes   SNES           snes;
73145595abSSatish Balay   PetscReal      h, norm, sum, umin, noise;
7479f36460SBarry Smith   PetscScalar    hs, dot;
75d2dddef7SLois Curfman McInnes   Vec            w, U, F;
76e6ed2b1fSLois Curfman McInnes   MPI_Comm       comm;
7777431f27SBarry Smith   PetscInt       iter;
785f80ce2aSJacob Faibussowitsch   PetscErrorCode (*eval_fct)(SNES, Vec, Vec);
79d2dddef7SLois Curfman McInnes 
8007a3eed9SLois Curfman McInnes   PetscFunctionBegin;
81d2dddef7SLois Curfman McInnes   /* We log matrix-free matrix-vector products separately, so that we can
82d2dddef7SLois Curfman McInnes      separate the performance monitoring from the cases that use conventional
83d2dddef7SLois Curfman McInnes      storage.  We may eventually modify event logging to associate events
84d2dddef7SLois Curfman McInnes      with particular objects, hence alleviating the more general problem. */
859566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MATMFFD_Mult, a, y, 0, 0));
86d2dddef7SLois Curfman McInnes 
879566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
889566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
89e6ed2b1fSLois Curfman McInnes   snes = ctx->snes;
90e6ed2b1fSLois Curfman McInnes   w    = ctx->w;
91e6ed2b1fSLois Curfman McInnes   umin = ctx->umin;
92e6ed2b1fSLois Curfman McInnes 
939566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &U));
94d2dddef7SLois Curfman McInnes   eval_fct = SNESComputeFunction;
959566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &F, NULL, NULL));
96295f7fbbSLois Curfman McInnes 
97d2dddef7SLois Curfman McInnes   /* Determine a "good" step size, h */
98d2dddef7SLois Curfman McInnes   if (ctx->need_h) {
99d2dddef7SLois Curfman McInnes     /* Use Jorge's method to compute h */
100d2dddef7SLois Curfman McInnes     if (ctx->jorge) {
1019566063dSJacob Faibussowitsch       PetscCall(SNESDiffParameterCompute_More(snes, ctx->data, U, a, &noise, &h));
102d2dddef7SLois Curfman McInnes 
103d2dddef7SLois Curfman McInnes       /* Use the Brown/Saad method to compute h */
104d2dddef7SLois Curfman McInnes     } else {
105295f7fbbSLois Curfman McInnes       /* Compute error if desired */
1069566063dSJacob Faibussowitsch       PetscCall(SNESGetIterationNumber(snes, &iter));
107145595abSSatish Balay       if ((ctx->need_err) || ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter - 1) % ctx->compute_err_freq)))) {
108d2dddef7SLois Curfman McInnes         /* Use Jorge's method to compute noise */
1099566063dSJacob Faibussowitsch         PetscCall(SNESDiffParameterCompute_More(snes, ctx->data, U, a, &noise, &h));
1101aa26658SKarl Rupp 
111369cc0aeSBarry Smith         ctx->error_rel = PetscSqrtReal(noise);
1121aa26658SKarl Rupp 
1139566063dSJacob 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));
1141aa26658SKarl Rupp 
115295f7fbbSLois Curfman McInnes         ctx->compute_err_iter = iter;
1163b2a6d2fSBarry Smith         ctx->need_err         = PETSC_FALSE;
117d2dddef7SLois Curfman McInnes       }
118e6ed2b1fSLois Curfman McInnes 
1199566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(U, a, &dot));
1209566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(a, NORM_1, &sum));
1219566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(a, NORM_2, &norm));
1229566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(U, a, &dot));
1239566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(a, NORM_1, &sum));
1249566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(a, NORM_2, &norm));
125e6ed2b1fSLois Curfman McInnes 
126d2dddef7SLois Curfman McInnes       /* Safeguard for step sizes too small */
1271aa26658SKarl Rupp       if (sum == 0.0) {
1281aa26658SKarl Rupp         dot  = 1.0;
1291aa26658SKarl Rupp         norm = 1.0;
1301aa26658SKarl Rupp       } else if (PetscAbsScalar(dot) < umin * sum && PetscRealPart(dot) >= 0.0) dot = umin * sum;
131145595abSSatish Balay       else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin * sum) dot = -umin * sum;
132145595abSSatish Balay       h = PetscRealPart(ctx->error_rel * dot / (norm * norm));
133d2dddef7SLois Curfman McInnes     }
1341aa26658SKarl Rupp   } else h = ctx->h;
1351aa26658SKarl Rupp 
1369566063dSJacob Faibussowitsch   if (!ctx->jorge || !ctx->need_h) PetscCall(PetscInfo(snes, "h = %g\n", (double)h));
137d2dddef7SLois Curfman McInnes 
138d2dddef7SLois Curfman McInnes   /* Evaluate function at F(u + ha) */
139db0b4b35SLois Curfman McInnes   hs = h;
1409566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(w, hs, a, U));
1419566063dSJacob Faibussowitsch   PetscCall(eval_fct(snes, w, y));
1429566063dSJacob Faibussowitsch   PetscCall(VecAXPY(y, -1.0, F));
1439566063dSJacob Faibussowitsch   PetscCall(VecScale(y, 1.0 / hs));
1449566063dSJacob Faibussowitsch   if (mat->nullsp) PetscCall(MatNullSpaceRemove(mat->nullsp, y));
145d2dddef7SLois Curfman McInnes 
1469566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MATMFFD_Mult, a, y, 0, 0));
14707a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
148d2dddef7SLois Curfman McInnes }
149d2dddef7SLois Curfman McInnes 
150d2dddef7SLois Curfman McInnes /*@C
15163dd3a1aSKris Buschelman    SNESMatrixFreeCreate2 - Creates a matrix-free matrix
152d2dddef7SLois Curfman McInnes    context for use with a SNES solver.  This matrix can be used as
153d2dddef7SLois Curfman McInnes    the Jacobian argument for the routine SNESSetJacobian().
154d2dddef7SLois Curfman McInnes 
155d2dddef7SLois Curfman McInnes    Input Parameters:
156a2b725a8SWilliam Gropp +  snes - the SNES context
157a2b725a8SWilliam Gropp -  x - vector where SNES solution is to be stored.
158d2dddef7SLois Curfman McInnes 
159d2dddef7SLois Curfman McInnes    Output Parameter:
160d2dddef7SLois Curfman McInnes .  J - the matrix-free matrix
161d2dddef7SLois Curfman McInnes 
16290c26489SBarry Smith    Level: advanced
16390c26489SBarry Smith 
164d2dddef7SLois Curfman McInnes    Notes:
165d2dddef7SLois Curfman McInnes    The matrix-free matrix context merely contains the function pointers
166d2dddef7SLois Curfman McInnes    and work space for performing finite difference approximations of
167d2dddef7SLois Curfman McInnes    Jacobian-vector products, J(u)*a, via
168d2dddef7SLois Curfman McInnes 
169d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h,
170d2dddef7SLois Curfman McInnes $   where by default
171d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
172d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
173d2dddef7SLois Curfman McInnes $   where
174d2dddef7SLois Curfman McInnes $        error_rel = square root of relative error in
175d2dddef7SLois Curfman McInnes $                    function evaluation
176d2dddef7SLois Curfman McInnes $        umin = minimum iterate parameter
177d2dddef7SLois Curfman McInnes $   Alternatively, the differencing parameter, h, can be set using
178e6ed2b1fSLois Curfman McInnes $   Jorge's nifty new strategy if one specifies the option
179e6ed2b1fSLois Curfman McInnes $          -snes_mf_jorge
180d2dddef7SLois Curfman McInnes 
1813ec795f1SBarry Smith    The user can set these parameters via MatMFFDSetFunctionError().
182a7f22e61SSatish Balay    See Users-Manual: ch_snes for details
183d2dddef7SLois Curfman McInnes 
184d2dddef7SLois Curfman McInnes    The user should call MatDestroy() when finished with the matrix-free
185d2dddef7SLois Curfman McInnes    matrix context.
186d2dddef7SLois Curfman McInnes 
187d2dddef7SLois Curfman McInnes    Options Database Keys:
188d2dddef7SLois Curfman McInnes $  -snes_mf_err <error_rel>
1891e2ae407SBarry Smith $  -snes_mf_umin <umin>
190295f7fbbSLois Curfman McInnes $  -snes_mf_compute_err
191295f7fbbSLois Curfman McInnes $  -snes_mf_freq_err <freq>
192e6ed2b1fSLois Curfman McInnes $  -snes_mf_jorge
193d2dddef7SLois Curfman McInnes 
194db781477SPatrick Sanan .seealso: `MatDestroy()`, `MatMFFDSetFunctionError()`
195d2dddef7SLois Curfman McInnes @*/
1969371c9d4SSatish Balay PetscErrorCode SNESDefaultMatrixFreeCreate2(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;
2049566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes, &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 
2579566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)*J, (PetscObject)mfctx->w));
2589566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)snes, (PetscObject)*J));
25907a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
260d2dddef7SLois Curfman McInnes }
261d2dddef7SLois Curfman McInnes 
26290c26489SBarry Smith /*@C
263758f395bSBarry Smith    SNESDefaultMatrixFreeSetParameters2 - Sets the parameters for the approximation of
264d2dddef7SLois Curfman McInnes    matrix-vector products using finite differences.
265d2dddef7SLois Curfman McInnes 
266d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h where
267d2dddef7SLois Curfman McInnes 
268d2dddef7SLois Curfman McInnes    either the user sets h directly here, or this parameter is computed via
269d2dddef7SLois Curfman McInnes 
270d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
271d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
272d2dddef7SLois Curfman McInnes $
273d2dddef7SLois Curfman McInnes 
274d2dddef7SLois Curfman McInnes    Input Parameters:
275758f395bSBarry Smith +  mat - the matrix
276d2dddef7SLois Curfman McInnes .  error_rel - relative error (should be set to the square root of
277d2dddef7SLois Curfman McInnes      the relative error in the function evaluations)
278d2dddef7SLois Curfman McInnes .  umin - minimum allowable u-value
279758f395bSBarry Smith -  h - differencing parameter
280d2dddef7SLois Curfman McInnes 
28190c26489SBarry Smith    Level: advanced
28290c26489SBarry Smith 
283d2dddef7SLois Curfman McInnes    Notes:
284d2dddef7SLois Curfman McInnes    If the user sets the parameter h directly, then this value will be used
285d2dddef7SLois Curfman McInnes    instead of the default computation indicated above.
286d2dddef7SLois Curfman McInnes 
287db781477SPatrick Sanan .seealso: `MatCreateSNESMF()`
288d2dddef7SLois Curfman McInnes @*/
2899371c9d4SSatish Balay PetscErrorCode SNESDefaultMatrixFreeSetParameters2(Mat mat, PetscReal error, PetscReal umin, PetscReal h) {
290d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
291d2dddef7SLois Curfman McInnes 
29207a3eed9SLois Curfman McInnes   PetscFunctionBegin;
2939566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
294d2dddef7SLois Curfman McInnes   if (ctx) {
295d2dddef7SLois Curfman McInnes     if (error != PETSC_DEFAULT) ctx->error_rel = error;
296d2dddef7SLois Curfman McInnes     if (umin != PETSC_DEFAULT) ctx->umin = umin;
297d2dddef7SLois Curfman McInnes     if (h != PETSC_DEFAULT) {
298d2dddef7SLois Curfman McInnes       ctx->h      = h;
299a7cc72afSBarry Smith       ctx->need_h = PETSC_FALSE;
300d2dddef7SLois Curfman McInnes     }
301d2dddef7SLois Curfman McInnes   }
30207a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
303d2dddef7SLois Curfman McInnes }
304d2dddef7SLois Curfman McInnes 
3059371c9d4SSatish Balay PetscErrorCode SNESUnSetMatrixFreeParameter(SNES snes) {
306d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
307d2dddef7SLois Curfman McInnes   Mat            mat;
308d2dddef7SLois Curfman McInnes 
30907a3eed9SLois Curfman McInnes   PetscFunctionBegin;
3109566063dSJacob Faibussowitsch   PetscCall(SNESGetJacobian(snes, &mat, NULL, NULL, NULL));
3119566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
312a7cc72afSBarry Smith   if (ctx) ctx->need_h = PETSC_TRUE;
31307a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
314d2dddef7SLois Curfman McInnes }
315