xref: /petsc/src/mat/impls/mffd/wp.c (revision 009bbdc485cd9ad46be9940d3549e2dde9cdc322)
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) {ierr =  PetscViewerASCIIPrintf(viewer,"    Computes normU\n");CHKERRQ(ierr);}
111     else                    {ierr =  PetscViewerASCIIPrintf(viewer,"    Does not compute normU\n");CHKERRQ(ierr);}
112   }
113   PetscFunctionReturn(0);
114 }
115 
116 #undef __FUNCT__
117 #define __FUNCT__ "MatMFFDSetFromOptions_WP"
118 /*
119    MatMFFDSetFromOptions_WP - Looks in the options database for
120      any options appropriate for this method
121 
122   Input Parameter:
123 .  ctx - the matrix free context
124 
125 */
126 static PetscErrorCode MatMFFDSetFromOptions_WP(MatMFFD ctx)
127 {
128   PetscErrorCode ierr;
129   MatMFFD_WP     *hctx = (MatMFFD_WP*)ctx->hctx;
130 
131   PetscFunctionBegin;
132   ierr = PetscOptionsHead("Walker-Pernice options");CHKERRQ(ierr);
133     ierr = PetscOptionsBool("-mat_mffd_compute_normu","Compute the norm of u","MatMFFDWPSetComputeNormU", hctx->computenormU,&hctx->computenormU,0);CHKERRQ(ierr);
134   ierr = PetscOptionsTail();CHKERRQ(ierr);
135   PetscFunctionReturn(0);
136 }
137 
138 #undef __FUNCT__
139 #define __FUNCT__ "MatMFFDDestroy_WP"
140 /*
141    MatMFFDDestroy_WP - Frees the space allocated by
142        MatCreateMFFD_WP().
143 
144   Input Parameter:
145 .  ctx - the matrix free context
146 
147    Notes: does not free the ctx, that is handled by the calling routine
148 
149 */
150 static PetscErrorCode MatMFFDDestroy_WP(MatMFFD ctx)
151 {
152   PetscErrorCode ierr;
153 
154   PetscFunctionBegin;
155   ierr = PetscFree(ctx->hctx);CHKERRQ(ierr);
156   PetscFunctionReturn(0);
157 }
158 
159 EXTERN_C_BEGIN
160 #undef __FUNCT__
161 #define __FUNCT__ "MatMFFDWPSetComputeNormU_P"
162 PetscErrorCode  MatMFFDWPSetComputeNormU_P(Mat mat,PetscBool  flag)
163 {
164   MatMFFD     ctx = (MatMFFD)mat->data;
165   MatMFFD_WP  *hctx = (MatMFFD_WP*)ctx->hctx;
166 
167   PetscFunctionBegin;
168   hctx->computenormU = flag;
169   PetscFunctionReturn(0);
170 }
171 EXTERN_C_END
172 
173 #undef __FUNCT__
174 #define __FUNCT__ "MatMFFDWPSetComputeNormU"
175 /*@
176     MatMFFDWPSetComputeNormU - Sets whether it computes the ||U|| used by the WP
177              PETSc routine for computing h. With any Krylov solver this need only
178              be computed during the first iteration and kept for later.
179 
180   Input Parameters:
181 +   A - the matrix created with MatCreateSNESMF()
182 -   flag - PETSC_TRUE causes it to compute ||U||, PETSC_FALSE uses the previous value
183 
184   Options Database Key:
185 .   -mat_mffd_compute_normu <true,false> - true by default, false can save calculations but you
186               must be sure that ||U|| has not changed in the mean time.
187 
188   Level: advanced
189 
190   Notes:
191    See the manual page for MATMFFD_WP for a complete description of the
192    algorithm used to compute h.
193 
194 .seealso: MatMFFDSetFunctionError(), MatCreateSNESMF()
195 
196 @*/
197 PetscErrorCode  MatMFFDWPSetComputeNormU(Mat A,PetscBool  flag)
198 {
199   PetscErrorCode ierr;
200 
201   PetscFunctionBegin;
202   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
203   ierr = PetscTryMethod(A,"MatMFFDWPSetComputeNormU_C",(Mat,PetscBool),(A,flag));CHKERRQ(ierr);
204   PetscFunctionReturn(0);
205 }
206 
207 EXTERN_C_BEGIN
208 #undef __FUNCT__
209 #define __FUNCT__ "MatCreateMFFD_WP"
210 /*
211      MatCreateMFFD_WP - Standard PETSc code for
212    computing h with matrix-free finite differences.
213 
214    Input Parameter:
215 .  ctx - the matrix free context created by MatCreateMFFD()
216 
217 */
218 PetscErrorCode  MatCreateMFFD_WP(MatMFFD ctx)
219 {
220   PetscErrorCode ierr;
221   MatMFFD_WP     *hctx;
222 
223   PetscFunctionBegin;
224   /* allocate my own private data structure */
225   ierr               = PetscNewLog(ctx,MatMFFD_WP,&hctx);CHKERRQ(ierr);
226   ctx->hctx          = (void*)hctx;
227   hctx->computenormU = PETSC_FALSE;
228 
229   /* set the functions I am providing */
230   ctx->ops->compute        = MatMFFDCompute_WP;
231   ctx->ops->destroy        = MatMFFDDestroy_WP;
232   ctx->ops->view           = MatMFFDView_WP;
233   ctx->ops->setfromoptions = MatMFFDSetFromOptions_WP;
234 
235   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatMFFDWPSetComputeNormU_C","MatMFFDWPSetComputeNormU_P",MatMFFDWPSetComputeNormU_P);CHKERRQ(ierr);
236   PetscFunctionReturn(0);
237 }
238 EXTERN_C_END
239 
240 
241 
242