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