xref: /petsc/src/mat/impls/mffd/wp.c (revision 95452b02e12c0ee11232c7ff2b24b568a8e07e43)
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 
24*95452b02SPatrick Sanan    Notes:
25*95452b02SPatrick Sanan     Requires no global collectives when used with GMRES
26e884886fSBarry Smith 
27e884886fSBarry Smith    Formula used:
28e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
29e884886fSBarry Smith 
30e884886fSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatCreateSNESMF(), MATMFFD_DS
31e884886fSBarry Smith 
32e884886fSBarry Smith M*/
33e884886fSBarry Smith 
34e884886fSBarry Smith /*
35e884886fSBarry Smith     This include file defines the data structure  MatMFFD that
36e884886fSBarry Smith    includes information about the computation of h. It is shared by
37e884886fSBarry Smith    all implementations that people provide.
38e884886fSBarry Smith 
39e884886fSBarry Smith    See snesmfjdef.c for  a full set of comments on the routines below.
40e884886fSBarry Smith */
41af0996ceSBarry Smith #include <petsc/private/matimpl.h>
42c6db04a5SJed Brown #include <../src/mat/impls/mffd/mffdimpl.h>   /*I  "petscmat.h"   I*/
43e884886fSBarry Smith 
44e884886fSBarry Smith typedef struct {
45e884886fSBarry Smith   PetscReal normUfact;                    /* previous sqrt(1.0 + || U ||) */
46ace3abfcSBarry Smith   PetscBool computenormU;
47e884886fSBarry Smith } MatMFFD_WP;
48e884886fSBarry Smith 
49e884886fSBarry Smith /*
50e884886fSBarry Smith      MatMFFDCompute_WP - Standard PETSc code for
51e884886fSBarry Smith    computing h with matrix-free finite differences.
52e884886fSBarry Smith 
53e884886fSBarry Smith   Input Parameters:
54e884886fSBarry Smith +   ctx - the matrix free context
55e884886fSBarry Smith .   U - the location at which you want the Jacobian
56e884886fSBarry Smith -   a - the direction you want the derivative
57e884886fSBarry Smith 
58e884886fSBarry Smith   Output Parameter:
59e884886fSBarry Smith .   h - the scale computed
60e884886fSBarry Smith 
61e884886fSBarry Smith */
62ace3abfcSBarry Smith static PetscErrorCode MatMFFDCompute_WP(MatMFFD ctx,Vec U,Vec a,PetscScalar *h,PetscBool  *zeroa)
63e884886fSBarry Smith {
64e884886fSBarry Smith   MatMFFD_WP     *hctx = (MatMFFD_WP*)ctx->hctx;
6551a79602SBarry Smith   PetscReal      normU,norma;
66e884886fSBarry Smith   PetscErrorCode ierr;
67e884886fSBarry Smith 
68e884886fSBarry Smith   PetscFunctionBegin;
69e884886fSBarry Smith   if (!(ctx->count % ctx->recomputeperiod)) {
7051a79602SBarry Smith     if (hctx->computenormU || !ctx->ncurrenth) {
71e884886fSBarry Smith       ierr            = VecNorm(U,NORM_2,&normU);CHKERRQ(ierr);
728f1a2a5eSBarry Smith       hctx->normUfact = PetscSqrtReal(1.0+normU);
73e884886fSBarry Smith     }
7451a79602SBarry Smith     ierr = VecNorm(a,NORM_2,&norma);CHKERRQ(ierr);
75e884886fSBarry Smith     if (norma == 0.0) {
76e884886fSBarry Smith       *zeroa = PETSC_TRUE;
77e884886fSBarry Smith       PetscFunctionReturn(0);
78e884886fSBarry Smith     }
79e884886fSBarry Smith     *zeroa = PETSC_FALSE;
80e884886fSBarry Smith     *h     = ctx->error_rel*hctx->normUfact/norma;
81e884886fSBarry Smith   } else {
82e884886fSBarry Smith     *h = ctx->currenth;
83e884886fSBarry Smith   }
84e884886fSBarry Smith   PetscFunctionReturn(0);
85e884886fSBarry Smith }
86e884886fSBarry Smith 
87e884886fSBarry Smith /*
88e884886fSBarry Smith    MatMFFDView_WP - Prints information about this particular
89e884886fSBarry Smith      method for computing h. Note that this does not print the general
90e884886fSBarry Smith      information about the matrix free, that is printed by the calling
91e884886fSBarry Smith      routine.
92e884886fSBarry Smith 
93e884886fSBarry Smith   Input Parameters:
94e884886fSBarry Smith +   ctx - the matrix free context
95e884886fSBarry Smith -   viewer - the PETSc viewer
96e884886fSBarry Smith 
97e884886fSBarry Smith */
98e884886fSBarry Smith static PetscErrorCode MatMFFDView_WP(MatMFFD ctx,PetscViewer viewer)
99e884886fSBarry Smith {
100e884886fSBarry Smith   MatMFFD_WP     *hctx = (MatMFFD_WP*)ctx->hctx;
101e884886fSBarry Smith   PetscErrorCode ierr;
102ace3abfcSBarry Smith   PetscBool      iascii;
103e884886fSBarry Smith 
104e884886fSBarry Smith   PetscFunctionBegin;
105251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
106e884886fSBarry Smith   if (iascii) {
1072205254eSKarl Rupp     if (hctx->computenormU) {
1082205254eSKarl Rupp       ierr =  PetscViewerASCIIPrintf(viewer,"    Computes normU\n");CHKERRQ(ierr);
1092205254eSKarl Rupp     } else {
1102205254eSKarl Rupp       ierr =  PetscViewerASCIIPrintf(viewer,"    Does not compute normU\n");CHKERRQ(ierr);
1112205254eSKarl Rupp     }
11211aeaf0aSBarry Smith   }
113e884886fSBarry Smith   PetscFunctionReturn(0);
114e884886fSBarry Smith }
115e884886fSBarry Smith 
116e884886fSBarry Smith /*
117e884886fSBarry Smith    MatMFFDSetFromOptions_WP - Looks in the options database for
118e884886fSBarry Smith      any options appropriate for this method
119e884886fSBarry Smith 
120e884886fSBarry Smith   Input Parameter:
121e884886fSBarry Smith .  ctx - the matrix free context
122e884886fSBarry Smith 
123e884886fSBarry Smith */
1244416b707SBarry Smith static PetscErrorCode MatMFFDSetFromOptions_WP(PetscOptionItems *PetscOptionsObject,MatMFFD ctx)
125e884886fSBarry Smith {
126e884886fSBarry Smith   PetscErrorCode ierr;
127e884886fSBarry Smith   MatMFFD_WP     *hctx = (MatMFFD_WP*)ctx->hctx;
128e884886fSBarry Smith 
129e884886fSBarry Smith   PetscFunctionBegin;
130e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Walker-Pernice options");CHKERRQ(ierr);
1318afaa268SBarry Smith   ierr = PetscOptionsBool("-mat_mffd_compute_normu","Compute the norm of u","MatMFFDWPSetComputeNormU", hctx->computenormU,&hctx->computenormU,NULL);CHKERRQ(ierr);
132e884886fSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
133e884886fSBarry Smith   PetscFunctionReturn(0);
134e884886fSBarry Smith }
135e884886fSBarry Smith 
136e884886fSBarry Smith /*
137e884886fSBarry Smith    MatMFFDDestroy_WP - Frees the space allocated by
1381d0fab5eSBarry Smith        MatCreateMFFD_WP().
139e884886fSBarry Smith 
140e884886fSBarry Smith   Input Parameter:
141e884886fSBarry Smith .  ctx - the matrix free context
142e884886fSBarry Smith 
143*95452b02SPatrick Sanan    Notes:
144*95452b02SPatrick Sanan     does not free the ctx, that is handled by the calling routine
145e884886fSBarry Smith 
146e884886fSBarry Smith */
147e884886fSBarry Smith static PetscErrorCode MatMFFDDestroy_WP(MatMFFD ctx)
148e884886fSBarry Smith {
149e884886fSBarry Smith   PetscErrorCode ierr;
1506e111a19SKarl Rupp 
151e884886fSBarry Smith   PetscFunctionBegin;
152e884886fSBarry Smith   ierr = PetscFree(ctx->hctx);CHKERRQ(ierr);
153e884886fSBarry Smith   PetscFunctionReturn(0);
154e884886fSBarry Smith }
155e884886fSBarry Smith 
1567087cfbeSBarry Smith PetscErrorCode  MatMFFDWPSetComputeNormU_P(Mat mat,PetscBool flag)
157e884886fSBarry Smith {
158e884886fSBarry Smith   MatMFFD    ctx   = (MatMFFD)mat->data;
159e884886fSBarry Smith   MatMFFD_WP *hctx = (MatMFFD_WP*)ctx->hctx;
160e884886fSBarry Smith 
161e884886fSBarry Smith   PetscFunctionBegin;
162e884886fSBarry Smith   hctx->computenormU = flag;
163e884886fSBarry Smith   PetscFunctionReturn(0);
164e884886fSBarry Smith }
165e884886fSBarry Smith 
166e884886fSBarry Smith /*@
167e884886fSBarry Smith     MatMFFDWPSetComputeNormU - Sets whether it computes the ||U|| used by the WP
168e884886fSBarry Smith              PETSc routine for computing h. With any Krylov solver this need only
169e884886fSBarry Smith              be computed during the first iteration and kept for later.
170e884886fSBarry Smith 
171e884886fSBarry Smith   Input Parameters:
172e884886fSBarry Smith +   A - the matrix created with MatCreateSNESMF()
173e884886fSBarry Smith -   flag - PETSC_TRUE causes it to compute ||U||, PETSC_FALSE uses the previous value
174e884886fSBarry Smith 
175e884886fSBarry Smith   Options Database Key:
176e884886fSBarry Smith .   -mat_mffd_compute_normu <true,false> - true by default, false can save calculations but you
177e884886fSBarry Smith               must be sure that ||U|| has not changed in the mean time.
178e884886fSBarry Smith 
179e884886fSBarry Smith   Level: advanced
180e884886fSBarry Smith 
181e884886fSBarry Smith   Notes:
182e884886fSBarry Smith    See the manual page for MATMFFD_WP for a complete description of the
183e884886fSBarry Smith    algorithm used to compute h.
184e884886fSBarry Smith 
185e884886fSBarry Smith .seealso: MatMFFDSetFunctionError(), MatCreateSNESMF()
186e884886fSBarry Smith 
187e884886fSBarry Smith @*/
1887087cfbeSBarry Smith PetscErrorCode  MatMFFDWPSetComputeNormU(Mat A,PetscBool flag)
189e884886fSBarry Smith {
1904ac538c5SBarry Smith   PetscErrorCode ierr;
191e884886fSBarry Smith 
192e884886fSBarry Smith   PetscFunctionBegin;
1930700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1944ac538c5SBarry Smith   ierr = PetscTryMethod(A,"MatMFFDWPSetComputeNormU_C",(Mat,PetscBool),(A,flag));CHKERRQ(ierr);
195e884886fSBarry Smith   PetscFunctionReturn(0);
196e884886fSBarry Smith }
197e884886fSBarry Smith 
198e884886fSBarry Smith /*
1991d0fab5eSBarry Smith      MatCreateMFFD_WP - Standard PETSc code for
200e884886fSBarry Smith    computing h with matrix-free finite differences.
201e884886fSBarry Smith 
202e884886fSBarry Smith    Input Parameter:
2031d0fab5eSBarry Smith .  ctx - the matrix free context created by MatCreateMFFD()
204e884886fSBarry Smith 
205e884886fSBarry Smith */
2068cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD_WP(MatMFFD ctx)
207e884886fSBarry Smith {
208e884886fSBarry Smith   PetscErrorCode ierr;
209e884886fSBarry Smith   MatMFFD_WP     *hctx;
210e884886fSBarry Smith 
211e884886fSBarry Smith   PetscFunctionBegin;
212e884886fSBarry Smith   /* allocate my own private data structure */
213b00a9115SJed Brown   ierr               = PetscNewLog(ctx,&hctx);CHKERRQ(ierr);
214e884886fSBarry Smith   ctx->hctx          = (void*)hctx;
215e884886fSBarry Smith   hctx->computenormU = PETSC_FALSE;
216e884886fSBarry Smith 
217e884886fSBarry Smith   /* set the functions I am providing */
218e884886fSBarry Smith   ctx->ops->compute        = MatMFFDCompute_WP;
219e884886fSBarry Smith   ctx->ops->destroy        = MatMFFDDestroy_WP;
220e884886fSBarry Smith   ctx->ops->view           = MatMFFDView_WP;
221e884886fSBarry Smith   ctx->ops->setfromoptions = MatMFFDSetFromOptions_WP;
222e884886fSBarry Smith 
223bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ctx->mat,"MatMFFDWPSetComputeNormU_C",MatMFFDWPSetComputeNormU_P);CHKERRQ(ierr);
224e884886fSBarry Smith   PetscFunctionReturn(0);
225e884886fSBarry Smith }
226e884886fSBarry Smith 
227e884886fSBarry Smith 
228e884886fSBarry Smith 
229