xref: /petsc/src/mat/impls/mffd/mffddef.c (revision bcda9346efad4e5ba2d553af84eb238771ba1e25)
1e884886fSBarry Smith /*
2e884886fSBarry Smith   Implements the DS PETSc approach for computing the h
3e884886fSBarry Smith   parameter used with the finite difference based matrix-free
4e884886fSBarry Smith   Jacobian-vector products.
5e884886fSBarry Smith 
6e884886fSBarry Smith   To make your own: clone this file and modify for your needs.
7e884886fSBarry Smith 
8e884886fSBarry Smith   Mandatory functions:
9e884886fSBarry Smith   -------------------
10e884886fSBarry Smith       MatMFFDCompute_  - for a given point and direction computes h
11e884886fSBarry Smith 
121d0fab5eSBarry Smith       MatCreateMFFD _ - fills in the MatMFFD data structure
13e884886fSBarry Smith                            for this particular implementation
14e884886fSBarry Smith 
15e884886fSBarry Smith    Optional functions:
16e884886fSBarry Smith    -------------------
17e884886fSBarry Smith       MatMFFDView_ - prints information about the parameters being used.
18e884886fSBarry Smith                        This is called when SNESView() or -snes_view is used.
19e884886fSBarry Smith 
20e884886fSBarry Smith       MatMFFDSetFromOptions_ - checks the options database for options that
21e884886fSBarry Smith                                apply to this method.
22e884886fSBarry Smith 
23e884886fSBarry Smith       MatMFFDDestroy_ - frees any space allocated by the routines above
24e884886fSBarry Smith 
25e884886fSBarry Smith */
26e884886fSBarry Smith 
27e884886fSBarry Smith /*
28e884886fSBarry Smith     This include file defines the data structure  MatMFFD that
29e884886fSBarry Smith    includes information about the computation of h. It is shared by
30e884886fSBarry Smith    all implementations that people provide
31e884886fSBarry Smith */
32af0996ceSBarry Smith #include <petsc/private/matimpl.h>
33c6db04a5SJed Brown #include <../src/mat/impls/mffd/mffdimpl.h> /*I  "petscmat.h"   I*/
34e884886fSBarry Smith 
35e884886fSBarry Smith /*
36e884886fSBarry Smith       The  method has one parameter that is used to
37e884886fSBarry Smith    "cutoff" very small values. This is stored in a data structure
38e884886fSBarry Smith    that is only visible to this file. If your method has no parameters
39e884886fSBarry Smith    it can omit this, if it has several simply reorganize the data structure.
40e884886fSBarry Smith    The data structure is "hung-off" the MatMFFD data structure in
41e884886fSBarry Smith    the void *hctx; field.
42e884886fSBarry Smith */
43e884886fSBarry Smith typedef struct {
44e884886fSBarry Smith   PetscReal umin; /* minimum allowable u'a value relative to |u|_1 */
45e884886fSBarry Smith } MatMFFD_DS;
46e884886fSBarry Smith 
MatMFFDCompute_DS(MatMFFD ctx,Vec U,Vec a,PetscScalar * h,PetscBool * zeroa)47d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDCompute_DS(MatMFFD ctx, Vec U, Vec a, PetscScalar *h, PetscBool *zeroa)
48d71ae5a4SJacob Faibussowitsch {
49e884886fSBarry Smith   MatMFFD_DS *hctx = (MatMFFD_DS *)ctx->hctx;
50e884886fSBarry Smith   PetscReal   nrm, sum, umin = hctx->umin;
51e884886fSBarry Smith   PetscScalar dot;
52e884886fSBarry Smith 
53e884886fSBarry Smith   PetscFunctionBegin;
54e884886fSBarry Smith   if (!(ctx->count % ctx->recomputeperiod)) {
55e884886fSBarry Smith     /*
56e884886fSBarry Smith      This algorithm requires 2 norms and 1 inner product. Rather than
57e884886fSBarry Smith      use directly the VecNorm() and VecDot() routines (and thus have
58e884886fSBarry Smith      three separate collective operations, we use the VecxxxBegin/End() routines
59e884886fSBarry Smith     */
609566063dSJacob Faibussowitsch     PetscCall(VecDotBegin(U, a, &dot));
619566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(a, NORM_1, &sum));
629566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(a, NORM_2, &nrm));
639566063dSJacob Faibussowitsch     PetscCall(VecDotEnd(U, a, &dot));
649566063dSJacob Faibussowitsch     PetscCall(VecNormEnd(a, NORM_1, &sum));
659566063dSJacob Faibussowitsch     PetscCall(VecNormEnd(a, NORM_2, &nrm));
66e884886fSBarry Smith 
67e884886fSBarry Smith     if (nrm == 0.0) {
68e884886fSBarry Smith       *zeroa = PETSC_TRUE;
693ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
70e884886fSBarry Smith     }
71e884886fSBarry Smith     *zeroa = PETSC_FALSE;
72e884886fSBarry Smith 
73e884886fSBarry Smith     /*
74e884886fSBarry Smith       Safeguard for step sizes that are "too small"
75e884886fSBarry Smith     */
76e884886fSBarry Smith     if (PetscAbsScalar(dot) < umin * sum && PetscRealPart(dot) >= 0.0) dot = umin * sum;
77e884886fSBarry Smith     else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin * sum) dot = -umin * sum;
78883424caSPierre Jolivet     *h = dot * ctx->error_rel / (nrm * nrm);
79aed4548fSBarry 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);
80e884886fSBarry Smith   } else {
81e884886fSBarry Smith     *h = ctx->currenth;
82e884886fSBarry Smith   }
83e884886fSBarry Smith   ctx->count++;
843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
85e884886fSBarry Smith }
86e884886fSBarry Smith 
MatMFFDView_DS(MatMFFD ctx,PetscViewer viewer)87d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDView_DS(MatMFFD ctx, PetscViewer viewer)
88d71ae5a4SJacob Faibussowitsch {
89e884886fSBarry Smith   MatMFFD_DS *hctx = (MatMFFD_DS *)ctx->hctx;
90*9f196a02SMartin Diehl   PetscBool   isascii;
91e884886fSBarry Smith 
92e884886fSBarry Smith   PetscFunctionBegin;
93e884886fSBarry Smith   /*
94e884886fSBarry Smith      Currently this only handles the ascii file viewers, others
95e884886fSBarry Smith      could be added, but for this type of object other viewers
96e884886fSBarry Smith      make less sense
97e884886fSBarry Smith   */
98*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
99*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "    umin=%g (minimum iterate parameter)\n", (double)hctx->umin));
1003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
101e884886fSBarry Smith }
102e884886fSBarry Smith 
MatMFFDSetFromOptions_DS(MatMFFD ctx,PetscOptionItems PetscOptionsObject)103ce78bad3SBarry Smith static PetscErrorCode MatMFFDSetFromOptions_DS(MatMFFD ctx, PetscOptionItems PetscOptionsObject)
104d71ae5a4SJacob Faibussowitsch {
105e884886fSBarry Smith   MatMFFD_DS *hctx = (MatMFFD_DS *)ctx->hctx;
106e884886fSBarry Smith 
107e884886fSBarry Smith   PetscFunctionBegin;
10801c1178eSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Finite difference matrix-free parameters");
1099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_mffd_umin", "umin", "MatMFFDDSSetUmin", hctx->umin, &hctx->umin, NULL));
110d0609cedSBarry Smith   PetscOptionsHeadEnd();
1113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
112e884886fSBarry Smith }
113e884886fSBarry Smith 
MatMFFDDestroy_DS(MatMFFD ctx)114d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDDestroy_DS(MatMFFD ctx)
115d71ae5a4SJacob Faibussowitsch {
116e884886fSBarry Smith   PetscFunctionBegin;
1179566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx->hctx));
1183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
119e884886fSBarry Smith }
120e884886fSBarry Smith 
121e884886fSBarry Smith /*
122e884886fSBarry Smith    The following two routines use the PetscObjectCompose() and PetscObjectQuery()
123e884886fSBarry Smith    mechanism to allow the user to change the Umin parameter used in this method.
124e884886fSBarry Smith */
MatMFFDDSSetUmin_DS(Mat mat,PetscReal umin)12566976f2fSJacob Faibussowitsch static PetscErrorCode MatMFFDDSSetUmin_DS(Mat mat, PetscReal umin)
126d71ae5a4SJacob Faibussowitsch {
127440c47e6SSatish Balay   MatMFFD     ctx = NULL;
128e884886fSBarry Smith   MatMFFD_DS *hctx;
129e884886fSBarry Smith 
130e884886fSBarry Smith   PetscFunctionBegin;
1319566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &ctx));
13228b400f6SJacob Faibussowitsch   PetscCheck(ctx, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MatMFFDDSSetUmin() attached to non-shell matrix");
133e884886fSBarry Smith   hctx       = (MatMFFD_DS *)ctx->hctx;
134e884886fSBarry Smith   hctx->umin = umin;
1353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
136e884886fSBarry Smith }
137e884886fSBarry Smith 
138e884886fSBarry Smith /*@
139e884886fSBarry Smith   MatMFFDDSSetUmin - Sets the "umin" parameter used by the
140e884886fSBarry Smith   PETSc routine for computing the differencing parameter, h, which is used
14111a5261eSBarry Smith   for matrix-free Jacobian-vector products for a `MATMFFD` matrix.
142e884886fSBarry Smith 
143e884886fSBarry Smith   Input Parameters:
14411a5261eSBarry Smith + A    - the `MATMFFD` matrix
145e884886fSBarry Smith - umin - the parameter
146e884886fSBarry Smith 
147e884886fSBarry Smith   Level: advanced
148e884886fSBarry Smith 
14911a5261eSBarry Smith   Note:
15011a5261eSBarry Smith   See the manual page for `MatCreateSNESMF()` for a complete description of the
151e884886fSBarry Smith   algorithm used to compute h.
152e884886fSBarry Smith 
15311a5261eSBarry Smith .seealso: `MATMFFD`, `MatMFFDSetFunctionError()`, `MatCreateSNESMF()`
154e884886fSBarry Smith @*/
MatMFFDDSSetUmin(Mat A,PetscReal umin)155d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDDSSetUmin(Mat A, PetscReal umin)
156d71ae5a4SJacob Faibussowitsch {
157e884886fSBarry Smith   PetscFunctionBegin;
1580700a824SBarry Smith   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
159cac4c232SBarry Smith   PetscTryMethod(A, "MatMFFDDSSetUmin_C", (Mat, PetscReal), (A, umin));
1603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
161e884886fSBarry Smith }
162e884886fSBarry Smith 
163e884886fSBarry Smith /*MC
16411a5261eSBarry Smith      MATMFFD_DS - algorithm for compute the "h" used in the finite difference matrix-free matrix vector product, `MatMult()`.
165e884886fSBarry Smith 
166e884886fSBarry Smith    Options Database Keys:
16711a5261eSBarry Smith .  -mat_mffd_umin <umin> - see `MatMFFDDSSetUmin()`
168e884886fSBarry Smith 
169e884886fSBarry Smith    Level: intermediate
170e884886fSBarry Smith 
17195452b02SPatrick Sanan    Notes:
17295452b02SPatrick Sanan     Requires 2 norms and 1 inner product, but they are computed together
17311a5261eSBarry Smith        so only one parallel collective operation is needed. See `MATMFFD_WP` for a method
17411a5261eSBarry Smith        (with `KSPGMRES`) that requires NO collective operations.
175e884886fSBarry Smith 
176e884886fSBarry Smith    Formula used:
177e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
178e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
179e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
180e884886fSBarry Smith   where
181e884886fSBarry Smith      error_rel = square root of relative error in function evaluation
182e884886fSBarry Smith      umin = minimum iterate parameter
183e884886fSBarry Smith 
1841d27aa22SBarry Smith   Method taken from {cite}`dennis:83`
185e884886fSBarry Smith 
18611a5261eSBarry Smith .seealso: `MATMFFD`, `MATMFFD_WP`, `MatCreateMFFD()`, `MatCreateSNESMF()`, `MATMFFD_WP`, `MatMFFDDSSetUmin()`
187e884886fSBarry Smith M*/
MatCreateMFFD_DS(MatMFFD ctx)188d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreateMFFD_DS(MatMFFD ctx)
189d71ae5a4SJacob Faibussowitsch {
190e884886fSBarry Smith   MatMFFD_DS *hctx;
191e884886fSBarry Smith 
192e884886fSBarry Smith   PetscFunctionBegin;
193e884886fSBarry Smith   /* allocate my own private data structure */
1944dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&hctx));
195e884886fSBarry Smith   ctx->hctx = (void *)hctx;
196e884886fSBarry Smith   /* set a default for my parameter */
197e884886fSBarry Smith   hctx->umin = 1.e-6;
198e884886fSBarry Smith 
199e884886fSBarry Smith   /* set the functions I am providing */
200e884886fSBarry Smith   ctx->ops->compute        = MatMFFDCompute_DS;
201e884886fSBarry Smith   ctx->ops->destroy        = MatMFFDDestroy_DS;
202e884886fSBarry Smith   ctx->ops->view           = MatMFFDView_DS;
203e884886fSBarry Smith   ctx->ops->setfromoptions = MatMFFDSetFromOptions_DS;
204e884886fSBarry Smith 
2059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ctx->mat, "MatMFFDDSSetUmin_C", MatMFFDDSSetUmin_DS));
2063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
207e884886fSBarry Smith }
208