1e884886fSBarry Smith /*MC
211a5261eSBarry Smith MATMFFD_WP - Implements an approach for computing the differencing parameter
31d27aa22SBarry Smith h used with the finite difference based matrix-free Jacobian. From Walker-Pernice {cite}`pw98`
4e884886fSBarry Smith
51d27aa22SBarry Smith $$
61d27aa22SBarry Smith h = error_rel * sqrt(1 + ||u||) / ||a||
71d27aa22SBarry Smith $$
8e884886fSBarry Smith
911a5261eSBarry Smith Options Database Key:
1011a5261eSBarry Smith . -mat_mffd_compute_normu -Compute the norm of u every time see `MatMFFDWPSetComputeNormU()`
11e884886fSBarry Smith
12e884886fSBarry Smith Level: intermediate
13e884886fSBarry Smith
1495452b02SPatrick Sanan Notes:
151d27aa22SBarry Smith $ || U || $ does not change between linear iterations so is reused
1611a5261eSBarry Smith
171d27aa22SBarry Smith In `KSPGMRES` $ || a || == 1 $ and so does not need to ever be computed except at restart
181cb3f120SPierre Jolivet when it is recomputed. Thus requires no global collectives when used with `KSPGMRES`
19e884886fSBarry Smith
2011a5261eSBarry Smith .seealso: `MATMFFD`, `MATMFFD_DS`, `MatCreateMFFD()`, `MatCreateSNESMF()`, `MATMFFD_DS`
21e884886fSBarry Smith M*/
22e884886fSBarry Smith
23e884886fSBarry Smith /*
24e884886fSBarry Smith This include file defines the data structure MatMFFD that
25e884886fSBarry Smith includes information about the computation of h. It is shared by
26e884886fSBarry Smith all implementations that people provide.
27e884886fSBarry Smith
28e884886fSBarry Smith See snesmfjdef.c for a full set of comments on the routines below.
29e884886fSBarry Smith */
30af0996ceSBarry Smith #include <petsc/private/matimpl.h>
31c6db04a5SJed Brown #include <../src/mat/impls/mffd/mffdimpl.h> /*I "petscmat.h" I*/
32e884886fSBarry Smith
33e884886fSBarry Smith typedef struct {
34e884886fSBarry Smith PetscReal normUfact; /* previous sqrt(1.0 + || U ||) */
35ace3abfcSBarry Smith PetscBool computenormU;
36e884886fSBarry Smith } MatMFFD_WP;
37e884886fSBarry Smith
38e884886fSBarry Smith /*
3911a5261eSBarry Smith MatMFFDCompute_WP - code for
40e884886fSBarry Smith computing h with matrix-free finite differences.
41e884886fSBarry Smith
42e884886fSBarry Smith Input Parameters:
4301c1178eSBarry Smith + ctx - the matrix-free context
44e884886fSBarry Smith . U - the location at which you want the Jacobian
45e884886fSBarry Smith - a - the direction you want the derivative
46e884886fSBarry Smith
47e884886fSBarry Smith Output Parameter:
48e884886fSBarry Smith . h - the scale computed
49e884886fSBarry Smith */
MatMFFDCompute_WP(MatMFFD ctx,Vec U,Vec a,PetscScalar * h,PetscBool * zeroa)50d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDCompute_WP(MatMFFD ctx, Vec U, Vec a, PetscScalar *h, PetscBool *zeroa)
51d71ae5a4SJacob Faibussowitsch {
52e884886fSBarry Smith MatMFFD_WP *hctx = (MatMFFD_WP *)ctx->hctx;
5351a79602SBarry Smith PetscReal normU, norma;
54e884886fSBarry Smith
55e884886fSBarry Smith PetscFunctionBegin;
56e884886fSBarry Smith if (!(ctx->count % ctx->recomputeperiod)) {
5751a79602SBarry Smith if (hctx->computenormU || !ctx->ncurrenth) {
589566063dSJacob Faibussowitsch PetscCall(VecNorm(U, NORM_2, &normU));
598f1a2a5eSBarry Smith hctx->normUfact = PetscSqrtReal(1.0 + normU);
60e884886fSBarry Smith }
619566063dSJacob Faibussowitsch PetscCall(VecNorm(a, NORM_2, &norma));
62e884886fSBarry Smith if (norma == 0.0) {
63e884886fSBarry Smith *zeroa = PETSC_TRUE;
643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
65e884886fSBarry Smith }
66e884886fSBarry Smith *zeroa = PETSC_FALSE;
67e884886fSBarry Smith *h = ctx->error_rel * hctx->normUfact / norma;
68e884886fSBarry Smith } else {
69e884886fSBarry Smith *h = ctx->currenth;
70e884886fSBarry Smith }
713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
72e884886fSBarry Smith }
73e884886fSBarry Smith
74e884886fSBarry Smith /*
75e884886fSBarry Smith MatMFFDView_WP - Prints information about this particular
76e884886fSBarry Smith method for computing h. Note that this does not print the general
7701c1178eSBarry Smith information about the matrix-free, that is printed by the calling
78e884886fSBarry Smith routine.
79e884886fSBarry Smith
80e884886fSBarry Smith Input Parameters:
8101c1178eSBarry Smith + ctx - the matrix-free context
82e884886fSBarry Smith - viewer - the PETSc viewer
83e884886fSBarry Smith
84e884886fSBarry Smith */
MatMFFDView_WP(MatMFFD ctx,PetscViewer viewer)85d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDView_WP(MatMFFD ctx, PetscViewer viewer)
86d71ae5a4SJacob Faibussowitsch {
87e884886fSBarry Smith MatMFFD_WP *hctx = (MatMFFD_WP *)ctx->hctx;
88*9f196a02SMartin Diehl PetscBool isascii;
89e884886fSBarry Smith
90e884886fSBarry Smith PetscFunctionBegin;
91*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
92*9f196a02SMartin Diehl if (isascii) {
932205254eSKarl Rupp if (hctx->computenormU) {
949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Computes normU\n"));
952205254eSKarl Rupp } else {
969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Does not compute normU\n"));
972205254eSKarl Rupp }
9811aeaf0aSBarry Smith }
993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
100e884886fSBarry Smith }
101e884886fSBarry Smith
102e884886fSBarry Smith /*
103e884886fSBarry Smith MatMFFDSetFromOptions_WP - Looks in the options database for
104e884886fSBarry Smith any options appropriate for this method
105e884886fSBarry Smith
106e884886fSBarry Smith Input Parameter:
10701c1178eSBarry Smith . ctx - the matrix-free context
108e884886fSBarry Smith
109e884886fSBarry Smith */
MatMFFDSetFromOptions_WP(MatMFFD ctx,PetscOptionItems PetscOptionsObject)110ce78bad3SBarry Smith static PetscErrorCode MatMFFDSetFromOptions_WP(MatMFFD ctx, PetscOptionItems PetscOptionsObject)
111d71ae5a4SJacob Faibussowitsch {
112e884886fSBarry Smith MatMFFD_WP *hctx = (MatMFFD_WP *)ctx->hctx;
113e884886fSBarry Smith
114e884886fSBarry Smith PetscFunctionBegin;
115d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Walker-Pernice options");
1169566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_mffd_compute_normu", "Compute the norm of u", "MatMFFDWPSetComputeNormU", hctx->computenormU, &hctx->computenormU, NULL));
117d0609cedSBarry Smith PetscOptionsHeadEnd();
1183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
119e884886fSBarry Smith }
120e884886fSBarry Smith
MatMFFDDestroy_WP(MatMFFD ctx)121d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDDestroy_WP(MatMFFD ctx)
122d71ae5a4SJacob Faibussowitsch {
123e884886fSBarry Smith PetscFunctionBegin;
1242e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ctx->mat, "MatMFFDWPSetComputeNormU_C", NULL));
1259566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx->hctx));
1263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
127e884886fSBarry Smith }
128e884886fSBarry Smith
MatMFFDWPSetComputeNormU_P(Mat mat,PetscBool flag)12966976f2fSJacob Faibussowitsch static PetscErrorCode MatMFFDWPSetComputeNormU_P(Mat mat, PetscBool flag)
130d71ae5a4SJacob Faibussowitsch {
131e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data;
132e884886fSBarry Smith MatMFFD_WP *hctx = (MatMFFD_WP *)ctx->hctx;
133e884886fSBarry Smith
134e884886fSBarry Smith PetscFunctionBegin;
135e884886fSBarry Smith hctx->computenormU = flag;
1363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
137e884886fSBarry Smith }
138e884886fSBarry Smith
139e884886fSBarry Smith /*@
1401d27aa22SBarry Smith MatMFFDWPSetComputeNormU - Sets whether it computes the ||U|| used by the Walker-Pernice {cite}`pw98`
141e884886fSBarry Smith PETSc routine for computing h. With any Krylov solver this need only
142e884886fSBarry Smith be computed during the first iteration and kept for later.
143e884886fSBarry Smith
144e884886fSBarry Smith Input Parameters:
14511a5261eSBarry Smith + A - the `MATMFFD` matrix
1461d27aa22SBarry Smith - flag - `PETSC_TRUE` causes it to compute $||U||$, `PETSC_FALSE` uses the previous value
147e884886fSBarry Smith
148e884886fSBarry Smith Options Database Key:
149e884886fSBarry Smith . -mat_mffd_compute_normu <true,false> - true by default, false can save calculations but you
1501d27aa22SBarry Smith must be sure that $||U||$ has not changed in the mean time.
151e884886fSBarry Smith
152e884886fSBarry Smith Level: advanced
153e884886fSBarry Smith
15411a5261eSBarry Smith Note:
15511a5261eSBarry Smith See the manual page for `MATMFFD_WP` for a complete description of the
156e884886fSBarry Smith algorithm used to compute h.
157e884886fSBarry Smith
15811a5261eSBarry Smith .seealso: `MATMFFD_WP`, `MATMFFD`, `MatMFFDSetFunctionError()`, `MatCreateSNESMF()`
159e884886fSBarry Smith @*/
MatMFFDWPSetComputeNormU(Mat A,PetscBool flag)160d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDWPSetComputeNormU(Mat A, PetscBool flag)
161d71ae5a4SJacob Faibussowitsch {
162e884886fSBarry Smith PetscFunctionBegin;
1630700a824SBarry Smith PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
164cac4c232SBarry Smith PetscTryMethod(A, "MatMFFDWPSetComputeNormU_C", (Mat, PetscBool), (A, flag));
1653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
166e884886fSBarry Smith }
167e884886fSBarry Smith
168e884886fSBarry Smith /*
1691d0fab5eSBarry Smith MatCreateMFFD_WP - Standard PETSc code for
170e884886fSBarry Smith computing h with matrix-free finite differences.
171e884886fSBarry Smith
172e884886fSBarry Smith Input Parameter:
17301c1178eSBarry Smith . ctx - the matrix-free context created by MatCreateMFFD()
174e884886fSBarry Smith
175e884886fSBarry Smith */
MatCreateMFFD_WP(MatMFFD ctx)176d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreateMFFD_WP(MatMFFD ctx)
177d71ae5a4SJacob Faibussowitsch {
178e884886fSBarry Smith MatMFFD_WP *hctx;
179e884886fSBarry Smith
180e884886fSBarry Smith PetscFunctionBegin;
181e884886fSBarry Smith /* allocate my own private data structure */
1824dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&hctx));
183e884886fSBarry Smith ctx->hctx = (void *)hctx;
184e884886fSBarry Smith hctx->computenormU = PETSC_FALSE;
185e884886fSBarry Smith
186e884886fSBarry Smith /* set the functions I am providing */
187e884886fSBarry Smith ctx->ops->compute = MatMFFDCompute_WP;
188e884886fSBarry Smith ctx->ops->destroy = MatMFFDDestroy_WP;
189e884886fSBarry Smith ctx->ops->view = MatMFFDView_WP;
190e884886fSBarry Smith ctx->ops->setfromoptions = MatMFFDSetFromOptions_WP;
191e884886fSBarry Smith
1929566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ctx->mat, "MatMFFDWPSetComputeNormU_C", MatMFFDWPSetComputeNormU_P));
1933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
194e884886fSBarry Smith }
195