xref: /petsc/src/mat/impls/mffd/wp.c (revision 2e956fe4fc852fabc23b437482e1fb7b77fddb0d)
1e884886fSBarry Smith 
2e884886fSBarry Smith /*MC
3e884886fSBarry Smith      MATMFFD_WP - Implements an alternative approach for computing the differencing parameter
4e884886fSBarry Smith         h used with the finite difference based matrix-free Jacobian.  This code
5e884886fSBarry Smith         implements the strategy of M. Pernice and H. Walker:
6e884886fSBarry Smith 
7e884886fSBarry Smith       h = error_rel * sqrt(1 + ||U||) / ||a||
8e884886fSBarry Smith 
9e884886fSBarry Smith       Notes:
10e884886fSBarry Smith         1) || U || does not change between linear iterations so is reused
11e884886fSBarry Smith         2) In GMRES || a || == 1 and so does not need to ever be computed except at restart
12e884886fSBarry Smith            when it is recomputed.
13e884886fSBarry Smith 
14e884886fSBarry Smith       Reference:  M. Pernice and H. F. Walker, "NITSOL: A Newton Iterative
15e884886fSBarry Smith       Solver for Nonlinear Systems", SIAM J. Sci. Stat. Comput.", 1998,
16e884886fSBarry Smith       vol 19, pp. 302--318.
17e884886fSBarry Smith 
18e884886fSBarry Smith    Options Database Keys:
19e884886fSBarry Smith .   -mat_mffd_compute_normu -Compute the norm of u every time see MatMFFDWPSetComputeNormU()
20e884886fSBarry Smith 
21e884886fSBarry Smith    Level: intermediate
22e884886fSBarry Smith 
2395452b02SPatrick Sanan    Notes:
2495452b02SPatrick Sanan     Requires no global collectives when used with GMRES
25e884886fSBarry Smith 
26e884886fSBarry Smith    Formula used:
27e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
28e884886fSBarry Smith 
29db781477SPatrick Sanan .seealso: `MATMFFD`, `MatCreateMFFD()`, `MatCreateSNESMF()`, `MATMFFD_DS`
30e884886fSBarry Smith 
31e884886fSBarry Smith M*/
32e884886fSBarry Smith 
33e884886fSBarry Smith /*
34e884886fSBarry Smith     This include file defines the data structure  MatMFFD that
35e884886fSBarry Smith    includes information about the computation of h. It is shared by
36e884886fSBarry Smith    all implementations that people provide.
37e884886fSBarry Smith 
38e884886fSBarry Smith    See snesmfjdef.c for  a full set of comments on the routines below.
39e884886fSBarry Smith */
40af0996ceSBarry Smith #include <petsc/private/matimpl.h>
41c6db04a5SJed Brown #include <../src/mat/impls/mffd/mffdimpl.h>   /*I  "petscmat.h"   I*/
42e884886fSBarry Smith 
43e884886fSBarry Smith typedef struct {
44e884886fSBarry Smith   PetscReal normUfact;                    /* previous sqrt(1.0 + || U ||) */
45ace3abfcSBarry Smith   PetscBool computenormU;
46e884886fSBarry Smith } MatMFFD_WP;
47e884886fSBarry Smith 
48e884886fSBarry Smith /*
49e884886fSBarry Smith      MatMFFDCompute_WP - Standard PETSc code for
50e884886fSBarry Smith    computing h with matrix-free finite differences.
51e884886fSBarry Smith 
52e884886fSBarry Smith   Input Parameters:
53e884886fSBarry Smith +   ctx - the matrix free context
54e884886fSBarry Smith .   U - the location at which you want the Jacobian
55e884886fSBarry Smith -   a - the direction you want the derivative
56e884886fSBarry Smith 
57e884886fSBarry Smith   Output Parameter:
58e884886fSBarry Smith .   h - the scale computed
59e884886fSBarry Smith 
60e884886fSBarry Smith */
61ace3abfcSBarry Smith static PetscErrorCode MatMFFDCompute_WP(MatMFFD ctx,Vec U,Vec a,PetscScalar *h,PetscBool  *zeroa)
62e884886fSBarry Smith {
63e884886fSBarry Smith   MatMFFD_WP     *hctx = (MatMFFD_WP*)ctx->hctx;
6451a79602SBarry Smith   PetscReal      normU,norma;
65e884886fSBarry Smith 
66e884886fSBarry Smith   PetscFunctionBegin;
67e884886fSBarry Smith   if (!(ctx->count % ctx->recomputeperiod)) {
6851a79602SBarry Smith     if (hctx->computenormU || !ctx->ncurrenth) {
699566063dSJacob Faibussowitsch       PetscCall(VecNorm(U,NORM_2,&normU));
708f1a2a5eSBarry Smith       hctx->normUfact = PetscSqrtReal(1.0+normU);
71e884886fSBarry Smith     }
729566063dSJacob Faibussowitsch     PetscCall(VecNorm(a,NORM_2,&norma));
73e884886fSBarry Smith     if (norma == 0.0) {
74e884886fSBarry Smith       *zeroa = PETSC_TRUE;
75e884886fSBarry Smith       PetscFunctionReturn(0);
76e884886fSBarry Smith     }
77e884886fSBarry Smith     *zeroa = PETSC_FALSE;
78e884886fSBarry Smith     *h     = ctx->error_rel*hctx->normUfact/norma;
79e884886fSBarry Smith   } else {
80e884886fSBarry Smith     *h = ctx->currenth;
81e884886fSBarry Smith   }
82e884886fSBarry Smith   PetscFunctionReturn(0);
83e884886fSBarry Smith }
84e884886fSBarry Smith 
85e884886fSBarry Smith /*
86e884886fSBarry Smith    MatMFFDView_WP - Prints information about this particular
87e884886fSBarry Smith      method for computing h. Note that this does not print the general
88e884886fSBarry Smith      information about the matrix free, that is printed by the calling
89e884886fSBarry Smith      routine.
90e884886fSBarry Smith 
91e884886fSBarry Smith   Input Parameters:
92e884886fSBarry Smith +   ctx - the matrix free context
93e884886fSBarry Smith -   viewer - the PETSc viewer
94e884886fSBarry Smith 
95e884886fSBarry Smith */
96e884886fSBarry Smith static PetscErrorCode MatMFFDView_WP(MatMFFD ctx,PetscViewer viewer)
97e884886fSBarry Smith {
98e884886fSBarry Smith   MatMFFD_WP     *hctx = (MatMFFD_WP*)ctx->hctx;
99ace3abfcSBarry Smith   PetscBool      iascii;
100e884886fSBarry Smith 
101e884886fSBarry Smith   PetscFunctionBegin;
1029566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
103e884886fSBarry Smith   if (iascii) {
1042205254eSKarl Rupp     if (hctx->computenormU) {
1059566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"    Computes normU\n"));
1062205254eSKarl Rupp     } else {
1079566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"    Does not compute normU\n"));
1082205254eSKarl Rupp     }
10911aeaf0aSBarry Smith   }
110e884886fSBarry Smith   PetscFunctionReturn(0);
111e884886fSBarry Smith }
112e884886fSBarry Smith 
113e884886fSBarry Smith /*
114e884886fSBarry Smith    MatMFFDSetFromOptions_WP - Looks in the options database for
115e884886fSBarry Smith      any options appropriate for this method
116e884886fSBarry Smith 
117e884886fSBarry Smith   Input Parameter:
118e884886fSBarry Smith .  ctx - the matrix free context
119e884886fSBarry Smith 
120e884886fSBarry Smith */
1214416b707SBarry Smith static PetscErrorCode MatMFFDSetFromOptions_WP(PetscOptionItems *PetscOptionsObject,MatMFFD ctx)
122e884886fSBarry Smith {
123e884886fSBarry Smith   MatMFFD_WP     *hctx = (MatMFFD_WP*)ctx->hctx;
124e884886fSBarry Smith 
125e884886fSBarry Smith   PetscFunctionBegin;
126d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"Walker-Pernice options");
1279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_mffd_compute_normu","Compute the norm of u","MatMFFDWPSetComputeNormU", hctx->computenormU,&hctx->computenormU,NULL));
128d0609cedSBarry Smith   PetscOptionsHeadEnd();
129e884886fSBarry Smith   PetscFunctionReturn(0);
130e884886fSBarry Smith }
131e884886fSBarry Smith 
132e884886fSBarry Smith /*
133e884886fSBarry Smith    MatMFFDDestroy_WP - Frees the space allocated by
1341d0fab5eSBarry Smith        MatCreateMFFD_WP().
135e884886fSBarry Smith 
136e884886fSBarry Smith   Input Parameter:
137e884886fSBarry Smith .  ctx - the matrix free context
138e884886fSBarry Smith 
13995452b02SPatrick Sanan    Notes:
14095452b02SPatrick Sanan     does not free the ctx, that is handled by the calling routine
141e884886fSBarry Smith 
142e884886fSBarry Smith */
143e884886fSBarry Smith static PetscErrorCode MatMFFDDestroy_WP(MatMFFD ctx)
144e884886fSBarry Smith {
145e884886fSBarry Smith   PetscFunctionBegin;
146*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ctx->mat,"MatMFFDWPSetComputeNormU_C",NULL));
1479566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx->hctx));
148e884886fSBarry Smith   PetscFunctionReturn(0);
149e884886fSBarry Smith }
150e884886fSBarry Smith 
1517087cfbeSBarry Smith PetscErrorCode  MatMFFDWPSetComputeNormU_P(Mat mat,PetscBool flag)
152e884886fSBarry Smith {
153e884886fSBarry Smith   MatMFFD    ctx   = (MatMFFD)mat->data;
154e884886fSBarry Smith   MatMFFD_WP *hctx = (MatMFFD_WP*)ctx->hctx;
155e884886fSBarry Smith 
156e884886fSBarry Smith   PetscFunctionBegin;
157e884886fSBarry Smith   hctx->computenormU = flag;
158e884886fSBarry Smith   PetscFunctionReturn(0);
159e884886fSBarry Smith }
160e884886fSBarry Smith 
161e884886fSBarry Smith /*@
162e884886fSBarry Smith     MatMFFDWPSetComputeNormU - Sets whether it computes the ||U|| used by the WP
163e884886fSBarry Smith              PETSc routine for computing h. With any Krylov solver this need only
164e884886fSBarry Smith              be computed during the first iteration and kept for later.
165e884886fSBarry Smith 
166e884886fSBarry Smith   Input Parameters:
167e884886fSBarry Smith +   A - the matrix created with MatCreateSNESMF()
168e884886fSBarry Smith -   flag - PETSC_TRUE causes it to compute ||U||, PETSC_FALSE uses the previous value
169e884886fSBarry Smith 
170e884886fSBarry Smith   Options Database Key:
171e884886fSBarry Smith .   -mat_mffd_compute_normu <true,false> - true by default, false can save calculations but you
172e884886fSBarry Smith               must be sure that ||U|| has not changed in the mean time.
173e884886fSBarry Smith 
174e884886fSBarry Smith   Level: advanced
175e884886fSBarry Smith 
176e884886fSBarry Smith   Notes:
177e884886fSBarry Smith    See the manual page for MATMFFD_WP for a complete description of the
178e884886fSBarry Smith    algorithm used to compute h.
179e884886fSBarry Smith 
180db781477SPatrick Sanan .seealso: `MatMFFDSetFunctionError()`, `MatCreateSNESMF()`
181e884886fSBarry Smith 
182e884886fSBarry Smith @*/
1837087cfbeSBarry Smith PetscErrorCode  MatMFFDWPSetComputeNormU(Mat A,PetscBool flag)
184e884886fSBarry Smith {
185e884886fSBarry Smith   PetscFunctionBegin;
1860700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
187cac4c232SBarry Smith   PetscTryMethod(A,"MatMFFDWPSetComputeNormU_C",(Mat,PetscBool),(A,flag));
188e884886fSBarry Smith   PetscFunctionReturn(0);
189e884886fSBarry Smith }
190e884886fSBarry Smith 
191e884886fSBarry Smith /*
1921d0fab5eSBarry Smith      MatCreateMFFD_WP - Standard PETSc code for
193e884886fSBarry Smith    computing h with matrix-free finite differences.
194e884886fSBarry Smith 
195e884886fSBarry Smith    Input Parameter:
1961d0fab5eSBarry Smith .  ctx - the matrix free context created by MatCreateMFFD()
197e884886fSBarry Smith 
198e884886fSBarry Smith */
1998cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD_WP(MatMFFD ctx)
200e884886fSBarry Smith {
201e884886fSBarry Smith   MatMFFD_WP     *hctx;
202e884886fSBarry Smith 
203e884886fSBarry Smith   PetscFunctionBegin;
204e884886fSBarry Smith   /* allocate my own private data structure */
2059566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ctx,&hctx));
206e884886fSBarry Smith   ctx->hctx          = (void*)hctx;
207e884886fSBarry Smith   hctx->computenormU = PETSC_FALSE;
208e884886fSBarry Smith 
209e884886fSBarry Smith   /* set the functions I am providing */
210e884886fSBarry Smith   ctx->ops->compute        = MatMFFDCompute_WP;
211e884886fSBarry Smith   ctx->ops->destroy        = MatMFFDDestroy_WP;
212e884886fSBarry Smith   ctx->ops->view           = MatMFFDView_WP;
213e884886fSBarry Smith   ctx->ops->setfromoptions = MatMFFDSetFromOptions_WP;
214e884886fSBarry Smith 
2159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ctx->mat,"MatMFFDWPSetComputeNormU_C",MatMFFDWPSetComputeNormU_P));
216e884886fSBarry Smith   PetscFunctionReturn(0);
217e884886fSBarry Smith }
218