xref: /petsc/src/mat/impls/mffd/mffddef.c (revision e6e75211d226c622f451867f53ce5d558649ff4f)
1 
2 /*
3   Implements the DS PETSc approach for computing the h
4   parameter used with the finite difference based matrix-free
5   Jacobian-vector products.
6 
7   To make your own: clone this file and modify for your needs.
8 
9   Mandatory functions:
10   -------------------
11       MatMFFDCompute_  - for a given point and direction computes h
12 
13       MatCreateMFFD _ - fills in the MatMFFD data structure
14                            for this particular implementation
15 
16 
17    Optional functions:
18    -------------------
19       MatMFFDView_ - prints information about the parameters being used.
20                        This is called when SNESView() or -snes_view is used.
21 
22       MatMFFDSetFromOptions_ - checks the options database for options that
23                                apply to this method.
24 
25       MatMFFDDestroy_ - frees any space allocated by the routines above
26 
27 */
28 
29 /*
30     This include file defines the data structure  MatMFFD that
31    includes information about the computation of h. It is shared by
32    all implementations that people provide
33 */
34 #include <petsc/private/matimpl.h>
35 #include <../src/mat/impls/mffd/mffdimpl.h>   /*I  "petscmat.h"   I*/
36 
37 /*
38       The  method has one parameter that is used to
39    "cutoff" very small values. This is stored in a data structure
40    that is only visible to this file. If your method has no parameters
41    it can omit this, if it has several simply reorganize the data structure.
42    The data structure is "hung-off" the MatMFFD data structure in
43    the void *hctx; field.
44 */
45 typedef struct {
46   PetscReal umin;          /* minimum allowable u'a value relative to |u|_1 */
47 } MatMFFD_DS;
48 
49 #undef __FUNCT__
50 #define __FUNCT__ "MatMFFDCompute_DS"
51 /*
52    MatMFFDCompute_DS - Standard PETSc code for computing the
53    differencing paramter (h) for use with matrix-free finite differences.
54 
55    Input Parameters:
56 +  ctx - the matrix free context
57 .  U - the location at which you want the Jacobian
58 -  a - the direction you want the derivative
59 
60 
61    Output Parameter:
62 .  h - the scale computed
63 
64 */
65 static PetscErrorCode MatMFFDCompute_DS(MatMFFD ctx,Vec U,Vec a,PetscScalar *h,PetscBool  *zeroa)
66 {
67   MatMFFD_DS     *hctx = (MatMFFD_DS*)ctx->hctx;
68   PetscReal      nrm,sum,umin = hctx->umin;
69   PetscScalar    dot;
70   PetscErrorCode ierr;
71 
72   PetscFunctionBegin;
73   if (!(ctx->count % ctx->recomputeperiod)) {
74     /*
75      This algorithm requires 2 norms and 1 inner product. Rather than
76      use directly the VecNorm() and VecDot() routines (and thus have
77      three separate collective operations, we use the VecxxxBegin/End() routines
78     */
79     ierr = VecDotBegin(U,a,&dot);CHKERRQ(ierr);
80     ierr = VecNormBegin(a,NORM_1,&sum);CHKERRQ(ierr);
81     ierr = VecNormBegin(a,NORM_2,&nrm);CHKERRQ(ierr);
82     ierr = VecDotEnd(U,a,&dot);CHKERRQ(ierr);
83     ierr = VecNormEnd(a,NORM_1,&sum);CHKERRQ(ierr);
84     ierr = VecNormEnd(a,NORM_2,&nrm);CHKERRQ(ierr);
85 
86     if (nrm == 0.0) {
87       *zeroa = PETSC_TRUE;
88       PetscFunctionReturn(0);
89     }
90     *zeroa = PETSC_FALSE;
91 
92     /*
93       Safeguard for step sizes that are "too small"
94     */
95     if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum;
96     else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum;
97     *h = ctx->error_rel*dot/(nrm*nrm);
98   } else {
99     *h = ctx->currenth;
100   }
101   if (*h != *h) SETERRQ3(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);
102   ctx->count++;
103   PetscFunctionReturn(0);
104 }
105 
106 #undef __FUNCT__
107 #define __FUNCT__ "MatMFFDView_DS"
108 /*
109    MatMFFDView_DS - Prints information about this particular
110    method for computing h. Note that this does not print the general
111    information about the matrix-free method, as such info is printed
112    by the calling routine.
113 
114    Input Parameters:
115 +  ctx - the matrix free context
116 -  viewer - the PETSc viewer
117 */
118 static PetscErrorCode MatMFFDView_DS(MatMFFD ctx,PetscViewer viewer)
119 {
120   MatMFFD_DS     *hctx = (MatMFFD_DS*)ctx->hctx;
121   PetscErrorCode ierr;
122   PetscBool      iascii;
123 
124   PetscFunctionBegin;
125   /*
126      Currently this only handles the ascii file viewers, others
127      could be added, but for this type of object other viewers
128      make less sense
129   */
130   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
131   if (iascii) {
132     ierr = PetscViewerASCIIPrintf(viewer,"    umin=%g (minimum iterate parameter)\n",(double)hctx->umin);CHKERRQ(ierr);
133   }
134   PetscFunctionReturn(0);
135 }
136 
137 #undef __FUNCT__
138 #define __FUNCT__ "MatMFFDSetFromOptions_DS"
139 /*
140    MatMFFDSetFromOptions_DS - Looks in the options database for
141    any options appropriate for this method.
142 
143    Input Parameter:
144 .  ctx - the matrix free context
145 
146 */
147 static PetscErrorCode MatMFFDSetFromOptions_DS(PetscOptions *PetscOptionsObject,MatMFFD ctx)
148 {
149   PetscErrorCode ierr;
150   MatMFFD_DS     *hctx = (MatMFFD_DS*)ctx->hctx;
151 
152   PetscFunctionBegin;
153   ierr = PetscOptionsHead(PetscOptionsObject,"Finite difference matrix free parameters");CHKERRQ(ierr);
154   ierr = PetscOptionsReal("-mat_mffd_umin","umin","MatMFFDDSSetUmin",hctx->umin,&hctx->umin,0);CHKERRQ(ierr);
155   ierr = PetscOptionsTail();CHKERRQ(ierr);
156   PetscFunctionReturn(0);
157 }
158 
159 #undef __FUNCT__
160 #define __FUNCT__ "MatMFFDDestroy_DS"
161 /*
162    MatMFFDDestroy_DS - Frees the space allocated by
163    MatCreateMFFD_DS().
164 
165    Input Parameter:
166 .  ctx - the matrix free context
167 
168    Notes:
169    Does not free the ctx, that is handled by the calling routine
170 */
171 static PetscErrorCode MatMFFDDestroy_DS(MatMFFD ctx)
172 {
173   PetscErrorCode ierr;
174 
175   PetscFunctionBegin;
176   ierr = PetscFree(ctx->hctx);CHKERRQ(ierr);
177   PetscFunctionReturn(0);
178 }
179 
180 #undef __FUNCT__
181 #define __FUNCT__ "MatMFFDDSSetUmin_DS"
182 /*
183    The following two routines use the PetscObjectCompose() and PetscObjectQuery()
184    mechanism to allow the user to change the Umin parameter used in this method.
185 */
186 PetscErrorCode MatMFFDDSSetUmin_DS(Mat mat,PetscReal umin)
187 {
188   MatMFFD    ctx = (MatMFFD)mat->data;
189   MatMFFD_DS *hctx;
190 
191   PetscFunctionBegin;
192   if (!ctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MatMFFDDSSetUmin() attached to non-shell matrix");
193   hctx       = (MatMFFD_DS*)ctx->hctx;
194   hctx->umin = umin;
195   PetscFunctionReturn(0);
196 }
197 
198 #undef __FUNCT__
199 #define __FUNCT__ "MatMFFDDSSetUmin"
200 /*@
201     MatMFFDDSSetUmin - Sets the "umin" parameter used by the
202     PETSc routine for computing the differencing parameter, h, which is used
203     for matrix-free Jacobian-vector products.
204 
205    Input Parameters:
206 +  A - the matrix created with MatCreateSNESMF()
207 -  umin - the parameter
208 
209    Level: advanced
210 
211    Notes:
212    See the manual page for MatCreateSNESMF() for a complete description of the
213    algorithm used to compute h.
214 
215 .seealso: MatMFFDSetFunctionError(), MatCreateSNESMF()
216 
217 @*/
218 PetscErrorCode  MatMFFDDSSetUmin(Mat A,PetscReal umin)
219 {
220   PetscErrorCode ierr;
221 
222   PetscFunctionBegin;
223   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
224   ierr = PetscTryMethod(A,"MatMFFDDSSetUmin_C",(Mat,PetscReal),(A,umin));CHKERRQ(ierr);
225   PetscFunctionReturn(0);
226 }
227 
228 /*MC
229      MATMFFD_DS - the code for compute the "h" used in the finite difference
230             matrix-free matrix vector product.  This code
231         implements the strategy in Dennis and Schnabel, "Numerical Methods for Unconstrained
232         Optimization and Nonlinear Equations".
233 
234    Options Database Keys:
235 .  -mat_mffd_umin <umin> see MatMFFDDSSetUmin()
236 
237    Level: intermediate
238 
239    Notes: Requires 2 norms and 1 inner product, but they are computed together
240        so only one parallel collective operation is needed. See MATMFFD_WP for a method
241        (with GMRES) that requires NO collective operations.
242 
243    Formula used:
244      F'(u)*a = [F(u+h*a) - F(u)]/h where
245      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
246        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
247  where
248      error_rel = square root of relative error in function evaluation
249      umin = minimum iterate parameter
250 
251 .seealso: MATMFFD, MatCreateMFFD(), MatCreateSNESMF(), MATMFFD_WP, MatMFFDDSSetUmin()
252 
253 M*/
254 #undef __FUNCT__
255 #define __FUNCT__ "MatCreateMFFD_DS"
256 PETSC_EXTERN PetscErrorCode MatCreateMFFD_DS(MatMFFD ctx)
257 {
258   MatMFFD_DS     *hctx;
259   PetscErrorCode ierr;
260 
261   PetscFunctionBegin;
262   /* allocate my own private data structure */
263   ierr      = PetscNewLog(ctx,&hctx);CHKERRQ(ierr);
264   ctx->hctx = (void*)hctx;
265   /* set a default for my parameter */
266   hctx->umin = 1.e-6;
267 
268   /* set the functions I am providing */
269   ctx->ops->compute        = MatMFFDCompute_DS;
270   ctx->ops->destroy        = MatMFFDDestroy_DS;
271   ctx->ops->view           = MatMFFDView_DS;
272   ctx->ops->setfromoptions = MatMFFDSetFromOptions_DS;
273 
274   ierr = PetscObjectComposeFunction((PetscObject)ctx->mat,"MatMFFDDSSetUmin_C",MatMFFDDSSetUmin_DS);CHKERRQ(ierr);
275   PetscFunctionReturn(0);
276 }
277 
278 
279 
280 
281 
282 
283 
284