xref: /petsc/src/mat/impls/mffd/wp.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1e884886fSBarry Smith 
2e884886fSBarry Smith /*MC
311a5261eSBarry Smith      MATMFFD_WP - Implements an approach for computing the differencing parameter
411a5261eSBarry Smith         h used with the finite difference based matrix-free Jacobian.
5e884886fSBarry Smith 
6e884886fSBarry Smith       h = error_rel * sqrt(1 + ||U||) / ||a||
7e884886fSBarry Smith 
811a5261eSBarry Smith    Options Database Key:
911a5261eSBarry Smith .   -mat_mffd_compute_normu -Compute the norm of u every time see `MatMFFDWPSetComputeNormU()`
10e884886fSBarry Smith 
11e884886fSBarry Smith    Level: intermediate
12e884886fSBarry Smith 
1395452b02SPatrick Sanan    Notes:
1411a5261eSBarry Smith    || U || does not change between linear iterations so is reused
1511a5261eSBarry Smith 
1611a5261eSBarry Smith    In `KSPGMRES` || a || == 1 and so does not need to ever be computed except at restart
1711a5261eSBarry Smith     when it is recomputed.  Thus equires no global collectives when used with `KSPGMRES`
18e884886fSBarry Smith 
19e884886fSBarry Smith    Formula used:
20e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
21e884886fSBarry Smith 
2211a5261eSBarry Smith    Reference:
2311a5261eSBarry Smith .  * -  M. Pernice and H. F. Walker, "NITSOL: A Newton Iterative
2411a5261eSBarry Smith       Solver for Nonlinear Systems", SIAM J. Sci. Stat. Comput.", 1998,
2511a5261eSBarry Smith       vol 19, pp. 302--318.
26e884886fSBarry Smith 
2711a5261eSBarry Smith .seealso: `MATMFFD`, `MATMFFD_DS`, `MatCreateMFFD()`, `MatCreateSNESMF()`, `MATMFFD_DS`
28e884886fSBarry Smith M*/
29e884886fSBarry Smith 
30e884886fSBarry Smith /*
31e884886fSBarry Smith     This include file defines the data structure  MatMFFD that
32e884886fSBarry Smith    includes information about the computation of h. It is shared by
33e884886fSBarry Smith    all implementations that people provide.
34e884886fSBarry Smith 
35e884886fSBarry Smith    See snesmfjdef.c for  a full set of comments on the routines below.
36e884886fSBarry Smith */
37af0996ceSBarry Smith #include <petsc/private/matimpl.h>
38c6db04a5SJed Brown #include <../src/mat/impls/mffd/mffdimpl.h> /*I  "petscmat.h"   I*/
39e884886fSBarry Smith 
40e884886fSBarry Smith typedef struct {
41e884886fSBarry Smith   PetscReal normUfact; /* previous sqrt(1.0 + || U ||) */
42ace3abfcSBarry Smith   PetscBool computenormU;
43e884886fSBarry Smith } MatMFFD_WP;
44e884886fSBarry Smith 
45e884886fSBarry Smith /*
4611a5261eSBarry Smith      MatMFFDCompute_WP - code for
47e884886fSBarry Smith    computing h with matrix-free finite differences.
48e884886fSBarry Smith 
49e884886fSBarry Smith   Input Parameters:
50e884886fSBarry Smith +   ctx - the matrix free context
51e884886fSBarry Smith .   U - the location at which you want the Jacobian
52e884886fSBarry Smith -   a - the direction you want the derivative
53e884886fSBarry Smith 
54e884886fSBarry Smith   Output Parameter:
55e884886fSBarry Smith .   h - the scale computed
56e884886fSBarry Smith 
57e884886fSBarry Smith */
58*d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDCompute_WP(MatMFFD ctx, Vec U, Vec a, PetscScalar *h, PetscBool *zeroa)
59*d71ae5a4SJacob Faibussowitsch {
60e884886fSBarry Smith   MatMFFD_WP *hctx = (MatMFFD_WP *)ctx->hctx;
6151a79602SBarry Smith   PetscReal   normU, norma;
62e884886fSBarry Smith 
63e884886fSBarry Smith   PetscFunctionBegin;
64e884886fSBarry Smith   if (!(ctx->count % ctx->recomputeperiod)) {
6551a79602SBarry Smith     if (hctx->computenormU || !ctx->ncurrenth) {
669566063dSJacob Faibussowitsch       PetscCall(VecNorm(U, NORM_2, &normU));
678f1a2a5eSBarry Smith       hctx->normUfact = PetscSqrtReal(1.0 + normU);
68e884886fSBarry Smith     }
699566063dSJacob Faibussowitsch     PetscCall(VecNorm(a, NORM_2, &norma));
70e884886fSBarry Smith     if (norma == 0.0) {
71e884886fSBarry Smith       *zeroa = PETSC_TRUE;
72e884886fSBarry Smith       PetscFunctionReturn(0);
73e884886fSBarry Smith     }
74e884886fSBarry Smith     *zeroa = PETSC_FALSE;
75e884886fSBarry Smith     *h     = ctx->error_rel * hctx->normUfact / norma;
76e884886fSBarry Smith   } else {
77e884886fSBarry Smith     *h = ctx->currenth;
78e884886fSBarry Smith   }
79e884886fSBarry Smith   PetscFunctionReturn(0);
80e884886fSBarry Smith }
81e884886fSBarry Smith 
82e884886fSBarry Smith /*
83e884886fSBarry Smith    MatMFFDView_WP - Prints information about this particular
84e884886fSBarry Smith      method for computing h. Note that this does not print the general
85e884886fSBarry Smith      information about the matrix free, that is printed by the calling
86e884886fSBarry Smith      routine.
87e884886fSBarry Smith 
88e884886fSBarry Smith   Input Parameters:
89e884886fSBarry Smith +   ctx - the matrix free context
90e884886fSBarry Smith -   viewer - the PETSc viewer
91e884886fSBarry Smith 
92e884886fSBarry Smith */
93*d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDView_WP(MatMFFD ctx, PetscViewer viewer)
94*d71ae5a4SJacob Faibussowitsch {
95e884886fSBarry Smith   MatMFFD_WP *hctx = (MatMFFD_WP *)ctx->hctx;
96ace3abfcSBarry Smith   PetscBool   iascii;
97e884886fSBarry Smith 
98e884886fSBarry Smith   PetscFunctionBegin;
999566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
100e884886fSBarry Smith   if (iascii) {
1012205254eSKarl Rupp     if (hctx->computenormU) {
1029566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Computes normU\n"));
1032205254eSKarl Rupp     } else {
1049566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Does not compute normU\n"));
1052205254eSKarl Rupp     }
10611aeaf0aSBarry Smith   }
107e884886fSBarry Smith   PetscFunctionReturn(0);
108e884886fSBarry Smith }
109e884886fSBarry Smith 
110e884886fSBarry Smith /*
111e884886fSBarry Smith    MatMFFDSetFromOptions_WP - Looks in the options database for
112e884886fSBarry Smith      any options appropriate for this method
113e884886fSBarry Smith 
114e884886fSBarry Smith   Input Parameter:
115e884886fSBarry Smith .  ctx - the matrix free context
116e884886fSBarry Smith 
117e884886fSBarry Smith */
118*d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFromOptions_WP(MatMFFD ctx, PetscOptionItems *PetscOptionsObject)
119*d71ae5a4SJacob Faibussowitsch {
120e884886fSBarry Smith   MatMFFD_WP *hctx = (MatMFFD_WP *)ctx->hctx;
121e884886fSBarry Smith 
122e884886fSBarry Smith   PetscFunctionBegin;
123d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Walker-Pernice options");
1249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_mffd_compute_normu", "Compute the norm of u", "MatMFFDWPSetComputeNormU", hctx->computenormU, &hctx->computenormU, NULL));
125d0609cedSBarry Smith   PetscOptionsHeadEnd();
126e884886fSBarry Smith   PetscFunctionReturn(0);
127e884886fSBarry Smith }
128e884886fSBarry Smith 
129e884886fSBarry Smith /*
130e884886fSBarry Smith    MatMFFDDestroy_WP - Frees the space allocated by
1311d0fab5eSBarry Smith        MatCreateMFFD_WP().
132e884886fSBarry Smith 
133e884886fSBarry Smith   Input Parameter:
134e884886fSBarry Smith .  ctx - the matrix free context
135e884886fSBarry Smith 
13695452b02SPatrick Sanan    Notes:
13795452b02SPatrick Sanan     does not free the ctx, that is handled by the calling routine
138e884886fSBarry Smith 
139e884886fSBarry Smith */
140*d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDDestroy_WP(MatMFFD ctx)
141*d71ae5a4SJacob Faibussowitsch {
142e884886fSBarry Smith   PetscFunctionBegin;
1432e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ctx->mat, "MatMFFDWPSetComputeNormU_C", NULL));
1449566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx->hctx));
145e884886fSBarry Smith   PetscFunctionReturn(0);
146e884886fSBarry Smith }
147e884886fSBarry Smith 
148*d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDWPSetComputeNormU_P(Mat mat, PetscBool flag)
149*d71ae5a4SJacob Faibussowitsch {
150e884886fSBarry Smith   MatMFFD     ctx  = (MatMFFD)mat->data;
151e884886fSBarry Smith   MatMFFD_WP *hctx = (MatMFFD_WP *)ctx->hctx;
152e884886fSBarry Smith 
153e884886fSBarry Smith   PetscFunctionBegin;
154e884886fSBarry Smith   hctx->computenormU = flag;
155e884886fSBarry Smith   PetscFunctionReturn(0);
156e884886fSBarry Smith }
157e884886fSBarry Smith 
158e884886fSBarry Smith /*@
15911a5261eSBarry Smith     MatMFFDWPSetComputeNormU - Sets whether it computes the ||U|| used by the Walker-Pernice
160e884886fSBarry Smith              PETSc routine for computing h. With any Krylov solver this need only
161e884886fSBarry Smith              be computed during the first iteration and kept for later.
162e884886fSBarry Smith 
163e884886fSBarry Smith   Input Parameters:
16411a5261eSBarry Smith +   A - the `MATMFFD` matrix
16511a5261eSBarry Smith -   flag - `PETSC_TRUE` causes it to compute ||U||, `PETSC_FALSE` uses the previous value
166e884886fSBarry Smith 
167e884886fSBarry Smith   Options Database Key:
168e884886fSBarry Smith .   -mat_mffd_compute_normu <true,false> - true by default, false can save calculations but you
169e884886fSBarry Smith               must be sure that ||U|| has not changed in the mean time.
170e884886fSBarry Smith 
171e884886fSBarry Smith   Level: advanced
172e884886fSBarry Smith 
17311a5261eSBarry Smith   Note:
17411a5261eSBarry Smith    See the manual page for `MATMFFD_WP` for a complete description of the
175e884886fSBarry Smith    algorithm used to compute h.
176e884886fSBarry Smith 
17711a5261eSBarry Smith .seealso: `MATMFFD_WP`, `MATMFFD`, `MatMFFDSetFunctionError()`, `MatCreateSNESMF()`
178e884886fSBarry Smith @*/
179*d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDWPSetComputeNormU(Mat A, PetscBool flag)
180*d71ae5a4SJacob Faibussowitsch {
181e884886fSBarry Smith   PetscFunctionBegin;
1820700a824SBarry Smith   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
183cac4c232SBarry Smith   PetscTryMethod(A, "MatMFFDWPSetComputeNormU_C", (Mat, PetscBool), (A, flag));
184e884886fSBarry Smith   PetscFunctionReturn(0);
185e884886fSBarry Smith }
186e884886fSBarry Smith 
187e884886fSBarry Smith /*
1881d0fab5eSBarry Smith      MatCreateMFFD_WP - Standard PETSc code for
189e884886fSBarry Smith    computing h with matrix-free finite differences.
190e884886fSBarry Smith 
191e884886fSBarry Smith    Input Parameter:
1921d0fab5eSBarry Smith .  ctx - the matrix free context created by MatCreateMFFD()
193e884886fSBarry Smith 
194e884886fSBarry Smith */
195*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreateMFFD_WP(MatMFFD ctx)
196*d71ae5a4SJacob Faibussowitsch {
197e884886fSBarry Smith   MatMFFD_WP *hctx;
198e884886fSBarry Smith 
199e884886fSBarry Smith   PetscFunctionBegin;
200e884886fSBarry Smith   /* allocate my own private data structure */
2014dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&hctx));
202e884886fSBarry Smith   ctx->hctx          = (void *)hctx;
203e884886fSBarry Smith   hctx->computenormU = PETSC_FALSE;
204e884886fSBarry Smith 
205e884886fSBarry Smith   /* set the functions I am providing */
206e884886fSBarry Smith   ctx->ops->compute        = MatMFFDCompute_WP;
207e884886fSBarry Smith   ctx->ops->destroy        = MatMFFDDestroy_WP;
208e884886fSBarry Smith   ctx->ops->view           = MatMFFDView_WP;
209e884886fSBarry Smith   ctx->ops->setfromoptions = MatMFFDSetFromOptions_WP;
210e884886fSBarry Smith 
2119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ctx->mat, "MatMFFDWPSetComputeNormU_C", MatMFFDWPSetComputeNormU_P));
212e884886fSBarry Smith   PetscFunctionReturn(0);
213e884886fSBarry Smith }
214