xref: /petsc/src/mat/impls/mffd/mffddef.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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 */
61ace3abfcSBarry Smith static PetscErrorCode MatMFFDCompute_DS(MatMFFD ctx,Vec U,Vec a,PetscScalar *h,PetscBool  *zeroa)
62e884886fSBarry Smith {
63e884886fSBarry Smith   MatMFFD_DS     *hctx = (MatMFFD_DS*)ctx->hctx;
64e884886fSBarry Smith   PetscReal      nrm,sum,umin = hctx->umin;
65e884886fSBarry Smith   PetscScalar    dot;
66e884886fSBarry Smith 
67e884886fSBarry Smith   PetscFunctionBegin;
68e884886fSBarry Smith   if (!(ctx->count % ctx->recomputeperiod)) {
69e884886fSBarry Smith     /*
70e884886fSBarry Smith      This algorithm requires 2 norms and 1 inner product. Rather than
71e884886fSBarry Smith      use directly the VecNorm() and VecDot() routines (and thus have
72e884886fSBarry Smith      three separate collective operations, we use the VecxxxBegin/End() routines
73e884886fSBarry Smith     */
745f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDotBegin(U,a,&dot));
755f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNormBegin(a,NORM_1,&sum));
765f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNormBegin(a,NORM_2,&nrm));
775f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDotEnd(U,a,&dot));
785f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNormEnd(a,NORM_1,&sum));
795f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNormEnd(a,NORM_2,&nrm));
80e884886fSBarry Smith 
81e884886fSBarry Smith     if (nrm == 0.0) {
82e884886fSBarry Smith       *zeroa = PETSC_TRUE;
83e884886fSBarry Smith       PetscFunctionReturn(0);
84e884886fSBarry Smith     }
85e884886fSBarry Smith     *zeroa = PETSC_FALSE;
86e884886fSBarry Smith 
87e884886fSBarry Smith     /*
88e884886fSBarry Smith       Safeguard for step sizes that are "too small"
89e884886fSBarry Smith     */
90e884886fSBarry Smith     if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum;
91e884886fSBarry Smith     else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum;
92e884886fSBarry Smith     *h = ctx->error_rel*dot/(nrm*nrm);
932c71b3e2SJacob Faibussowitsch     PetscCheckFalse(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);
94e884886fSBarry Smith   } else {
95e884886fSBarry Smith     *h = ctx->currenth;
96e884886fSBarry Smith   }
97e884886fSBarry Smith   ctx->count++;
98e884886fSBarry Smith   PetscFunctionReturn(0);
99e884886fSBarry Smith }
100e884886fSBarry Smith 
101e884886fSBarry Smith /*
102e884886fSBarry Smith    MatMFFDView_DS - Prints information about this particular
103e884886fSBarry Smith    method for computing h. Note that this does not print the general
104e884886fSBarry Smith    information about the matrix-free method, as such info is printed
105e884886fSBarry Smith    by the calling routine.
106e884886fSBarry Smith 
107e884886fSBarry Smith    Input Parameters:
108e884886fSBarry Smith +  ctx - the matrix free context
109e884886fSBarry Smith -  viewer - the PETSc viewer
110e884886fSBarry Smith */
111e884886fSBarry Smith static PetscErrorCode MatMFFDView_DS(MatMFFD ctx,PetscViewer viewer)
112e884886fSBarry Smith {
113e884886fSBarry Smith   MatMFFD_DS     *hctx = (MatMFFD_DS*)ctx->hctx;
114ace3abfcSBarry Smith   PetscBool      iascii;
115e884886fSBarry Smith 
116e884886fSBarry Smith   PetscFunctionBegin;
117e884886fSBarry Smith   /*
118e884886fSBarry Smith      Currently this only handles the ascii file viewers, others
119e884886fSBarry Smith      could be added, but for this type of object other viewers
120e884886fSBarry Smith      make less sense
121e884886fSBarry Smith   */
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
123e884886fSBarry Smith   if (iascii) {
1245f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"    umin=%g (minimum iterate parameter)\n",(double)hctx->umin));
125e884886fSBarry Smith   }
126e884886fSBarry Smith   PetscFunctionReturn(0);
127e884886fSBarry Smith }
128e884886fSBarry Smith 
129e884886fSBarry Smith /*
130e884886fSBarry Smith    MatMFFDSetFromOptions_DS - Looks in the options database for
131e884886fSBarry Smith    any options appropriate for this method.
132e884886fSBarry Smith 
133e884886fSBarry Smith    Input Parameter:
134e884886fSBarry Smith .  ctx - the matrix free context
135e884886fSBarry Smith 
136e884886fSBarry Smith */
1374416b707SBarry Smith static PetscErrorCode MatMFFDSetFromOptions_DS(PetscOptionItems *PetscOptionsObject,MatMFFD ctx)
138e884886fSBarry Smith {
139e884886fSBarry Smith   MatMFFD_DS     *hctx = (MatMFFD_DS*)ctx->hctx;
140e884886fSBarry Smith 
141e884886fSBarry Smith   PetscFunctionBegin;
1425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHead(PetscOptionsObject,"Finite difference matrix free parameters"));
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-mat_mffd_umin","umin","MatMFFDDSSetUmin",hctx->umin,&hctx->umin,NULL));
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsTail());
145e884886fSBarry Smith   PetscFunctionReturn(0);
146e884886fSBarry Smith }
147e884886fSBarry Smith 
148e884886fSBarry Smith /*
149e884886fSBarry Smith    MatMFFDDestroy_DS - Frees the space allocated by
1501d0fab5eSBarry Smith    MatCreateMFFD_DS().
151e884886fSBarry Smith 
152e884886fSBarry Smith    Input Parameter:
153e884886fSBarry Smith .  ctx - the matrix free context
154e884886fSBarry Smith 
155e884886fSBarry Smith    Notes:
156e884886fSBarry Smith    Does not free the ctx, that is handled by the calling routine
157e884886fSBarry Smith */
158e884886fSBarry Smith static PetscErrorCode MatMFFDDestroy_DS(MatMFFD ctx)
159e884886fSBarry Smith {
160e884886fSBarry Smith   PetscFunctionBegin;
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(ctx->hctx));
162e884886fSBarry Smith   PetscFunctionReturn(0);
163e884886fSBarry Smith }
164e884886fSBarry Smith 
165e884886fSBarry Smith /*
166e884886fSBarry Smith    The following two routines use the PetscObjectCompose() and PetscObjectQuery()
167e884886fSBarry Smith    mechanism to allow the user to change the Umin parameter used in this method.
168e884886fSBarry Smith */
169304ddee0SJed Brown PetscErrorCode MatMFFDDSSetUmin_DS(Mat mat,PetscReal umin)
170e884886fSBarry Smith {
171440c47e6SSatish Balay   MatMFFD        ctx=NULL;
172e884886fSBarry Smith   MatMFFD_DS     *hctx;
173e884886fSBarry Smith 
174e884886fSBarry Smith   PetscFunctionBegin;
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(mat,&ctx));
176*28b400f6SJacob Faibussowitsch   PetscCheck(ctx,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MatMFFDDSSetUmin() attached to non-shell matrix");
177e884886fSBarry Smith   hctx       = (MatMFFD_DS*)ctx->hctx;
178e884886fSBarry Smith   hctx->umin = umin;
179e884886fSBarry Smith   PetscFunctionReturn(0);
180e884886fSBarry Smith }
181e884886fSBarry Smith 
182e884886fSBarry Smith /*@
183e884886fSBarry Smith     MatMFFDDSSetUmin - Sets the "umin" parameter used by the
184e884886fSBarry Smith     PETSc routine for computing the differencing parameter, h, which is used
185e884886fSBarry Smith     for matrix-free Jacobian-vector products.
186e884886fSBarry Smith 
187e884886fSBarry Smith    Input Parameters:
188e884886fSBarry Smith +  A - the matrix created with MatCreateSNESMF()
189e884886fSBarry Smith -  umin - the parameter
190e884886fSBarry Smith 
191e884886fSBarry Smith    Level: advanced
192e884886fSBarry Smith 
193e884886fSBarry Smith    Notes:
194e884886fSBarry Smith    See the manual page for MatCreateSNESMF() for a complete description of the
195e884886fSBarry Smith    algorithm used to compute h.
196e884886fSBarry Smith 
197e884886fSBarry Smith .seealso: MatMFFDSetFunctionError(), MatCreateSNESMF()
198e884886fSBarry Smith 
199e884886fSBarry Smith @*/
2007087cfbeSBarry Smith PetscErrorCode  MatMFFDDSSetUmin(Mat A,PetscReal umin)
201e884886fSBarry Smith {
202e884886fSBarry Smith   PetscFunctionBegin;
2030700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
2045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTryMethod(A,"MatMFFDDSSetUmin_C",(Mat,PetscReal),(A,umin)));
205e884886fSBarry Smith   PetscFunctionReturn(0);
206e884886fSBarry Smith }
207e884886fSBarry Smith 
208e884886fSBarry Smith /*MC
209e884886fSBarry Smith      MATMFFD_DS - the code for compute the "h" used in the finite difference
210e884886fSBarry Smith             matrix-free matrix vector product.  This code
211e884886fSBarry Smith         implements the strategy in Dennis and Schnabel, "Numerical Methods for Unconstrained
212e884886fSBarry Smith         Optimization and Nonlinear Equations".
213e884886fSBarry Smith 
214e884886fSBarry Smith    Options Database Keys:
21567b8a455SSatish Balay .  -mat_mffd_umin <umin> - see MatMFFDDSSetUmin()
216e884886fSBarry Smith 
217e884886fSBarry Smith    Level: intermediate
218e884886fSBarry Smith 
21995452b02SPatrick Sanan    Notes:
22095452b02SPatrick Sanan     Requires 2 norms and 1 inner product, but they are computed together
221e884886fSBarry Smith        so only one parallel collective operation is needed. See MATMFFD_WP for a method
222e884886fSBarry Smith        (with GMRES) that requires NO collective operations.
223e884886fSBarry Smith 
224e884886fSBarry Smith    Formula used:
225e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
226e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
227e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
228e884886fSBarry Smith  where
229e884886fSBarry Smith      error_rel = square root of relative error in function evaluation
230e884886fSBarry Smith      umin = minimum iterate parameter
231e884886fSBarry Smith 
232e884886fSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatCreateSNESMF(), MATMFFD_WP, MatMFFDDSSetUmin()
233e884886fSBarry Smith 
234e884886fSBarry Smith M*/
2358cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD_DS(MatMFFD ctx)
236e884886fSBarry Smith {
237e884886fSBarry Smith   MatMFFD_DS     *hctx;
238e884886fSBarry Smith 
239e884886fSBarry Smith   PetscFunctionBegin;
240e884886fSBarry Smith   /* allocate my own private data structure */
2415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNewLog(ctx,&hctx));
242e884886fSBarry Smith   ctx->hctx = (void*)hctx;
243e884886fSBarry Smith   /* set a default for my parameter */
244e884886fSBarry Smith   hctx->umin = 1.e-6;
245e884886fSBarry Smith 
246e884886fSBarry Smith   /* set the functions I am providing */
247e884886fSBarry Smith   ctx->ops->compute        = MatMFFDCompute_DS;
248e884886fSBarry Smith   ctx->ops->destroy        = MatMFFDDestroy_DS;
249e884886fSBarry Smith   ctx->ops->view           = MatMFFDView_DS;
250e884886fSBarry Smith   ctx->ops->setfromoptions = MatMFFDSetFromOptions_DS;
251e884886fSBarry Smith 
2525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)ctx->mat,"MatMFFDDSSetUmin_C",MatMFFDDSSetUmin_DS));
253e884886fSBarry Smith   PetscFunctionReturn(0);
254e884886fSBarry Smith }
255