xref: /petsc/src/mat/impls/mffd/wp.c (revision e884886f040fd140b63a22c479f657cf835a1983)
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 "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   PetscTruth computenorma,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,PetscTruth *zeroa)
65 {
66   MatMFFD_WP    *hctx = (MatMFFD_WP*)ctx->hctx;
67   PetscReal      normU,norma = 1.0;
68   PetscErrorCode ierr;
69 
70   PetscFunctionBegin;
71 
72   if (!(ctx->count % ctx->recomputeperiod)) {
73     if (hctx->computenorma && (hctx->computenormU || !ctx->ncurrenth)) {
74       ierr = VecNormBegin(U,NORM_2,&normU);CHKERRQ(ierr);
75       ierr = VecNormBegin(a,NORM_2,&norma);CHKERRQ(ierr);
76       ierr = VecNormEnd(U,NORM_2,&normU);CHKERRQ(ierr);
77       ierr = VecNormEnd(a,NORM_2,&norma);CHKERRQ(ierr);
78       hctx->normUfact = sqrt(1.0+normU);
79     } else if (hctx->computenormU || !ctx->ncurrenth) {
80       ierr = VecNorm(U,NORM_2,&normU);CHKERRQ(ierr);
81       hctx->normUfact = sqrt(1.0+normU);
82     } else if (hctx->computenorma) {
83       ierr = VecNorm(a,NORM_2,&norma);CHKERRQ(ierr);
84     }
85     if (norma == 0.0) {
86       *zeroa = PETSC_TRUE;
87       PetscFunctionReturn(0);
88     }
89     *zeroa = PETSC_FALSE;
90     *h = ctx->error_rel*hctx->normUfact/norma;
91   } else {
92     *h = ctx->currenth;
93   }
94   PetscFunctionReturn(0);
95 }
96 
97 #undef __FUNCT__
98 #define __FUNCT__ "MatMFFDView_WP"
99 /*
100    MatMFFDView_WP - Prints information about this particular
101      method for computing h. Note that this does not print the general
102      information about the matrix free, that is printed by the calling
103      routine.
104 
105   Input Parameters:
106 +   ctx - the matrix free context
107 -   viewer - the PETSc viewer
108 
109 */
110 static PetscErrorCode MatMFFDView_WP(MatMFFD ctx,PetscViewer viewer)
111 {
112   MatMFFD_WP     *hctx = (MatMFFD_WP *)ctx->hctx;
113   PetscErrorCode ierr;
114   PetscTruth     iascii;
115 
116   PetscFunctionBegin;
117   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
118   if (iascii) {
119     if (hctx->computenorma){ierr = PetscViewerASCIIPrintf(viewer,"    Computes normA\n");CHKERRQ(ierr);}
120     else                   {ierr =  PetscViewerASCIIPrintf(viewer,"    Does not compute normA\n");CHKERRQ(ierr);}
121     if (hctx->computenormU){ierr =  PetscViewerASCIIPrintf(viewer,"    Computes normU\n");CHKERRQ(ierr);}
122     else                   {ierr =  PetscViewerASCIIPrintf(viewer,"    Does not compute normU\n");CHKERRQ(ierr);}
123   } else {
124     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES matrix-free WP",((PetscObject)viewer)->type_name);
125   }
126   PetscFunctionReturn(0);
127 }
128 
129 #undef __FUNCT__
130 #define __FUNCT__ "MatMFFDSetFromOptions_WP"
131 /*
132    MatMFFDSetFromOptions_WP - Looks in the options database for
133      any options appropriate for this method
134 
135   Input Parameter:
136 .  ctx - the matrix free context
137 
138 */
139 static PetscErrorCode MatMFFDSetFromOptions_WP(MatMFFD ctx)
140 {
141   PetscErrorCode ierr;
142   MatMFFD_WP     *hctx = (MatMFFD_WP*)ctx->hctx;
143 
144   PetscFunctionBegin;
145   ierr = PetscOptionsHead("Walker-Pernice options");CHKERRQ(ierr);
146     ierr = PetscOptionsTruth("-mat_mffd_compute_normu","Compute the norm of u","MatMFFDWPSetComputeNormU",
147                           hctx->computenorma,&hctx->computenorma,0);CHKERRQ(ierr);
148   ierr = PetscOptionsTail();CHKERRQ(ierr);
149   PetscFunctionReturn(0);
150 }
151 
152 #undef __FUNCT__
153 #define __FUNCT__ "MatMFFDDestroy_WP"
154 /*
155    MatMFFDDestroy_WP - Frees the space allocated by
156        MatMFFDCreate_WP().
157 
158   Input Parameter:
159 .  ctx - the matrix free context
160 
161    Notes: does not free the ctx, that is handled by the calling routine
162 
163 */
164 static PetscErrorCode MatMFFDDestroy_WP(MatMFFD ctx)
165 {
166   PetscErrorCode ierr;
167   PetscFunctionBegin;
168   ierr = PetscFree(ctx->hctx);CHKERRQ(ierr);
169   PetscFunctionReturn(0);
170 }
171 
172 EXTERN_C_BEGIN
173 #undef __FUNCT__
174 #define __FUNCT__ "MatMFFDWPSetComputeNormU_P"
175 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU_P(Mat mat,PetscTruth flag)
176 {
177   MatMFFD     ctx = (MatMFFD)mat->data;
178   MatMFFD_WP  *hctx = (MatMFFD_WP*)ctx->hctx;
179 
180   PetscFunctionBegin;
181   hctx->computenormU = flag;
182   PetscFunctionReturn(0);
183 }
184 EXTERN_C_END
185 
186 #undef __FUNCT__
187 #define __FUNCT__ "MatMFFDWPSetComputeNormU"
188 /*@
189     MatMFFDWPSetComputeNormU - Sets whether it computes the ||U|| used by the WP
190              PETSc routine for computing h. With any Krylov solver this need only
191              be computed during the first iteration and kept for later.
192 
193   Input Parameters:
194 +   A - the matrix created with MatCreateSNESMF()
195 -   flag - PETSC_TRUE causes it to compute ||U||, PETSC_FALSE uses the previous value
196 
197   Options Database Key:
198 .   -mat_mffd_compute_normu <true,false> - true by default, false can save calculations but you
199               must be sure that ||U|| has not changed in the mean time.
200 
201   Level: advanced
202 
203   Notes:
204    See the manual page for MATMFFD_WP for a complete description of the
205    algorithm used to compute h.
206 
207 .seealso: MatMFFDSetFunctionError(), MatCreateSNESMF()
208 
209 @*/
210 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat A,PetscTruth flag)
211 {
212   PetscErrorCode ierr,(*f)(Mat,PetscTruth);
213 
214   PetscFunctionBegin;
215   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
216   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMFFDWPSetComputeNormU_C",(void (**)(void))&f);CHKERRQ(ierr);
217   if (f) {
218     ierr = (*f)(A,flag);CHKERRQ(ierr);
219   }
220   PetscFunctionReturn(0);
221 }
222 
223 EXTERN_C_BEGIN
224 #undef __FUNCT__
225 #define __FUNCT__ "MatMFFDCreate_WP"
226 /*
227      MatMFFDCreate_WP - Standard PETSc code for
228    computing h with matrix-free finite differences.
229 
230    Input Parameter:
231 .  ctx - the matrix free context created by MatMFFDCreate()
232 
233 */
234 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCreate_WP(MatMFFD ctx)
235 {
236   PetscErrorCode ierr;
237   MatMFFD_WP     *hctx;
238 
239   PetscFunctionBegin;
240 
241   /* allocate my own private data structure */
242   ierr               = PetscNew(MatMFFD_WP,&hctx);CHKERRQ(ierr);
243   ctx->hctx          = (void*)hctx;
244   hctx->computenormU = PETSC_FALSE;
245   hctx->computenorma = PETSC_TRUE;
246 
247   /* set the functions I am providing */
248   ctx->ops->compute        = MatMFFDCompute_WP;
249   ctx->ops->destroy        = MatMFFDDestroy_WP;
250   ctx->ops->view           = MatMFFDView_WP;
251   ctx->ops->setfromoptions = MatMFFDSetFromOptions_WP;
252 
253   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatMFFDWPSetComputeNormU_C",
254                             "MatMFFDWPSetComputeNormU_P",
255                              MatMFFDWPSetComputeNormU_P);CHKERRQ(ierr);
256 
257   PetscFunctionReturn(0);
258 }
259 EXTERN_C_END
260 
261 
262 
263