xref: /petsc/src/snes/interface/noise/snesmfj2.c (revision acfcf0e5ba359b58e6a8a7af3f239cd7334278e8)
163dd3a1aSKris Buschelman #define PETSCSNES_DLL
2d2dddef7SLois Curfman McInnes 
37c4f633dSBarry Smith #include "private/snesimpl.h"   /*I  "petscsnes.h"   I*/
46370053bSBarry Smith /* matimpl.h is needed only for logging of matrix operation */
57c4f633dSBarry Smith #include "private/matimpl.h"
6d2dddef7SLois Curfman McInnes 
7dfbe8321SBarry Smith EXTERN PetscErrorCode DiffParameterCreate_More(SNES,Vec,void**);
8dfbe8321SBarry Smith EXTERN PetscErrorCode DiffParameterCompute_More(SNES,void*,Vec,Vec,PetscReal*,PetscReal*);
9dfbe8321SBarry Smith EXTERN PetscErrorCode 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 */
1474637425SBarry Smith   MatNullSpace sp;               /* null space context */
15145595abSSatish Balay   PetscReal    error_rel;        /* square root of relative error in computing function */
16145595abSSatish Balay   PetscReal    umin;             /* minimum allowable u'a value relative to |u|_1 */
17ace3abfcSBarry Smith   PetscBool    jorge;            /* flag indicating use of Jorge's method for determining
18d2dddef7SLois Curfman McInnes                                    the differencing parameter */
19145595abSSatish Balay   PetscReal    h;                /* differencing parameter */
20ace3abfcSBarry Smith   PetscBool    need_h;           /* flag indicating whether we must compute h */
21ace3abfcSBarry Smith   PetscBool    need_err;         /* flag indicating whether we must currently compute error_rel */
22ace3abfcSBarry Smith   PetscBool    compute_err;      /* flag indicating whether we must ever compute error_rel */
23a7cc72afSBarry Smith   PetscInt     compute_err_iter; /* last iter where we've computer error_rel */
24a7cc72afSBarry Smith   PetscInt     compute_err_freq; /* frequency of computing error_rel */
25d2dddef7SLois Curfman McInnes   void         *data;            /* implementation-specific data */
26d2dddef7SLois Curfman McInnes } MFCtx_Private;
27d2dddef7SLois Curfman McInnes 
284a2ae208SSatish Balay #undef __FUNCT__
294a2ae208SSatish Balay #define __FUNCT__ "SNESMatrixFreeDestroy2_Private" /* ADIC Ignore */
30dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeDestroy2_Private(Mat mat)
31d2dddef7SLois Curfman McInnes {
32dfbe8321SBarry Smith   PetscErrorCode 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);
3874637425SBarry Smith   if (ctx->sp) {ierr = MatNullSpaceDestroy(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 
444a2ae208SSatish Balay #undef __FUNCT__
454a2ae208SSatish Balay #define __FUNCT__ "SNESMatrixFreeView2_Private" /* ADIC Ignore */
46d2dddef7SLois Curfman McInnes /*
47d2dddef7SLois Curfman McInnes    SNESMatrixFreeView2_Private - Views matrix-free parameters.
48d2dddef7SLois Curfman McInnes  */
49dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeView2_Private(Mat J,PetscViewer viewer)
50d2dddef7SLois Curfman McInnes {
51dfbe8321SBarry Smith   PetscErrorCode ierr;
52d2dddef7SLois Curfman McInnes   MFCtx_Private  *ctx;
53ace3abfcSBarry Smith   PetscBool      iascii;
54d2dddef7SLois Curfman McInnes 
5507a3eed9SLois Curfman McInnes   PetscFunctionBegin;
56d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(J,(void **)&ctx);CHKERRQ(ierr);
572692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
5832077d6dSBarry Smith   if (iascii) {
59b0a32e0cSBarry Smith      ierr = PetscViewerASCIIPrintf(viewer,"  SNES matrix-free approximation:\n");CHKERRQ(ierr);
60758f395bSBarry Smith      if (ctx->jorge) {
61b0a32e0cSBarry Smith        ierr = PetscViewerASCIIPrintf(viewer,"    using Jorge's method of determining differencing parameter\n");CHKERRQ(ierr);
62758f395bSBarry Smith      }
63a83599f4SBarry Smith      ierr = PetscViewerASCIIPrintf(viewer,"    err=%G (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr);
64a83599f4SBarry Smith      ierr = PetscViewerASCIIPrintf(viewer,"    umin=%G (minimum iterate parameter)\n",ctx->umin);CHKERRQ(ierr);
65758f395bSBarry Smith      if (ctx->compute_err) {
6677431f27SBarry Smith        ierr = PetscViewerASCIIPrintf(viewer,"    freq_err=%D (frequency for computing err)\n",ctx->compute_err_freq);CHKERRQ(ierr);
67d2dddef7SLois Curfman McInnes      }
68758f395bSBarry Smith   } else {
69e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"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 
744a2ae208SSatish Balay #undef __FUNCT__
754a2ae208SSatish Balay #define __FUNCT__ "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 */
84dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y)
85d2dddef7SLois Curfman McInnes {
86d2dddef7SLois Curfman McInnes   MFCtx_Private  *ctx;
87d2dddef7SLois Curfman McInnes   SNES           snes;
88145595abSSatish Balay   PetscReal      h,norm,sum,umin,noise;
8979f36460SBarry Smith   PetscScalar    hs,dot;
90d2dddef7SLois Curfman McInnes   Vec            w,U,F;
91dfbe8321SBarry Smith   PetscErrorCode ierr,(*eval_fct)(SNES,Vec,Vec);
92e6ed2b1fSLois Curfman McInnes   MPI_Comm       comm;
9377431f27SBarry Smith   PetscInt       iter;
94d2dddef7SLois Curfman McInnes 
9507a3eed9SLois Curfman McInnes   PetscFunctionBegin;
9607a3eed9SLois Curfman McInnes 
97d2dddef7SLois Curfman McInnes   /* We log matrix-free matrix-vector products separately, so that we can
98d2dddef7SLois Curfman McInnes      separate the performance monitoring from the cases that use conventional
99d2dddef7SLois Curfman McInnes      storage.  We may eventually modify event logging to associate events
100d2dddef7SLois Curfman McInnes      with particular objects, hence alleviating the more general problem. */
1013ec795f1SBarry Smith   ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
102d2dddef7SLois Curfman McInnes 
1032d0c0e3bSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
104e6ed2b1fSLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
105e6ed2b1fSLois Curfman McInnes   snes = ctx->snes;
106e6ed2b1fSLois Curfman McInnes   w    = ctx->w;
107e6ed2b1fSLois Curfman McInnes   umin = ctx->umin;
108e6ed2b1fSLois Curfman McInnes 
109d2dddef7SLois Curfman McInnes   ierr = SNESGetSolution(snes,&U);CHKERRQ(ierr);
110d2dddef7SLois Curfman McInnes   eval_fct = SNESComputeFunction;
11174637425SBarry Smith   ierr = SNESGetFunction(snes,&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
112295f7fbbSLois Curfman McInnes 
113d2dddef7SLois Curfman McInnes   /* Determine a "good" step size, h */
114d2dddef7SLois Curfman McInnes   if (ctx->need_h) {
115d2dddef7SLois Curfman McInnes 
116d2dddef7SLois Curfman McInnes     /* Use Jorge's method to compute h */
117d2dddef7SLois Curfman McInnes     if (ctx->jorge) {
118d2dddef7SLois Curfman McInnes       ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr);
119d2dddef7SLois Curfman McInnes 
120d2dddef7SLois Curfman McInnes     /* Use the Brown/Saad method to compute h */
121d2dddef7SLois Curfman McInnes     } else {
122295f7fbbSLois Curfman McInnes       /* Compute error if desired */
123295f7fbbSLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
124145595abSSatish Balay       if ((ctx->need_err) || ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter-1)%ctx->compute_err_freq)))) {
125d2dddef7SLois Curfman McInnes         /* Use Jorge's method to compute noise */
126d2dddef7SLois Curfman McInnes         ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr);
127d2dddef7SLois Curfman McInnes         ctx->error_rel = sqrt(noise);
128ae15b995SBarry Smith         ierr = PetscInfo3(snes,"Using Jorge's noise: noise=%G, sqrt(noise)=%G, h_more=%G\n",noise,ctx->error_rel,h);CHKERRQ(ierr);
129295f7fbbSLois Curfman McInnes         ctx->compute_err_iter = iter;
1303b2a6d2fSBarry Smith         ctx->need_err = PETSC_FALSE;
131d2dddef7SLois Curfman McInnes       }
132e6ed2b1fSLois Curfman McInnes 
133691598edSBarry Smith       ierr = VecDotBegin(U,a,&dot);CHKERRQ(ierr);
134691598edSBarry Smith       ierr = VecNormBegin(a,NORM_1,&sum);CHKERRQ(ierr);
135691598edSBarry Smith       ierr = VecNormBegin(a,NORM_2,&norm);CHKERRQ(ierr);
136691598edSBarry Smith       ierr = VecDotEnd(U,a,&dot);CHKERRQ(ierr);
137691598edSBarry Smith       ierr = VecNormEnd(a,NORM_1,&sum);CHKERRQ(ierr);
138691598edSBarry Smith       ierr = VecNormEnd(a,NORM_2,&norm);CHKERRQ(ierr);
139e6ed2b1fSLois Curfman McInnes 
140d2dddef7SLois Curfman McInnes 
141d2dddef7SLois Curfman McInnes       /* Safeguard for step sizes too small */
142d2dddef7SLois Curfman McInnes       if (sum == 0.0) {dot = 1.0; norm = 1.0;}
143aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
144145595abSSatish Balay       else if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum;
145145595abSSatish Balay       else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum;
146d2dddef7SLois Curfman McInnes #else
147d2dddef7SLois Curfman McInnes       else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
148d2dddef7SLois Curfman McInnes       else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
149d2dddef7SLois Curfman McInnes #endif
150145595abSSatish Balay       h = PetscRealPart(ctx->error_rel*dot/(norm*norm));
151d2dddef7SLois Curfman McInnes     }
152d2dddef7SLois Curfman McInnes   } else {
153d2dddef7SLois Curfman McInnes     h = ctx->h;
154d2dddef7SLois Curfman McInnes   }
155ae15b995SBarry Smith   if (!ctx->jorge || !ctx->need_h) {ierr = PetscInfo1(snes,"h = %G\n",h);CHKERRQ(ierr);}
156d2dddef7SLois Curfman McInnes 
157d2dddef7SLois Curfman McInnes   /* Evaluate function at F(u + ha) */
158db0b4b35SLois Curfman McInnes   hs = h;
1592dcb1b2aSMatthew Knepley   ierr = VecWAXPY(w,hs,a,U);CHKERRQ(ierr);
160d2dddef7SLois Curfman McInnes   ierr = eval_fct(snes,w,y);CHKERRQ(ierr);
16179f36460SBarry Smith   ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr);
16279f36460SBarry Smith   ierr = VecScale(y,1.0/hs);CHKERRQ(ierr);
16374637425SBarry Smith   if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);CHKERRQ(ierr);}
164d2dddef7SLois Curfman McInnes 
1653ec795f1SBarry Smith   ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
16607a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
167d2dddef7SLois Curfman McInnes }
168d2dddef7SLois Curfman McInnes 
1694a2ae208SSatish Balay #undef __FUNCT__
17063dd3a1aSKris Buschelman #define __FUNCT__ "SNESMatrixFreeCreate2"
171d2dddef7SLois Curfman McInnes /*@C
17263dd3a1aSKris Buschelman    SNESMatrixFreeCreate2 - Creates a matrix-free matrix
173d2dddef7SLois Curfman McInnes    context for use with a SNES solver.  This matrix can be used as
174d2dddef7SLois Curfman McInnes    the Jacobian argument for the routine SNESSetJacobian().
175d2dddef7SLois Curfman McInnes 
176d2dddef7SLois Curfman McInnes    Input Parameters:
177d2dddef7SLois Curfman McInnes .  snes - the SNES context
178d2dddef7SLois Curfman McInnes .  x - vector where SNES solution is to be stored.
179d2dddef7SLois Curfman McInnes 
180d2dddef7SLois Curfman McInnes    Output Parameter:
181d2dddef7SLois Curfman McInnes .  J - the matrix-free matrix
182d2dddef7SLois Curfman McInnes 
18390c26489SBarry Smith    Level: advanced
18490c26489SBarry Smith 
185d2dddef7SLois Curfman McInnes    Notes:
186d2dddef7SLois Curfman McInnes    The matrix-free matrix context merely contains the function pointers
187d2dddef7SLois Curfman McInnes    and work space for performing finite difference approximations of
188d2dddef7SLois Curfman McInnes    Jacobian-vector products, J(u)*a, via
189d2dddef7SLois Curfman McInnes 
190d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h,
191d2dddef7SLois Curfman McInnes $   where by default
192d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
193d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
194d2dddef7SLois Curfman McInnes $   where
195d2dddef7SLois Curfman McInnes $        error_rel = square root of relative error in
196d2dddef7SLois Curfman McInnes $                    function evaluation
197d2dddef7SLois Curfman McInnes $        umin = minimum iterate parameter
198d2dddef7SLois Curfman McInnes $   Alternatively, the differencing parameter, h, can be set using
199e6ed2b1fSLois Curfman McInnes $   Jorge's nifty new strategy if one specifies the option
200e6ed2b1fSLois Curfman McInnes $          -snes_mf_jorge
201d2dddef7SLois Curfman McInnes 
2023ec795f1SBarry Smith    The user can set these parameters via MatMFFDSetFunctionError().
2030598bfebSBarry Smith    See the <A href="../../docs/manual.pdf#nameddest=ch_snes">SNES chapter of the users manual</A>.
204d2dddef7SLois Curfman McInnes 
205d2dddef7SLois Curfman McInnes    The user should call MatDestroy() when finished with the matrix-free
206d2dddef7SLois Curfman McInnes    matrix context.
207d2dddef7SLois Curfman McInnes 
208d2dddef7SLois Curfman McInnes    Options Database Keys:
209d2dddef7SLois Curfman McInnes $  -snes_mf_err <error_rel>
210d2dddef7SLois Curfman McInnes $  -snes_mf_unim <umin>
211295f7fbbSLois Curfman McInnes $  -snes_mf_compute_err
212295f7fbbSLois Curfman McInnes $  -snes_mf_freq_err <freq>
213e6ed2b1fSLois Curfman McInnes $  -snes_mf_jorge
214d2dddef7SLois Curfman McInnes 
215d2dddef7SLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix
216d2dddef7SLois Curfman McInnes 
2173ec795f1SBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError()
218d2dddef7SLois Curfman McInnes @*/
21963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES snes,Vec x,Mat *J)
220d2dddef7SLois Curfman McInnes {
221d2dddef7SLois Curfman McInnes   MPI_Comm       comm;
222d2dddef7SLois Curfman McInnes   MFCtx_Private  *mfctx;
2236849ba73SBarry Smith   PetscErrorCode ierr;
224a7cc72afSBarry Smith   PetscInt       n,nloc;
225ace3abfcSBarry Smith   PetscBool      flg;
226d2dddef7SLois Curfman McInnes   char           p[64];
227d2dddef7SLois Curfman McInnes 
22807a3eed9SLois Curfman McInnes   PetscFunctionBegin;
22938f2d2fdSLisandro Dalcin   ierr = PetscNewLog(snes,MFCtx_Private,&mfctx);CHKERRQ(ierr);
230d2dddef7SLois Curfman McInnes   mfctx->sp   = 0;
231d2dddef7SLois Curfman McInnes   mfctx->snes = snes;
23277d8c4bbSBarry Smith   mfctx->error_rel        = PETSC_SQRT_MACHINE_EPSILON;
233d2dddef7SLois Curfman McInnes   mfctx->umin             = 1.e-6;
234d2dddef7SLois Curfman McInnes   mfctx->h                = 0.0;
235a7cc72afSBarry Smith   mfctx->need_h           = PETSC_TRUE;
236a7cc72afSBarry Smith   mfctx->need_err         = PETSC_FALSE;
237a7cc72afSBarry Smith   mfctx->compute_err      = PETSC_FALSE;
238295f7fbbSLois Curfman McInnes   mfctx->compute_err_freq = 0;
239295f7fbbSLois Curfman McInnes   mfctx->compute_err_iter = -1;
24090d69ab7SBarry Smith   mfctx->compute_err      = PETSC_FALSE;
24190d69ab7SBarry Smith   mfctx->jorge            = PETSC_FALSE;
2427adad957SLisandro Dalcin   ierr = PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_mf_err",&mfctx->error_rel,PETSC_NULL);CHKERRQ(ierr);
2437adad957SLisandro Dalcin   ierr = PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_mf_umin",&mfctx->umin,PETSC_NULL);CHKERRQ(ierr);
244*acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_mf_jorge",&mfctx->jorge,PETSC_NULL);CHKERRQ(ierr);
245*acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_mf_compute_err",&mfctx->compute_err,PETSC_NULL);CHKERRQ(ierr);
2467adad957SLisandro Dalcin   ierr = PetscOptionsGetInt(((PetscObject)snes)->prefix,"-snes_mf_freq_err",&mfctx->compute_err_freq,&flg);CHKERRQ(ierr);
247295f7fbbSLois Curfman McInnes   if (flg) {
248295f7fbbSLois Curfman McInnes     if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0;
2493b2a6d2fSBarry Smith     mfctx->compute_err = PETSC_TRUE;
250295f7fbbSLois Curfman McInnes   }
2513b2a6d2fSBarry Smith   if (mfctx->compute_err) mfctx->need_err = PETSC_TRUE;
252295f7fbbSLois Curfman McInnes   if (mfctx->jorge || mfctx->compute_err) {
253d2dddef7SLois Curfman McInnes     ierr = DiffParameterCreate_More(snes,x,&mfctx->data);CHKERRQ(ierr);
254d2dddef7SLois Curfman McInnes   } else mfctx->data = 0;
255d2dddef7SLois Curfman McInnes 
256b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr);
257549d3d68SSatish Balay   ierr = PetscStrcpy(p,"-");CHKERRQ(ierr);
2587adad957SLisandro Dalcin   if (((PetscObject)snes)->prefix) PetscStrcat(p,((PetscObject)snes)->prefix);
259d2dddef7SLois Curfman McInnes   if (flg) {
2607adad957SLisandro Dalcin     ierr = PetscPrintf(((PetscObject)snes)->comm," Matrix-free Options (via SNES):\n");CHKERRQ(ierr);
2617adad957SLisandro Dalcin     ierr = PetscPrintf(((PetscObject)snes)->comm,"   %ssnes_mf_err <err>: set sqrt of relative error in function (default %G)\n",p,mfctx->error_rel);CHKERRQ(ierr);
2627adad957SLisandro Dalcin     ierr = PetscPrintf(((PetscObject)snes)->comm,"   %ssnes_mf_umin <umin>: see users manual (default %G)\n",p,mfctx->umin);CHKERRQ(ierr);
2637adad957SLisandro Dalcin     ierr = PetscPrintf(((PetscObject)snes)->comm,"   %ssnes_mf_jorge: use Jorge More's method\n",p);CHKERRQ(ierr);
2647adad957SLisandro Dalcin     ierr = PetscPrintf(((PetscObject)snes)->comm,"   %ssnes_mf_compute_err: compute sqrt or relative error in function\n",p);CHKERRQ(ierr);
2657adad957SLisandro Dalcin     ierr = PetscPrintf(((PetscObject)snes)->comm,"   %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)\n",p);CHKERRQ(ierr);
2667adad957SLisandro Dalcin     ierr = PetscPrintf(((PetscObject)snes)->comm,"   %ssnes_mf_noise_file <file>: set file for printing noise info\n",p);CHKERRQ(ierr);
267d2dddef7SLois Curfman McInnes   }
268d2dddef7SLois Curfman McInnes   ierr = VecDuplicate(x,&mfctx->w);CHKERRQ(ierr);
269d2dddef7SLois Curfman McInnes   ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr);
270d2dddef7SLois Curfman McInnes   ierr = VecGetSize(x,&n);CHKERRQ(ierr);
271d2dddef7SLois Curfman McInnes   ierr = VecGetLocalSize(x,&nloc);CHKERRQ(ierr);
272f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,J);CHKERRQ(ierr);
273f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*J,nloc,n,n,n);CHKERRQ(ierr);
274f204ca49SKris Buschelman   ierr = MatSetType(*J,MATSHELL);CHKERRQ(ierr);
275f204ca49SKris Buschelman   ierr = MatShellSetContext(*J,mfctx);CHKERRQ(ierr);
276c134de8dSSatish Balay   ierr = MatShellSetOperation(*J,MATOP_MULT,(void(*)(void))SNESMatrixFreeMult2_Private);CHKERRQ(ierr);
277c134de8dSSatish Balay   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void(*)(void))SNESMatrixFreeDestroy2_Private);CHKERRQ(ierr);
278c134de8dSSatish Balay   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void(*)(void))SNESMatrixFreeView2_Private);CHKERRQ(ierr);
279d2dddef7SLois Curfman McInnes 
28052e6d16bSBarry Smith   ierr = PetscLogObjectParent(*J,mfctx->w);CHKERRQ(ierr);
28152e6d16bSBarry Smith   ierr = PetscLogObjectParent(snes,*J);CHKERRQ(ierr);
28207a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
283d2dddef7SLois Curfman McInnes }
284d2dddef7SLois Curfman McInnes 
2854a2ae208SSatish Balay #undef __FUNCT__
2864a2ae208SSatish Balay #define __FUNCT__ "SNESDefaultMatrixFreeSetParameters2"
28790c26489SBarry Smith /*@C
288758f395bSBarry Smith    SNESDefaultMatrixFreeSetParameters2 - Sets the parameters for the approximation of
289d2dddef7SLois Curfman McInnes    matrix-vector products using finite differences.
290d2dddef7SLois Curfman McInnes 
291d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h where
292d2dddef7SLois Curfman McInnes 
293d2dddef7SLois Curfman McInnes    either the user sets h directly here, or this parameter is computed via
294d2dddef7SLois Curfman McInnes 
295d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
296d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
297d2dddef7SLois Curfman McInnes $
298d2dddef7SLois Curfman McInnes 
299d2dddef7SLois Curfman McInnes    Input Parameters:
300758f395bSBarry Smith +  mat - the matrix
301d2dddef7SLois Curfman McInnes .  error_rel - relative error (should be set to the square root of
302d2dddef7SLois Curfman McInnes      the relative error in the function evaluations)
303d2dddef7SLois Curfman McInnes .  umin - minimum allowable u-value
304758f395bSBarry Smith -  h - differencing parameter
305d2dddef7SLois Curfman McInnes 
30690c26489SBarry Smith    Level: advanced
30790c26489SBarry Smith 
308d2dddef7SLois Curfman McInnes    Notes:
309d2dddef7SLois Curfman McInnes    If the user sets the parameter h directly, then this value will be used
310d2dddef7SLois Curfman McInnes    instead of the default computation indicated above.
311d2dddef7SLois Curfman McInnes 
312d2dddef7SLois Curfman McInnes .keywords: SNES, matrix-free, parameters
313d2dddef7SLois Curfman McInnes 
3145a655dc6SBarry Smith .seealso: MatCreateSNESMF()
315d2dddef7SLois Curfman McInnes @*/
31663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeSetParameters2(Mat mat,PetscReal error,PetscReal umin,PetscReal h)
317d2dddef7SLois Curfman McInnes {
318d2dddef7SLois Curfman McInnes   MFCtx_Private  *ctx;
319dfbe8321SBarry Smith   PetscErrorCode ierr;
320d2dddef7SLois Curfman McInnes 
32107a3eed9SLois Curfman McInnes   PetscFunctionBegin;
322d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
323d2dddef7SLois Curfman McInnes   if (ctx) {
324d2dddef7SLois Curfman McInnes     if (error != PETSC_DEFAULT) ctx->error_rel = error;
325d2dddef7SLois Curfman McInnes     if (umin  != PETSC_DEFAULT) ctx->umin = umin;
326d2dddef7SLois Curfman McInnes     if (h     != PETSC_DEFAULT) {
327d2dddef7SLois Curfman McInnes       ctx->h = h;
328a7cc72afSBarry Smith       ctx->need_h = PETSC_FALSE;
329d2dddef7SLois Curfman McInnes     }
330d2dddef7SLois Curfman McInnes   }
33107a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
332d2dddef7SLois Curfman McInnes }
333d2dddef7SLois Curfman McInnes 
33463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESUnSetMatrixFreeParameter(SNES snes)
335d2dddef7SLois Curfman McInnes {
336d2dddef7SLois Curfman McInnes   MFCtx_Private  *ctx;
337dfbe8321SBarry Smith   PetscErrorCode ierr;
338d2dddef7SLois Curfman McInnes   Mat            mat;
339d2dddef7SLois Curfman McInnes 
34007a3eed9SLois Curfman McInnes   PetscFunctionBegin;
34174637425SBarry Smith   ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
342d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
343a7cc72afSBarry Smith   if (ctx) ctx->need_h = PETSC_TRUE;
34407a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
345d2dddef7SLois Curfman McInnes }
346d2dddef7SLois Curfman McInnes 
347