xref: /petsc/src/mat/impls/mffd/wp.c (revision e884886f040fd140b63a22c479f657cf835a1983)
1*e884886fSBarry Smith #define PETSCMAT_DLL
2*e884886fSBarry Smith 
3*e884886fSBarry Smith /*MC
4*e884886fSBarry Smith      MATMFFD_WP - Implements an alternative approach for computing the differencing parameter
5*e884886fSBarry Smith         h used with the finite difference based matrix-free Jacobian.  This code
6*e884886fSBarry Smith         implements the strategy of M. Pernice and H. Walker:
7*e884886fSBarry Smith 
8*e884886fSBarry Smith       h = error_rel * sqrt(1 + ||U||) / ||a||
9*e884886fSBarry Smith 
10*e884886fSBarry Smith       Notes:
11*e884886fSBarry Smith         1) || U || does not change between linear iterations so is reused
12*e884886fSBarry Smith         2) In GMRES || a || == 1 and so does not need to ever be computed except at restart
13*e884886fSBarry Smith            when it is recomputed.
14*e884886fSBarry Smith 
15*e884886fSBarry Smith       Reference:  M. Pernice and H. F. Walker, "NITSOL: A Newton Iterative
16*e884886fSBarry Smith       Solver for Nonlinear Systems", SIAM J. Sci. Stat. Comput.", 1998,
17*e884886fSBarry Smith       vol 19, pp. 302--318.
18*e884886fSBarry Smith 
19*e884886fSBarry Smith    Options Database Keys:
20*e884886fSBarry Smith .   -mat_mffd_compute_normu -Compute the norm of u everytime see MatMFFDWPSetComputeNormU()
21*e884886fSBarry Smith 
22*e884886fSBarry Smith 
23*e884886fSBarry Smith    Level: intermediate
24*e884886fSBarry Smith 
25*e884886fSBarry Smith    Notes: Requires no global collectives when used with GMRES
26*e884886fSBarry Smith 
27*e884886fSBarry Smith    Formula used:
28*e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
29*e884886fSBarry Smith 
30*e884886fSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatCreateSNESMF(), MATMFFD_DS
31*e884886fSBarry Smith 
32*e884886fSBarry Smith M*/
33*e884886fSBarry Smith 
34*e884886fSBarry Smith /*
35*e884886fSBarry Smith     This include file defines the data structure  MatMFFD that
36*e884886fSBarry Smith    includes information about the computation of h. It is shared by
37*e884886fSBarry Smith    all implementations that people provide.
38*e884886fSBarry Smith 
39*e884886fSBarry Smith    See snesmfjdef.c for  a full set of comments on the routines below.
40*e884886fSBarry Smith */
41*e884886fSBarry Smith #include "include/private/matimpl.h"
42*e884886fSBarry Smith #include "src/mat/impls/mffd/mffdimpl.h"   /*I  "petscmat.h"   I*/
43*e884886fSBarry Smith 
44*e884886fSBarry Smith typedef struct {
45*e884886fSBarry Smith   PetscReal  normUfact;                   /* previous sqrt(1.0 + || U ||) */
46*e884886fSBarry Smith   PetscTruth computenorma,computenormU;
47*e884886fSBarry Smith } MatMFFD_WP;
48*e884886fSBarry Smith 
49*e884886fSBarry Smith #undef __FUNCT__
50*e884886fSBarry Smith #define __FUNCT__ "MatMFFDCompute_WP"
51*e884886fSBarry Smith /*
52*e884886fSBarry Smith      MatMFFDCompute_WP - Standard PETSc code for
53*e884886fSBarry Smith    computing h with matrix-free finite differences.
54*e884886fSBarry Smith 
55*e884886fSBarry Smith   Input Parameters:
56*e884886fSBarry Smith +   ctx - the matrix free context
57*e884886fSBarry Smith .   U - the location at which you want the Jacobian
58*e884886fSBarry Smith -   a - the direction you want the derivative
59*e884886fSBarry Smith 
60*e884886fSBarry Smith   Output Parameter:
61*e884886fSBarry Smith .   h - the scale computed
62*e884886fSBarry Smith 
63*e884886fSBarry Smith */
64*e884886fSBarry Smith static PetscErrorCode MatMFFDCompute_WP(MatMFFD ctx,Vec U,Vec a,PetscScalar *h,PetscTruth *zeroa)
65*e884886fSBarry Smith {
66*e884886fSBarry Smith   MatMFFD_WP    *hctx = (MatMFFD_WP*)ctx->hctx;
67*e884886fSBarry Smith   PetscReal      normU,norma = 1.0;
68*e884886fSBarry Smith   PetscErrorCode ierr;
69*e884886fSBarry Smith 
70*e884886fSBarry Smith   PetscFunctionBegin;
71*e884886fSBarry Smith 
72*e884886fSBarry Smith   if (!(ctx->count % ctx->recomputeperiod)) {
73*e884886fSBarry Smith     if (hctx->computenorma && (hctx->computenormU || !ctx->ncurrenth)) {
74*e884886fSBarry Smith       ierr = VecNormBegin(U,NORM_2,&normU);CHKERRQ(ierr);
75*e884886fSBarry Smith       ierr = VecNormBegin(a,NORM_2,&norma);CHKERRQ(ierr);
76*e884886fSBarry Smith       ierr = VecNormEnd(U,NORM_2,&normU);CHKERRQ(ierr);
77*e884886fSBarry Smith       ierr = VecNormEnd(a,NORM_2,&norma);CHKERRQ(ierr);
78*e884886fSBarry Smith       hctx->normUfact = sqrt(1.0+normU);
79*e884886fSBarry Smith     } else if (hctx->computenormU || !ctx->ncurrenth) {
80*e884886fSBarry Smith       ierr = VecNorm(U,NORM_2,&normU);CHKERRQ(ierr);
81*e884886fSBarry Smith       hctx->normUfact = sqrt(1.0+normU);
82*e884886fSBarry Smith     } else if (hctx->computenorma) {
83*e884886fSBarry Smith       ierr = VecNorm(a,NORM_2,&norma);CHKERRQ(ierr);
84*e884886fSBarry Smith     }
85*e884886fSBarry Smith     if (norma == 0.0) {
86*e884886fSBarry Smith       *zeroa = PETSC_TRUE;
87*e884886fSBarry Smith       PetscFunctionReturn(0);
88*e884886fSBarry Smith     }
89*e884886fSBarry Smith     *zeroa = PETSC_FALSE;
90*e884886fSBarry Smith     *h = ctx->error_rel*hctx->normUfact/norma;
91*e884886fSBarry Smith   } else {
92*e884886fSBarry Smith     *h = ctx->currenth;
93*e884886fSBarry Smith   }
94*e884886fSBarry Smith   PetscFunctionReturn(0);
95*e884886fSBarry Smith }
96*e884886fSBarry Smith 
97*e884886fSBarry Smith #undef __FUNCT__
98*e884886fSBarry Smith #define __FUNCT__ "MatMFFDView_WP"
99*e884886fSBarry Smith /*
100*e884886fSBarry Smith    MatMFFDView_WP - Prints information about this particular
101*e884886fSBarry Smith      method for computing h. Note that this does not print the general
102*e884886fSBarry Smith      information about the matrix free, that is printed by the calling
103*e884886fSBarry Smith      routine.
104*e884886fSBarry Smith 
105*e884886fSBarry Smith   Input Parameters:
106*e884886fSBarry Smith +   ctx - the matrix free context
107*e884886fSBarry Smith -   viewer - the PETSc viewer
108*e884886fSBarry Smith 
109*e884886fSBarry Smith */
110*e884886fSBarry Smith static PetscErrorCode MatMFFDView_WP(MatMFFD ctx,PetscViewer viewer)
111*e884886fSBarry Smith {
112*e884886fSBarry Smith   MatMFFD_WP     *hctx = (MatMFFD_WP *)ctx->hctx;
113*e884886fSBarry Smith   PetscErrorCode ierr;
114*e884886fSBarry Smith   PetscTruth     iascii;
115*e884886fSBarry Smith 
116*e884886fSBarry Smith   PetscFunctionBegin;
117*e884886fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
118*e884886fSBarry Smith   if (iascii) {
119*e884886fSBarry Smith     if (hctx->computenorma){ierr = PetscViewerASCIIPrintf(viewer,"    Computes normA\n");CHKERRQ(ierr);}
120*e884886fSBarry Smith     else                   {ierr =  PetscViewerASCIIPrintf(viewer,"    Does not compute normA\n");CHKERRQ(ierr);}
121*e884886fSBarry Smith     if (hctx->computenormU){ierr =  PetscViewerASCIIPrintf(viewer,"    Computes normU\n");CHKERRQ(ierr);}
122*e884886fSBarry Smith     else                   {ierr =  PetscViewerASCIIPrintf(viewer,"    Does not compute normU\n");CHKERRQ(ierr);}
123*e884886fSBarry Smith   } else {
124*e884886fSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES matrix-free WP",((PetscObject)viewer)->type_name);
125*e884886fSBarry Smith   }
126*e884886fSBarry Smith   PetscFunctionReturn(0);
127*e884886fSBarry Smith }
128*e884886fSBarry Smith 
129*e884886fSBarry Smith #undef __FUNCT__
130*e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFromOptions_WP"
131*e884886fSBarry Smith /*
132*e884886fSBarry Smith    MatMFFDSetFromOptions_WP - Looks in the options database for
133*e884886fSBarry Smith      any options appropriate for this method
134*e884886fSBarry Smith 
135*e884886fSBarry Smith   Input Parameter:
136*e884886fSBarry Smith .  ctx - the matrix free context
137*e884886fSBarry Smith 
138*e884886fSBarry Smith */
139*e884886fSBarry Smith static PetscErrorCode MatMFFDSetFromOptions_WP(MatMFFD ctx)
140*e884886fSBarry Smith {
141*e884886fSBarry Smith   PetscErrorCode ierr;
142*e884886fSBarry Smith   MatMFFD_WP     *hctx = (MatMFFD_WP*)ctx->hctx;
143*e884886fSBarry Smith 
144*e884886fSBarry Smith   PetscFunctionBegin;
145*e884886fSBarry Smith   ierr = PetscOptionsHead("Walker-Pernice options");CHKERRQ(ierr);
146*e884886fSBarry Smith     ierr = PetscOptionsTruth("-mat_mffd_compute_normu","Compute the norm of u","MatMFFDWPSetComputeNormU",
147*e884886fSBarry Smith                           hctx->computenorma,&hctx->computenorma,0);CHKERRQ(ierr);
148*e884886fSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
149*e884886fSBarry Smith   PetscFunctionReturn(0);
150*e884886fSBarry Smith }
151*e884886fSBarry Smith 
152*e884886fSBarry Smith #undef __FUNCT__
153*e884886fSBarry Smith #define __FUNCT__ "MatMFFDDestroy_WP"
154*e884886fSBarry Smith /*
155*e884886fSBarry Smith    MatMFFDDestroy_WP - Frees the space allocated by
156*e884886fSBarry Smith        MatMFFDCreate_WP().
157*e884886fSBarry Smith 
158*e884886fSBarry Smith   Input Parameter:
159*e884886fSBarry Smith .  ctx - the matrix free context
160*e884886fSBarry Smith 
161*e884886fSBarry Smith    Notes: does not free the ctx, that is handled by the calling routine
162*e884886fSBarry Smith 
163*e884886fSBarry Smith */
164*e884886fSBarry Smith static PetscErrorCode MatMFFDDestroy_WP(MatMFFD ctx)
165*e884886fSBarry Smith {
166*e884886fSBarry Smith   PetscErrorCode ierr;
167*e884886fSBarry Smith   PetscFunctionBegin;
168*e884886fSBarry Smith   ierr = PetscFree(ctx->hctx);CHKERRQ(ierr);
169*e884886fSBarry Smith   PetscFunctionReturn(0);
170*e884886fSBarry Smith }
171*e884886fSBarry Smith 
172*e884886fSBarry Smith EXTERN_C_BEGIN
173*e884886fSBarry Smith #undef __FUNCT__
174*e884886fSBarry Smith #define __FUNCT__ "MatMFFDWPSetComputeNormU_P"
175*e884886fSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU_P(Mat mat,PetscTruth flag)
176*e884886fSBarry Smith {
177*e884886fSBarry Smith   MatMFFD     ctx = (MatMFFD)mat->data;
178*e884886fSBarry Smith   MatMFFD_WP  *hctx = (MatMFFD_WP*)ctx->hctx;
179*e884886fSBarry Smith 
180*e884886fSBarry Smith   PetscFunctionBegin;
181*e884886fSBarry Smith   hctx->computenormU = flag;
182*e884886fSBarry Smith   PetscFunctionReturn(0);
183*e884886fSBarry Smith }
184*e884886fSBarry Smith EXTERN_C_END
185*e884886fSBarry Smith 
186*e884886fSBarry Smith #undef __FUNCT__
187*e884886fSBarry Smith #define __FUNCT__ "MatMFFDWPSetComputeNormU"
188*e884886fSBarry Smith /*@
189*e884886fSBarry Smith     MatMFFDWPSetComputeNormU - Sets whether it computes the ||U|| used by the WP
190*e884886fSBarry Smith              PETSc routine for computing h. With any Krylov solver this need only
191*e884886fSBarry Smith              be computed during the first iteration and kept for later.
192*e884886fSBarry Smith 
193*e884886fSBarry Smith   Input Parameters:
194*e884886fSBarry Smith +   A - the matrix created with MatCreateSNESMF()
195*e884886fSBarry Smith -   flag - PETSC_TRUE causes it to compute ||U||, PETSC_FALSE uses the previous value
196*e884886fSBarry Smith 
197*e884886fSBarry Smith   Options Database Key:
198*e884886fSBarry Smith .   -mat_mffd_compute_normu <true,false> - true by default, false can save calculations but you
199*e884886fSBarry Smith               must be sure that ||U|| has not changed in the mean time.
200*e884886fSBarry Smith 
201*e884886fSBarry Smith   Level: advanced
202*e884886fSBarry Smith 
203*e884886fSBarry Smith   Notes:
204*e884886fSBarry Smith    See the manual page for MATMFFD_WP for a complete description of the
205*e884886fSBarry Smith    algorithm used to compute h.
206*e884886fSBarry Smith 
207*e884886fSBarry Smith .seealso: MatMFFDSetFunctionError(), MatCreateSNESMF()
208*e884886fSBarry Smith 
209*e884886fSBarry Smith @*/
210*e884886fSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat A,PetscTruth flag)
211*e884886fSBarry Smith {
212*e884886fSBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscTruth);
213*e884886fSBarry Smith 
214*e884886fSBarry Smith   PetscFunctionBegin;
215*e884886fSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
216*e884886fSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMFFDWPSetComputeNormU_C",(void (**)(void))&f);CHKERRQ(ierr);
217*e884886fSBarry Smith   if (f) {
218*e884886fSBarry Smith     ierr = (*f)(A,flag);CHKERRQ(ierr);
219*e884886fSBarry Smith   }
220*e884886fSBarry Smith   PetscFunctionReturn(0);
221*e884886fSBarry Smith }
222*e884886fSBarry Smith 
223*e884886fSBarry Smith EXTERN_C_BEGIN
224*e884886fSBarry Smith #undef __FUNCT__
225*e884886fSBarry Smith #define __FUNCT__ "MatMFFDCreate_WP"
226*e884886fSBarry Smith /*
227*e884886fSBarry Smith      MatMFFDCreate_WP - Standard PETSc code for
228*e884886fSBarry Smith    computing h with matrix-free finite differences.
229*e884886fSBarry Smith 
230*e884886fSBarry Smith    Input Parameter:
231*e884886fSBarry Smith .  ctx - the matrix free context created by MatMFFDCreate()
232*e884886fSBarry Smith 
233*e884886fSBarry Smith */
234*e884886fSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCreate_WP(MatMFFD ctx)
235*e884886fSBarry Smith {
236*e884886fSBarry Smith   PetscErrorCode ierr;
237*e884886fSBarry Smith   MatMFFD_WP     *hctx;
238*e884886fSBarry Smith 
239*e884886fSBarry Smith   PetscFunctionBegin;
240*e884886fSBarry Smith 
241*e884886fSBarry Smith   /* allocate my own private data structure */
242*e884886fSBarry Smith   ierr               = PetscNew(MatMFFD_WP,&hctx);CHKERRQ(ierr);
243*e884886fSBarry Smith   ctx->hctx          = (void*)hctx;
244*e884886fSBarry Smith   hctx->computenormU = PETSC_FALSE;
245*e884886fSBarry Smith   hctx->computenorma = PETSC_TRUE;
246*e884886fSBarry Smith 
247*e884886fSBarry Smith   /* set the functions I am providing */
248*e884886fSBarry Smith   ctx->ops->compute        = MatMFFDCompute_WP;
249*e884886fSBarry Smith   ctx->ops->destroy        = MatMFFDDestroy_WP;
250*e884886fSBarry Smith   ctx->ops->view           = MatMFFDView_WP;
251*e884886fSBarry Smith   ctx->ops->setfromoptions = MatMFFDSetFromOptions_WP;
252*e884886fSBarry Smith 
253*e884886fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatMFFDWPSetComputeNormU_C",
254*e884886fSBarry Smith                             "MatMFFDWPSetComputeNormU_P",
255*e884886fSBarry Smith                              MatMFFDWPSetComputeNormU_P);CHKERRQ(ierr);
256*e884886fSBarry Smith 
257*e884886fSBarry Smith   PetscFunctionReturn(0);
258*e884886fSBarry Smith }
259*e884886fSBarry Smith EXTERN_C_END
260*e884886fSBarry Smith 
261*e884886fSBarry Smith 
262*e884886fSBarry Smith 
263