xref: /petsc/src/mat/impls/mffd/wp.c (revision 4e278199b78715991f5c71ebbd945c1489263e6c)
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    Level: intermediate
22 
23    Notes:
24     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 /*
49      MatMFFDCompute_WP - Standard PETSc code for
50    computing h with matrix-free finite differences.
51 
52   Input Parameters:
53 +   ctx - the matrix free context
54 .   U - the location at which you want the Jacobian
55 -   a - the direction you want the derivative
56 
57   Output Parameter:
58 .   h - the scale computed
59 
60 */
61 static PetscErrorCode MatMFFDCompute_WP(MatMFFD ctx,Vec U,Vec a,PetscScalar *h,PetscBool  *zeroa)
62 {
63   MatMFFD_WP     *hctx = (MatMFFD_WP*)ctx->hctx;
64   PetscReal      normU,norma;
65   PetscErrorCode ierr;
66 
67   PetscFunctionBegin;
68   if (!(ctx->count % ctx->recomputeperiod)) {
69     if (hctx->computenormU || !ctx->ncurrenth) {
70       ierr            = VecNorm(U,NORM_2,&normU);CHKERRQ(ierr);
71       hctx->normUfact = PetscSqrtReal(1.0+normU);
72     }
73     ierr = VecNorm(a,NORM_2,&norma);CHKERRQ(ierr);
74     if (norma == 0.0) {
75       *zeroa = PETSC_TRUE;
76       PetscFunctionReturn(0);
77     }
78     *zeroa = PETSC_FALSE;
79     *h     = ctx->error_rel*hctx->normUfact/norma;
80   } else {
81     *h = ctx->currenth;
82   }
83   PetscFunctionReturn(0);
84 }
85 
86 /*
87    MatMFFDView_WP - Prints information about this particular
88      method for computing h. Note that this does not print the general
89      information about the matrix free, that is printed by the calling
90      routine.
91 
92   Input Parameters:
93 +   ctx - the matrix free context
94 -   viewer - the PETSc viewer
95 
96 */
97 static PetscErrorCode MatMFFDView_WP(MatMFFD ctx,PetscViewer viewer)
98 {
99   MatMFFD_WP     *hctx = (MatMFFD_WP*)ctx->hctx;
100   PetscErrorCode ierr;
101   PetscBool      iascii;
102 
103   PetscFunctionBegin;
104   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
105   if (iascii) {
106     if (hctx->computenormU) {
107       ierr =  PetscViewerASCIIPrintf(viewer,"    Computes normU\n");CHKERRQ(ierr);
108     } else {
109       ierr =  PetscViewerASCIIPrintf(viewer,"    Does not compute normU\n");CHKERRQ(ierr);
110     }
111   }
112   PetscFunctionReturn(0);
113 }
114 
115 /*
116    MatMFFDSetFromOptions_WP - Looks in the options database for
117      any options appropriate for this method
118 
119   Input Parameter:
120 .  ctx - the matrix free context
121 
122 */
123 static PetscErrorCode MatMFFDSetFromOptions_WP(PetscOptionItems *PetscOptionsObject,MatMFFD ctx)
124 {
125   PetscErrorCode ierr;
126   MatMFFD_WP     *hctx = (MatMFFD_WP*)ctx->hctx;
127 
128   PetscFunctionBegin;
129   ierr = PetscOptionsHead(PetscOptionsObject,"Walker-Pernice options");CHKERRQ(ierr);
130   ierr = PetscOptionsBool("-mat_mffd_compute_normu","Compute the norm of u","MatMFFDWPSetComputeNormU", hctx->computenormU,&hctx->computenormU,NULL);CHKERRQ(ierr);
131   ierr = PetscOptionsTail();CHKERRQ(ierr);
132   PetscFunctionReturn(0);
133 }
134 
135 /*
136    MatMFFDDestroy_WP - Frees the space allocated by
137        MatCreateMFFD_WP().
138 
139   Input Parameter:
140 .  ctx - the matrix free context
141 
142    Notes:
143     does not free the ctx, that is handled by the calling routine
144 
145 */
146 static PetscErrorCode MatMFFDDestroy_WP(MatMFFD ctx)
147 {
148   PetscErrorCode ierr;
149 
150   PetscFunctionBegin;
151   ierr = PetscFree(ctx->hctx);CHKERRQ(ierr);
152   PetscFunctionReturn(0);
153 }
154 
155 PetscErrorCode  MatMFFDWPSetComputeNormU_P(Mat mat,PetscBool flag)
156 {
157   MatMFFD    ctx   = (MatMFFD)mat->data;
158   MatMFFD_WP *hctx = (MatMFFD_WP*)ctx->hctx;
159 
160   PetscFunctionBegin;
161   hctx->computenormU = flag;
162   PetscFunctionReturn(0);
163 }
164 
165 /*@
166     MatMFFDWPSetComputeNormU - Sets whether it computes the ||U|| used by the WP
167              PETSc routine for computing h. With any Krylov solver this need only
168              be computed during the first iteration and kept for later.
169 
170   Input Parameters:
171 +   A - the matrix created with MatCreateSNESMF()
172 -   flag - PETSC_TRUE causes it to compute ||U||, PETSC_FALSE uses the previous value
173 
174   Options Database Key:
175 .   -mat_mffd_compute_normu <true,false> - true by default, false can save calculations but you
176               must be sure that ||U|| has not changed in the mean time.
177 
178   Level: advanced
179 
180   Notes:
181    See the manual page for MATMFFD_WP for a complete description of the
182    algorithm used to compute h.
183 
184 .seealso: MatMFFDSetFunctionError(), MatCreateSNESMF()
185 
186 @*/
187 PetscErrorCode  MatMFFDWPSetComputeNormU(Mat A,PetscBool flag)
188 {
189   PetscErrorCode ierr;
190 
191   PetscFunctionBegin;
192   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
193   ierr = PetscTryMethod(A,"MatMFFDWPSetComputeNormU_C",(Mat,PetscBool),(A,flag));CHKERRQ(ierr);
194   PetscFunctionReturn(0);
195 }
196 
197 /*
198      MatCreateMFFD_WP - Standard PETSc code for
199    computing h with matrix-free finite differences.
200 
201    Input Parameter:
202 .  ctx - the matrix free context created by MatCreateMFFD()
203 
204 */
205 PETSC_EXTERN PetscErrorCode MatCreateMFFD_WP(MatMFFD ctx)
206 {
207   PetscErrorCode ierr;
208   MatMFFD_WP     *hctx;
209 
210   PetscFunctionBegin;
211   /* allocate my own private data structure */
212   ierr               = PetscNewLog(ctx,&hctx);CHKERRQ(ierr);
213   ctx->hctx          = (void*)hctx;
214   hctx->computenormU = PETSC_FALSE;
215 
216   /* set the functions I am providing */
217   ctx->ops->compute        = MatMFFDCompute_WP;
218   ctx->ops->destroy        = MatMFFDDestroy_WP;
219   ctx->ops->view           = MatMFFDView_WP;
220   ctx->ops->setfromoptions = MatMFFDSetFromOptions_WP;
221 
222   ierr = PetscObjectComposeFunction((PetscObject)ctx->mat,"MatMFFDWPSetComputeNormU_C",MatMFFDWPSetComputeNormU_P);CHKERRQ(ierr);
223   PetscFunctionReturn(0);
224 }
225 
226