xref: /petsc/src/mat/impls/mffd/mffddef.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1e884886fSBarry Smith 
2e884886fSBarry Smith /*
3e884886fSBarry Smith   Implements the DS PETSc approach for computing the h
4e884886fSBarry Smith   parameter used with the finite difference based matrix-free
5e884886fSBarry Smith   Jacobian-vector products.
6e884886fSBarry Smith 
7e884886fSBarry Smith   To make your own: clone this file and modify for your needs.
8e884886fSBarry Smith 
9e884886fSBarry Smith   Mandatory functions:
10e884886fSBarry Smith   -------------------
11e884886fSBarry Smith       MatMFFDCompute_  - for a given point and direction computes h
12e884886fSBarry Smith 
131d0fab5eSBarry Smith       MatCreateMFFD _ - fills in the MatMFFD data structure
14e884886fSBarry Smith                            for this particular implementation
15e884886fSBarry Smith 
16e884886fSBarry Smith    Optional functions:
17e884886fSBarry Smith    -------------------
18e884886fSBarry Smith       MatMFFDView_ - prints information about the parameters being used.
19e884886fSBarry Smith                        This is called when SNESView() or -snes_view is used.
20e884886fSBarry Smith 
21e884886fSBarry Smith       MatMFFDSetFromOptions_ - checks the options database for options that
22e884886fSBarry Smith                                apply to this method.
23e884886fSBarry Smith 
24e884886fSBarry Smith       MatMFFDDestroy_ - frees any space allocated by the routines above
25e884886fSBarry Smith 
26e884886fSBarry Smith */
27e884886fSBarry Smith 
28e884886fSBarry Smith /*
29e884886fSBarry Smith     This include file defines the data structure  MatMFFD that
30e884886fSBarry Smith    includes information about the computation of h. It is shared by
31e884886fSBarry Smith    all implementations that people provide
32e884886fSBarry Smith */
33af0996ceSBarry Smith #include <petsc/private/matimpl.h>
34c6db04a5SJed Brown #include <../src/mat/impls/mffd/mffdimpl.h> /*I  "petscmat.h"   I*/
35e884886fSBarry Smith 
36e884886fSBarry Smith /*
37e884886fSBarry Smith       The  method has one parameter that is used to
38e884886fSBarry Smith    "cutoff" very small values. This is stored in a data structure
39e884886fSBarry Smith    that is only visible to this file. If your method has no parameters
40e884886fSBarry Smith    it can omit this, if it has several simply reorganize the data structure.
41e884886fSBarry Smith    The data structure is "hung-off" the MatMFFD data structure in
42e884886fSBarry Smith    the void *hctx; field.
43e884886fSBarry Smith */
44e884886fSBarry Smith typedef struct {
45e884886fSBarry Smith   PetscReal umin; /* minimum allowable u'a value relative to |u|_1 */
46e884886fSBarry Smith } MatMFFD_DS;
47e884886fSBarry Smith 
48e884886fSBarry Smith /*
49e884886fSBarry Smith    MatMFFDCompute_DS - Standard PETSc code for computing the
50fd292e60Sprj-    differencing parameter (h) for use with matrix-free finite differences.
51e884886fSBarry Smith 
52e884886fSBarry Smith    Input Parameters:
53e884886fSBarry Smith +  ctx - the matrix free context
54e884886fSBarry Smith .  U - the location at which you want the Jacobian
55e884886fSBarry Smith -  a - the direction you want the derivative
56e884886fSBarry Smith 
57e884886fSBarry Smith    Output Parameter:
58e884886fSBarry Smith .  h - the scale computed
59e884886fSBarry Smith 
60e884886fSBarry Smith */
619371c9d4SSatish Balay static PetscErrorCode MatMFFDCompute_DS(MatMFFD ctx, Vec U, Vec a, PetscScalar *h, PetscBool *zeroa) {
62e884886fSBarry Smith   MatMFFD_DS *hctx = (MatMFFD_DS *)ctx->hctx;
63e884886fSBarry Smith   PetscReal   nrm, sum, umin = hctx->umin;
64e884886fSBarry Smith   PetscScalar dot;
65e884886fSBarry Smith 
66e884886fSBarry Smith   PetscFunctionBegin;
67e884886fSBarry Smith   if (!(ctx->count % ctx->recomputeperiod)) {
68e884886fSBarry Smith     /*
69e884886fSBarry Smith      This algorithm requires 2 norms and 1 inner product. Rather than
70e884886fSBarry Smith      use directly the VecNorm() and VecDot() routines (and thus have
71e884886fSBarry Smith      three separate collective operations, we use the VecxxxBegin/End() routines
72e884886fSBarry Smith     */
739566063dSJacob Faibussowitsch     PetscCall(VecDotBegin(U, a, &dot));
749566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(a, NORM_1, &sum));
759566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(a, NORM_2, &nrm));
769566063dSJacob Faibussowitsch     PetscCall(VecDotEnd(U, a, &dot));
779566063dSJacob Faibussowitsch     PetscCall(VecNormEnd(a, NORM_1, &sum));
789566063dSJacob Faibussowitsch     PetscCall(VecNormEnd(a, NORM_2, &nrm));
79e884886fSBarry Smith 
80e884886fSBarry Smith     if (nrm == 0.0) {
81e884886fSBarry Smith       *zeroa = PETSC_TRUE;
82e884886fSBarry Smith       PetscFunctionReturn(0);
83e884886fSBarry Smith     }
84e884886fSBarry Smith     *zeroa = PETSC_FALSE;
85e884886fSBarry Smith 
86e884886fSBarry Smith     /*
87e884886fSBarry Smith       Safeguard for step sizes that are "too small"
88e884886fSBarry Smith     */
89e884886fSBarry Smith     if (PetscAbsScalar(dot) < umin * sum && PetscRealPart(dot) >= 0.0) dot = umin * sum;
90e884886fSBarry Smith     else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin * sum) dot = -umin * sum;
91e884886fSBarry Smith     *h = ctx->error_rel * dot / (nrm * nrm);
92aed4548fSBarry Smith     PetscCheck(!PetscIsInfOrNanScalar(*h), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Differencing parameter is not a number sum = %g dot = %g norm = %g", (double)sum, (double)PetscRealPart(dot), (double)nrm);
93e884886fSBarry Smith   } else {
94e884886fSBarry Smith     *h = ctx->currenth;
95e884886fSBarry Smith   }
96e884886fSBarry Smith   ctx->count++;
97e884886fSBarry Smith   PetscFunctionReturn(0);
98e884886fSBarry Smith }
99e884886fSBarry Smith 
100e884886fSBarry Smith /*
101e884886fSBarry Smith    MatMFFDView_DS - Prints information about this particular
102e884886fSBarry Smith    method for computing h. Note that this does not print the general
103e884886fSBarry Smith    information about the matrix-free method, as such info is printed
104e884886fSBarry Smith    by the calling routine.
105e884886fSBarry Smith 
106e884886fSBarry Smith    Input Parameters:
107e884886fSBarry Smith +  ctx - the matrix free context
108e884886fSBarry Smith -  viewer - the PETSc viewer
109e884886fSBarry Smith */
1109371c9d4SSatish Balay static PetscErrorCode MatMFFDView_DS(MatMFFD ctx, PetscViewer viewer) {
111e884886fSBarry Smith   MatMFFD_DS *hctx = (MatMFFD_DS *)ctx->hctx;
112ace3abfcSBarry Smith   PetscBool   iascii;
113e884886fSBarry Smith 
114e884886fSBarry Smith   PetscFunctionBegin;
115e884886fSBarry Smith   /*
116e884886fSBarry Smith      Currently this only handles the ascii file viewers, others
117e884886fSBarry Smith      could be added, but for this type of object other viewers
118e884886fSBarry Smith      make less sense
119e884886fSBarry Smith   */
1209566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
12148a46eb9SPierre Jolivet   if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, "    umin=%g (minimum iterate parameter)\n", (double)hctx->umin));
122e884886fSBarry Smith   PetscFunctionReturn(0);
123e884886fSBarry Smith }
124e884886fSBarry Smith 
125e884886fSBarry Smith /*
126e884886fSBarry Smith    MatMFFDSetFromOptions_DS - Looks in the options database for
127e884886fSBarry Smith    any options appropriate for this method.
128e884886fSBarry Smith 
129e884886fSBarry Smith    Input Parameter:
130e884886fSBarry Smith .  ctx - the matrix free context
131e884886fSBarry Smith 
132e884886fSBarry Smith */
1339371c9d4SSatish Balay static PetscErrorCode MatMFFDSetFromOptions_DS(MatMFFD ctx, PetscOptionItems *PetscOptionsObject) {
134e884886fSBarry Smith   MatMFFD_DS *hctx = (MatMFFD_DS *)ctx->hctx;
135e884886fSBarry Smith 
136e884886fSBarry Smith   PetscFunctionBegin;
137d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Finite difference matrix free parameters");
1389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mffd_umin", "umin", "MatMFFDDSSetUmin", hctx->umin, &hctx->umin, NULL));
139d0609cedSBarry Smith   PetscOptionsHeadEnd();
140e884886fSBarry Smith   PetscFunctionReturn(0);
141e884886fSBarry Smith }
142e884886fSBarry Smith 
143e884886fSBarry Smith /*
144e884886fSBarry Smith    MatMFFDDestroy_DS - Frees the space allocated by
1451d0fab5eSBarry Smith    MatCreateMFFD_DS().
146e884886fSBarry Smith 
147e884886fSBarry Smith    Input Parameter:
148e884886fSBarry Smith .  ctx - the matrix free context
149e884886fSBarry Smith 
15011a5261eSBarry Smith    Note:
151e884886fSBarry Smith    Does not free the ctx, that is handled by the calling routine
152e884886fSBarry Smith */
1539371c9d4SSatish Balay static PetscErrorCode MatMFFDDestroy_DS(MatMFFD ctx) {
154e884886fSBarry Smith   PetscFunctionBegin;
1559566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx->hctx));
156e884886fSBarry Smith   PetscFunctionReturn(0);
157e884886fSBarry Smith }
158e884886fSBarry Smith 
159e884886fSBarry Smith /*
160e884886fSBarry Smith    The following two routines use the PetscObjectCompose() and PetscObjectQuery()
161e884886fSBarry Smith    mechanism to allow the user to change the Umin parameter used in this method.
162e884886fSBarry Smith */
1639371c9d4SSatish Balay PetscErrorCode MatMFFDDSSetUmin_DS(Mat mat, PetscReal umin) {
164440c47e6SSatish Balay   MatMFFD     ctx = NULL;
165e884886fSBarry Smith   MatMFFD_DS *hctx;
166e884886fSBarry Smith 
167e884886fSBarry Smith   PetscFunctionBegin;
1689566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
16928b400f6SJacob Faibussowitsch   PetscCheck(ctx, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MatMFFDDSSetUmin() attached to non-shell matrix");
170e884886fSBarry Smith   hctx       = (MatMFFD_DS *)ctx->hctx;
171e884886fSBarry Smith   hctx->umin = umin;
172e884886fSBarry Smith   PetscFunctionReturn(0);
173e884886fSBarry Smith }
174e884886fSBarry Smith 
175e884886fSBarry Smith /*@
176e884886fSBarry Smith     MatMFFDDSSetUmin - Sets the "umin" parameter used by the
177e884886fSBarry Smith     PETSc routine for computing the differencing parameter, h, which is used
17811a5261eSBarry Smith     for matrix-free Jacobian-vector products for a `MATMFFD` matrix.
179e884886fSBarry Smith 
180e884886fSBarry Smith    Input Parameters:
18111a5261eSBarry Smith +  A - the `MATMFFD` matrix
182e884886fSBarry Smith -  umin - the parameter
183e884886fSBarry Smith 
184e884886fSBarry Smith    Level: advanced
185e884886fSBarry Smith 
18611a5261eSBarry Smith    Note:
18711a5261eSBarry Smith    See the manual page for `MatCreateSNESMF()` for a complete description of the
188e884886fSBarry Smith    algorithm used to compute h.
189e884886fSBarry Smith 
19011a5261eSBarry Smith .seealso: `MATMFFD`, `MatMFFDSetFunctionError()`, `MatCreateSNESMF()`
191e884886fSBarry Smith @*/
1929371c9d4SSatish Balay PetscErrorCode MatMFFDDSSetUmin(Mat A, PetscReal umin) {
193e884886fSBarry Smith   PetscFunctionBegin;
1940700a824SBarry Smith   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
195cac4c232SBarry Smith   PetscTryMethod(A, "MatMFFDDSSetUmin_C", (Mat, PetscReal), (A, umin));
196e884886fSBarry Smith   PetscFunctionReturn(0);
197e884886fSBarry Smith }
198e884886fSBarry Smith 
199e884886fSBarry Smith /*MC
20011a5261eSBarry Smith      MATMFFD_DS - algorithm for compute the "h" used in the finite difference matrix-free matrix vector product, `MatMult()`.
201e884886fSBarry Smith 
202e884886fSBarry Smith    Options Database Keys:
20311a5261eSBarry Smith .  -mat_mffd_umin <umin> - see `MatMFFDDSSetUmin()`
204e884886fSBarry Smith 
205e884886fSBarry Smith    Level: intermediate
206e884886fSBarry Smith 
20795452b02SPatrick Sanan    Notes:
20895452b02SPatrick Sanan     Requires 2 norms and 1 inner product, but they are computed together
20911a5261eSBarry Smith        so only one parallel collective operation is needed. See `MATMFFD_WP` for a method
21011a5261eSBarry Smith        (with `KSPGMRES`) that requires NO collective operations.
211e884886fSBarry Smith 
212e884886fSBarry Smith    Formula used:
213e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
214e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
215e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
216e884886fSBarry Smith  where
217e884886fSBarry Smith      error_rel = square root of relative error in function evaluation
218e884886fSBarry Smith      umin = minimum iterate parameter
219e884886fSBarry Smith 
22011a5261eSBarry Smith   References:
22111a5261eSBarry Smith .  * -  Dennis and Schnabel, "Numerical Methods for Unconstrained Optimization and Nonlinear Equations"
222e884886fSBarry Smith 
22311a5261eSBarry Smith .seealso: `MATMFFD`, `MATMFFD_WP`, `MatCreateMFFD()`, `MatCreateSNESMF()`, `MATMFFD_WP`, `MatMFFDDSSetUmin()`
224e884886fSBarry Smith M*/
2259371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatCreateMFFD_DS(MatMFFD ctx) {
226e884886fSBarry Smith   MatMFFD_DS *hctx;
227e884886fSBarry Smith 
228e884886fSBarry Smith   PetscFunctionBegin;
229e884886fSBarry Smith   /* allocate my own private data structure */
230*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&hctx));
231e884886fSBarry Smith   ctx->hctx  = (void *)hctx;
232e884886fSBarry Smith   /* set a default for my parameter */
233e884886fSBarry Smith   hctx->umin = 1.e-6;
234e884886fSBarry Smith 
235e884886fSBarry Smith   /* set the functions I am providing */
236e884886fSBarry Smith   ctx->ops->compute        = MatMFFDCompute_DS;
237e884886fSBarry Smith   ctx->ops->destroy        = MatMFFDDestroy_DS;
238e884886fSBarry Smith   ctx->ops->view           = MatMFFDView_DS;
239e884886fSBarry Smith   ctx->ops->setfromoptions = MatMFFDSetFromOptions_DS;
240e884886fSBarry Smith 
2419566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ctx->mat, "MatMFFDDSSetUmin_C", MatMFFDDSSetUmin_DS));
242e884886fSBarry Smith   PetscFunctionReturn(0);
243e884886fSBarry Smith }
244