xref: /petsc/src/snes/interface/noise/snesmfj2.c (revision a2b725a8db0d6bf6cc2a1c6df7dd8029aadfff6e)
1d2dddef7SLois Curfman McInnes 
2af0996ceSBarry Smith #include <petsc/private/snesimpl.h>   /*I  "petscsnes.h"   I*/
36370053bSBarry Smith /* matimpl.h is needed only for logging of matrix operation */
4af0996ceSBarry Smith #include <petsc/private/matimpl.h>
5d2dddef7SLois Curfman McInnes 
625acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESUnSetMatrixFreeParameter(SNES);
725acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
825acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDefaultMatrixFreeSetParameters2(Mat,PetscReal,PetscReal,PetscReal);
925acbd8eSLisandro Dalcin 
1025acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterCreate_More(SNES,Vec,void**);
1125acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterCompute_More(SNES,void*,Vec,Vec,PetscReal*,PetscReal*);
1225acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterDestroy_More(void*);
13d2dddef7SLois Curfman McInnes 
14d2dddef7SLois Curfman McInnes typedef struct {  /* default context for matrix-free SNES */
15d2dddef7SLois Curfman McInnes   SNES         snes;             /* SNES context */
16d2dddef7SLois Curfman McInnes   Vec          w;                /* work vector */
1774637425SBarry Smith   MatNullSpace sp;               /* null space context */
18145595abSSatish Balay   PetscReal    error_rel;        /* square root of relative error in computing function */
19145595abSSatish Balay   PetscReal    umin;             /* minimum allowable u'a value relative to |u|_1 */
20f5af7f23SKarl Rupp   PetscBool    jorge;            /* flag indicating use of Jorge's method for determining the differencing parameter */
21145595abSSatish Balay   PetscReal    h;                /* differencing parameter */
22ace3abfcSBarry Smith   PetscBool    need_h;           /* flag indicating whether we must compute h */
23ace3abfcSBarry Smith   PetscBool    need_err;         /* flag indicating whether we must currently compute error_rel */
24ace3abfcSBarry Smith   PetscBool    compute_err;      /* flag indicating whether we must ever compute error_rel */
25a7cc72afSBarry Smith   PetscInt     compute_err_iter; /* last iter where we've computer error_rel */
26a7cc72afSBarry Smith   PetscInt     compute_err_freq; /* frequency of computing error_rel */
27d2dddef7SLois Curfman McInnes   void         *data;            /* implementation-specific data */
28d2dddef7SLois Curfman McInnes } MFCtx_Private;
29d2dddef7SLois Curfman McInnes 
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;
3622d28d08SBarry Smith   ierr = MatShellGetContext(mat,(void**)&ctx);CHKERRQ(ierr);
376bf464f9SBarry Smith   ierr = VecDestroy(&ctx->w);CHKERRQ(ierr);
386bf464f9SBarry Smith   ierr = MatNullSpaceDestroy(&ctx->sp);CHKERRQ(ierr);
398b5db460SBarry Smith   if (ctx->jorge || ctx->compute_err) {ierr = SNESDiffParameterDestroy_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 /*
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;
51ace3abfcSBarry Smith   PetscBool      iascii;
52d2dddef7SLois Curfman McInnes 
5307a3eed9SLois Curfman McInnes   PetscFunctionBegin;
54d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(J,(void**)&ctx);CHKERRQ(ierr);
55251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&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     }
6157622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    err=%g (relative error in function evaluation)\n",(double)ctx->error_rel);CHKERRQ(ierr);
6257622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    umin=%g (minimum iterate parameter)\n",(double)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   }
6707a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
68d2dddef7SLois Curfman McInnes }
69d2dddef7SLois Curfman McInnes 
70d2dddef7SLois Curfman McInnes /*
71d2dddef7SLois Curfman McInnes   SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector
72d2dddef7SLois Curfman McInnes   product, y = F'(u)*a:
73d2dddef7SLois Curfman McInnes         y = (F(u + ha) - F(u)) /h,
74d2dddef7SLois Curfman McInnes   where F = nonlinear function, as set by SNESSetFunction()
75d2dddef7SLois Curfman McInnes         u = current iterate
76d2dddef7SLois Curfman McInnes         h = difference interval
77d2dddef7SLois Curfman McInnes */
78dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y)
79d2dddef7SLois Curfman McInnes {
80d2dddef7SLois Curfman McInnes   MFCtx_Private  *ctx;
81d2dddef7SLois Curfman McInnes   SNES           snes;
82145595abSSatish Balay   PetscReal      h,norm,sum,umin,noise;
8379f36460SBarry Smith   PetscScalar    hs,dot;
84d2dddef7SLois Curfman McInnes   Vec            w,U,F;
85dfbe8321SBarry Smith   PetscErrorCode ierr,(*eval_fct)(SNES,Vec,Vec);
86e6ed2b1fSLois Curfman McInnes   MPI_Comm       comm;
8777431f27SBarry Smith   PetscInt       iter;
88d2dddef7SLois Curfman McInnes 
8907a3eed9SLois Curfman McInnes   PetscFunctionBegin;
90d2dddef7SLois Curfman McInnes   /* We log matrix-free matrix-vector products separately, so that we can
91d2dddef7SLois Curfman McInnes      separate the performance monitoring from the cases that use conventional
92d2dddef7SLois Curfman McInnes      storage.  We may eventually modify event logging to associate events
93d2dddef7SLois Curfman McInnes      with particular objects, hence alleviating the more general problem. */
943ec795f1SBarry Smith   ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
95d2dddef7SLois Curfman McInnes 
962d0c0e3bSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
97e6ed2b1fSLois Curfman McInnes   ierr = MatShellGetContext(mat,(void**)&ctx);CHKERRQ(ierr);
98e6ed2b1fSLois Curfman McInnes   snes = ctx->snes;
99e6ed2b1fSLois Curfman McInnes   w    = ctx->w;
100e6ed2b1fSLois Curfman McInnes   umin = ctx->umin;
101e6ed2b1fSLois Curfman McInnes 
102d2dddef7SLois Curfman McInnes   ierr     = SNESGetSolution(snes,&U);CHKERRQ(ierr);
103d2dddef7SLois Curfman McInnes   eval_fct = SNESComputeFunction;
1040298fd71SBarry Smith   ierr     = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr);
105295f7fbbSLois Curfman McInnes 
106d2dddef7SLois Curfman McInnes   /* Determine a "good" step size, h */
107d2dddef7SLois Curfman McInnes   if (ctx->need_h) {
108d2dddef7SLois Curfman McInnes 
109d2dddef7SLois Curfman McInnes     /* Use Jorge's method to compute h */
110d2dddef7SLois Curfman McInnes     if (ctx->jorge) {
1118b5db460SBarry Smith       ierr = SNESDiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr);
112d2dddef7SLois Curfman McInnes 
113d2dddef7SLois Curfman McInnes       /* Use the Brown/Saad method to compute h */
114d2dddef7SLois Curfman McInnes     } else {
115295f7fbbSLois Curfman McInnes       /* Compute error if desired */
116295f7fbbSLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
117145595abSSatish Balay       if ((ctx->need_err) || ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter-1)%ctx->compute_err_freq)))) {
118d2dddef7SLois Curfman McInnes         /* Use Jorge's method to compute noise */
1198b5db460SBarry Smith         ierr = SNESDiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr);
1201aa26658SKarl Rupp 
121369cc0aeSBarry Smith         ctx->error_rel = PetscSqrtReal(noise);
1221aa26658SKarl Rupp 
12357622a8eSBarry Smith         ierr = PetscInfo3(snes,"Using Jorge's noise: noise=%g, sqrt(noise)=%g, h_more=%g\n",(double)noise,(double)ctx->error_rel,(double)h);CHKERRQ(ierr);
1241aa26658SKarl Rupp 
125295f7fbbSLois Curfman McInnes         ctx->compute_err_iter = iter;
1263b2a6d2fSBarry Smith         ctx->need_err         = PETSC_FALSE;
127d2dddef7SLois Curfman McInnes       }
128e6ed2b1fSLois Curfman McInnes 
129691598edSBarry Smith       ierr = VecDotBegin(U,a,&dot);CHKERRQ(ierr);
130691598edSBarry Smith       ierr = VecNormBegin(a,NORM_1,&sum);CHKERRQ(ierr);
131691598edSBarry Smith       ierr = VecNormBegin(a,NORM_2,&norm);CHKERRQ(ierr);
132691598edSBarry Smith       ierr = VecDotEnd(U,a,&dot);CHKERRQ(ierr);
133691598edSBarry Smith       ierr = VecNormEnd(a,NORM_1,&sum);CHKERRQ(ierr);
134691598edSBarry Smith       ierr = VecNormEnd(a,NORM_2,&norm);CHKERRQ(ierr);
135e6ed2b1fSLois Curfman McInnes 
136d2dddef7SLois Curfman McInnes 
137d2dddef7SLois Curfman McInnes       /* Safeguard for step sizes too small */
1381aa26658SKarl Rupp       if (sum == 0.0) {
1391aa26658SKarl Rupp         dot  = 1.0;
1401aa26658SKarl Rupp         norm = 1.0;
1411aa26658SKarl Rupp       } else if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum;
142145595abSSatish Balay       else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum;
143145595abSSatish Balay       h = PetscRealPart(ctx->error_rel*dot/(norm*norm));
144d2dddef7SLois Curfman McInnes     }
1451aa26658SKarl Rupp   } else h = ctx->h;
1461aa26658SKarl Rupp 
14757622a8eSBarry Smith   if (!ctx->jorge || !ctx->need_h) {ierr = PetscInfo1(snes,"h = %g\n",(double)h);CHKERRQ(ierr);}
148d2dddef7SLois Curfman McInnes 
149d2dddef7SLois Curfman McInnes   /* Evaluate function at F(u + ha) */
150db0b4b35SLois Curfman McInnes   hs   = h;
1512dcb1b2aSMatthew Knepley   ierr = VecWAXPY(w,hs,a,U);CHKERRQ(ierr);
152d2dddef7SLois Curfman McInnes   ierr = eval_fct(snes,w,y);CHKERRQ(ierr);
15379f36460SBarry Smith   ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr);
15479f36460SBarry Smith   ierr = VecScale(y,1.0/hs);CHKERRQ(ierr);
15539601f4eSBarry Smith   if (mat->nullsp) {ierr = MatNullSpaceRemove(mat->nullsp,y);CHKERRQ(ierr);}
156d2dddef7SLois Curfman McInnes 
1573ec795f1SBarry Smith   ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
15807a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
159d2dddef7SLois Curfman McInnes }
160d2dddef7SLois Curfman McInnes 
161d2dddef7SLois Curfman McInnes /*@C
16263dd3a1aSKris Buschelman    SNESMatrixFreeCreate2 - Creates a matrix-free matrix
163d2dddef7SLois Curfman McInnes    context for use with a SNES solver.  This matrix can be used as
164d2dddef7SLois Curfman McInnes    the Jacobian argument for the routine SNESSetJacobian().
165d2dddef7SLois Curfman McInnes 
166d2dddef7SLois Curfman McInnes    Input Parameters:
167*a2b725a8SWilliam Gropp +  snes - the SNES context
168*a2b725a8SWilliam Gropp -  x - vector where SNES solution is to be stored.
169d2dddef7SLois Curfman McInnes 
170d2dddef7SLois Curfman McInnes    Output Parameter:
171d2dddef7SLois Curfman McInnes .  J - the matrix-free matrix
172d2dddef7SLois Curfman McInnes 
17390c26489SBarry Smith    Level: advanced
17490c26489SBarry Smith 
175d2dddef7SLois Curfman McInnes    Notes:
176d2dddef7SLois Curfman McInnes    The matrix-free matrix context merely contains the function pointers
177d2dddef7SLois Curfman McInnes    and work space for performing finite difference approximations of
178d2dddef7SLois Curfman McInnes    Jacobian-vector products, J(u)*a, via
179d2dddef7SLois Curfman McInnes 
180d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h,
181d2dddef7SLois Curfman McInnes $   where by default
182d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
183d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
184d2dddef7SLois Curfman McInnes $   where
185d2dddef7SLois Curfman McInnes $        error_rel = square root of relative error in
186d2dddef7SLois Curfman McInnes $                    function evaluation
187d2dddef7SLois Curfman McInnes $        umin = minimum iterate parameter
188d2dddef7SLois Curfman McInnes $   Alternatively, the differencing parameter, h, can be set using
189e6ed2b1fSLois Curfman McInnes $   Jorge's nifty new strategy if one specifies the option
190e6ed2b1fSLois Curfman McInnes $          -snes_mf_jorge
191d2dddef7SLois Curfman McInnes 
1923ec795f1SBarry Smith    The user can set these parameters via MatMFFDSetFunctionError().
193a7f22e61SSatish Balay    See Users-Manual: ch_snes for details
194d2dddef7SLois Curfman McInnes 
195d2dddef7SLois Curfman McInnes    The user should call MatDestroy() when finished with the matrix-free
196d2dddef7SLois Curfman McInnes    matrix context.
197d2dddef7SLois Curfman McInnes 
198d2dddef7SLois Curfman McInnes    Options Database Keys:
199d2dddef7SLois Curfman McInnes $  -snes_mf_err <error_rel>
200d2dddef7SLois Curfman McInnes $  -snes_mf_unim <umin>
201295f7fbbSLois Curfman McInnes $  -snes_mf_compute_err
202295f7fbbSLois Curfman McInnes $  -snes_mf_freq_err <freq>
203e6ed2b1fSLois Curfman McInnes $  -snes_mf_jorge
204d2dddef7SLois Curfman McInnes 
205d2dddef7SLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix
206d2dddef7SLois Curfman McInnes 
2073ec795f1SBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError()
208d2dddef7SLois Curfman McInnes @*/
2097087cfbeSBarry Smith PetscErrorCode  SNESDefaultMatrixFreeCreate2(SNES snes,Vec x,Mat *J)
210d2dddef7SLois Curfman McInnes {
211d2dddef7SLois Curfman McInnes   MPI_Comm       comm;
212d2dddef7SLois Curfman McInnes   MFCtx_Private  *mfctx;
2136849ba73SBarry Smith   PetscErrorCode ierr;
214a7cc72afSBarry Smith   PetscInt       n,nloc;
215ace3abfcSBarry Smith   PetscBool      flg;
216d2dddef7SLois Curfman McInnes   char           p[64];
217d2dddef7SLois Curfman McInnes 
21807a3eed9SLois Curfman McInnes   PetscFunctionBegin;
219b00a9115SJed Brown   ierr                    = PetscNewLog(snes,&mfctx);CHKERRQ(ierr);
220d2dddef7SLois Curfman McInnes   mfctx->sp               = 0;
221d2dddef7SLois Curfman McInnes   mfctx->snes             = snes;
22277d8c4bbSBarry Smith   mfctx->error_rel        = PETSC_SQRT_MACHINE_EPSILON;
223d2dddef7SLois Curfman McInnes   mfctx->umin             = 1.e-6;
224d2dddef7SLois Curfman McInnes   mfctx->h                = 0.0;
225a7cc72afSBarry Smith   mfctx->need_h           = PETSC_TRUE;
226a7cc72afSBarry Smith   mfctx->need_err         = PETSC_FALSE;
227a7cc72afSBarry Smith   mfctx->compute_err      = PETSC_FALSE;
228295f7fbbSLois Curfman McInnes   mfctx->compute_err_freq = 0;
229295f7fbbSLois Curfman McInnes   mfctx->compute_err_iter = -1;
23090d69ab7SBarry Smith   mfctx->compute_err      = PETSC_FALSE;
23190d69ab7SBarry Smith   mfctx->jorge            = PETSC_FALSE;
2321aa26658SKarl Rupp 
233c5929fdfSBarry Smith   ierr = PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_err",&mfctx->error_rel,NULL);CHKERRQ(ierr);
234c5929fdfSBarry Smith   ierr = PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_umin",&mfctx->umin,NULL);CHKERRQ(ierr);
235c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_jorge",&mfctx->jorge,NULL);CHKERRQ(ierr);
236c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_compute_err",&mfctx->compute_err,NULL);CHKERRQ(ierr);
237c5929fdfSBarry Smith   ierr = PetscOptionsGetInt(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_freq_err",&mfctx->compute_err_freq,&flg);CHKERRQ(ierr);
238295f7fbbSLois Curfman McInnes   if (flg) {
239295f7fbbSLois Curfman McInnes     if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0;
2403b2a6d2fSBarry Smith     mfctx->compute_err = PETSC_TRUE;
241295f7fbbSLois Curfman McInnes   }
2423b2a6d2fSBarry Smith   if (mfctx->compute_err) mfctx->need_err = PETSC_TRUE;
243295f7fbbSLois Curfman McInnes   if (mfctx->jorge || mfctx->compute_err) {
2448b5db460SBarry Smith     ierr = SNESDiffParameterCreate_More(snes,x,&mfctx->data);CHKERRQ(ierr);
245d2dddef7SLois Curfman McInnes   } else mfctx->data = 0;
246d2dddef7SLois Curfman McInnes 
2472d747510SLisandro Dalcin   ierr = PetscOptionsHasHelp(((PetscObject)snes)->options,&flg);CHKERRQ(ierr);
248a126751eSBarry Smith   ierr = PetscStrncpy(p,"-",sizeof(p));CHKERRQ(ierr);
249a126751eSBarry Smith   if (((PetscObject)snes)->prefix) {ierr = PetscStrlcat(p,((PetscObject)snes)->prefix,sizeof(p));CHKERRQ(ierr);}
250d2dddef7SLois Curfman McInnes   if (flg) {
251ce94432eSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)snes)," Matrix-free Options (via SNES):\n");CHKERRQ(ierr);
25257622a8eSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"   %ssnes_mf_err <err>: set sqrt of relative error in function (default %g)\n",p,(double)mfctx->error_rel);CHKERRQ(ierr);
25357622a8eSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"   %ssnes_mf_umin <umin>: see users manual (default %g)\n",p,(double)mfctx->umin);CHKERRQ(ierr);
254ce94432eSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"   %ssnes_mf_jorge: use Jorge More's method\n",p);CHKERRQ(ierr);
255ce94432eSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"   %ssnes_mf_compute_err: compute sqrt or relative error in function\n",p);CHKERRQ(ierr);
256ce94432eSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"   %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)\n",p);CHKERRQ(ierr);
257ce94432eSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"   %ssnes_mf_noise_file <file>: set file for printing noise info\n",p);CHKERRQ(ierr);
258d2dddef7SLois Curfman McInnes   }
259d2dddef7SLois Curfman McInnes   ierr = VecDuplicate(x,&mfctx->w);CHKERRQ(ierr);
260d2dddef7SLois Curfman McInnes   ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr);
261d2dddef7SLois Curfman McInnes   ierr = VecGetSize(x,&n);CHKERRQ(ierr);
262d2dddef7SLois Curfman McInnes   ierr = VecGetLocalSize(x,&nloc);CHKERRQ(ierr);
263f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,J);CHKERRQ(ierr);
264f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*J,nloc,n,n,n);CHKERRQ(ierr);
265f204ca49SKris Buschelman   ierr = MatSetType(*J,MATSHELL);CHKERRQ(ierr);
266f204ca49SKris Buschelman   ierr = MatShellSetContext(*J,mfctx);CHKERRQ(ierr);
267c134de8dSSatish Balay   ierr = MatShellSetOperation(*J,MATOP_MULT,(void (*)(void))SNESMatrixFreeMult2_Private);CHKERRQ(ierr);
268c134de8dSSatish Balay   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void (*)(void))SNESMatrixFreeDestroy2_Private);CHKERRQ(ierr);
269c134de8dSSatish Balay   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void (*)(void))SNESMatrixFreeView2_Private);CHKERRQ(ierr);
2707831f9fcSJed Brown   ierr = MatSetUp(*J);CHKERRQ(ierr);
271d2dddef7SLois Curfman McInnes 
2723bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)*J,(PetscObject)mfctx->w);CHKERRQ(ierr);
2733bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)snes,(PetscObject)*J);CHKERRQ(ierr);
27407a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
275d2dddef7SLois Curfman McInnes }
276d2dddef7SLois Curfman McInnes 
27790c26489SBarry Smith /*@C
278758f395bSBarry Smith    SNESDefaultMatrixFreeSetParameters2 - Sets the parameters for the approximation of
279d2dddef7SLois Curfman McInnes    matrix-vector products using finite differences.
280d2dddef7SLois Curfman McInnes 
281d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h where
282d2dddef7SLois Curfman McInnes 
283d2dddef7SLois Curfman McInnes    either the user sets h directly here, or this parameter is computed via
284d2dddef7SLois Curfman McInnes 
285d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
286d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
287d2dddef7SLois Curfman McInnes $
288d2dddef7SLois Curfman McInnes 
289d2dddef7SLois Curfman McInnes    Input Parameters:
290758f395bSBarry Smith +  mat - the matrix
291d2dddef7SLois Curfman McInnes .  error_rel - relative error (should be set to the square root of
292d2dddef7SLois Curfman McInnes      the relative error in the function evaluations)
293d2dddef7SLois Curfman McInnes .  umin - minimum allowable u-value
294758f395bSBarry Smith -  h - differencing parameter
295d2dddef7SLois Curfman McInnes 
29690c26489SBarry Smith    Level: advanced
29790c26489SBarry Smith 
298d2dddef7SLois Curfman McInnes    Notes:
299d2dddef7SLois Curfman McInnes    If the user sets the parameter h directly, then this value will be used
300d2dddef7SLois Curfman McInnes    instead of the default computation indicated above.
301d2dddef7SLois Curfman McInnes 
302d2dddef7SLois Curfman McInnes .keywords: SNES, matrix-free, parameters
303d2dddef7SLois Curfman McInnes 
3045a655dc6SBarry Smith .seealso: MatCreateSNESMF()
305d2dddef7SLois Curfman McInnes @*/
3067087cfbeSBarry Smith PetscErrorCode  SNESDefaultMatrixFreeSetParameters2(Mat mat,PetscReal error,PetscReal umin,PetscReal h)
307d2dddef7SLois Curfman McInnes {
308d2dddef7SLois Curfman McInnes   MFCtx_Private  *ctx;
309dfbe8321SBarry Smith   PetscErrorCode ierr;
310d2dddef7SLois Curfman McInnes 
31107a3eed9SLois Curfman McInnes   PetscFunctionBegin;
312d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void**)&ctx);CHKERRQ(ierr);
313d2dddef7SLois Curfman McInnes   if (ctx) {
314d2dddef7SLois Curfman McInnes     if (error != PETSC_DEFAULT) ctx->error_rel = error;
315d2dddef7SLois Curfman McInnes     if (umin  != PETSC_DEFAULT) ctx->umin = umin;
316d2dddef7SLois Curfman McInnes     if (h     != PETSC_DEFAULT) {
317d2dddef7SLois Curfman McInnes       ctx->h      = h;
318a7cc72afSBarry Smith       ctx->need_h = PETSC_FALSE;
319d2dddef7SLois Curfman McInnes     }
320d2dddef7SLois Curfman McInnes   }
32107a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
322d2dddef7SLois Curfman McInnes }
323d2dddef7SLois Curfman McInnes 
3247087cfbeSBarry Smith PetscErrorCode  SNESUnSetMatrixFreeParameter(SNES snes)
325d2dddef7SLois Curfman McInnes {
326d2dddef7SLois Curfman McInnes   MFCtx_Private  *ctx;
327dfbe8321SBarry Smith   PetscErrorCode ierr;
328d2dddef7SLois Curfman McInnes   Mat            mat;
329d2dddef7SLois Curfman McInnes 
33007a3eed9SLois Curfman McInnes   PetscFunctionBegin;
3310298fd71SBarry Smith   ierr = SNESGetJacobian(snes,&mat,NULL,NULL,NULL);CHKERRQ(ierr);
332d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void**)&ctx);CHKERRQ(ierr);
333a7cc72afSBarry Smith   if (ctx) ctx->need_h = PETSC_TRUE;
33407a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
335d2dddef7SLois Curfman McInnes }
336d2dddef7SLois Curfman McInnes 
337