xref: /petsc/src/snes/interface/noise/snesmfj2.c (revision e6ed2b1f1d6cf9347963e2d08df9328b0218e4f5)
1d2dddef7SLois Curfman McInnes #ifndef lint
2*e6ed2b1fSLois Curfman McInnes static char vcid[] = "$Id: snesmfj2.c,v 1.1 1997/07/04 03:35:30 curfman Exp curfman $";
3d2dddef7SLois Curfman McInnes #endif
4d2dddef7SLois Curfman McInnes 
5d2dddef7SLois Curfman McInnes #include "src/snes/snesimpl.h"   /*I  "snes.h"   I*/
6d2dddef7SLois Curfman McInnes #include "pinclude/pviewer.h"
7d2dddef7SLois Curfman McInnes #include <math.h>
8d2dddef7SLois Curfman McInnes 
9d2dddef7SLois Curfman McInnes extern int DiffParameterCreate_More(SNES,Vec,void**);
10d2dddef7SLois Curfman McInnes extern int DiffParameterCompute_More(SNES,void*,Vec,Vec,double*,double*);
11d2dddef7SLois Curfman McInnes extern int DiffParameterDestroy_More(void*);
12d2dddef7SLois Curfman McInnes 
13d2dddef7SLois Curfman McInnes typedef struct {  /* default context for matrix-free SNES */
14d2dddef7SLois Curfman McInnes   SNES        snes;      /* SNES context */
15d2dddef7SLois Curfman McInnes   Vec         w;         /* work vector */
16d2dddef7SLois Curfman McInnes   PCNullSpace sp;        /* null space context */
17d2dddef7SLois Curfman McInnes   double      error_rel; /* square root of relative error in computing function */
18d2dddef7SLois Curfman McInnes   double      umin;      /* minimum allowable u'a value relative to |u|_1 */
19d2dddef7SLois Curfman McInnes   int         jorge;     /* flag indicating use of Jorge's method for determining
20d2dddef7SLois Curfman McInnes                             the differencing parameter */
21d2dddef7SLois Curfman McInnes   Scalar      h;         /* differencing parameter */
22d2dddef7SLois Curfman McInnes   int         need_h;    /* flag indicating whether we must compute h */
23d2dddef7SLois Curfman McInnes   int         need_err;  /* flag indicating whether we must compute error_rel */
24d2dddef7SLois Curfman McInnes   void        *data;     /* implementation-specific data */
25d2dddef7SLois Curfman McInnes } MFCtx_Private;
26d2dddef7SLois Curfman McInnes 
27d2dddef7SLois Curfman McInnes #undef __FUNC__
28d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeDestroy2_Private" /* ADIC Ignore */
29d2dddef7SLois Curfman McInnes int SNESMatrixFreeDestroy2_Private(PetscObject obj)
30d2dddef7SLois Curfman McInnes {
31d2dddef7SLois Curfman McInnes   int           ierr;
32d2dddef7SLois Curfman McInnes   Mat           mat = (Mat) obj;
33d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
34d2dddef7SLois Curfman McInnes 
35d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx);
36d2dddef7SLois Curfman McInnes   ierr = VecDestroy(ctx->w); CHKERRQ(ierr);
37d2dddef7SLois Curfman McInnes   if (ctx->sp) {ierr = PCNullSpaceDestroy(ctx->sp); CHKERRQ(ierr);}
38d2dddef7SLois Curfman McInnes   if (ctx->jorge || ctx->need_err) {ierr = DiffParameterDestroy_More(ctx->data); CHKERRQ(ierr);}
39d2dddef7SLois Curfman McInnes   PetscFree(ctx);
40d2dddef7SLois Curfman McInnes   return 0;
41d2dddef7SLois Curfman McInnes }
42d2dddef7SLois Curfman McInnes 
43d2dddef7SLois Curfman McInnes #undef __FUNC__
44d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeView2_Private" /* ADIC Ignore */
45d2dddef7SLois Curfman McInnes /*
46d2dddef7SLois Curfman McInnes    SNESMatrixFreeView2_Private - Views matrix-free parameters.
47d2dddef7SLois Curfman McInnes  */
48d2dddef7SLois Curfman McInnes int SNESMatrixFreeView2_Private(Mat J,Viewer viewer)
49d2dddef7SLois Curfman McInnes {
50d2dddef7SLois Curfman McInnes   int           ierr;
51d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
52d2dddef7SLois Curfman McInnes   MPI_Comm      comm;
53d2dddef7SLois Curfman McInnes   FILE          *fd;
54d2dddef7SLois Curfman McInnes   ViewerType    vtype;
55d2dddef7SLois Curfman McInnes 
56d2dddef7SLois Curfman McInnes   PetscObjectGetComm((PetscObject)J,&comm);
57d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
58d2dddef7SLois Curfman McInnes   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
59d2dddef7SLois Curfman McInnes   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
60d2dddef7SLois Curfman McInnes   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
61d2dddef7SLois Curfman McInnes      PetscFPrintf(comm,fd,"  SNES matrix-free approximation:\n");
62d2dddef7SLois Curfman McInnes      if (ctx->jorge)
63d2dddef7SLois Curfman McInnes        PetscFPrintf(comm,fd,"    using Jorge's method of determining differencing parameter\n");
64d2dddef7SLois Curfman McInnes      PetscFPrintf(comm,fd,"    err=%g (relative error in function evaluation)\n",ctx->error_rel);
65d2dddef7SLois Curfman McInnes      PetscFPrintf(comm,fd,"    umin=%g (minimum iterate parameter)\n",ctx->umin);
66d2dddef7SLois Curfman McInnes   }
67d2dddef7SLois Curfman McInnes   return 0;
68d2dddef7SLois Curfman McInnes }
69d2dddef7SLois Curfman McInnes 
70*e6ed2b1fSLois Curfman McInnes extern int VecDot_Seq(Vec,Vec,Scalar *);
71*e6ed2b1fSLois Curfman McInnes extern int VecNorm_Seq(Vec,NormType,double *);
72*e6ed2b1fSLois Curfman McInnes 
73d2dddef7SLois Curfman McInnes #undef __FUNC__
74d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeMult2_Private"
75d2dddef7SLois Curfman McInnes /*
76d2dddef7SLois Curfman McInnes   SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector
77d2dddef7SLois Curfman McInnes   product, y = F'(u)*a:
78d2dddef7SLois Curfman McInnes         y = ( F(u + ha) - F(u) ) /h,
79d2dddef7SLois Curfman McInnes   where F = nonlinear function, as set by SNESSetFunction()
80d2dddef7SLois Curfman McInnes         u = current iterate
81d2dddef7SLois Curfman McInnes         h = difference interval
82d2dddef7SLois Curfman McInnes */
83d2dddef7SLois Curfman McInnes int SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y)
84d2dddef7SLois Curfman McInnes {
85d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
86d2dddef7SLois Curfman McInnes   SNES          snes;
87*e6ed2b1fSLois Curfman McInnes   double        norm, sum, umin, noise, ovalues[3],values[3];
88d2dddef7SLois Curfman McInnes   Scalar        h, dot, mone = -1.0;
89d2dddef7SLois Curfman McInnes   Vec           w,U,F;
90d2dddef7SLois Curfman McInnes   int           ierr, (*eval_fct)(SNES,Vec,Vec);
91*e6ed2b1fSLois Curfman McInnes   MPI_Comm      comm;
92d2dddef7SLois Curfman McInnes 
93d2dddef7SLois Curfman McInnes 
94d2dddef7SLois Curfman McInnes   /* We log matrix-free matrix-vector products separately, so that we can
95d2dddef7SLois Curfman McInnes      separate the performance monitoring from the cases that use conventional
96d2dddef7SLois Curfman McInnes      storage.  We may eventually modify event logging to associate events
97d2dddef7SLois Curfman McInnes      with particular objects, hence alleviating the more general problem. */
98d2dddef7SLois Curfman McInnes   PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0);
99d2dddef7SLois Curfman McInnes 
100*e6ed2b1fSLois Curfman McInnes   PetscObjectGetComm((PetscObject)mat,&comm);
101*e6ed2b1fSLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
102*e6ed2b1fSLois Curfman McInnes   snes = ctx->snes;
103*e6ed2b1fSLois Curfman McInnes   w    = ctx->w;
104*e6ed2b1fSLois Curfman McInnes   umin = ctx->umin;
105*e6ed2b1fSLois Curfman McInnes 
106d2dddef7SLois Curfman McInnes   ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr);
107d2dddef7SLois Curfman McInnes   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
108d2dddef7SLois Curfman McInnes     eval_fct = SNESComputeFunction;
109d2dddef7SLois Curfman McInnes     ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr);
110d2dddef7SLois Curfman McInnes   }
111d2dddef7SLois Curfman McInnes   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
112d2dddef7SLois Curfman McInnes     eval_fct = SNESComputeGradient;
113d2dddef7SLois Curfman McInnes     ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr);
114d2dddef7SLois Curfman McInnes   }
115d2dddef7SLois Curfman McInnes   else SETERRQ(1,0,"Invalid method class");
116d2dddef7SLois Curfman McInnes 
117d2dddef7SLois Curfman McInnes   /* Determine a "good" step size, h */
118d2dddef7SLois Curfman McInnes   if (ctx->need_h) {
119d2dddef7SLois Curfman McInnes 
120d2dddef7SLois Curfman McInnes     /* Use Jorge's method to compute h */
121d2dddef7SLois Curfman McInnes     if (ctx->jorge) {
122d2dddef7SLois Curfman McInnes       ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h); CHKERRQ(ierr);
123d2dddef7SLois Curfman McInnes 
124d2dddef7SLois Curfman McInnes     /* Use the Brown/Saad method to compute h */
125d2dddef7SLois Curfman McInnes     } else {
126d2dddef7SLois Curfman McInnes       if (ctx->need_err) {
127d2dddef7SLois Curfman McInnes         /* Use Jorge's method to compute noise */
128d2dddef7SLois Curfman McInnes         ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h); CHKERRQ(ierr);
129d2dddef7SLois Curfman McInnes         ctx->error_rel = sqrt(noise);
130d2dddef7SLois Curfman McInnes         PLogInfo(snes,"SNESMatrixFreeMult2_Private: Using Jorge's noise: noise=%g, sqrt(noise)=%g, h_more=%g\n",
131d2dddef7SLois Curfman McInnes             noise,ctx->error_rel,h);
132d2dddef7SLois Curfman McInnes 
133d2dddef7SLois Curfman McInnes         ctx->need_err = 0;
134d2dddef7SLois Curfman McInnes       }
135*e6ed2b1fSLois Curfman McInnes 
136*e6ed2b1fSLois Curfman McInnes       /*
137d2dddef7SLois Curfman McInnes       ierr = VecDot(U,a,&dot); CHKERRQ(ierr);
138d2dddef7SLois Curfman McInnes       ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr);
139d2dddef7SLois Curfman McInnes       ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr);
140*e6ed2b1fSLois Curfman McInnes       */
141*e6ed2b1fSLois Curfman McInnes 
142*e6ed2b1fSLois Curfman McInnes      /*
143*e6ed2b1fSLois Curfman McInnes         Call the Seq Vector routines and then do a single reduction
144*e6ed2b1fSLois Curfman McInnes         to reduce the number of communications required
145*e6ed2b1fSLois Curfman McInnes      */
146*e6ed2b1fSLois Curfman McInnes 
147*e6ed2b1fSLois Curfman McInnes #if !defined(PETSC_COMPLEX)
148*e6ed2b1fSLois Curfman McInnes      PLogEventBegin(VEC_Dot,U,a,0,0);
149*e6ed2b1fSLois Curfman McInnes      ierr = VecDot_Seq(U,a,ovalues); CHKERRQ(ierr);
150*e6ed2b1fSLois Curfman McInnes      PLogEventEnd(VEC_Dot,U,a,0,0);
151*e6ed2b1fSLois Curfman McInnes      PLogEventBegin(VEC_Norm,a,0,0,0);
152*e6ed2b1fSLois Curfman McInnes      ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr);
153*e6ed2b1fSLois Curfman McInnes      ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr);
154*e6ed2b1fSLois Curfman McInnes      ovalues[2] = ovalues[2]*ovalues[2];
155*e6ed2b1fSLois Curfman McInnes      MPI_Allreduce(ovalues,values,3,MPI_DOUBLE,MPI_SUM,comm );
156*e6ed2b1fSLois Curfman McInnes      dot = values[0]; sum = values[1]; norm = sqrt(values[2]);
157*e6ed2b1fSLois Curfman McInnes      PLogEventEnd(VEC_Norm,a,0,0,0);
158*e6ed2b1fSLois Curfman McInnes #else
159*e6ed2b1fSLois Curfman McInnes      {
160*e6ed2b1fSLois Curfman McInnes        Scalar cvalues[3],covalues[3];
161*e6ed2b1fSLois Curfman McInnes 
162*e6ed2b1fSLois Curfman McInnes        PLogEventBegin(VEC_Dot,U,a,0,0);
163*e6ed2b1fSLois Curfman McInnes        ierr = VecDot_Seq(U,a,covalues); CHKERRQ(ierr);
164*e6ed2b1fSLois Curfman McInnes        PLogEventEnd(VEC_Dot,U,a,0,0);
165*e6ed2b1fSLois Curfman McInnes        PLogEventBegin(VEC_Norm,a,0,0,0);
166*e6ed2b1fSLois Curfman McInnes        ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr);
167*e6ed2b1fSLois Curfman McInnes        ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr);
168*e6ed2b1fSLois Curfman McInnes        covalues[1] = ovalues[1];
169*e6ed2b1fSLois Curfman McInnes        covalues[2] = ovalues[2]*ovalues[2];
170*e6ed2b1fSLois Curfman McInnes        MPI_Allreduce(covalues,cvalues,6,MPI_DOUBLE,MPI_SUM,comm );
171*e6ed2b1fSLois Curfman McInnes        dot = cvalues[0]; sum = PetscReal(cvalues[1]); norm = sqrt(PetscReal(cvalues[2]));
172*e6ed2b1fSLois Curfman McInnes        PLogEventEnd(VEC_Norm,a,0,0,0);
173*e6ed2b1fSLois Curfman McInnes      }
174*e6ed2b1fSLois Curfman McInnes #endif
175d2dddef7SLois Curfman McInnes 
176d2dddef7SLois Curfman McInnes       /* Safeguard for step sizes too small */
177d2dddef7SLois Curfman McInnes       if (sum == 0.0) {dot = 1.0; norm = 1.0;}
178d2dddef7SLois Curfman McInnes #if defined(PETSC_COMPLEX)
179d2dddef7SLois Curfman McInnes       else if (abs(dot) < umin*sum && real(dot) >= 0.0) dot = umin*sum;
180d2dddef7SLois Curfman McInnes       else if (abs(dot) < 0.0 && real(dot) > -umin*sum) dot = -umin*sum;
181d2dddef7SLois Curfman McInnes #else
182d2dddef7SLois Curfman McInnes       else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
183d2dddef7SLois Curfman McInnes       else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
184d2dddef7SLois Curfman McInnes #endif
185d2dddef7SLois Curfman McInnes       h = ctx->error_rel*dot/(norm*norm);
186d2dddef7SLois Curfman McInnes     }
187d2dddef7SLois Curfman McInnes   } else {
188d2dddef7SLois Curfman McInnes     h = ctx->h;
189d2dddef7SLois Curfman McInnes   }
190d2dddef7SLois Curfman McInnes   if (!ctx->jorge || !ctx->need_h) PLogInfo(snes,"SNESMatrixFreeMult2_Private: h = %g\n",h);
191d2dddef7SLois Curfman McInnes 
192d2dddef7SLois Curfman McInnes   /* Evaluate function at F(u + ha) */
193d2dddef7SLois Curfman McInnes   ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr);
194d2dddef7SLois Curfman McInnes   ierr = eval_fct(snes,w,y); CHKERRQ(ierr);
195d2dddef7SLois Curfman McInnes   ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr);
196d2dddef7SLois Curfman McInnes   h = 1.0/h;
197d2dddef7SLois Curfman McInnes   ierr = VecScale(&h,y); CHKERRQ(ierr);
198d2dddef7SLois Curfman McInnes   if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);}
199d2dddef7SLois Curfman McInnes 
200d2dddef7SLois Curfman McInnes   PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0);
201d2dddef7SLois Curfman McInnes   return 0;
202d2dddef7SLois Curfman McInnes }
203d2dddef7SLois Curfman McInnes 
204d2dddef7SLois Curfman McInnes #undef __FUNC__
205d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESDefaultMatrixFreeMatCreate2"
206d2dddef7SLois Curfman McInnes /*@C
207d2dddef7SLois Curfman McInnes    SNESMatrixFreeMatCreate2 - Creates a matrix-free matrix
208d2dddef7SLois Curfman McInnes    context for use with a SNES solver.  This matrix can be used as
209d2dddef7SLois Curfman McInnes    the Jacobian argument for the routine SNESSetJacobian().
210d2dddef7SLois Curfman McInnes 
211d2dddef7SLois Curfman McInnes    Input Parameters:
212d2dddef7SLois Curfman McInnes .  snes - the SNES context
213d2dddef7SLois Curfman McInnes .  x - vector where SNES solution is to be stored.
214d2dddef7SLois Curfman McInnes 
215d2dddef7SLois Curfman McInnes    Output Parameter:
216d2dddef7SLois Curfman McInnes .  J - the matrix-free matrix
217d2dddef7SLois Curfman McInnes 
218d2dddef7SLois Curfman McInnes    Notes:
219d2dddef7SLois Curfman McInnes    The matrix-free matrix context merely contains the function pointers
220d2dddef7SLois Curfman McInnes    and work space for performing finite difference approximations of
221d2dddef7SLois Curfman McInnes    Jacobian-vector products, J(u)*a, via
222d2dddef7SLois Curfman McInnes 
223d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h,
224d2dddef7SLois Curfman McInnes $   where by default
225d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
226d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
227d2dddef7SLois Curfman McInnes $   where
228d2dddef7SLois Curfman McInnes $        error_rel = square root of relative error in
229d2dddef7SLois Curfman McInnes $                    function evaluation
230d2dddef7SLois Curfman McInnes $        umin = minimum iterate parameter
231d2dddef7SLois Curfman McInnes $   Alternatively, the differencing parameter, h, can be set using
232*e6ed2b1fSLois Curfman McInnes $   Jorge's nifty new strategy if one specifies the option
233*e6ed2b1fSLois Curfman McInnes $          -snes_mf_jorge
234d2dddef7SLois Curfman McInnes 
235d2dddef7SLois Curfman McInnes    The user can set these parameters via SNESSetMatrixFreeParameters().
236d2dddef7SLois Curfman McInnes    See the nonlinear solvers chapter of the users manual for details.
237d2dddef7SLois Curfman McInnes 
238d2dddef7SLois Curfman McInnes    The user should call MatDestroy() when finished with the matrix-free
239d2dddef7SLois Curfman McInnes    matrix context.
240d2dddef7SLois Curfman McInnes 
241d2dddef7SLois Curfman McInnes    Options Database Keys:
242*e6ed2b1fSLois Curfman McInnes $  -snes_mf_compute_err
243d2dddef7SLois Curfman McInnes $  -snes_mf_err <error_rel>
244d2dddef7SLois Curfman McInnes $  -snes_mf_unim <umin>
245*e6ed2b1fSLois Curfman McInnes $  -snes_mf_jorge
246d2dddef7SLois Curfman McInnes 
247d2dddef7SLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix
248d2dddef7SLois Curfman McInnes 
249d2dddef7SLois Curfman McInnes .seealso: MatDestroy(), SNESSetMatrixFreeParameters()
250d2dddef7SLois Curfman McInnes @*/
251d2dddef7SLois Curfman McInnes int SNESMatrixFreeMatCreate2(SNES snes,Vec x, Mat *J)
252d2dddef7SLois Curfman McInnes {
253d2dddef7SLois Curfman McInnes   MPI_Comm      comm;
254d2dddef7SLois Curfman McInnes   MFCtx_Private *mfctx;
255d2dddef7SLois Curfman McInnes   int           n, nloc, ierr, flg;
256d2dddef7SLois Curfman McInnes   char          p[64];
257d2dddef7SLois Curfman McInnes 
258d2dddef7SLois Curfman McInnes   mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx);
259d2dddef7SLois Curfman McInnes   PLogObjectMemory(snes,sizeof(MFCtx_Private));
260d2dddef7SLois Curfman McInnes   mfctx->sp   = 0;
261d2dddef7SLois Curfman McInnes   mfctx->snes = snes;
262d2dddef7SLois Curfman McInnes   mfctx->error_rel = 1.e-8; /* assumes double precision */
263d2dddef7SLois Curfman McInnes   mfctx->umin      = 1.e-6;
264d2dddef7SLois Curfman McInnes   mfctx->h         = 0.0;
265d2dddef7SLois Curfman McInnes   mfctx->need_h    = 1;
266d2dddef7SLois Curfman McInnes   mfctx->need_err  = 0;
267d2dddef7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr);
268d2dddef7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr);
269*e6ed2b1fSLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf_jorge",&mfctx->jorge); CHKERRQ(ierr);
270*e6ed2b1fSLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf_compute_err",&mfctx->need_err); CHKERRQ(ierr);
271d2dddef7SLois Curfman McInnes   if (mfctx->jorge || mfctx->need_err) {
272d2dddef7SLois Curfman McInnes     ierr = DiffParameterCreate_More(snes,x,&mfctx->data); CHKERRQ(ierr);
273d2dddef7SLois Curfman McInnes   } else mfctx->data = 0;
274d2dddef7SLois Curfman McInnes 
275d2dddef7SLois Curfman McInnes   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
276d2dddef7SLois Curfman McInnes   PetscStrcpy(p,"-");
277d2dddef7SLois Curfman McInnes   if (snes->prefix) PetscStrcat(p,snes->prefix);
278d2dddef7SLois Curfman McInnes   if (flg) {
279*e6ed2b1fSLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_jorge: use Jorge More's method\n",p);
280*e6ed2b1fSLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_compute_err: compute sqrt rel error in function\n",p);
281d2dddef7SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel);
282d2dddef7SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_umin <umin>: see users manual (default %g)\n",p,mfctx->umin);
283d2dddef7SLois Curfman McInnes   }
284d2dddef7SLois Curfman McInnes   ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr);
285d2dddef7SLois Curfman McInnes   ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr);
286d2dddef7SLois Curfman McInnes   ierr = VecGetSize(x,&n); CHKERRQ(ierr);
287d2dddef7SLois Curfman McInnes   ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr);
288d2dddef7SLois Curfman McInnes   ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr);
289d2dddef7SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult2_Private);CHKERRQ(ierr);
290d2dddef7SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy2_Private);CHKERRQ(ierr);
291d2dddef7SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView2_Private); CHKERRQ(ierr);
292d2dddef7SLois Curfman McInnes 
293d2dddef7SLois Curfman McInnes   PLogObjectParent(*J,mfctx->w);
294d2dddef7SLois Curfman McInnes   PLogObjectParent(snes,*J);
295d2dddef7SLois Curfman McInnes   return 0;
296d2dddef7SLois Curfman McInnes }
297d2dddef7SLois Curfman McInnes 
298d2dddef7SLois Curfman McInnes #undef __FUNC__
299d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESSetMatrixFreeParameters2"
300d2dddef7SLois Curfman McInnes /*@
301d2dddef7SLois Curfman McInnes    SNESSetMatrixFreeParameters2 - Sets the parameters for the approximation of
302d2dddef7SLois Curfman McInnes    matrix-vector products using finite differences.
303d2dddef7SLois Curfman McInnes 
304d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h where
305d2dddef7SLois Curfman McInnes 
306d2dddef7SLois Curfman McInnes    either the user sets h directly here, or this parameter is computed via
307d2dddef7SLois Curfman McInnes 
308d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
309d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
310d2dddef7SLois Curfman McInnes $
311d2dddef7SLois Curfman McInnes 
312d2dddef7SLois Curfman McInnes    Input Parameters:
313d2dddef7SLois Curfman McInnes .  snes - the SNES context
314d2dddef7SLois Curfman McInnes .  error_rel - relative error (should be set to the square root of
315d2dddef7SLois Curfman McInnes      the relative error in the function evaluations)
316d2dddef7SLois Curfman McInnes .  umin - minimum allowable u-value
317d2dddef7SLois Curfman McInnes .  h - differencing parameter
318d2dddef7SLois Curfman McInnes 
319d2dddef7SLois Curfman McInnes    Notes:
320d2dddef7SLois Curfman McInnes    If the user sets the parameter h directly, then this value will be used
321d2dddef7SLois Curfman McInnes    instead of the default computation indicated above.
322d2dddef7SLois Curfman McInnes 
323d2dddef7SLois Curfman McInnes .keywords: SNES, matrix-free, parameters
324d2dddef7SLois Curfman McInnes 
325d2dddef7SLois Curfman McInnes .seealso: SNESDefaultMatrixFreeMatCreate()
326d2dddef7SLois Curfman McInnes @*/
327d2dddef7SLois Curfman McInnes int SNESSetMatrixFreeParameters2(SNES snes,double error,double umin,double h)
328d2dddef7SLois Curfman McInnes {
329d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
330d2dddef7SLois Curfman McInnes   int           ierr;
331d2dddef7SLois Curfman McInnes   Mat           mat;
332d2dddef7SLois Curfman McInnes 
333d2dddef7SLois Curfman McInnes   ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr);
334d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
335d2dddef7SLois Curfman McInnes   if (ctx) {
336d2dddef7SLois Curfman McInnes     if (error != PETSC_DEFAULT) ctx->error_rel = error;
337d2dddef7SLois Curfman McInnes     if (umin  != PETSC_DEFAULT) ctx->umin = umin;
338d2dddef7SLois Curfman McInnes     if (h     != PETSC_DEFAULT) {
339d2dddef7SLois Curfman McInnes       ctx->h = h;
340d2dddef7SLois Curfman McInnes       ctx->need_h = 0;
341d2dddef7SLois Curfman McInnes     }
342d2dddef7SLois Curfman McInnes   }
343d2dddef7SLois Curfman McInnes   return 0;
344d2dddef7SLois Curfman McInnes }
345d2dddef7SLois Curfman McInnes 
346d2dddef7SLois Curfman McInnes int SNESUnSetMatrixFreeParameter(SNES snes)
347d2dddef7SLois Curfman McInnes {
348d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
349d2dddef7SLois Curfman McInnes   int           ierr;
350d2dddef7SLois Curfman McInnes   Mat           mat;
351d2dddef7SLois Curfman McInnes 
352d2dddef7SLois Curfman McInnes   ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr);
353d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
354d2dddef7SLois Curfman McInnes   if (ctx) ctx->need_h = 1;
355d2dddef7SLois Curfman McInnes   return 0;
356d2dddef7SLois Curfman McInnes }
357d2dddef7SLois Curfman McInnes 
358