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