xref: /petsc/src/snes/interface/noise/snesmfj2.c (revision a83599f43e878c398a63bdbdbc1ebde0ef34fa4d)
163dd3a1aSKris Buschelman #define PETSCSNES_DLL
2d2dddef7SLois Curfman McInnes 
30a835dfdSSatish Balay #include "src/snes/snesimpl.h"   /*I  "petscsnes.h"   I*/
4d2dddef7SLois Curfman McInnes 
5dfbe8321SBarry Smith EXTERN PetscErrorCode DiffParameterCreate_More(SNES,Vec,void**);
6dfbe8321SBarry Smith EXTERN PetscErrorCode DiffParameterCompute_More(SNES,void*,Vec,Vec,PetscReal*,PetscReal*);
7dfbe8321SBarry Smith EXTERN PetscErrorCode DiffParameterDestroy_More(void*);
8d2dddef7SLois Curfman McInnes 
9d2dddef7SLois Curfman McInnes typedef struct {  /* default context for matrix-free SNES */
10d2dddef7SLois Curfman McInnes   SNES         snes;             /* SNES context */
11d2dddef7SLois Curfman McInnes   Vec          w;                /* work vector */
1274637425SBarry Smith   MatNullSpace sp;               /* null space context */
13145595abSSatish Balay   PetscReal    error_rel;        /* square root of relative error in computing function */
14145595abSSatish Balay   PetscReal    umin;             /* minimum allowable u'a value relative to |u|_1 */
15a7cc72afSBarry Smith   PetscTruth   jorge;            /* flag indicating use of Jorge's method for determining
16d2dddef7SLois Curfman McInnes                                    the differencing parameter */
17145595abSSatish Balay   PetscReal    h;                /* differencing parameter */
18a7cc72afSBarry Smith   PetscTruth   need_h;           /* flag indicating whether we must compute h */
19a7cc72afSBarry Smith   PetscTruth   need_err;         /* flag indicating whether we must currently compute error_rel */
20a7cc72afSBarry Smith   PetscTruth   compute_err;      /* flag indicating whether we must ever compute error_rel */
21a7cc72afSBarry Smith   PetscInt     compute_err_iter; /* last iter where we've computer error_rel */
22a7cc72afSBarry Smith   PetscInt     compute_err_freq; /* frequency of computing error_rel */
23d2dddef7SLois Curfman McInnes   void         *data;            /* implementation-specific data */
24d2dddef7SLois Curfman McInnes } MFCtx_Private;
25d2dddef7SLois Curfman McInnes 
264a2ae208SSatish Balay #undef __FUNCT__
274a2ae208SSatish Balay #define __FUNCT__ "SNESMatrixFreeDestroy2_Private" /* ADIC Ignore */
28dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeDestroy2_Private(Mat mat)
29d2dddef7SLois Curfman McInnes {
30dfbe8321SBarry Smith   PetscErrorCode ierr;
31d2dddef7SLois Curfman McInnes   MFCtx_Private  *ctx;
32d2dddef7SLois Curfman McInnes 
3307a3eed9SLois Curfman McInnes   PetscFunctionBegin;
34d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx);
35d2dddef7SLois Curfman McInnes   ierr = VecDestroy(ctx->w);CHKERRQ(ierr);
3674637425SBarry Smith   if (ctx->sp) {ierr = MatNullSpaceDestroy(ctx->sp);CHKERRQ(ierr);}
37295f7fbbSLois Curfman McInnes   if (ctx->jorge || ctx->compute_err) {ierr = DiffParameterDestroy_More(ctx->data);CHKERRQ(ierr);}
38606d414cSSatish Balay   ierr = PetscFree(ctx);CHKERRQ(ierr);
3907a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
40d2dddef7SLois Curfman McInnes }
41d2dddef7SLois Curfman McInnes 
424a2ae208SSatish Balay #undef __FUNCT__
434a2ae208SSatish Balay #define __FUNCT__ "SNESMatrixFreeView2_Private" /* ADIC Ignore */
44d2dddef7SLois Curfman McInnes /*
45d2dddef7SLois Curfman McInnes    SNESMatrixFreeView2_Private - Views matrix-free parameters.
46d2dddef7SLois Curfman McInnes  */
47dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeView2_Private(Mat J,PetscViewer viewer)
48d2dddef7SLois Curfman McInnes {
49dfbe8321SBarry Smith   PetscErrorCode ierr;
50d2dddef7SLois Curfman McInnes   MFCtx_Private  *ctx;
5132077d6dSBarry Smith   PetscTruth     iascii;
52d2dddef7SLois Curfman McInnes 
5307a3eed9SLois Curfman McInnes   PetscFunctionBegin;
54d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(J,(void **)&ctx);CHKERRQ(ierr);
5532077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
5632077d6dSBarry Smith   if (iascii) {
57b0a32e0cSBarry Smith      ierr = PetscViewerASCIIPrintf(viewer,"  SNES matrix-free approximation:\n");CHKERRQ(ierr);
58758f395bSBarry Smith      if (ctx->jorge) {
59b0a32e0cSBarry Smith        ierr = PetscViewerASCIIPrintf(viewer,"    using Jorge's method of determining differencing parameter\n");CHKERRQ(ierr);
60758f395bSBarry Smith      }
61*a83599f4SBarry Smith      ierr = PetscViewerASCIIPrintf(viewer,"    err=%G (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr);
62*a83599f4SBarry Smith      ierr = PetscViewerASCIIPrintf(viewer,"    umin=%G (minimum iterate parameter)\n",ctx->umin);CHKERRQ(ierr);
63758f395bSBarry Smith      if (ctx->compute_err) {
6477431f27SBarry Smith        ierr = PetscViewerASCIIPrintf(viewer,"    freq_err=%D (frequency for computing err)\n",ctx->compute_err_freq);CHKERRQ(ierr);
65d2dddef7SLois Curfman McInnes      }
66758f395bSBarry Smith   } else {
6779a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by SNES matrix free Jorge",((PetscObject)viewer)->type_name);
68758f395bSBarry Smith   }
6907a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
70d2dddef7SLois Curfman McInnes }
71d2dddef7SLois Curfman McInnes 
724a2ae208SSatish Balay #undef __FUNCT__
734a2ae208SSatish Balay #define __FUNCT__ "SNESMatrixFreeMult2_Private"
74d2dddef7SLois Curfman McInnes /*
75d2dddef7SLois Curfman McInnes   SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector
76d2dddef7SLois Curfman McInnes   product, y = F'(u)*a:
77d2dddef7SLois Curfman McInnes         y = (F(u + ha) - F(u)) /h,
78d2dddef7SLois Curfman McInnes   where F = nonlinear function, as set by SNESSetFunction()
79d2dddef7SLois Curfman McInnes         u = current iterate
80d2dddef7SLois Curfman McInnes         h = difference interval
81d2dddef7SLois Curfman McInnes */
82dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y)
83d2dddef7SLois Curfman McInnes {
84d2dddef7SLois Curfman McInnes   MFCtx_Private  *ctx;
85d2dddef7SLois Curfman McInnes   SNES           snes;
86145595abSSatish Balay   PetscReal      h,norm,sum,umin,noise;
8779f36460SBarry Smith   PetscScalar    hs,dot;
88d2dddef7SLois Curfman McInnes   Vec            w,U,F;
89dfbe8321SBarry Smith   PetscErrorCode ierr,(*eval_fct)(SNES,Vec,Vec);
90e6ed2b1fSLois Curfman McInnes   MPI_Comm       comm;
9177431f27SBarry Smith   PetscInt       iter;
92d2dddef7SLois Curfman McInnes 
9307a3eed9SLois Curfman McInnes   PetscFunctionBegin;
9407a3eed9SLois Curfman McInnes 
95d2dddef7SLois Curfman McInnes   /* We log matrix-free matrix-vector products separately, so that we can
96d2dddef7SLois Curfman McInnes      separate the performance monitoring from the cases that use conventional
97d2dddef7SLois Curfman McInnes      storage.  We may eventually modify event logging to associate events
98d2dddef7SLois Curfman McInnes      with particular objects, hence alleviating the more general problem. */
9946129b97SKris Buschelman   ierr = PetscLogEventBegin(MATSNESMF_Mult,a,y,0,0);CHKERRQ(ierr);
100d2dddef7SLois Curfman McInnes 
1012d0c0e3bSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
102e6ed2b1fSLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
103e6ed2b1fSLois Curfman McInnes   snes = ctx->snes;
104e6ed2b1fSLois Curfman McInnes   w    = ctx->w;
105e6ed2b1fSLois Curfman McInnes   umin = ctx->umin;
106e6ed2b1fSLois Curfman McInnes 
107d2dddef7SLois Curfman McInnes   ierr = SNESGetSolution(snes,&U);CHKERRQ(ierr);
108d2dddef7SLois Curfman McInnes   eval_fct = SNESComputeFunction;
10974637425SBarry Smith   ierr = SNESGetFunction(snes,&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
110295f7fbbSLois Curfman McInnes 
111d2dddef7SLois Curfman McInnes   /* Determine a "good" step size, h */
112d2dddef7SLois Curfman McInnes   if (ctx->need_h) {
113d2dddef7SLois Curfman McInnes 
114d2dddef7SLois Curfman McInnes     /* Use Jorge's method to compute h */
115d2dddef7SLois Curfman McInnes     if (ctx->jorge) {
116d2dddef7SLois Curfman McInnes       ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr);
117d2dddef7SLois Curfman McInnes 
118d2dddef7SLois Curfman McInnes     /* Use the Brown/Saad method to compute h */
119d2dddef7SLois Curfman McInnes     } else {
120295f7fbbSLois Curfman McInnes       /* Compute error if desired */
121295f7fbbSLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
122145595abSSatish Balay       if ((ctx->need_err) || ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter-1)%ctx->compute_err_freq)))) {
123d2dddef7SLois Curfman McInnes         /* Use Jorge's method to compute noise */
124d2dddef7SLois Curfman McInnes         ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr);
125d2dddef7SLois Curfman McInnes         ctx->error_rel = sqrt(noise);
126*a83599f4SBarry Smith         ierr = PetscVerboseInfo((snes,"SNESMatrixFreeMult2_Private: Using Jorge's noise: noise=%G, sqrt(noise)=%G, h_more=%G\n",noise,ctx->error_rel,h));CHKERRQ(ierr);
127295f7fbbSLois Curfman McInnes         ctx->compute_err_iter = iter;
1283b2a6d2fSBarry Smith         ctx->need_err = PETSC_FALSE;
129d2dddef7SLois Curfman McInnes       }
130e6ed2b1fSLois Curfman McInnes 
131691598edSBarry Smith       ierr = VecDotBegin(U,a,&dot);CHKERRQ(ierr);
132691598edSBarry Smith       ierr = VecNormBegin(a,NORM_1,&sum);CHKERRQ(ierr);
133691598edSBarry Smith       ierr = VecNormBegin(a,NORM_2,&norm);CHKERRQ(ierr);
134691598edSBarry Smith       ierr = VecDotEnd(U,a,&dot);CHKERRQ(ierr);
135691598edSBarry Smith       ierr = VecNormEnd(a,NORM_1,&sum);CHKERRQ(ierr);
136691598edSBarry Smith       ierr = VecNormEnd(a,NORM_2,&norm);CHKERRQ(ierr);
137e6ed2b1fSLois Curfman McInnes 
138d2dddef7SLois Curfman McInnes 
139d2dddef7SLois Curfman McInnes       /* Safeguard for step sizes too small */
140d2dddef7SLois Curfman McInnes       if (sum == 0.0) {dot = 1.0; norm = 1.0;}
141aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
142145595abSSatish Balay       else if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum;
143145595abSSatish Balay       else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum;
144d2dddef7SLois Curfman McInnes #else
145d2dddef7SLois Curfman McInnes       else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
146d2dddef7SLois Curfman McInnes       else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
147d2dddef7SLois Curfman McInnes #endif
148145595abSSatish Balay       h = PetscRealPart(ctx->error_rel*dot/(norm*norm));
149d2dddef7SLois Curfman McInnes     }
150d2dddef7SLois Curfman McInnes   } else {
151d2dddef7SLois Curfman McInnes     h = ctx->h;
152d2dddef7SLois Curfman McInnes   }
153*a83599f4SBarry Smith   if (!ctx->jorge || !ctx->need_h) {ierr = PetscVerboseInfo((snes,"SNESMatrixFreeMult2_Private: h = %G\n",h));CHKERRQ(ierr);}
154d2dddef7SLois Curfman McInnes 
155d2dddef7SLois Curfman McInnes   /* Evaluate function at F(u + ha) */
156db0b4b35SLois Curfman McInnes   hs = h;
1572dcb1b2aSMatthew Knepley   ierr = VecWAXPY(w,hs,a,U);CHKERRQ(ierr);
158d2dddef7SLois Curfman McInnes   ierr = eval_fct(snes,w,y);CHKERRQ(ierr);
15979f36460SBarry Smith   ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr);
16079f36460SBarry Smith   ierr = VecScale(y,1.0/hs);CHKERRQ(ierr);
16174637425SBarry Smith   if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);CHKERRQ(ierr);}
162d2dddef7SLois Curfman McInnes 
16346129b97SKris Buschelman   ierr = PetscLogEventEnd(MATSNESMF_Mult,a,y,0,0);CHKERRQ(ierr);
16407a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
165d2dddef7SLois Curfman McInnes }
166d2dddef7SLois Curfman McInnes 
1674a2ae208SSatish Balay #undef __FUNCT__
16863dd3a1aSKris Buschelman #define __FUNCT__ "SNESMatrixFreeCreate2"
169d2dddef7SLois Curfman McInnes /*@C
17063dd3a1aSKris Buschelman    SNESMatrixFreeCreate2 - Creates a matrix-free matrix
171d2dddef7SLois Curfman McInnes    context for use with a SNES solver.  This matrix can be used as
172d2dddef7SLois Curfman McInnes    the Jacobian argument for the routine SNESSetJacobian().
173d2dddef7SLois Curfman McInnes 
174d2dddef7SLois Curfman McInnes    Input Parameters:
175d2dddef7SLois Curfman McInnes .  snes - the SNES context
176d2dddef7SLois Curfman McInnes .  x - vector where SNES solution is to be stored.
177d2dddef7SLois Curfman McInnes 
178d2dddef7SLois Curfman McInnes    Output Parameter:
179d2dddef7SLois Curfman McInnes .  J - the matrix-free matrix
180d2dddef7SLois Curfman McInnes 
18190c26489SBarry Smith    Level: advanced
18290c26489SBarry Smith 
183d2dddef7SLois Curfman McInnes    Notes:
184d2dddef7SLois Curfman McInnes    The matrix-free matrix context merely contains the function pointers
185d2dddef7SLois Curfman McInnes    and work space for performing finite difference approximations of
186d2dddef7SLois Curfman McInnes    Jacobian-vector products, J(u)*a, via
187d2dddef7SLois Curfman McInnes 
188d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h,
189d2dddef7SLois Curfman McInnes $   where by default
190d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
191d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
192d2dddef7SLois Curfman McInnes $   where
193d2dddef7SLois Curfman McInnes $        error_rel = square root of relative error in
194d2dddef7SLois Curfman McInnes $                    function evaluation
195d2dddef7SLois Curfman McInnes $        umin = minimum iterate parameter
196d2dddef7SLois Curfman McInnes $   Alternatively, the differencing parameter, h, can be set using
197e6ed2b1fSLois Curfman McInnes $   Jorge's nifty new strategy if one specifies the option
198e6ed2b1fSLois Curfman McInnes $          -snes_mf_jorge
199d2dddef7SLois Curfman McInnes 
2005a655dc6SBarry Smith    The user can set these parameters via MatSNESMFSetFunctionError().
201d2dddef7SLois Curfman McInnes    See the nonlinear solvers chapter of the users manual for details.
202d2dddef7SLois Curfman McInnes 
203d2dddef7SLois Curfman McInnes    The user should call MatDestroy() when finished with the matrix-free
204d2dddef7SLois Curfman McInnes    matrix context.
205d2dddef7SLois Curfman McInnes 
206d2dddef7SLois Curfman McInnes    Options Database Keys:
207d2dddef7SLois Curfman McInnes $  -snes_mf_err <error_rel>
208d2dddef7SLois Curfman McInnes $  -snes_mf_unim <umin>
209295f7fbbSLois Curfman McInnes $  -snes_mf_compute_err
210295f7fbbSLois Curfman McInnes $  -snes_mf_freq_err <freq>
211e6ed2b1fSLois Curfman McInnes $  -snes_mf_jorge
212d2dddef7SLois Curfman McInnes 
213d2dddef7SLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix
214d2dddef7SLois Curfman McInnes 
2155a655dc6SBarry Smith .seealso: MatDestroy(), MatSNESMFSetFunctionError()
216d2dddef7SLois Curfman McInnes @*/
21763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES snes,Vec x,Mat *J)
218d2dddef7SLois Curfman McInnes {
219d2dddef7SLois Curfman McInnes   MPI_Comm       comm;
220d2dddef7SLois Curfman McInnes   MFCtx_Private  *mfctx;
2216849ba73SBarry Smith   PetscErrorCode ierr;
222a7cc72afSBarry Smith   PetscInt       n,nloc;
223145595abSSatish Balay   PetscTruth     flg;
224d2dddef7SLois Curfman McInnes   char           p[64];
225d2dddef7SLois Curfman McInnes 
22607a3eed9SLois Curfman McInnes   PetscFunctionBegin;
22782502324SSatish Balay   ierr = PetscNew(MFCtx_Private,&mfctx);CHKERRQ(ierr);
22852e6d16bSBarry Smith   ierr = PetscLogObjectMemory(snes,sizeof(MFCtx_Private));CHKERRQ(ierr);
229d2dddef7SLois Curfman McInnes   mfctx->sp   = 0;
230d2dddef7SLois Curfman McInnes   mfctx->snes = snes;
23177d8c4bbSBarry Smith   mfctx->error_rel        = PETSC_SQRT_MACHINE_EPSILON;
232d2dddef7SLois Curfman McInnes   mfctx->umin             = 1.e-6;
233d2dddef7SLois Curfman McInnes   mfctx->h                = 0.0;
234a7cc72afSBarry Smith   mfctx->need_h           = PETSC_TRUE;
235a7cc72afSBarry Smith   mfctx->need_err         = PETSC_FALSE;
236a7cc72afSBarry Smith   mfctx->compute_err      = PETSC_FALSE;
237295f7fbbSLois Curfman McInnes   mfctx->compute_err_freq = 0;
238295f7fbbSLois Curfman McInnes   mfctx->compute_err_iter = -1;
23987828ca2SBarry Smith   ierr = PetscOptionsGetReal(snes->prefix,"-snes_mf_err",&mfctx->error_rel,PETSC_NULL);CHKERRQ(ierr);
24087828ca2SBarry Smith   ierr = PetscOptionsGetReal(snes->prefix,"-snes_mf_umin",&mfctx->umin,PETSC_NULL);CHKERRQ(ierr);
241a7cc72afSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_jorge",&mfctx->jorge);CHKERRQ(ierr);
242a7cc72afSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_compute_err",&mfctx->compute_err);CHKERRQ(ierr);
243b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(snes->prefix,"-snes_mf_freq_err",&mfctx->compute_err_freq,&flg);CHKERRQ(ierr);
244295f7fbbSLois Curfman McInnes   if (flg) {
245295f7fbbSLois Curfman McInnes     if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0;
2463b2a6d2fSBarry Smith     mfctx->compute_err = PETSC_TRUE;
247295f7fbbSLois Curfman McInnes   }
2483b2a6d2fSBarry Smith   if (mfctx->compute_err) mfctx->need_err = PETSC_TRUE;
249295f7fbbSLois Curfman McInnes   if (mfctx->jorge || mfctx->compute_err) {
250d2dddef7SLois Curfman McInnes     ierr = DiffParameterCreate_More(snes,x,&mfctx->data);CHKERRQ(ierr);
251d2dddef7SLois Curfman McInnes   } else mfctx->data = 0;
252d2dddef7SLois Curfman McInnes 
253b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr);
254549d3d68SSatish Balay   ierr = PetscStrcpy(p,"-");CHKERRQ(ierr);
255d2dddef7SLois Curfman McInnes   if (snes->prefix) PetscStrcat(p,snes->prefix);
256d2dddef7SLois Curfman McInnes   if (flg) {
257d132466eSBarry Smith     ierr = PetscPrintf(snes->comm," Matrix-free Options (via SNES):\n");CHKERRQ(ierr);
258*a83599f4SBarry 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);
259*a83599f4SBarry Smith     ierr = PetscPrintf(snes->comm,"   %ssnes_mf_umin <umin>: see users manual (default %G)\n",p,mfctx->umin);CHKERRQ(ierr);
260d132466eSBarry Smith     ierr = PetscPrintf(snes->comm,"   %ssnes_mf_jorge: use Jorge More's method\n",p);CHKERRQ(ierr);
261d132466eSBarry Smith     ierr = PetscPrintf(snes->comm,"   %ssnes_mf_compute_err: compute sqrt or relative error in function\n",p);CHKERRQ(ierr);
262d132466eSBarry Smith     ierr = PetscPrintf(snes->comm,"   %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)\n",p);CHKERRQ(ierr);
263d132466eSBarry Smith     ierr = PetscPrintf(snes->comm,"   %ssnes_mf_noise_file <file>: set file for printing noise info\n",p);CHKERRQ(ierr);
264d2dddef7SLois Curfman McInnes   }
265d2dddef7SLois Curfman McInnes   ierr = VecDuplicate(x,&mfctx->w);CHKERRQ(ierr);
266d2dddef7SLois Curfman McInnes   ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr);
267d2dddef7SLois Curfman McInnes   ierr = VecGetSize(x,&n);CHKERRQ(ierr);
268d2dddef7SLois Curfman McInnes   ierr = VecGetLocalSize(x,&nloc);CHKERRQ(ierr);
269f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,J);CHKERRQ(ierr);
270f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*J,nloc,n,n,n);CHKERRQ(ierr);
271f204ca49SKris Buschelman   ierr = MatSetType(*J,MATSHELL);CHKERRQ(ierr);
272f204ca49SKris Buschelman   ierr = MatShellSetContext(*J,mfctx);CHKERRQ(ierr);
273c134de8dSSatish Balay   ierr = MatShellSetOperation(*J,MATOP_MULT,(void(*)(void))SNESMatrixFreeMult2_Private);CHKERRQ(ierr);
274c134de8dSSatish Balay   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void(*)(void))SNESMatrixFreeDestroy2_Private);CHKERRQ(ierr);
275c134de8dSSatish Balay   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void(*)(void))SNESMatrixFreeView2_Private);CHKERRQ(ierr);
276d2dddef7SLois Curfman McInnes 
27752e6d16bSBarry Smith   ierr = PetscLogObjectParent(*J,mfctx->w);CHKERRQ(ierr);
27852e6d16bSBarry Smith   ierr = PetscLogObjectParent(snes,*J);CHKERRQ(ierr);
27907a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
280d2dddef7SLois Curfman McInnes }
281d2dddef7SLois Curfman McInnes 
2824a2ae208SSatish Balay #undef __FUNCT__
2834a2ae208SSatish Balay #define __FUNCT__ "SNESDefaultMatrixFreeSetParameters2"
28490c26489SBarry Smith /*@C
285758f395bSBarry Smith    SNESDefaultMatrixFreeSetParameters2 - Sets the parameters for the approximation of
286d2dddef7SLois Curfman McInnes    matrix-vector products using finite differences.
287d2dddef7SLois Curfman McInnes 
288d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h where
289d2dddef7SLois Curfman McInnes 
290d2dddef7SLois Curfman McInnes    either the user sets h directly here, or this parameter is computed via
291d2dddef7SLois Curfman McInnes 
292d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
293d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
294d2dddef7SLois Curfman McInnes $
295d2dddef7SLois Curfman McInnes 
296d2dddef7SLois Curfman McInnes    Input Parameters:
297758f395bSBarry Smith +  mat - the matrix
298d2dddef7SLois Curfman McInnes .  error_rel - relative error (should be set to the square root of
299d2dddef7SLois Curfman McInnes      the relative error in the function evaluations)
300d2dddef7SLois Curfman McInnes .  umin - minimum allowable u-value
301758f395bSBarry Smith -  h - differencing parameter
302d2dddef7SLois Curfman McInnes 
30390c26489SBarry Smith    Level: advanced
30490c26489SBarry Smith 
305d2dddef7SLois Curfman McInnes    Notes:
306d2dddef7SLois Curfman McInnes    If the user sets the parameter h directly, then this value will be used
307d2dddef7SLois Curfman McInnes    instead of the default computation indicated above.
308d2dddef7SLois Curfman McInnes 
309d2dddef7SLois Curfman McInnes .keywords: SNES, matrix-free, parameters
310d2dddef7SLois Curfman McInnes 
3115a655dc6SBarry Smith .seealso: MatCreateSNESMF()
312d2dddef7SLois Curfman McInnes @*/
31363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeSetParameters2(Mat mat,PetscReal error,PetscReal umin,PetscReal h)
314d2dddef7SLois Curfman McInnes {
315d2dddef7SLois Curfman McInnes   MFCtx_Private  *ctx;
316dfbe8321SBarry Smith   PetscErrorCode ierr;
317d2dddef7SLois Curfman McInnes 
31807a3eed9SLois Curfman McInnes   PetscFunctionBegin;
319d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
320d2dddef7SLois Curfman McInnes   if (ctx) {
321d2dddef7SLois Curfman McInnes     if (error != PETSC_DEFAULT) ctx->error_rel = error;
322d2dddef7SLois Curfman McInnes     if (umin  != PETSC_DEFAULT) ctx->umin = umin;
323d2dddef7SLois Curfman McInnes     if (h     != PETSC_DEFAULT) {
324d2dddef7SLois Curfman McInnes       ctx->h = h;
325a7cc72afSBarry Smith       ctx->need_h = PETSC_FALSE;
326d2dddef7SLois Curfman McInnes     }
327d2dddef7SLois Curfman McInnes   }
32807a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
329d2dddef7SLois Curfman McInnes }
330d2dddef7SLois Curfman McInnes 
33163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESUnSetMatrixFreeParameter(SNES snes)
332d2dddef7SLois Curfman McInnes {
333d2dddef7SLois Curfman McInnes   MFCtx_Private  *ctx;
334dfbe8321SBarry Smith   PetscErrorCode ierr;
335d2dddef7SLois Curfman McInnes   Mat            mat;
336d2dddef7SLois Curfman McInnes 
33707a3eed9SLois Curfman McInnes   PetscFunctionBegin;
33874637425SBarry Smith   ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
339d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
340a7cc72afSBarry Smith   if (ctx) ctx->need_h = PETSC_TRUE;
34107a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
342d2dddef7SLois Curfman McInnes }
343d2dddef7SLois Curfman McInnes 
344