xref: /petsc/src/snes/interface/noise/snesmfj2.c (revision 549d3d68a6ae470532d58d544870024f02ff2d7c)
1d834d2e4SLois Curfman McInnes #ifdef PETSC_RCS_HEADER
2*549d3d68SSatish Balay static char vcid[] = "$Id: snesmfj2.c,v 1.15 1999/04/19 22:15:33 bsmith Exp balay $";
3d2dddef7SLois Curfman McInnes #endif
4d2dddef7SLois Curfman McInnes 
5d2dddef7SLois Curfman McInnes #include "src/snes/snesimpl.h"   /*I  "snes.h"   I*/
6d2dddef7SLois Curfman McInnes 
7d2dddef7SLois Curfman McInnes extern int DiffParameterCreate_More(SNES,Vec,void**);
8d2dddef7SLois Curfman McInnes extern int DiffParameterCompute_More(SNES,void*,Vec,Vec,double*,double*);
9d2dddef7SLois Curfman McInnes extern int DiffParameterDestroy_More(void*);
10d2dddef7SLois Curfman McInnes 
11d2dddef7SLois Curfman McInnes typedef struct {  /* default context for matrix-free SNES */
12d2dddef7SLois Curfman McInnes   SNES        snes;             /* SNES context */
13d2dddef7SLois Curfman McInnes   Vec         w;                /* work vector */
14d2dddef7SLois Curfman McInnes   PCNullSpace sp;               /* null space context */
15d2dddef7SLois Curfman McInnes   double      error_rel;        /* square root of relative error in computing function */
16d2dddef7SLois Curfman McInnes   double      umin;             /* minimum allowable u'a value relative to |u|_1 */
17d2dddef7SLois Curfman McInnes   int         jorge;            /* flag indicating use of Jorge's method for determining
18d2dddef7SLois Curfman McInnes                                    the differencing parameter */
19db0b4b35SLois Curfman McInnes   double      h;                /* differencing parameter */
20d2dddef7SLois Curfman McInnes   int         need_h;           /* flag indicating whether we must compute h */
21295f7fbbSLois Curfman McInnes   int         need_err;         /* flag indicating whether we must currently compute error_rel */
22295f7fbbSLois Curfman McInnes   int         compute_err;      /* flag indicating whether we must ever compute error_rel */
23295f7fbbSLois Curfman McInnes   int         compute_err_iter; /* last iter where we've computer error_rel */
24295f7fbbSLois Curfman McInnes   int         compute_err_freq; /* frequency of computing error_rel */
25d2dddef7SLois Curfman McInnes   void        *data;            /* implementation-specific data */
26d2dddef7SLois Curfman McInnes } MFCtx_Private;
27d2dddef7SLois Curfman McInnes 
28d2dddef7SLois Curfman McInnes #undef __FUNC__
29d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeDestroy2_Private" /* ADIC Ignore */
30f1d4f8a3SBarry Smith int SNESMatrixFreeDestroy2_Private(Mat mat)
31d2dddef7SLois Curfman McInnes {
32d2dddef7SLois Curfman McInnes   int           ierr;
33d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
34d2dddef7SLois Curfman McInnes 
3507a3eed9SLois Curfman McInnes   PetscFunctionBegin;
36d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx);
37d2dddef7SLois Curfman McInnes   ierr = VecDestroy(ctx->w);CHKERRQ(ierr);
38d2dddef7SLois Curfman McInnes   if (ctx->sp) {ierr = PCNullSpaceDestroy(ctx->sp);CHKERRQ(ierr);}
39295f7fbbSLois Curfman McInnes   if (ctx->jorge || ctx->compute_err) {ierr = DiffParameterDestroy_More(ctx->data);CHKERRQ(ierr);}
40d2dddef7SLois Curfman McInnes   PetscFree(ctx);
4107a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
42d2dddef7SLois Curfman McInnes }
43d2dddef7SLois Curfman McInnes 
44d2dddef7SLois Curfman McInnes #undef __FUNC__
45d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeView2_Private" /* ADIC Ignore */
46d2dddef7SLois Curfman McInnes /*
47d2dddef7SLois Curfman McInnes    SNESMatrixFreeView2_Private - Views matrix-free parameters.
48d2dddef7SLois Curfman McInnes  */
49d2dddef7SLois Curfman McInnes int SNESMatrixFreeView2_Private(Mat J,Viewer viewer)
50d2dddef7SLois Curfman McInnes {
51d2dddef7SLois Curfman McInnes   int           ierr;
52d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
53d2dddef7SLois Curfman McInnes   MPI_Comm      comm;
54d2dddef7SLois Curfman McInnes   FILE          *fd;
55d2dddef7SLois Curfman McInnes   ViewerType    vtype;
56d2dddef7SLois Curfman McInnes 
5707a3eed9SLois Curfman McInnes   PetscFunctionBegin;
582d0c0e3bSBarry Smith   ierr = PetscObjectGetComm((PetscObject)J,&comm);CHKERRQ(ierr);
59d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(J,(void **)&ctx);CHKERRQ(ierr);
60d2dddef7SLois Curfman McInnes   ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr);
61d2dddef7SLois Curfman McInnes   ierr = ViewerASCIIGetPointer(viewer,&fd);CHKERRQ(ierr);
62758f395bSBarry Smith   if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
63d132466eSBarry Smith      ierr = PetscFPrintf(comm,fd,"  SNES matrix-free approximation:\n");CHKERRQ(ierr);
64758f395bSBarry Smith      if (ctx->jorge) {
65d132466eSBarry Smith        ierr = PetscFPrintf(comm,fd,"    using Jorge's method of determining differencing parameter\n");CHKERRQ(ierr);
66758f395bSBarry Smith      }
67d132466eSBarry Smith      ierr = PetscFPrintf(comm,fd,"    err=%g (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr);
68d132466eSBarry Smith      ierr = PetscFPrintf(comm,fd,"    umin=%g (minimum iterate parameter)\n",ctx->umin);CHKERRQ(ierr);
69758f395bSBarry Smith      if (ctx->compute_err) {
70d132466eSBarry Smith        ierr = PetscFPrintf(comm,fd,"    freq_err=%d (frequency for computing err)\n",ctx->compute_err_freq);CHKERRQ(ierr);
71d2dddef7SLois Curfman McInnes      }
72758f395bSBarry Smith   } else {
73758f395bSBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
74758f395bSBarry Smith   }
7507a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
76d2dddef7SLois Curfman McInnes }
77d2dddef7SLois Curfman McInnes 
78d2dddef7SLois Curfman McInnes #undef __FUNC__
79d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeMult2_Private"
80d2dddef7SLois Curfman McInnes /*
81d2dddef7SLois Curfman McInnes   SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector
82d2dddef7SLois Curfman McInnes   product, y = F'(u)*a:
83d2dddef7SLois Curfman McInnes         y = ( F(u + ha) - F(u) ) /h,
84d2dddef7SLois Curfman McInnes   where F = nonlinear function, as set by SNESSetFunction()
85d2dddef7SLois Curfman McInnes         u = current iterate
86d2dddef7SLois Curfman McInnes         h = difference interval
87d2dddef7SLois Curfman McInnes */
88d2dddef7SLois Curfman McInnes int SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y)
89d2dddef7SLois Curfman McInnes {
90d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
91d2dddef7SLois Curfman McInnes   SNES          snes;
92691598edSBarry Smith   double        h, norm, sum, umin, noise;
93db0b4b35SLois Curfman McInnes   Scalar        hs, dot, mone = -1.0;
94d2dddef7SLois Curfman McInnes   Vec           w,U,F;
95295f7fbbSLois Curfman McInnes   int           ierr, iter, (*eval_fct)(SNES,Vec,Vec);
96e6ed2b1fSLois Curfman McInnes   MPI_Comm      comm;
97d2dddef7SLois Curfman McInnes 
9807a3eed9SLois Curfman McInnes   PetscFunctionBegin;
9907a3eed9SLois Curfman McInnes 
100d2dddef7SLois Curfman McInnes   /* We log matrix-free matrix-vector products separately, so that we can
101d2dddef7SLois Curfman McInnes      separate the performance monitoring from the cases that use conventional
102d2dddef7SLois Curfman McInnes      storage.  We may eventually modify event logging to associate events
103d2dddef7SLois Curfman McInnes      with particular objects, hence alleviating the more general problem. */
104d2dddef7SLois Curfman McInnes   PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0);
105d2dddef7SLois Curfman McInnes 
1062d0c0e3bSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
107e6ed2b1fSLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
108e6ed2b1fSLois Curfman McInnes   snes = ctx->snes;
109e6ed2b1fSLois Curfman McInnes   w    = ctx->w;
110e6ed2b1fSLois Curfman McInnes   umin = ctx->umin;
111e6ed2b1fSLois Curfman McInnes 
112d2dddef7SLois Curfman McInnes   ierr = SNESGetSolution(snes,&U);CHKERRQ(ierr);
113d2dddef7SLois Curfman McInnes   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
114d2dddef7SLois Curfman McInnes     eval_fct = SNESComputeFunction;
115d2dddef7SLois Curfman McInnes     ierr = SNESGetFunction(snes,&F);CHKERRQ(ierr);
116d2dddef7SLois Curfman McInnes   }
117d2dddef7SLois Curfman McInnes   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
118d2dddef7SLois Curfman McInnes     eval_fct = SNESComputeGradient;
119d2dddef7SLois Curfman McInnes     ierr = SNESGetGradient(snes,&F);CHKERRQ(ierr);
120d2dddef7SLois Curfman McInnes   }
121d2dddef7SLois Curfman McInnes   else SETERRQ(1,0,"Invalid method class");
122d2dddef7SLois Curfman McInnes 
123295f7fbbSLois Curfman McInnes 
124d2dddef7SLois Curfman McInnes   /* Determine a "good" step size, h */
125d2dddef7SLois Curfman McInnes   if (ctx->need_h) {
126d2dddef7SLois Curfman McInnes 
127d2dddef7SLois Curfman McInnes     /* Use Jorge's method to compute h */
128d2dddef7SLois Curfman McInnes     if (ctx->jorge) {
129d2dddef7SLois Curfman McInnes       ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr);
130d2dddef7SLois Curfman McInnes 
131d2dddef7SLois Curfman McInnes     /* Use the Brown/Saad method to compute h */
132d2dddef7SLois Curfman McInnes     } else {
133295f7fbbSLois Curfman McInnes       /* Compute error if desired */
134295f7fbbSLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
135295f7fbbSLois Curfman McInnes       if ((ctx->need_err) ||
136295f7fbbSLois Curfman McInnes           ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter-1)%ctx->compute_err_freq)))) {
137d2dddef7SLois Curfman McInnes         /* Use Jorge's method to compute noise */
138d2dddef7SLois Curfman McInnes         ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr);
139d2dddef7SLois Curfman McInnes         ctx->error_rel = sqrt(noise);
140d2dddef7SLois Curfman McInnes         PLogInfo(snes,"SNESMatrixFreeMult2_Private: Using Jorge's noise: noise=%g, sqrt(noise)=%g, h_more=%g\n",
141d2dddef7SLois Curfman McInnes             noise,ctx->error_rel,h);
142295f7fbbSLois Curfman McInnes         ctx->compute_err_iter = iter;
143d2dddef7SLois Curfman McInnes         ctx->need_err = 0;
144d2dddef7SLois Curfman McInnes       }
145e6ed2b1fSLois Curfman McInnes 
146691598edSBarry Smith       ierr = VecDotBegin(U,a,&dot);CHKERRQ(ierr);
147691598edSBarry Smith       ierr = VecNormBegin(a,NORM_1,&sum);CHKERRQ(ierr);
148691598edSBarry Smith       ierr = VecNormBegin(a,NORM_2,&norm);CHKERRQ(ierr);
149691598edSBarry Smith       ierr = VecDotEnd(U,a,&dot);CHKERRQ(ierr);
150691598edSBarry Smith       ierr = VecNormEnd(a,NORM_1,&sum);CHKERRQ(ierr);
151691598edSBarry Smith       ierr = VecNormEnd(a,NORM_2,&norm);CHKERRQ(ierr);
152e6ed2b1fSLois Curfman McInnes 
153d2dddef7SLois Curfman McInnes 
154d2dddef7SLois Curfman McInnes       /* Safeguard for step sizes too small */
155d2dddef7SLois Curfman McInnes       if (sum == 0.0) {dot = 1.0; norm = 1.0;}
156db0b4b35SLois Curfman McInnes #if defined(USE_PETSC_COMPLEX)
157e20fef11SSatish Balay       else if (PetscAbsScalar(dot) < umin*sum && PetscReal(dot) >= 0.0) dot = umin*sum;
158e20fef11SSatish Balay       else if (PetscAbsScalar(dot) < 0.0 && PetscReal(dot) > -umin*sum) dot = -umin*sum;
159d2dddef7SLois Curfman McInnes #else
160d2dddef7SLois Curfman McInnes       else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
161d2dddef7SLois Curfman McInnes       else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
162d2dddef7SLois Curfman McInnes #endif
163db0b4b35SLois Curfman McInnes       h = PetscReal(ctx->error_rel*dot/(norm*norm));
164d2dddef7SLois Curfman McInnes     }
165d2dddef7SLois Curfman McInnes   } else {
166d2dddef7SLois Curfman McInnes     h = ctx->h;
167d2dddef7SLois Curfman McInnes   }
168d2dddef7SLois Curfman McInnes   if (!ctx->jorge || !ctx->need_h) PLogInfo(snes,"SNESMatrixFreeMult2_Private: h = %g\n",h);
169d2dddef7SLois Curfman McInnes 
170d2dddef7SLois Curfman McInnes   /* Evaluate function at F(u + ha) */
171db0b4b35SLois Curfman McInnes   hs = h;
172db0b4b35SLois Curfman McInnes   ierr = VecWAXPY(&hs,a,U,w);CHKERRQ(ierr);
173d2dddef7SLois Curfman McInnes   ierr = eval_fct(snes,w,y);CHKERRQ(ierr);
174d2dddef7SLois Curfman McInnes   ierr = VecAXPY(&mone,F,y);CHKERRQ(ierr);
175db0b4b35SLois Curfman McInnes   hs = 1.0/hs;
176db0b4b35SLois Curfman McInnes   ierr = VecScale(&hs,y);CHKERRQ(ierr);
177d2dddef7SLois Curfman McInnes   if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y);CHKERRQ(ierr);}
178d2dddef7SLois Curfman McInnes 
179d2dddef7SLois Curfman McInnes   PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0);
18007a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
181d2dddef7SLois Curfman McInnes }
182d2dddef7SLois Curfman McInnes 
183d2dddef7SLois Curfman McInnes #undef __FUNC__
1841a7c28d5SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeMatCreate2"
185d2dddef7SLois Curfman McInnes /*@C
186d2dddef7SLois Curfman McInnes    SNESMatrixFreeMatCreate2 - Creates a matrix-free matrix
187d2dddef7SLois Curfman McInnes    context for use with a SNES solver.  This matrix can be used as
188d2dddef7SLois Curfman McInnes    the Jacobian argument for the routine SNESSetJacobian().
189d2dddef7SLois Curfman McInnes 
190d2dddef7SLois Curfman McInnes    Input Parameters:
191d2dddef7SLois Curfman McInnes .  snes - the SNES context
192d2dddef7SLois Curfman McInnes .  x - vector where SNES solution is to be stored.
193d2dddef7SLois Curfman McInnes 
194d2dddef7SLois Curfman McInnes    Output Parameter:
195d2dddef7SLois Curfman McInnes .  J - the matrix-free matrix
196d2dddef7SLois Curfman McInnes 
197d2dddef7SLois Curfman McInnes    Notes:
198d2dddef7SLois Curfman McInnes    The matrix-free matrix context merely contains the function pointers
199d2dddef7SLois Curfman McInnes    and work space for performing finite difference approximations of
200d2dddef7SLois Curfman McInnes    Jacobian-vector products, J(u)*a, via
201d2dddef7SLois Curfman McInnes 
202d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h,
203d2dddef7SLois Curfman McInnes $   where by default
204d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
205d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
206d2dddef7SLois Curfman McInnes $   where
207d2dddef7SLois Curfman McInnes $        error_rel = square root of relative error in
208d2dddef7SLois Curfman McInnes $                    function evaluation
209d2dddef7SLois Curfman McInnes $        umin = minimum iterate parameter
210d2dddef7SLois Curfman McInnes $   Alternatively, the differencing parameter, h, can be set using
211e6ed2b1fSLois Curfman McInnes $   Jorge's nifty new strategy if one specifies the option
212e6ed2b1fSLois Curfman McInnes $          -snes_mf_jorge
213d2dddef7SLois Curfman McInnes 
2145a655dc6SBarry Smith    The user can set these parameters via MatSNESMFSetFunctionError().
215d2dddef7SLois Curfman McInnes    See the nonlinear solvers chapter of the users manual for details.
216d2dddef7SLois Curfman McInnes 
217d2dddef7SLois Curfman McInnes    The user should call MatDestroy() when finished with the matrix-free
218d2dddef7SLois Curfman McInnes    matrix context.
219d2dddef7SLois Curfman McInnes 
220d2dddef7SLois Curfman McInnes    Options Database Keys:
221d2dddef7SLois Curfman McInnes $  -snes_mf_err <error_rel>
222d2dddef7SLois Curfman McInnes $  -snes_mf_unim <umin>
223295f7fbbSLois Curfman McInnes $  -snes_mf_compute_err
224295f7fbbSLois Curfman McInnes $  -snes_mf_freq_err <freq>
225e6ed2b1fSLois Curfman McInnes $  -snes_mf_jorge
226d2dddef7SLois Curfman McInnes 
227d2dddef7SLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix
228d2dddef7SLois Curfman McInnes 
2295a655dc6SBarry Smith .seealso: MatDestroy(), MatSNESMFSetFunctionError()
230d2dddef7SLois Curfman McInnes @*/
231758f395bSBarry Smith int SNESDefaultMatrixFreeCreate2(SNES snes,Vec x, Mat *J)
232d2dddef7SLois Curfman McInnes {
233d2dddef7SLois Curfman McInnes   MPI_Comm      comm;
234d2dddef7SLois Curfman McInnes   MFCtx_Private *mfctx;
235d2dddef7SLois Curfman McInnes   int           n, nloc, ierr, flg;
236d2dddef7SLois Curfman McInnes   char          p[64];
237d2dddef7SLois Curfman McInnes 
23807a3eed9SLois Curfman McInnes   PetscFunctionBegin;
239d2dddef7SLois Curfman McInnes   mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private));CHKPTRQ(mfctx);
240*549d3d68SSatish Balay   ierr  = PetscMemzero(mfctx,sizeof(MFCtx_Private));CHKERRQ(ierr);
241d2dddef7SLois Curfman McInnes   PLogObjectMemory(snes,sizeof(MFCtx_Private));
242d2dddef7SLois Curfman McInnes   mfctx->sp   = 0;
243d2dddef7SLois Curfman McInnes   mfctx->snes = snes;
244d2dddef7SLois Curfman McInnes   mfctx->error_rel        = 1.e-8; /* assumes double precision */
245d2dddef7SLois Curfman McInnes   mfctx->umin             = 1.e-6;
246d2dddef7SLois Curfman McInnes   mfctx->h                = 0.0;
247d2dddef7SLois Curfman McInnes   mfctx->need_h           = 1;
248d2dddef7SLois Curfman McInnes   mfctx->need_err         = 0;
249295f7fbbSLois Curfman McInnes   mfctx->compute_err      = 0;
250295f7fbbSLois Curfman McInnes   mfctx->compute_err_freq = 0;
251295f7fbbSLois Curfman McInnes   mfctx->compute_err_iter = -1;
252d2dddef7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg);CHKERRQ(ierr);
253d2dddef7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg);CHKERRQ(ierr);
254e6ed2b1fSLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf_jorge",&mfctx->jorge);CHKERRQ(ierr);
255295f7fbbSLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf_compute_err",&mfctx->compute_err);CHKERRQ(ierr);
256295f7fbbSLois Curfman McInnes   ierr = OptionsGetInt(snes->prefix,"-snes_mf_freq_err",&mfctx->compute_err_freq,&flg);CHKERRQ(ierr);
257295f7fbbSLois Curfman McInnes   if (flg) {
258295f7fbbSLois Curfman McInnes     if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0;
259295f7fbbSLois Curfman McInnes     mfctx->compute_err = 1;
260295f7fbbSLois Curfman McInnes   }
261295f7fbbSLois Curfman McInnes   if (mfctx->compute_err == 1) mfctx->need_err = 1;
262295f7fbbSLois Curfman McInnes   if (mfctx->jorge || mfctx->compute_err) {
263d2dddef7SLois Curfman McInnes     ierr = DiffParameterCreate_More(snes,x,&mfctx->data);CHKERRQ(ierr);
264d2dddef7SLois Curfman McInnes   } else mfctx->data = 0;
265d2dddef7SLois Curfman McInnes 
266d2dddef7SLois Curfman McInnes   ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr);
267*549d3d68SSatish Balay   ierr = PetscStrcpy(p,"-");CHKERRQ(ierr);
268d2dddef7SLois Curfman McInnes   if (snes->prefix) PetscStrcat(p,snes->prefix);
269d2dddef7SLois Curfman McInnes   if (flg) {
270d132466eSBarry Smith     ierr = PetscPrintf(snes->comm," Matrix-free Options (via SNES):\n");CHKERRQ(ierr);
271d132466eSBarry Smith     ierr = PetscPrintf(snes->comm,"   %ssnes_mf_err <err>: set sqrt of relative error in function (default %g)\n",p,mfctx->error_rel);CHKERRQ(ierr);
272d132466eSBarry Smith     ierr = PetscPrintf(snes->comm,"   %ssnes_mf_umin <umin>: see users manual (default %g)\n",p,mfctx->umin);CHKERRQ(ierr);
273d132466eSBarry Smith     ierr = PetscPrintf(snes->comm,"   %ssnes_mf_jorge: use Jorge More's method\n",p);CHKERRQ(ierr);
274d132466eSBarry Smith     ierr = PetscPrintf(snes->comm,"   %ssnes_mf_compute_err: compute sqrt or relative error in function\n",p);CHKERRQ(ierr);
275d132466eSBarry Smith     ierr = PetscPrintf(snes->comm,"   %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)\n",p);CHKERRQ(ierr);
276d132466eSBarry Smith     ierr = PetscPrintf(snes->comm,"   %ssnes_mf_noise_file <file>: set file for printing noise info\n",p);CHKERRQ(ierr);
277d2dddef7SLois Curfman McInnes   }
278d2dddef7SLois Curfman McInnes   ierr = VecDuplicate(x,&mfctx->w);CHKERRQ(ierr);
279d2dddef7SLois Curfman McInnes   ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr);
280d2dddef7SLois Curfman McInnes   ierr = VecGetSize(x,&n);CHKERRQ(ierr);
281d2dddef7SLois Curfman McInnes   ierr = VecGetLocalSize(x,&nloc);CHKERRQ(ierr);
282d2dddef7SLois Curfman McInnes   ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J);CHKERRQ(ierr);
283d2dddef7SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult2_Private);CHKERRQ(ierr);
284d2dddef7SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy2_Private);CHKERRQ(ierr);
285d2dddef7SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView2_Private);CHKERRQ(ierr);
286d2dddef7SLois Curfman McInnes 
287d2dddef7SLois Curfman McInnes   PLogObjectParent(*J,mfctx->w);
288d2dddef7SLois Curfman McInnes   PLogObjectParent(snes,*J);
28907a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
290d2dddef7SLois Curfman McInnes }
291d2dddef7SLois Curfman McInnes 
292d2dddef7SLois Curfman McInnes #undef __FUNC__
293758f395bSBarry Smith #define __FUNC__ "SNESDefaultMatrixFreeSetParameters2"
294d2dddef7SLois Curfman McInnes /*@
295758f395bSBarry Smith    SNESDefaultMatrixFreeSetParameters2 - Sets the parameters for the approximation of
296d2dddef7SLois Curfman McInnes    matrix-vector products using finite differences.
297d2dddef7SLois Curfman McInnes 
298d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h where
299d2dddef7SLois Curfman McInnes 
300d2dddef7SLois Curfman McInnes    either the user sets h directly here, or this parameter is computed via
301d2dddef7SLois Curfman McInnes 
302d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
303d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
304d2dddef7SLois Curfman McInnes $
305d2dddef7SLois Curfman McInnes 
306d2dddef7SLois Curfman McInnes    Input Parameters:
307758f395bSBarry Smith +  mat - the matrix
308d2dddef7SLois Curfman McInnes .  error_rel - relative error (should be set to the square root of
309d2dddef7SLois Curfman McInnes      the relative error in the function evaluations)
310d2dddef7SLois Curfman McInnes .  umin - minimum allowable u-value
311758f395bSBarry Smith -  h - differencing parameter
312d2dddef7SLois Curfman McInnes 
313d2dddef7SLois Curfman McInnes    Notes:
314d2dddef7SLois Curfman McInnes    If the user sets the parameter h directly, then this value will be used
315d2dddef7SLois Curfman McInnes    instead of the default computation indicated above.
316d2dddef7SLois Curfman McInnes 
317d2dddef7SLois Curfman McInnes .keywords: SNES, matrix-free, parameters
318d2dddef7SLois Curfman McInnes 
3195a655dc6SBarry Smith .seealso: MatCreateSNESMF()
320d2dddef7SLois Curfman McInnes @*/
321758f395bSBarry Smith int SNESDefaultMatrixFreeSetParameters2(Mat mat,double error,double umin,double h)
322d2dddef7SLois Curfman McInnes {
323d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
324d2dddef7SLois Curfman McInnes   int           ierr;
325d2dddef7SLois Curfman McInnes 
32607a3eed9SLois Curfman McInnes   PetscFunctionBegin;
327d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
328d2dddef7SLois Curfman McInnes   if (ctx) {
329d2dddef7SLois Curfman McInnes     if (error != PETSC_DEFAULT) ctx->error_rel = error;
330d2dddef7SLois Curfman McInnes     if (umin  != PETSC_DEFAULT) ctx->umin = umin;
331d2dddef7SLois Curfman McInnes     if (h     != PETSC_DEFAULT) {
332d2dddef7SLois Curfman McInnes       ctx->h = h;
333d2dddef7SLois Curfman McInnes       ctx->need_h = 0;
334d2dddef7SLois Curfman McInnes     }
335d2dddef7SLois Curfman McInnes   }
33607a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
337d2dddef7SLois Curfman McInnes }
338d2dddef7SLois Curfman McInnes 
339d2dddef7SLois Curfman McInnes int SNESUnSetMatrixFreeParameter(SNES snes)
340d2dddef7SLois Curfman McInnes {
341d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
342d2dddef7SLois Curfman McInnes   int           ierr;
343d2dddef7SLois Curfman McInnes   Mat           mat;
344d2dddef7SLois Curfman McInnes 
34507a3eed9SLois Curfman McInnes   PetscFunctionBegin;
346d2dddef7SLois Curfman McInnes   ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
347d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
348d2dddef7SLois Curfman McInnes   if (ctx) ctx->need_h = 1;
34907a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
350d2dddef7SLois Curfman McInnes }
351d2dddef7SLois Curfman McInnes 
352