xref: /petsc/src/mat/impls/mffd/wp.c (revision 0700a8246d308f50502909ba325e6169d3ee27eb)
1e884886fSBarry Smith #define PETSCMAT_DLL
2e884886fSBarry Smith 
3e884886fSBarry Smith /*MC
4e884886fSBarry Smith      MATMFFD_WP - Implements an alternative approach for computing the differencing parameter
5e884886fSBarry Smith         h used with the finite difference based matrix-free Jacobian.  This code
6e884886fSBarry Smith         implements the strategy of M. Pernice and H. Walker:
7e884886fSBarry Smith 
8e884886fSBarry Smith       h = error_rel * sqrt(1 + ||U||) / ||a||
9e884886fSBarry Smith 
10e884886fSBarry Smith       Notes:
11e884886fSBarry Smith         1) || U || does not change between linear iterations so is reused
12e884886fSBarry Smith         2) In GMRES || a || == 1 and so does not need to ever be computed except at restart
13e884886fSBarry Smith            when it is recomputed.
14e884886fSBarry Smith 
15e884886fSBarry Smith       Reference:  M. Pernice and H. F. Walker, "NITSOL: A Newton Iterative
16e884886fSBarry Smith       Solver for Nonlinear Systems", SIAM J. Sci. Stat. Comput.", 1998,
17e884886fSBarry Smith       vol 19, pp. 302--318.
18e884886fSBarry Smith 
19e884886fSBarry Smith    Options Database Keys:
20e884886fSBarry Smith .   -mat_mffd_compute_normu -Compute the norm of u everytime see MatMFFDWPSetComputeNormU()
21e884886fSBarry Smith 
22e884886fSBarry Smith 
23e884886fSBarry Smith    Level: intermediate
24e884886fSBarry Smith 
25e884886fSBarry Smith    Notes: 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 */
417c4f633dSBarry Smith #include "private/matimpl.h"
427c4f633dSBarry Smith #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 ||) */
4651a79602SBarry Smith   PetscTruth computenormU;
47e884886fSBarry Smith } MatMFFD_WP;
48e884886fSBarry Smith 
49e884886fSBarry Smith #undef __FUNCT__
50e884886fSBarry Smith #define __FUNCT__ "MatMFFDCompute_WP"
51e884886fSBarry Smith /*
52e884886fSBarry Smith      MatMFFDCompute_WP - Standard PETSc code for
53e884886fSBarry Smith    computing h with matrix-free finite differences.
54e884886fSBarry Smith 
55e884886fSBarry Smith   Input Parameters:
56e884886fSBarry Smith +   ctx - the matrix free context
57e884886fSBarry Smith .   U - the location at which you want the Jacobian
58e884886fSBarry Smith -   a - the direction you want the derivative
59e884886fSBarry Smith 
60e884886fSBarry Smith   Output Parameter:
61e884886fSBarry Smith .   h - the scale computed
62e884886fSBarry Smith 
63e884886fSBarry Smith */
64e884886fSBarry Smith static PetscErrorCode MatMFFDCompute_WP(MatMFFD ctx,Vec U,Vec a,PetscScalar *h,PetscTruth *zeroa)
65e884886fSBarry Smith {
66e884886fSBarry Smith   MatMFFD_WP    *hctx = (MatMFFD_WP*)ctx->hctx;
6751a79602SBarry Smith   PetscReal      normU,norma;
68e884886fSBarry Smith   PetscErrorCode ierr;
69e884886fSBarry Smith 
70e884886fSBarry Smith   PetscFunctionBegin;
71e884886fSBarry Smith   if (!(ctx->count % ctx->recomputeperiod)) {
7251a79602SBarry Smith     if (hctx->computenormU || !ctx->ncurrenth) {
73e884886fSBarry Smith       ierr = VecNorm(U,NORM_2,&normU);CHKERRQ(ierr);
74e884886fSBarry Smith       hctx->normUfact = sqrt(1.0+normU);
75e884886fSBarry Smith     }
7651a79602SBarry Smith     ierr = VecNorm(a,NORM_2,&norma);CHKERRQ(ierr);
77e884886fSBarry Smith     if (norma == 0.0) {
78e884886fSBarry Smith       *zeroa = PETSC_TRUE;
79e884886fSBarry Smith       PetscFunctionReturn(0);
80e884886fSBarry Smith     }
81e884886fSBarry Smith     *zeroa = PETSC_FALSE;
82e884886fSBarry Smith     *h = ctx->error_rel*hctx->normUfact/norma;
83e884886fSBarry Smith   } else {
84e884886fSBarry Smith     *h = ctx->currenth;
85e884886fSBarry Smith   }
86e884886fSBarry Smith   PetscFunctionReturn(0);
87e884886fSBarry Smith }
88e884886fSBarry Smith 
89e884886fSBarry Smith #undef __FUNCT__
90e884886fSBarry Smith #define __FUNCT__ "MatMFFDView_WP"
91e884886fSBarry Smith /*
92e884886fSBarry Smith    MatMFFDView_WP - Prints information about this particular
93e884886fSBarry Smith      method for computing h. Note that this does not print the general
94e884886fSBarry Smith      information about the matrix free, that is printed by the calling
95e884886fSBarry Smith      routine.
96e884886fSBarry Smith 
97e884886fSBarry Smith   Input Parameters:
98e884886fSBarry Smith +   ctx - the matrix free context
99e884886fSBarry Smith -   viewer - the PETSc viewer
100e884886fSBarry Smith 
101e884886fSBarry Smith */
102e884886fSBarry Smith static PetscErrorCode MatMFFDView_WP(MatMFFD ctx,PetscViewer viewer)
103e884886fSBarry Smith {
104e884886fSBarry Smith   MatMFFD_WP     *hctx = (MatMFFD_WP *)ctx->hctx;
105e884886fSBarry Smith   PetscErrorCode ierr;
106e884886fSBarry Smith   PetscTruth     iascii;
107e884886fSBarry Smith 
108e884886fSBarry Smith   PetscFunctionBegin;
109e884886fSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
110e884886fSBarry Smith   if (iascii) {
111e884886fSBarry Smith     if (hctx->computenormU){ierr =  PetscViewerASCIIPrintf(viewer,"    Computes normU\n");CHKERRQ(ierr);}
112e884886fSBarry Smith     else                   {ierr =  PetscViewerASCIIPrintf(viewer,"    Does not compute normU\n");CHKERRQ(ierr);}
113e884886fSBarry Smith   } else {
114e884886fSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES matrix-free WP",((PetscObject)viewer)->type_name);
115e884886fSBarry 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);
136e884886fSBarry Smith     ierr = PetscOptionsTruth("-mat_mffd_compute_normu","Compute the norm of u","MatMFFDWPSetComputeNormU",
13751a79602SBarry Smith                           hctx->computenormU,&hctx->computenormU,0);CHKERRQ(ierr);
138e884886fSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
139e884886fSBarry Smith   PetscFunctionReturn(0);
140e884886fSBarry Smith }
141e884886fSBarry Smith 
142e884886fSBarry Smith #undef __FUNCT__
143e884886fSBarry Smith #define __FUNCT__ "MatMFFDDestroy_WP"
144e884886fSBarry Smith /*
145e884886fSBarry Smith    MatMFFDDestroy_WP - Frees the space allocated by
146e884886fSBarry Smith        MatMFFDCreate_WP().
147e884886fSBarry Smith 
148e884886fSBarry Smith   Input Parameter:
149e884886fSBarry Smith .  ctx - the matrix free context
150e884886fSBarry Smith 
151e884886fSBarry Smith    Notes: does not free the ctx, that is handled by the calling routine
152e884886fSBarry Smith 
153e884886fSBarry Smith */
154e884886fSBarry Smith static PetscErrorCode MatMFFDDestroy_WP(MatMFFD ctx)
155e884886fSBarry Smith {
156e884886fSBarry Smith   PetscErrorCode ierr;
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"
165e884886fSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU_P(Mat mat,PetscTruth 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 @*/
200e884886fSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat A,PetscTruth flag)
201e884886fSBarry Smith {
202e884886fSBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscTruth);
203e884886fSBarry Smith 
204e884886fSBarry Smith   PetscFunctionBegin;
205*0700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
206e884886fSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMFFDWPSetComputeNormU_C",(void (**)(void))&f);CHKERRQ(ierr);
207e884886fSBarry Smith   if (f) {
208e884886fSBarry Smith     ierr = (*f)(A,flag);CHKERRQ(ierr);
209e884886fSBarry Smith   }
210e884886fSBarry Smith   PetscFunctionReturn(0);
211e884886fSBarry Smith }
212e884886fSBarry Smith 
213e884886fSBarry Smith EXTERN_C_BEGIN
214e884886fSBarry Smith #undef __FUNCT__
215e884886fSBarry Smith #define __FUNCT__ "MatMFFDCreate_WP"
216e884886fSBarry Smith /*
217e884886fSBarry Smith      MatMFFDCreate_WP - Standard PETSc code for
218e884886fSBarry Smith    computing h with matrix-free finite differences.
219e884886fSBarry Smith 
220e884886fSBarry Smith    Input Parameter:
221e884886fSBarry Smith .  ctx - the matrix free context created by MatMFFDCreate()
222e884886fSBarry Smith 
223e884886fSBarry Smith */
224e884886fSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCreate_WP(MatMFFD ctx)
225e884886fSBarry Smith {
226e884886fSBarry Smith   PetscErrorCode ierr;
227e884886fSBarry Smith   MatMFFD_WP     *hctx;
228e884886fSBarry Smith 
229e884886fSBarry Smith   PetscFunctionBegin;
230e884886fSBarry Smith 
231e884886fSBarry Smith   /* allocate my own private data structure */
23238f2d2fdSLisandro Dalcin   ierr               = PetscNewLog(ctx,MatMFFD_WP,&hctx);CHKERRQ(ierr);
233e884886fSBarry Smith   ctx->hctx          = (void*)hctx;
234e884886fSBarry Smith   hctx->computenormU = PETSC_FALSE;
235e884886fSBarry Smith 
236e884886fSBarry Smith   /* set the functions I am providing */
237e884886fSBarry Smith   ctx->ops->compute        = MatMFFDCompute_WP;
238e884886fSBarry Smith   ctx->ops->destroy        = MatMFFDDestroy_WP;
239e884886fSBarry Smith   ctx->ops->view           = MatMFFDView_WP;
240e884886fSBarry Smith   ctx->ops->setfromoptions = MatMFFDSetFromOptions_WP;
241e884886fSBarry Smith 
242e884886fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatMFFDWPSetComputeNormU_C",
243e884886fSBarry Smith                             "MatMFFDWPSetComputeNormU_P",
244e884886fSBarry Smith                              MatMFFDWPSetComputeNormU_P);CHKERRQ(ierr);
245e884886fSBarry Smith 
246e884886fSBarry Smith   PetscFunctionReturn(0);
247e884886fSBarry Smith }
248e884886fSBarry Smith EXTERN_C_END
249e884886fSBarry Smith 
250e884886fSBarry Smith 
251e884886fSBarry Smith 
252