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