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