xref: /petsc/src/snes/interface/noise/snesmfj2.c (revision 6831982abb6453c2f3e25efecb5d0951d333371e)
1d834d2e4SLois Curfman McInnes #ifdef PETSC_RCS_HEADER
2*6831982aSBarry Smith static char vcid[] = "$Id: snesmfj2.c,v 1.21 1999/10/13 20:38:26 bsmith Exp bsmith $";
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);}
40606d414cSSatish Balay   ierr = PetscFree(ctx);CHKERRQ(ierr);
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;
53*6831982aSBarry Smith   PetscTruth    isascii;
54d2dddef7SLois Curfman McInnes 
5507a3eed9SLois Curfman McInnes   PetscFunctionBegin;
56d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(J,(void **)&ctx);CHKERRQ(ierr);
57*6831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
580f5bd95cSBarry Smith   if (isascii) {
59*6831982aSBarry Smith      ierr = ViewerASCIIPrintf(viewer,"  SNES matrix-free approximation:\n");CHKERRQ(ierr);
60758f395bSBarry Smith      if (ctx->jorge) {
61*6831982aSBarry Smith        ierr = ViewerASCIIPrintf(viewer,"    using Jorge's method of determining differencing parameter\n");CHKERRQ(ierr);
62758f395bSBarry Smith      }
63*6831982aSBarry Smith      ierr = ViewerASCIIPrintf(viewer,"    err=%g (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr);
64*6831982aSBarry Smith      ierr = ViewerASCIIPrintf(viewer,"    umin=%g (minimum iterate parameter)\n",ctx->umin);CHKERRQ(ierr);
65758f395bSBarry Smith      if (ctx->compute_err) {
66*6831982aSBarry Smith        ierr = ViewerASCIIPrintf(viewer,"    freq_err=%d (frequency for computing err)\n",ctx->compute_err_freq);CHKERRQ(ierr);
67d2dddef7SLois Curfman McInnes      }
68758f395bSBarry Smith   } else {
690f5bd95cSBarry Smith     SETERRQ1(1,1,"Viewer type %s not supported by SNES matrix free Jorge",((PetscObject)viewer)->type_name);
70758f395bSBarry Smith   }
7107a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
72d2dddef7SLois Curfman McInnes }
73d2dddef7SLois Curfman McInnes 
74d2dddef7SLois Curfman McInnes #undef __FUNC__
75d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeMult2_Private"
76d2dddef7SLois Curfman McInnes /*
77d2dddef7SLois Curfman McInnes   SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector
78d2dddef7SLois Curfman McInnes   product, y = F'(u)*a:
79d2dddef7SLois Curfman McInnes         y = ( F(u + ha) - F(u) ) /h,
80d2dddef7SLois Curfman McInnes   where F = nonlinear function, as set by SNESSetFunction()
81d2dddef7SLois Curfman McInnes         u = current iterate
82d2dddef7SLois Curfman McInnes         h = difference interval
83d2dddef7SLois Curfman McInnes */
84d2dddef7SLois Curfman McInnes int SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y)
85d2dddef7SLois Curfman McInnes {
86d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
87d2dddef7SLois Curfman McInnes   SNES          snes;
88691598edSBarry Smith   double        h, norm, sum, umin, noise;
89db0b4b35SLois Curfman McInnes   Scalar        hs, dot, mone = -1.0;
90d2dddef7SLois Curfman McInnes   Vec           w,U,F;
91295f7fbbSLois Curfman McInnes   int           ierr, iter, (*eval_fct)(SNES,Vec,Vec);
92e6ed2b1fSLois Curfman McInnes   MPI_Comm      comm;
93d2dddef7SLois Curfman McInnes 
9407a3eed9SLois Curfman McInnes   PetscFunctionBegin;
9507a3eed9SLois Curfman McInnes 
96d2dddef7SLois Curfman McInnes   /* We log matrix-free matrix-vector products separately, so that we can
97d2dddef7SLois Curfman McInnes      separate the performance monitoring from the cases that use conventional
98d2dddef7SLois Curfman McInnes      storage.  We may eventually modify event logging to associate events
99d2dddef7SLois Curfman McInnes      with particular objects, hence alleviating the more general problem. */
100d2dddef7SLois Curfman McInnes   PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0);
101d2dddef7SLois Curfman McInnes 
1022d0c0e3bSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
103e6ed2b1fSLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
104e6ed2b1fSLois Curfman McInnes   snes = ctx->snes;
105e6ed2b1fSLois Curfman McInnes   w    = ctx->w;
106e6ed2b1fSLois Curfman McInnes   umin = ctx->umin;
107e6ed2b1fSLois Curfman McInnes 
108d2dddef7SLois Curfman McInnes   ierr = SNESGetSolution(snes,&U);CHKERRQ(ierr);
109d2dddef7SLois Curfman McInnes   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
110d2dddef7SLois Curfman McInnes     eval_fct = SNESComputeFunction;
111184914b5SBarry Smith     ierr = SNESGetFunction(snes,&F,PETSC_NULL);CHKERRQ(ierr);
112d2dddef7SLois Curfman McInnes   }
113d2dddef7SLois Curfman McInnes   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
114d2dddef7SLois Curfman McInnes     eval_fct = SNESComputeGradient;
115184914b5SBarry Smith     ierr = SNESGetGradient(snes,&F,PETSC_NULL);CHKERRQ(ierr);
116d2dddef7SLois Curfman McInnes   }
117d2dddef7SLois Curfman McInnes   else SETERRQ(1,0,"Invalid method class");
118d2dddef7SLois Curfman McInnes 
119295f7fbbSLois Curfman McInnes 
120d2dddef7SLois Curfman McInnes   /* Determine a "good" step size, h */
121d2dddef7SLois Curfman McInnes   if (ctx->need_h) {
122d2dddef7SLois Curfman McInnes 
123d2dddef7SLois Curfman McInnes     /* Use Jorge's method to compute h */
124d2dddef7SLois Curfman McInnes     if (ctx->jorge) {
125d2dddef7SLois Curfman McInnes       ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr);
126d2dddef7SLois Curfman McInnes 
127d2dddef7SLois Curfman McInnes     /* Use the Brown/Saad method to compute h */
128d2dddef7SLois Curfman McInnes     } else {
129295f7fbbSLois Curfman McInnes       /* Compute error if desired */
130295f7fbbSLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
131295f7fbbSLois Curfman McInnes       if ((ctx->need_err) ||
132295f7fbbSLois Curfman McInnes           ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter-1)%ctx->compute_err_freq)))) {
133d2dddef7SLois Curfman McInnes         /* Use Jorge's method to compute noise */
134d2dddef7SLois Curfman McInnes         ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr);
135d2dddef7SLois Curfman McInnes         ctx->error_rel = sqrt(noise);
136d2dddef7SLois Curfman McInnes         PLogInfo(snes,"SNESMatrixFreeMult2_Private: Using Jorge's noise: noise=%g, sqrt(noise)=%g, h_more=%g\n",
137d2dddef7SLois Curfman McInnes             noise,ctx->error_rel,h);
138295f7fbbSLois Curfman McInnes         ctx->compute_err_iter = iter;
139d2dddef7SLois Curfman McInnes         ctx->need_err = 0;
140d2dddef7SLois Curfman McInnes       }
141e6ed2b1fSLois Curfman McInnes 
142691598edSBarry Smith       ierr = VecDotBegin(U,a,&dot);CHKERRQ(ierr);
143691598edSBarry Smith       ierr = VecNormBegin(a,NORM_1,&sum);CHKERRQ(ierr);
144691598edSBarry Smith       ierr = VecNormBegin(a,NORM_2,&norm);CHKERRQ(ierr);
145691598edSBarry Smith       ierr = VecDotEnd(U,a,&dot);CHKERRQ(ierr);
146691598edSBarry Smith       ierr = VecNormEnd(a,NORM_1,&sum);CHKERRQ(ierr);
147691598edSBarry Smith       ierr = VecNormEnd(a,NORM_2,&norm);CHKERRQ(ierr);
148e6ed2b1fSLois Curfman McInnes 
149d2dddef7SLois Curfman McInnes 
150d2dddef7SLois Curfman McInnes       /* Safeguard for step sizes too small */
151d2dddef7SLois Curfman McInnes       if (sum == 0.0) {dot = 1.0; norm = 1.0;}
152aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
153e20fef11SSatish Balay       else if (PetscAbsScalar(dot) < umin*sum && PetscReal(dot) >= 0.0) dot = umin*sum;
154e20fef11SSatish Balay       else if (PetscAbsScalar(dot) < 0.0 && PetscReal(dot) > -umin*sum) dot = -umin*sum;
155d2dddef7SLois Curfman McInnes #else
156d2dddef7SLois Curfman McInnes       else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
157d2dddef7SLois Curfman McInnes       else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
158d2dddef7SLois Curfman McInnes #endif
159db0b4b35SLois Curfman McInnes       h = PetscReal(ctx->error_rel*dot/(norm*norm));
160d2dddef7SLois Curfman McInnes     }
161d2dddef7SLois Curfman McInnes   } else {
162d2dddef7SLois Curfman McInnes     h = ctx->h;
163d2dddef7SLois Curfman McInnes   }
164d2dddef7SLois Curfman McInnes   if (!ctx->jorge || !ctx->need_h) PLogInfo(snes,"SNESMatrixFreeMult2_Private: h = %g\n",h);
165d2dddef7SLois Curfman McInnes 
166d2dddef7SLois Curfman McInnes   /* Evaluate function at F(u + ha) */
167db0b4b35SLois Curfman McInnes   hs = h;
168db0b4b35SLois Curfman McInnes   ierr = VecWAXPY(&hs,a,U,w);CHKERRQ(ierr);
169d2dddef7SLois Curfman McInnes   ierr = eval_fct(snes,w,y);CHKERRQ(ierr);
170d2dddef7SLois Curfman McInnes   ierr = VecAXPY(&mone,F,y);CHKERRQ(ierr);
171db0b4b35SLois Curfman McInnes   hs = 1.0/hs;
172db0b4b35SLois Curfman McInnes   ierr = VecScale(&hs,y);CHKERRQ(ierr);
173d2dddef7SLois Curfman McInnes   if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y);CHKERRQ(ierr);}
174d2dddef7SLois Curfman McInnes 
175d2dddef7SLois Curfman McInnes   PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0);
17607a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
177d2dddef7SLois Curfman McInnes }
178d2dddef7SLois Curfman McInnes 
179d2dddef7SLois Curfman McInnes #undef __FUNC__
1801a7c28d5SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeMatCreate2"
181d2dddef7SLois Curfman McInnes /*@C
182d2dddef7SLois Curfman McInnes    SNESMatrixFreeMatCreate2 - Creates a matrix-free matrix
183d2dddef7SLois Curfman McInnes    context for use with a SNES solver.  This matrix can be used as
184d2dddef7SLois Curfman McInnes    the Jacobian argument for the routine SNESSetJacobian().
185d2dddef7SLois Curfman McInnes 
186d2dddef7SLois Curfman McInnes    Input Parameters:
187d2dddef7SLois Curfman McInnes .  snes - the SNES context
188d2dddef7SLois Curfman McInnes .  x - vector where SNES solution is to be stored.
189d2dddef7SLois Curfman McInnes 
190d2dddef7SLois Curfman McInnes    Output Parameter:
191d2dddef7SLois Curfman McInnes .  J - the matrix-free matrix
192d2dddef7SLois Curfman McInnes 
193d2dddef7SLois Curfman McInnes    Notes:
194d2dddef7SLois Curfman McInnes    The matrix-free matrix context merely contains the function pointers
195d2dddef7SLois Curfman McInnes    and work space for performing finite difference approximations of
196d2dddef7SLois Curfman McInnes    Jacobian-vector products, J(u)*a, via
197d2dddef7SLois Curfman McInnes 
198d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h,
199d2dddef7SLois Curfman McInnes $   where by default
200d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
201d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
202d2dddef7SLois Curfman McInnes $   where
203d2dddef7SLois Curfman McInnes $        error_rel = square root of relative error in
204d2dddef7SLois Curfman McInnes $                    function evaluation
205d2dddef7SLois Curfman McInnes $        umin = minimum iterate parameter
206d2dddef7SLois Curfman McInnes $   Alternatively, the differencing parameter, h, can be set using
207e6ed2b1fSLois Curfman McInnes $   Jorge's nifty new strategy if one specifies the option
208e6ed2b1fSLois Curfman McInnes $          -snes_mf_jorge
209d2dddef7SLois Curfman McInnes 
2105a655dc6SBarry Smith    The user can set these parameters via MatSNESMFSetFunctionError().
211d2dddef7SLois Curfman McInnes    See the nonlinear solvers chapter of the users manual for details.
212d2dddef7SLois Curfman McInnes 
213d2dddef7SLois Curfman McInnes    The user should call MatDestroy() when finished with the matrix-free
214d2dddef7SLois Curfman McInnes    matrix context.
215d2dddef7SLois Curfman McInnes 
216d2dddef7SLois Curfman McInnes    Options Database Keys:
217d2dddef7SLois Curfman McInnes $  -snes_mf_err <error_rel>
218d2dddef7SLois Curfman McInnes $  -snes_mf_unim <umin>
219295f7fbbSLois Curfman McInnes $  -snes_mf_compute_err
220295f7fbbSLois Curfman McInnes $  -snes_mf_freq_err <freq>
221e6ed2b1fSLois Curfman McInnes $  -snes_mf_jorge
222d2dddef7SLois Curfman McInnes 
223d2dddef7SLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix
224d2dddef7SLois Curfman McInnes 
2255a655dc6SBarry Smith .seealso: MatDestroy(), MatSNESMFSetFunctionError()
226d2dddef7SLois Curfman McInnes @*/
227758f395bSBarry Smith int SNESDefaultMatrixFreeCreate2(SNES snes,Vec x, Mat *J)
228d2dddef7SLois Curfman McInnes {
229d2dddef7SLois Curfman McInnes   MPI_Comm      comm;
230d2dddef7SLois Curfman McInnes   MFCtx_Private *mfctx;
231d2dddef7SLois Curfman McInnes   int           n, nloc, ierr, flg;
232d2dddef7SLois Curfman McInnes   char          p[64];
233d2dddef7SLois Curfman McInnes 
23407a3eed9SLois Curfman McInnes   PetscFunctionBegin;
235d2dddef7SLois Curfman McInnes   mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private));CHKPTRQ(mfctx);
236549d3d68SSatish Balay   ierr  = PetscMemzero(mfctx,sizeof(MFCtx_Private));CHKERRQ(ierr);
237d2dddef7SLois Curfman McInnes   PLogObjectMemory(snes,sizeof(MFCtx_Private));
238d2dddef7SLois Curfman McInnes   mfctx->sp   = 0;
239d2dddef7SLois Curfman McInnes   mfctx->snes = snes;
240d2dddef7SLois Curfman McInnes   mfctx->error_rel        = 1.e-8; /* assumes double precision */
241d2dddef7SLois Curfman McInnes   mfctx->umin             = 1.e-6;
242d2dddef7SLois Curfman McInnes   mfctx->h                = 0.0;
243d2dddef7SLois Curfman McInnes   mfctx->need_h           = 1;
244d2dddef7SLois Curfman McInnes   mfctx->need_err         = 0;
245295f7fbbSLois Curfman McInnes   mfctx->compute_err      = 0;
246295f7fbbSLois Curfman McInnes   mfctx->compute_err_freq = 0;
247295f7fbbSLois Curfman McInnes   mfctx->compute_err_iter = -1;
248d2dddef7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg);CHKERRQ(ierr);
249d2dddef7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg);CHKERRQ(ierr);
250e6ed2b1fSLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf_jorge",&mfctx->jorge);CHKERRQ(ierr);
251295f7fbbSLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf_compute_err",&mfctx->compute_err);CHKERRQ(ierr);
252295f7fbbSLois Curfman McInnes   ierr = OptionsGetInt(snes->prefix,"-snes_mf_freq_err",&mfctx->compute_err_freq,&flg);CHKERRQ(ierr);
253295f7fbbSLois Curfman McInnes   if (flg) {
254295f7fbbSLois Curfman McInnes     if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0;
255295f7fbbSLois Curfman McInnes     mfctx->compute_err = 1;
256295f7fbbSLois Curfman McInnes   }
257295f7fbbSLois Curfman McInnes   if (mfctx->compute_err == 1) mfctx->need_err = 1;
258295f7fbbSLois Curfman McInnes   if (mfctx->jorge || mfctx->compute_err) {
259d2dddef7SLois Curfman McInnes     ierr = DiffParameterCreate_More(snes,x,&mfctx->data);CHKERRQ(ierr);
260d2dddef7SLois Curfman McInnes   } else mfctx->data = 0;
261d2dddef7SLois Curfman McInnes 
262d2dddef7SLois Curfman McInnes   ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr);
263549d3d68SSatish Balay   ierr = PetscStrcpy(p,"-");CHKERRQ(ierr);
264d2dddef7SLois Curfman McInnes   if (snes->prefix) PetscStrcat(p,snes->prefix);
265d2dddef7SLois Curfman McInnes   if (flg) {
266d132466eSBarry Smith     ierr = PetscPrintf(snes->comm," Matrix-free Options (via SNES):\n");CHKERRQ(ierr);
267d132466eSBarry 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);
268d132466eSBarry Smith     ierr = PetscPrintf(snes->comm,"   %ssnes_mf_umin <umin>: see users manual (default %g)\n",p,mfctx->umin);CHKERRQ(ierr);
269d132466eSBarry Smith     ierr = PetscPrintf(snes->comm,"   %ssnes_mf_jorge: use Jorge More's method\n",p);CHKERRQ(ierr);
270d132466eSBarry Smith     ierr = PetscPrintf(snes->comm,"   %ssnes_mf_compute_err: compute sqrt or relative error in function\n",p);CHKERRQ(ierr);
271d132466eSBarry Smith     ierr = PetscPrintf(snes->comm,"   %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)\n",p);CHKERRQ(ierr);
272d132466eSBarry Smith     ierr = PetscPrintf(snes->comm,"   %ssnes_mf_noise_file <file>: set file for printing noise info\n",p);CHKERRQ(ierr);
273d2dddef7SLois Curfman McInnes   }
274d2dddef7SLois Curfman McInnes   ierr = VecDuplicate(x,&mfctx->w);CHKERRQ(ierr);
275d2dddef7SLois Curfman McInnes   ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr);
276d2dddef7SLois Curfman McInnes   ierr = VecGetSize(x,&n);CHKERRQ(ierr);
277d2dddef7SLois Curfman McInnes   ierr = VecGetLocalSize(x,&nloc);CHKERRQ(ierr);
278d2dddef7SLois Curfman McInnes   ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J);CHKERRQ(ierr);
279d2dddef7SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult2_Private);CHKERRQ(ierr);
280d2dddef7SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy2_Private);CHKERRQ(ierr);
281d2dddef7SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView2_Private);CHKERRQ(ierr);
282d2dddef7SLois Curfman McInnes 
283d2dddef7SLois Curfman McInnes   PLogObjectParent(*J,mfctx->w);
284d2dddef7SLois Curfman McInnes   PLogObjectParent(snes,*J);
28507a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
286d2dddef7SLois Curfman McInnes }
287d2dddef7SLois Curfman McInnes 
288d2dddef7SLois Curfman McInnes #undef __FUNC__
289758f395bSBarry Smith #define __FUNC__ "SNESDefaultMatrixFreeSetParameters2"
290d2dddef7SLois Curfman McInnes /*@
291758f395bSBarry Smith    SNESDefaultMatrixFreeSetParameters2 - Sets the parameters for the approximation of
292d2dddef7SLois Curfman McInnes    matrix-vector products using finite differences.
293d2dddef7SLois Curfman McInnes 
294d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h where
295d2dddef7SLois Curfman McInnes 
296d2dddef7SLois Curfman McInnes    either the user sets h directly here, or this parameter is computed via
297d2dddef7SLois Curfman McInnes 
298d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
299d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
300d2dddef7SLois Curfman McInnes $
301d2dddef7SLois Curfman McInnes 
302d2dddef7SLois Curfman McInnes    Input Parameters:
303758f395bSBarry Smith +  mat - the matrix
304d2dddef7SLois Curfman McInnes .  error_rel - relative error (should be set to the square root of
305d2dddef7SLois Curfman McInnes      the relative error in the function evaluations)
306d2dddef7SLois Curfman McInnes .  umin - minimum allowable u-value
307758f395bSBarry Smith -  h - differencing parameter
308d2dddef7SLois Curfman McInnes 
309d2dddef7SLois Curfman McInnes    Notes:
310d2dddef7SLois Curfman McInnes    If the user sets the parameter h directly, then this value will be used
311d2dddef7SLois Curfman McInnes    instead of the default computation indicated above.
312d2dddef7SLois Curfman McInnes 
313d2dddef7SLois Curfman McInnes .keywords: SNES, matrix-free, parameters
314d2dddef7SLois Curfman McInnes 
3155a655dc6SBarry Smith .seealso: MatCreateSNESMF()
316d2dddef7SLois Curfman McInnes @*/
317758f395bSBarry Smith int SNESDefaultMatrixFreeSetParameters2(Mat mat,double error,double umin,double h)
318d2dddef7SLois Curfman McInnes {
319d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
320d2dddef7SLois Curfman McInnes   int           ierr;
321d2dddef7SLois Curfman McInnes 
32207a3eed9SLois Curfman McInnes   PetscFunctionBegin;
323d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
324d2dddef7SLois Curfman McInnes   if (ctx) {
325d2dddef7SLois Curfman McInnes     if (error != PETSC_DEFAULT) ctx->error_rel = error;
326d2dddef7SLois Curfman McInnes     if (umin  != PETSC_DEFAULT) ctx->umin = umin;
327d2dddef7SLois Curfman McInnes     if (h     != PETSC_DEFAULT) {
328d2dddef7SLois Curfman McInnes       ctx->h = h;
329d2dddef7SLois Curfman McInnes       ctx->need_h = 0;
330d2dddef7SLois Curfman McInnes     }
331d2dddef7SLois Curfman McInnes   }
33207a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
333d2dddef7SLois Curfman McInnes }
334d2dddef7SLois Curfman McInnes 
335d2dddef7SLois Curfman McInnes int SNESUnSetMatrixFreeParameter(SNES snes)
336d2dddef7SLois Curfman McInnes {
337d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
338d2dddef7SLois Curfman McInnes   int           ierr;
339d2dddef7SLois Curfman McInnes   Mat           mat;
340d2dddef7SLois Curfman McInnes 
34107a3eed9SLois Curfman McInnes   PetscFunctionBegin;
342d2dddef7SLois Curfman McInnes   ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
343d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
344d2dddef7SLois Curfman McInnes   if (ctx) ctx->need_h = 1;
34507a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
346d2dddef7SLois Curfman McInnes }
347d2dddef7SLois Curfman McInnes 
348