xref: /petsc/src/mat/impls/mffd/wp.c (revision 9895aa37ac365bac650f6bd8bf977519f7222510)
1 
2 /*MC
3      MATMFFD_WP - Implements an alternative approach for computing the differencing parameter
4         h used with the finite difference based matrix-free Jacobian.  This code
5         implements the strategy of M. Pernice and H. Walker:
6 
7       h = error_rel * sqrt(1 + ||U||) / ||a||
8 
9       Notes:
10         1) || U || does not change between linear iterations so is reused
11         2) In GMRES || a || == 1 and so does not need to ever be computed except at restart
12            when it is recomputed.
13 
14       Reference:  M. Pernice and H. F. Walker, "NITSOL: A Newton Iterative
15       Solver for Nonlinear Systems", SIAM J. Sci. Stat. Comput.", 1998,
16       vol 19, pp. 302--318.
17 
18    Options Database Keys:
19 .   -mat_mffd_compute_normu -Compute the norm of u everytime see MatMFFDWPSetComputeNormU()
20 
21 
22    Level: intermediate
23 
24    Notes: Requires no global collectives when used with GMRES
25 
26    Formula used:
27      F'(u)*a = [F(u+h*a) - F(u)]/h where
28 
29 .seealso: MATMFFD, MatCreateMFFD(), MatCreateSNESMF(), MATMFFD_DS
30 
31 M*/
32 
33 /*
34     This include file defines the data structure  MatMFFD that
35    includes information about the computation of h. It is shared by
36    all implementations that people provide.
37 
38    See snesmfjdef.c for  a full set of comments on the routines below.
39 */
40 #include <petsc-private/matimpl.h>
41 #include <../src/mat/impls/mffd/mffdimpl.h>   /*I  "petscmat.h"   I*/
42 
43 typedef struct {
44   PetscReal normUfact;                    /* previous sqrt(1.0 + || U ||) */
45   PetscBool computenormU;
46 } MatMFFD_WP;
47 
48 #undef __FUNCT__
49 #define __FUNCT__ "MatMFFDCompute_WP"
50 /*
51      MatMFFDCompute_WP - Standard PETSc code for
52    computing h with matrix-free finite differences.
53 
54   Input Parameters:
55 +   ctx - the matrix free context
56 .   U - the location at which you want the Jacobian
57 -   a - the direction you want the derivative
58 
59   Output Parameter:
60 .   h - the scale computed
61 
62 */
63 static PetscErrorCode MatMFFDCompute_WP(MatMFFD ctx,Vec U,Vec a,PetscScalar *h,PetscBool  *zeroa)
64 {
65   MatMFFD_WP     *hctx = (MatMFFD_WP*)ctx->hctx;
66   PetscReal      normU,norma;
67   PetscErrorCode ierr;
68 
69   PetscFunctionBegin;
70   if (!(ctx->count % ctx->recomputeperiod)) {
71     if (hctx->computenormU || !ctx->ncurrenth) {
72       ierr            = VecNorm(U,NORM_2,&normU);CHKERRQ(ierr);
73       hctx->normUfact = PetscSqrtReal(1.0+normU);
74     }
75     ierr = VecNorm(a,NORM_2,&norma);CHKERRQ(ierr);
76     if (norma == 0.0) {
77       *zeroa = PETSC_TRUE;
78       PetscFunctionReturn(0);
79     }
80     *zeroa = PETSC_FALSE;
81     *h     = ctx->error_rel*hctx->normUfact/norma;
82   } else {
83     *h = ctx->currenth;
84   }
85   PetscFunctionReturn(0);
86 }
87 
88 #undef __FUNCT__
89 #define __FUNCT__ "MatMFFDView_WP"
90 /*
91    MatMFFDView_WP - Prints information about this particular
92      method for computing h. Note that this does not print the general
93      information about the matrix free, that is printed by the calling
94      routine.
95 
96   Input Parameters:
97 +   ctx - the matrix free context
98 -   viewer - the PETSc viewer
99 
100 */
101 static PetscErrorCode MatMFFDView_WP(MatMFFD ctx,PetscViewer viewer)
102 {
103   MatMFFD_WP     *hctx = (MatMFFD_WP*)ctx->hctx;
104   PetscErrorCode ierr;
105   PetscBool      iascii;
106 
107   PetscFunctionBegin;
108   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
109   if (iascii) {
110     if (hctx->computenormU) {
111       ierr =  PetscViewerASCIIPrintf(viewer,"    Computes normU\n");CHKERRQ(ierr);
112     } else {
113       ierr =  PetscViewerASCIIPrintf(viewer,"    Does not compute normU\n");CHKERRQ(ierr);
114     }
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", hctx->computenormU,&hctx->computenormU,0);CHKERRQ(ierr);
137   ierr = PetscOptionsTail();CHKERRQ(ierr);
138   PetscFunctionReturn(0);
139 }
140 
141 #undef __FUNCT__
142 #define __FUNCT__ "MatMFFDDestroy_WP"
143 /*
144    MatMFFDDestroy_WP - Frees the space allocated by
145        MatCreateMFFD_WP().
146 
147   Input Parameter:
148 .  ctx - the matrix free context
149 
150    Notes: does not free the ctx, that is handled by the calling routine
151 
152 */
153 static PetscErrorCode MatMFFDDestroy_WP(MatMFFD ctx)
154 {
155   PetscErrorCode ierr;
156 
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  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  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  MatCreateMFFD_WP(MatMFFD ctx)
222 {
223   PetscErrorCode ierr;
224   MatMFFD_WP     *hctx;
225 
226   PetscFunctionBegin;
227   /* allocate my own private data structure */
228   ierr               = PetscNewLog(ctx,MatMFFD_WP,&hctx);CHKERRQ(ierr);
229   ctx->hctx          = (void*)hctx;
230   hctx->computenormU = PETSC_FALSE;
231 
232   /* set the functions I am providing */
233   ctx->ops->compute        = MatMFFDCompute_WP;
234   ctx->ops->destroy        = MatMFFDDestroy_WP;
235   ctx->ops->view           = MatMFFDView_WP;
236   ctx->ops->setfromoptions = MatMFFDSetFromOptions_WP;
237 
238   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatMFFDWPSetComputeNormU_C","MatMFFDWPSetComputeNormU_P",MatMFFDWPSetComputeNormU_P);CHKERRQ(ierr);
239   PetscFunctionReturn(0);
240 }
241 EXTERN_C_END
242 
243 
244 
245