1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
26370053bSBarry Smith /* matimpl.h is needed only for logging of matrix operation */
3af0996ceSBarry Smith #include <petsc/private/matimpl.h>
4d2dddef7SLois Curfman McInnes
525acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESUnSetMatrixFreeParameter(SNES);
625acbd8eSLisandro Dalcin
725acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterCreate_More(SNES, Vec, void **);
825acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterCompute_More(SNES, void *, Vec, Vec, PetscReal *, PetscReal *);
925acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterDestroy_More(void *);
10d2dddef7SLois Curfman McInnes
11d2dddef7SLois Curfman McInnes typedef struct { /* default context for matrix-free SNES */
12d2dddef7SLois Curfman McInnes SNES snes; /* SNES context */
13d2dddef7SLois Curfman McInnes Vec w; /* work vector */
1474637425SBarry Smith MatNullSpace sp; /* null space context */
15145595abSSatish Balay PetscReal error_rel; /* square root of relative error in computing function */
16145595abSSatish Balay PetscReal umin; /* minimum allowable u'a value relative to |u|_1 */
17f5af7f23SKarl Rupp PetscBool jorge; /* flag indicating use of Jorge's method for determining the differencing parameter */
18145595abSSatish Balay PetscReal h; /* differencing parameter */
19ace3abfcSBarry Smith PetscBool need_h; /* flag indicating whether we must compute h */
20ace3abfcSBarry Smith PetscBool need_err; /* flag indicating whether we must currently compute error_rel */
21ace3abfcSBarry Smith PetscBool compute_err; /* flag indicating whether we must ever compute error_rel */
22a7cc72afSBarry Smith PetscInt compute_err_iter; /* last iter where we've computer error_rel */
23a7cc72afSBarry Smith PetscInt compute_err_freq; /* frequency of computing error_rel */
24d2dddef7SLois Curfman McInnes void *data; /* implementation-specific data */
25d2dddef7SLois Curfman McInnes } MFCtx_Private;
26d2dddef7SLois Curfman McInnes
SNESMatrixFreeDestroy2_Private(Mat mat)2766976f2fSJacob Faibussowitsch static PetscErrorCode SNESMatrixFreeDestroy2_Private(Mat mat)
28d71ae5a4SJacob Faibussowitsch {
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));
373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
38d2dddef7SLois Curfman McInnes }
39d2dddef7SLois Curfman McInnes
40d2dddef7SLois Curfman McInnes /*
41d2dddef7SLois Curfman McInnes SNESMatrixFreeView2_Private - Views matrix-free parameters.
42d2dddef7SLois Curfman McInnes */
SNESMatrixFreeView2_Private(Mat J,PetscViewer viewer)4366976f2fSJacob Faibussowitsch static PetscErrorCode SNESMatrixFreeView2_Private(Mat J, PetscViewer viewer)
44d71ae5a4SJacob Faibussowitsch {
45d2dddef7SLois Curfman McInnes MFCtx_Private *ctx;
469f196a02SMartin Diehl PetscBool isascii;
47d2dddef7SLois Curfman McInnes
4807a3eed9SLois Curfman McInnes PetscFunctionBegin;
499566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J, &ctx));
509f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
519f196a02SMartin Diehl if (isascii) {
529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " SNES matrix-free approximation:\n"));
5348a46eb9SPierre Jolivet if (ctx->jorge) PetscCall(PetscViewerASCIIPrintf(viewer, " using Jorge's method of determining differencing parameter\n"));
549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " err=%g (relative error in function evaluation)\n", (double)ctx->error_rel));
559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " umin=%g (minimum iterate parameter)\n", (double)ctx->umin));
5648a46eb9SPierre Jolivet if (ctx->compute_err) PetscCall(PetscViewerASCIIPrintf(viewer, " freq_err=%" PetscInt_FMT " (frequency for computing err)\n", ctx->compute_err_freq));
57758f395bSBarry Smith }
583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
59d2dddef7SLois Curfman McInnes }
60d2dddef7SLois Curfman McInnes
61d2dddef7SLois Curfman McInnes /*
62d2dddef7SLois Curfman McInnes SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector
63d2dddef7SLois Curfman McInnes product, y = F'(u)*a:
64d2dddef7SLois Curfman McInnes y = (F(u + ha) - F(u)) /h,
65d2dddef7SLois Curfman McInnes where F = nonlinear function, as set by SNESSetFunction()
66d2dddef7SLois Curfman McInnes u = current iterate
67d2dddef7SLois Curfman McInnes h = difference interval
68d2dddef7SLois Curfman McInnes */
SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y)6966976f2fSJacob Faibussowitsch static PetscErrorCode SNESMatrixFreeMult2_Private(Mat mat, Vec a, Vec y)
70d71ae5a4SJacob Faibussowitsch {
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));
107f4f49eeaSPierre Jolivet 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));
1473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
148d2dddef7SLois Curfman McInnes }
149d2dddef7SLois Curfman McInnes
1505d83a8b1SBarry Smith /*@
151f6dfbefdSBarry Smith MatCreateSNESMFMore - Creates a matrix-free matrix
152f6dfbefdSBarry 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
153f6dfbefdSBarry Smith the Jacobian argument for the routine `SNESSetJacobian()`.
154d2dddef7SLois Curfman McInnes
155d2dddef7SLois Curfman McInnes Input Parameters:
156f6dfbefdSBarry Smith + snes - the `SNES` context
157f6dfbefdSBarry Smith - 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
162f6dfbefdSBarry Smith Options Database Keys:
163f6dfbefdSBarry Smith + -snes_mf_err <error_rel> - see `MatCreateSNESMF()`
164f6dfbefdSBarry Smith . -snes_mf_umin <umin> - see `MatCreateSNESMF()`
165f6dfbefdSBarry Smith . -snes_mf_compute_err - compute the square root or relative error in function
166f6dfbefdSBarry Smith . -snes_mf_freq_err <freq> - set the frequency to recompute the parameters
167f6dfbefdSBarry Smith - -snes_mf_jorge - use the method of Jorge More
168f6dfbefdSBarry Smith
16990c26489SBarry Smith Level: advanced
17090c26489SBarry Smith
171d2dddef7SLois Curfman McInnes Notes:
172f6dfbefdSBarry Smith This is an experimental approach, use `MatCreateSNESMF()`.
173f6dfbefdSBarry Smith
174d2dddef7SLois Curfman McInnes The matrix-free matrix context merely contains the function pointers
175d2dddef7SLois Curfman McInnes and work space for performing finite difference approximations of
176d2dddef7SLois Curfman McInnes Jacobian-vector products, J(u)*a, via
177d2dddef7SLois Curfman McInnes
17837fdd005SBarry Smith .vb
17937fdd005SBarry Smith J(u)*a = [J(u+h*a) - J(u)]/h,
18037fdd005SBarry Smith where by default
18137fdd005SBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1}
18237fdd005SBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise
18337fdd005SBarry Smith where
18437fdd005SBarry Smith error_rel = square root of relative error in function evaluation
18537fdd005SBarry Smith umin = minimum iterate parameter
18637fdd005SBarry Smith Alternatively, the differencing parameter, h, can be set using
18737fdd005SBarry Smith Jorge's nifty new strategy if one specifies the option
18837fdd005SBarry Smith -snes_mf_jorge
18937fdd005SBarry Smith .ve
190d2dddef7SLois Curfman McInnes
191f6dfbefdSBarry Smith The user can set these parameters via `MatMFFDSetFunctionError()`.
192d2dddef7SLois Curfman McInnes
193f6dfbefdSBarry Smith The user should call `MatDestroy()` when finished with the matrix-free
194d2dddef7SLois Curfman McInnes matrix context.
195d2dddef7SLois Curfman McInnes
196420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreateMF()`, `MatCreateMFFD()`, `MatDestroy()`, `MatMFFDSetFunctionError()`
197d2dddef7SLois Curfman McInnes @*/
MatCreateSNESMFMore(SNES snes,Vec x,Mat * J)198d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSNESMFMore(SNES snes, Vec x, Mat *J)
199d71ae5a4SJacob Faibussowitsch {
200d2dddef7SLois Curfman McInnes MPI_Comm comm;
201d2dddef7SLois Curfman McInnes MFCtx_Private *mfctx;
202a7cc72afSBarry Smith PetscInt n, nloc;
203ace3abfcSBarry Smith PetscBool flg;
204d2dddef7SLois Curfman McInnes char p[64];
205d2dddef7SLois Curfman McInnes
20607a3eed9SLois Curfman McInnes PetscFunctionBegin;
2074dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&mfctx));
2089e5d0892SLisandro Dalcin mfctx->sp = NULL;
209d2dddef7SLois Curfman McInnes mfctx->snes = snes;
21077d8c4bbSBarry Smith mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON;
211d2dddef7SLois Curfman McInnes mfctx->umin = 1.e-6;
212d2dddef7SLois Curfman McInnes mfctx->h = 0.0;
213a7cc72afSBarry Smith mfctx->need_h = PETSC_TRUE;
214a7cc72afSBarry Smith mfctx->need_err = PETSC_FALSE;
215a7cc72afSBarry Smith mfctx->compute_err = PETSC_FALSE;
216295f7fbbSLois Curfman McInnes mfctx->compute_err_freq = 0;
217295f7fbbSLois Curfman McInnes mfctx->compute_err_iter = -1;
21890d69ab7SBarry Smith mfctx->compute_err = PETSC_FALSE;
21990d69ab7SBarry Smith mfctx->jorge = PETSC_FALSE;
2201aa26658SKarl Rupp
2219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_mf_err", &mfctx->error_rel, NULL));
2229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_mf_umin", &mfctx->umin, NULL));
2239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_mf_jorge", &mfctx->jorge, NULL));
2249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_mf_compute_err", &mfctx->compute_err, NULL));
2259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_mf_freq_err", &mfctx->compute_err_freq, &flg));
226295f7fbbSLois Curfman McInnes if (flg) {
227295f7fbbSLois Curfman McInnes if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0;
2283b2a6d2fSBarry Smith mfctx->compute_err = PETSC_TRUE;
229295f7fbbSLois Curfman McInnes }
2303b2a6d2fSBarry Smith if (mfctx->compute_err) mfctx->need_err = PETSC_TRUE;
231ac530a7eSPierre Jolivet if (mfctx->jorge || mfctx->compute_err) PetscCall(SNESDiffParameterCreate_More(snes, x, &mfctx->data));
232ac530a7eSPierre Jolivet else mfctx->data = NULL;
233d2dddef7SLois Curfman McInnes
2349566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasHelp(((PetscObject)snes)->options, &flg));
2359566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(p, "-", sizeof(p)));
2369566063dSJacob Faibussowitsch if (((PetscObject)snes)->prefix) PetscCall(PetscStrlcat(p, ((PetscObject)snes)->prefix, sizeof(p)));
237d2dddef7SLois Curfman McInnes if (flg) {
2389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), " Matrix-free Options (via SNES):\n"));
2399566063dSJacob 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));
2409566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), " %ssnes_mf_umin <umin>: see users manual (default %g)\n", p, (double)mfctx->umin));
2419566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), " %ssnes_mf_jorge: use Jorge More's method\n", p));
2429566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), " %ssnes_mf_compute_err: compute sqrt or relative error in function\n", p));
2439566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), " %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)\n", p));
2449566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), " %ssnes_mf_noise_file <file>: set file for printing noise info\n", p));
245d2dddef7SLois Curfman McInnes }
2469566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &mfctx->w));
2479566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)x, &comm));
2489566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &n));
2499566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x, &nloc));
2509566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, J));
2519566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J, nloc, n, n, n));
2529566063dSJacob Faibussowitsch PetscCall(MatSetType(*J, MATSHELL));
2539566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(*J, mfctx));
254*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_MULT, (PetscErrorCodeFn *)SNESMatrixFreeMult2_Private));
255*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_DESTROY, (PetscErrorCodeFn *)SNESMatrixFreeDestroy2_Private));
256*57d50842SBarry Smith PetscCall(MatShellSetOperation(*J, MATOP_VIEW, (PetscErrorCodeFn *)SNESMatrixFreeView2_Private));
2579566063dSJacob Faibussowitsch PetscCall(MatSetUp(*J));
2583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
259d2dddef7SLois Curfman McInnes }
260d2dddef7SLois Curfman McInnes
261cc4c1da9SBarry Smith /*@
262f6dfbefdSBarry Smith MatSNESMFMoreSetParameters - Sets the parameters for the approximation of
263f6dfbefdSBarry Smith matrix-vector products using finite differences, see `MatCreateSNESMFMore()`
264d2dddef7SLois Curfman McInnes
265d2dddef7SLois Curfman McInnes Input Parameters:
266758f395bSBarry Smith + mat - the matrix
267e4094ef1SJacob Faibussowitsch . error - relative error (should be set to the square root of the relative error in the function evaluations)
268d2dddef7SLois Curfman McInnes . umin - minimum allowable u-value
269758f395bSBarry Smith - h - differencing parameter
270d2dddef7SLois Curfman McInnes
271f6dfbefdSBarry Smith Options Database Keys:
272f6dfbefdSBarry Smith + -snes_mf_err <error_rel> - see `MatCreateSNESMF()`
273f6dfbefdSBarry Smith . -snes_mf_umin <umin> - see `MatCreateSNESMF()`
274f6dfbefdSBarry Smith . -snes_mf_compute_err - compute the square root or relative error in function
275f6dfbefdSBarry Smith . -snes_mf_freq_err <freq> - set the frequency to recompute the parameters
276f6dfbefdSBarry Smith - -snes_mf_jorge - use the method of Jorge More
277f6dfbefdSBarry Smith
27890c26489SBarry Smith Level: advanced
27990c26489SBarry Smith
280f6dfbefdSBarry Smith Note:
281420bcc1bSBarry Smith If the user sets the parameter `h` directly, then this value will be used
282f6dfbefdSBarry Smith instead of the default computation as discussed in `MatCreateSNESMFMore()`
283d2dddef7SLois Curfman McInnes
284420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `MatCreateSNESMF()`, `MatCreateSNESMFMore()`
285d2dddef7SLois Curfman McInnes @*/
MatSNESMFMoreSetParameters(Mat mat,PetscReal error,PetscReal umin,PetscReal h)286d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSNESMFMoreSetParameters(Mat mat, PetscReal error, PetscReal umin, PetscReal h)
287d71ae5a4SJacob Faibussowitsch {
288d2dddef7SLois Curfman McInnes MFCtx_Private *ctx;
289d2dddef7SLois Curfman McInnes
29007a3eed9SLois Curfman McInnes PetscFunctionBegin;
2919566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx));
292d2dddef7SLois Curfman McInnes if (ctx) {
29313bcc0bdSJacob Faibussowitsch if (error != (PetscReal)PETSC_DEFAULT) ctx->error_rel = error;
29413bcc0bdSJacob Faibussowitsch if (umin != (PetscReal)PETSC_DEFAULT) ctx->umin = umin;
29513bcc0bdSJacob Faibussowitsch if (h != (PetscReal)PETSC_DEFAULT) {
296d2dddef7SLois Curfman McInnes ctx->h = h;
297a7cc72afSBarry Smith ctx->need_h = PETSC_FALSE;
298d2dddef7SLois Curfman McInnes }
299d2dddef7SLois Curfman McInnes }
3003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
301d2dddef7SLois Curfman McInnes }
302d2dddef7SLois Curfman McInnes
SNESUnSetMatrixFreeParameter(SNES snes)303d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESUnSetMatrixFreeParameter(SNES snes)
304d71ae5a4SJacob Faibussowitsch {
305d2dddef7SLois Curfman McInnes MFCtx_Private *ctx;
306d2dddef7SLois Curfman McInnes Mat mat;
307d2dddef7SLois Curfman McInnes
30807a3eed9SLois Curfman McInnes PetscFunctionBegin;
3099566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &mat, NULL, NULL, NULL));
3109566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx));
311a7cc72afSBarry Smith if (ctx) ctx->need_h = PETSC_TRUE;
3123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
313d2dddef7SLois Curfman McInnes }
314