xref: /petsc/src/mat/impls/mffd/wp.c (revision 09573ac72a50d3e7ecd55a2b7f0ef28450cd0a8b)
1 #define PETSCMAT_DLL
2 
3 /*MC
4      MATMFFD_WP - Implements an alternative approach for computing the differencing parameter
5         h used with the finite difference based matrix-free Jacobian.  This code
6         implements the strategy of M. Pernice and H. Walker:
7 
8       h = error_rel * sqrt(1 + ||U||) / ||a||
9 
10       Notes:
11         1) || U || does not change between linear iterations so is reused
12         2) In GMRES || a || == 1 and so does not need to ever be computed except at restart
13            when it is recomputed.
14 
15       Reference:  M. Pernice and H. F. Walker, "NITSOL: A Newton Iterative
16       Solver for Nonlinear Systems", SIAM J. Sci. Stat. Comput.", 1998,
17       vol 19, pp. 302--318.
18 
19    Options Database Keys:
20 .   -mat_mffd_compute_normu -Compute the norm of u everytime see MatMFFDWPSetComputeNormU()
21 
22 
23    Level: intermediate
24 
25    Notes: Requires no global collectives when used with GMRES
26 
27    Formula used:
28      F'(u)*a = [F(u+h*a) - F(u)]/h where
29 
30 .seealso: MATMFFD, MatCreateMFFD(), MatCreateSNESMF(), MATMFFD_DS
31 
32 M*/
33 
34 /*
35     This include file defines the data structure  MatMFFD that
36    includes information about the computation of h. It is shared by
37    all implementations that people provide.
38 
39    See snesmfjdef.c for  a full set of comments on the routines below.
40 */
41 #include "private/matimpl.h"
42 #include "../src/mat/impls/mffd/mffdimpl.h"   /*I  "petscmat.h"   I*/
43 
44 typedef struct {
45   PetscReal  normUfact;                   /* previous sqrt(1.0 + || U ||) */
46   PetscBool  computenormU;
47 } MatMFFD_WP;
48 
49 #undef __FUNCT__
50 #define __FUNCT__ "MatMFFDCompute_WP"
51 /*
52      MatMFFDCompute_WP - Standard PETSc code for
53    computing h 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   Output Parameter:
61 .   h - the scale computed
62 
63 */
64 static PetscErrorCode MatMFFDCompute_WP(MatMFFD ctx,Vec U,Vec a,PetscScalar *h,PetscBool  *zeroa)
65 {
66   MatMFFD_WP    *hctx = (MatMFFD_WP*)ctx->hctx;
67   PetscReal      normU,norma;
68   PetscErrorCode ierr;
69 
70   PetscFunctionBegin;
71   if (!(ctx->count % ctx->recomputeperiod)) {
72     if (hctx->computenormU || !ctx->ncurrenth) {
73       ierr = VecNorm(U,NORM_2,&normU);CHKERRQ(ierr);
74       hctx->normUfact = sqrt(1.0+normU);
75     }
76     ierr = VecNorm(a,NORM_2,&norma);CHKERRQ(ierr);
77     if (norma == 0.0) {
78       *zeroa = PETSC_TRUE;
79       PetscFunctionReturn(0);
80     }
81     *zeroa = PETSC_FALSE;
82     *h = ctx->error_rel*hctx->normUfact/norma;
83   } else {
84     *h = ctx->currenth;
85   }
86   PetscFunctionReturn(0);
87 }
88 
89 #undef __FUNCT__
90 #define __FUNCT__ "MatMFFDView_WP"
91 /*
92    MatMFFDView_WP - Prints information about this particular
93      method for computing h. Note that this does not print the general
94      information about the matrix free, that is printed by the calling
95      routine.
96 
97   Input Parameters:
98 +   ctx - the matrix free context
99 -   viewer - the PETSc viewer
100 
101 */
102 static PetscErrorCode MatMFFDView_WP(MatMFFD ctx,PetscViewer viewer)
103 {
104   MatMFFD_WP     *hctx = (MatMFFD_WP *)ctx->hctx;
105   PetscErrorCode ierr;
106   PetscBool      iascii;
107 
108   PetscFunctionBegin;
109   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
110   if (iascii) {
111     if (hctx->computenormU){ierr =  PetscViewerASCIIPrintf(viewer,"    Computes normU\n");CHKERRQ(ierr);}
112     else                   {ierr =  PetscViewerASCIIPrintf(viewer,"    Does not compute normU\n");CHKERRQ(ierr);}
113   } else {
114     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES matrix-free WP",((PetscObject)viewer)->type_name);
115   }
116   PetscFunctionReturn(0);
117 }
118 
119 #undef __FUNCT__
120 #define __FUNCT__ "MatMFFDSetFromOptions_WP"
121 /*
122    MatMFFDSetFromOptions_WP - Looks in the options database for
123      any options appropriate for this method
124 
125   Input Parameter:
126 .  ctx - the matrix free context
127 
128 */
129 static PetscErrorCode MatMFFDSetFromOptions_WP(MatMFFD ctx)
130 {
131   PetscErrorCode ierr;
132   MatMFFD_WP     *hctx = (MatMFFD_WP*)ctx->hctx;
133 
134   PetscFunctionBegin;
135   ierr = PetscOptionsHead("Walker-Pernice options");CHKERRQ(ierr);
136     ierr = PetscOptionsBool("-mat_mffd_compute_normu","Compute the norm of u","MatMFFDWPSetComputeNormU",
137                           hctx->computenormU,&hctx->computenormU,0);CHKERRQ(ierr);
138   ierr = PetscOptionsTail();CHKERRQ(ierr);
139   PetscFunctionReturn(0);
140 }
141 
142 #undef __FUNCT__
143 #define __FUNCT__ "MatMFFDDestroy_WP"
144 /*
145    MatMFFDDestroy_WP - Frees the space allocated by
146        MatCreateMFFD_WP().
147 
148   Input Parameter:
149 .  ctx - the matrix free context
150 
151    Notes: does not free the ctx, that is handled by the calling routine
152 
153 */
154 static PetscErrorCode MatMFFDDestroy_WP(MatMFFD ctx)
155 {
156   PetscErrorCode ierr;
157   PetscFunctionBegin;
158   ierr = PetscFree(ctx->hctx);CHKERRQ(ierr);
159   PetscFunctionReturn(0);
160 }
161 
162 EXTERN_C_BEGIN
163 #undef __FUNCT__
164 #define __FUNCT__ "MatMFFDWPSetComputeNormU_P"
165 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU_P(Mat mat,PetscBool  flag)
166 {
167   MatMFFD     ctx = (MatMFFD)mat->data;
168   MatMFFD_WP  *hctx = (MatMFFD_WP*)ctx->hctx;
169 
170   PetscFunctionBegin;
171   hctx->computenormU = flag;
172   PetscFunctionReturn(0);
173 }
174 EXTERN_C_END
175 
176 #undef __FUNCT__
177 #define __FUNCT__ "MatMFFDWPSetComputeNormU"
178 /*@
179     MatMFFDWPSetComputeNormU - Sets whether it computes the ||U|| used by the WP
180              PETSc routine for computing h. With any Krylov solver this need only
181              be computed during the first iteration and kept for later.
182 
183   Input Parameters:
184 +   A - the matrix created with MatCreateSNESMF()
185 -   flag - PETSC_TRUE causes it to compute ||U||, PETSC_FALSE uses the previous value
186 
187   Options Database Key:
188 .   -mat_mffd_compute_normu <true,false> - true by default, false can save calculations but you
189               must be sure that ||U|| has not changed in the mean time.
190 
191   Level: advanced
192 
193   Notes:
194    See the manual page for MATMFFD_WP for a complete description of the
195    algorithm used to compute h.
196 
197 .seealso: MatMFFDSetFunctionError(), MatCreateSNESMF()
198 
199 @*/
200 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat A,PetscBool  flag)
201 {
202   PetscErrorCode ierr;
203 
204   PetscFunctionBegin;
205   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
206   ierr = PetscTryMethod(A,"MatMFFDWPSetComputeNormU_C",(Mat,PetscBool),(A,flag));CHKERRQ(ierr);
207   PetscFunctionReturn(0);
208 }
209 
210 EXTERN_C_BEGIN
211 #undef __FUNCT__
212 #define __FUNCT__ "MatCreateMFFD_WP"
213 /*
214      MatCreateMFFD_WP - Standard PETSc code for
215    computing h with matrix-free finite differences.
216 
217    Input Parameter:
218 .  ctx - the matrix free context created by MatCreateMFFD()
219 
220 */
221 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD_WP(MatMFFD ctx)
222 {
223   PetscErrorCode ierr;
224   MatMFFD_WP     *hctx;
225 
226   PetscFunctionBegin;
227 
228   /* allocate my own private data structure */
229   ierr               = PetscNewLog(ctx,MatMFFD_WP,&hctx);CHKERRQ(ierr);
230   ctx->hctx          = (void*)hctx;
231   hctx->computenormU = PETSC_FALSE;
232 
233   /* set the functions I am providing */
234   ctx->ops->compute        = MatMFFDCompute_WP;
235   ctx->ops->destroy        = MatMFFDDestroy_WP;
236   ctx->ops->view           = MatMFFDView_WP;
237   ctx->ops->setfromoptions = MatMFFDSetFromOptions_WP;
238 
239   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatMFFDWPSetComputeNormU_C",
240                             "MatMFFDWPSetComputeNormU_P",
241                              MatMFFDWPSetComputeNormU_P);CHKERRQ(ierr);
242 
243   PetscFunctionReturn(0);
244 }
245 EXTERN_C_END
246 
247 
248 
249