xref: /petsc/src/mat/impls/mffd/wp.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
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 */
589371c9d4SSatish Balay static PetscErrorCode MatMFFDCompute_WP(MatMFFD ctx, Vec U, Vec a, PetscScalar *h, PetscBool *zeroa) {
59e884886fSBarry Smith   MatMFFD_WP *hctx = (MatMFFD_WP *)ctx->hctx;
6051a79602SBarry Smith   PetscReal   normU, norma;
61e884886fSBarry Smith 
62e884886fSBarry Smith   PetscFunctionBegin;
63e884886fSBarry Smith   if (!(ctx->count % ctx->recomputeperiod)) {
6451a79602SBarry Smith     if (hctx->computenormU || !ctx->ncurrenth) {
659566063dSJacob Faibussowitsch       PetscCall(VecNorm(U, NORM_2, &normU));
668f1a2a5eSBarry Smith       hctx->normUfact = PetscSqrtReal(1.0 + normU);
67e884886fSBarry Smith     }
689566063dSJacob Faibussowitsch     PetscCall(VecNorm(a, NORM_2, &norma));
69e884886fSBarry Smith     if (norma == 0.0) {
70e884886fSBarry Smith       *zeroa = PETSC_TRUE;
71e884886fSBarry Smith       PetscFunctionReturn(0);
72e884886fSBarry Smith     }
73e884886fSBarry Smith     *zeroa = PETSC_FALSE;
74e884886fSBarry Smith     *h     = ctx->error_rel * hctx->normUfact / norma;
75e884886fSBarry Smith   } else {
76e884886fSBarry Smith     *h = ctx->currenth;
77e884886fSBarry Smith   }
78e884886fSBarry Smith   PetscFunctionReturn(0);
79e884886fSBarry Smith }
80e884886fSBarry Smith 
81e884886fSBarry Smith /*
82e884886fSBarry Smith    MatMFFDView_WP - Prints information about this particular
83e884886fSBarry Smith      method for computing h. Note that this does not print the general
84e884886fSBarry Smith      information about the matrix free, that is printed by the calling
85e884886fSBarry Smith      routine.
86e884886fSBarry Smith 
87e884886fSBarry Smith   Input Parameters:
88e884886fSBarry Smith +   ctx - the matrix free context
89e884886fSBarry Smith -   viewer - the PETSc viewer
90e884886fSBarry Smith 
91e884886fSBarry Smith */
929371c9d4SSatish Balay static PetscErrorCode MatMFFDView_WP(MatMFFD ctx, PetscViewer viewer) {
93e884886fSBarry Smith   MatMFFD_WP *hctx = (MatMFFD_WP *)ctx->hctx;
94ace3abfcSBarry Smith   PetscBool   iascii;
95e884886fSBarry Smith 
96e884886fSBarry Smith   PetscFunctionBegin;
979566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
98e884886fSBarry Smith   if (iascii) {
992205254eSKarl Rupp     if (hctx->computenormU) {
1009566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Computes normU\n"));
1012205254eSKarl Rupp     } else {
1029566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Does not compute normU\n"));
1032205254eSKarl Rupp     }
10411aeaf0aSBarry Smith   }
105e884886fSBarry Smith   PetscFunctionReturn(0);
106e884886fSBarry Smith }
107e884886fSBarry Smith 
108e884886fSBarry Smith /*
109e884886fSBarry Smith    MatMFFDSetFromOptions_WP - Looks in the options database for
110e884886fSBarry Smith      any options appropriate for this method
111e884886fSBarry Smith 
112e884886fSBarry Smith   Input Parameter:
113e884886fSBarry Smith .  ctx - the matrix free context
114e884886fSBarry Smith 
115e884886fSBarry Smith */
1169371c9d4SSatish Balay static PetscErrorCode MatMFFDSetFromOptions_WP(MatMFFD ctx, PetscOptionItems *PetscOptionsObject) {
117e884886fSBarry Smith   MatMFFD_WP *hctx = (MatMFFD_WP *)ctx->hctx;
118e884886fSBarry Smith 
119e884886fSBarry Smith   PetscFunctionBegin;
120d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Walker-Pernice options");
1219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_mffd_compute_normu", "Compute the norm of u", "MatMFFDWPSetComputeNormU", hctx->computenormU, &hctx->computenormU, NULL));
122d0609cedSBarry Smith   PetscOptionsHeadEnd();
123e884886fSBarry Smith   PetscFunctionReturn(0);
124e884886fSBarry Smith }
125e884886fSBarry Smith 
126e884886fSBarry Smith /*
127e884886fSBarry Smith    MatMFFDDestroy_WP - Frees the space allocated by
1281d0fab5eSBarry Smith        MatCreateMFFD_WP().
129e884886fSBarry Smith 
130e884886fSBarry Smith   Input Parameter:
131e884886fSBarry Smith .  ctx - the matrix free context
132e884886fSBarry Smith 
13395452b02SPatrick Sanan    Notes:
13495452b02SPatrick Sanan     does not free the ctx, that is handled by the calling routine
135e884886fSBarry Smith 
136e884886fSBarry Smith */
1379371c9d4SSatish Balay static PetscErrorCode MatMFFDDestroy_WP(MatMFFD ctx) {
138e884886fSBarry Smith   PetscFunctionBegin;
1392e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ctx->mat, "MatMFFDWPSetComputeNormU_C", NULL));
1409566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx->hctx));
141e884886fSBarry Smith   PetscFunctionReturn(0);
142e884886fSBarry Smith }
143e884886fSBarry Smith 
1449371c9d4SSatish Balay PetscErrorCode MatMFFDWPSetComputeNormU_P(Mat mat, PetscBool flag) {
145e884886fSBarry Smith   MatMFFD     ctx  = (MatMFFD)mat->data;
146e884886fSBarry Smith   MatMFFD_WP *hctx = (MatMFFD_WP *)ctx->hctx;
147e884886fSBarry Smith 
148e884886fSBarry Smith   PetscFunctionBegin;
149e884886fSBarry Smith   hctx->computenormU = flag;
150e884886fSBarry Smith   PetscFunctionReturn(0);
151e884886fSBarry Smith }
152e884886fSBarry Smith 
153e884886fSBarry Smith /*@
15411a5261eSBarry Smith     MatMFFDWPSetComputeNormU - Sets whether it computes the ||U|| used by the Walker-Pernice
155e884886fSBarry Smith              PETSc routine for computing h. With any Krylov solver this need only
156e884886fSBarry Smith              be computed during the first iteration and kept for later.
157e884886fSBarry Smith 
158e884886fSBarry Smith   Input Parameters:
15911a5261eSBarry Smith +   A - the `MATMFFD` matrix
16011a5261eSBarry Smith -   flag - `PETSC_TRUE` causes it to compute ||U||, `PETSC_FALSE` uses the previous value
161e884886fSBarry Smith 
162e884886fSBarry Smith   Options Database Key:
163e884886fSBarry Smith .   -mat_mffd_compute_normu <true,false> - true by default, false can save calculations but you
164e884886fSBarry Smith               must be sure that ||U|| has not changed in the mean time.
165e884886fSBarry Smith 
166e884886fSBarry Smith   Level: advanced
167e884886fSBarry Smith 
16811a5261eSBarry Smith   Note:
16911a5261eSBarry Smith    See the manual page for `MATMFFD_WP` for a complete description of the
170e884886fSBarry Smith    algorithm used to compute h.
171e884886fSBarry Smith 
17211a5261eSBarry Smith .seealso: `MATMFFD_WP`, `MATMFFD`, `MatMFFDSetFunctionError()`, `MatCreateSNESMF()`
173e884886fSBarry Smith @*/
1749371c9d4SSatish Balay PetscErrorCode MatMFFDWPSetComputeNormU(Mat A, PetscBool flag) {
175e884886fSBarry Smith   PetscFunctionBegin;
1760700a824SBarry Smith   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
177cac4c232SBarry Smith   PetscTryMethod(A, "MatMFFDWPSetComputeNormU_C", (Mat, PetscBool), (A, flag));
178e884886fSBarry Smith   PetscFunctionReturn(0);
179e884886fSBarry Smith }
180e884886fSBarry Smith 
181e884886fSBarry Smith /*
1821d0fab5eSBarry Smith      MatCreateMFFD_WP - Standard PETSc code for
183e884886fSBarry Smith    computing h with matrix-free finite differences.
184e884886fSBarry Smith 
185e884886fSBarry Smith    Input Parameter:
1861d0fab5eSBarry Smith .  ctx - the matrix free context created by MatCreateMFFD()
187e884886fSBarry Smith 
188e884886fSBarry Smith */
1899371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatCreateMFFD_WP(MatMFFD ctx) {
190e884886fSBarry Smith   MatMFFD_WP *hctx;
191e884886fSBarry Smith 
192e884886fSBarry Smith   PetscFunctionBegin;
193e884886fSBarry Smith   /* allocate my own private data structure */
194*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&hctx));
195e884886fSBarry Smith   ctx->hctx          = (void *)hctx;
196e884886fSBarry Smith   hctx->computenormU = PETSC_FALSE;
197e884886fSBarry Smith 
198e884886fSBarry Smith   /* set the functions I am providing */
199e884886fSBarry Smith   ctx->ops->compute        = MatMFFDCompute_WP;
200e884886fSBarry Smith   ctx->ops->destroy        = MatMFFDDestroy_WP;
201e884886fSBarry Smith   ctx->ops->view           = MatMFFDView_WP;
202e884886fSBarry Smith   ctx->ops->setfromoptions = MatMFFDSetFromOptions_WP;
203e884886fSBarry Smith 
2049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ctx->mat, "MatMFFDWPSetComputeNormU_C", MatMFFDWPSetComputeNormU_P));
205e884886fSBarry Smith   PetscFunctionReturn(0);
206e884886fSBarry Smith }
207