xref: /petsc/src/snes/interface/noise/snesmfj2.c (revision 7831f9fc9edc069bf7ca487e05b8397f944b3ff8)
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 
68b5db460SBarry Smith extern PetscErrorCode SNESDiffParameterCreate_More(SNES,Vec,void**);
78b5db460SBarry Smith extern PetscErrorCode SNESDiffParameterCompute_More(SNES,void*,Vec,Vec,PetscReal*,PetscReal*);
88b5db460SBarry Smith extern PetscErrorCode SNESDiffParameterDestroy_More(void*);
9d2dddef7SLois Curfman McInnes 
10d2dddef7SLois Curfman McInnes typedef struct {  /* default context for matrix-free SNES */
11d2dddef7SLois Curfman McInnes   SNES         snes;             /* SNES context */
12d2dddef7SLois Curfman McInnes   Vec          w;                /* work vector */
1374637425SBarry Smith   MatNullSpace sp;               /* null space context */
14145595abSSatish Balay   PetscReal    error_rel;        /* square root of relative error in computing function */
15145595abSSatish Balay   PetscReal    umin;             /* minimum allowable u'a value relative to |u|_1 */
16f5af7f23SKarl Rupp   PetscBool    jorge;            /* flag indicating use of Jorge's method for determining the differencing parameter */
17145595abSSatish Balay   PetscReal    h;                /* differencing parameter */
18ace3abfcSBarry Smith   PetscBool    need_h;           /* flag indicating whether we must compute h */
19ace3abfcSBarry Smith   PetscBool    need_err;         /* flag indicating whether we must currently compute error_rel */
20ace3abfcSBarry Smith   PetscBool    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 
26dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeDestroy2_Private(Mat mat)
27d2dddef7SLois Curfman McInnes {
28dfbe8321SBarry Smith   PetscErrorCode ierr;
29d2dddef7SLois Curfman McInnes   MFCtx_Private  *ctx;
30d2dddef7SLois Curfman McInnes 
3107a3eed9SLois Curfman McInnes   PetscFunctionBegin;
3222d28d08SBarry Smith   ierr = MatShellGetContext(mat,(void**)&ctx);CHKERRQ(ierr);
336bf464f9SBarry Smith   ierr = VecDestroy(&ctx->w);CHKERRQ(ierr);
346bf464f9SBarry Smith   ierr = MatNullSpaceDestroy(&ctx->sp);CHKERRQ(ierr);
358b5db460SBarry Smith   if (ctx->jorge || ctx->compute_err) {ierr = SNESDiffParameterDestroy_More(ctx->data);CHKERRQ(ierr);}
36606d414cSSatish Balay   ierr = PetscFree(ctx);CHKERRQ(ierr);
3707a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
38d2dddef7SLois Curfman McInnes }
39d2dddef7SLois Curfman McInnes 
40d2dddef7SLois Curfman McInnes /*
41d2dddef7SLois Curfman McInnes    SNESMatrixFreeView2_Private - Views matrix-free parameters.
42d2dddef7SLois Curfman McInnes  */
43dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeView2_Private(Mat J,PetscViewer viewer)
44d2dddef7SLois Curfman McInnes {
45dfbe8321SBarry Smith   PetscErrorCode ierr;
46d2dddef7SLois Curfman McInnes   MFCtx_Private  *ctx;
47ace3abfcSBarry Smith   PetscBool      iascii;
48d2dddef7SLois Curfman McInnes 
4907a3eed9SLois Curfman McInnes   PetscFunctionBegin;
50d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(J,(void**)&ctx);CHKERRQ(ierr);
51251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
5232077d6dSBarry Smith   if (iascii) {
53b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  SNES matrix-free approximation:\n");CHKERRQ(ierr);
54758f395bSBarry Smith     if (ctx->jorge) {
55b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    using Jorge's method of determining differencing parameter\n");CHKERRQ(ierr);
56758f395bSBarry Smith     }
5757622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    err=%g (relative error in function evaluation)\n",(double)ctx->error_rel);CHKERRQ(ierr);
5857622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    umin=%g (minimum iterate parameter)\n",(double)ctx->umin);CHKERRQ(ierr);
59758f395bSBarry Smith     if (ctx->compute_err) {
6077431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    freq_err=%D (frequency for computing err)\n",ctx->compute_err_freq);CHKERRQ(ierr);
61d2dddef7SLois Curfman McInnes     }
62758f395bSBarry Smith   }
6307a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
64d2dddef7SLois Curfman McInnes }
65d2dddef7SLois Curfman McInnes 
66d2dddef7SLois Curfman McInnes /*
67d2dddef7SLois Curfman McInnes   SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector
68d2dddef7SLois Curfman McInnes   product, y = F'(u)*a:
69d2dddef7SLois Curfman McInnes         y = (F(u + ha) - F(u)) /h,
70d2dddef7SLois Curfman McInnes   where F = nonlinear function, as set by SNESSetFunction()
71d2dddef7SLois Curfman McInnes         u = current iterate
72d2dddef7SLois Curfman McInnes         h = difference interval
73d2dddef7SLois Curfman McInnes */
74dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y)
75d2dddef7SLois Curfman McInnes {
76d2dddef7SLois Curfman McInnes   MFCtx_Private  *ctx;
77d2dddef7SLois Curfman McInnes   SNES           snes;
78145595abSSatish Balay   PetscReal      h,norm,sum,umin,noise;
7979f36460SBarry Smith   PetscScalar    hs,dot;
80d2dddef7SLois Curfman McInnes   Vec            w,U,F;
81dfbe8321SBarry Smith   PetscErrorCode ierr,(*eval_fct)(SNES,Vec,Vec);
82e6ed2b1fSLois Curfman McInnes   MPI_Comm       comm;
8377431f27SBarry Smith   PetscInt       iter;
84d2dddef7SLois Curfman McInnes 
8507a3eed9SLois Curfman McInnes   PetscFunctionBegin;
86d2dddef7SLois Curfman McInnes   /* We log matrix-free matrix-vector products separately, so that we can
87d2dddef7SLois Curfman McInnes      separate the performance monitoring from the cases that use conventional
88d2dddef7SLois Curfman McInnes      storage.  We may eventually modify event logging to associate events
89d2dddef7SLois Curfman McInnes      with particular objects, hence alleviating the more general problem. */
903ec795f1SBarry Smith   ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
91d2dddef7SLois Curfman McInnes 
922d0c0e3bSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
93e6ed2b1fSLois Curfman McInnes   ierr = MatShellGetContext(mat,(void**)&ctx);CHKERRQ(ierr);
94e6ed2b1fSLois Curfman McInnes   snes = ctx->snes;
95e6ed2b1fSLois Curfman McInnes   w    = ctx->w;
96e6ed2b1fSLois Curfman McInnes   umin = ctx->umin;
97e6ed2b1fSLois Curfman McInnes 
98d2dddef7SLois Curfman McInnes   ierr     = SNESGetSolution(snes,&U);CHKERRQ(ierr);
99d2dddef7SLois Curfman McInnes   eval_fct = SNESComputeFunction;
1000298fd71SBarry Smith   ierr     = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr);
101295f7fbbSLois Curfman McInnes 
102d2dddef7SLois Curfman McInnes   /* Determine a "good" step size, h */
103d2dddef7SLois Curfman McInnes   if (ctx->need_h) {
104d2dddef7SLois Curfman McInnes 
105d2dddef7SLois Curfman McInnes     /* Use Jorge's method to compute h */
106d2dddef7SLois Curfman McInnes     if (ctx->jorge) {
1078b5db460SBarry Smith       ierr = SNESDiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr);
108d2dddef7SLois Curfman McInnes 
109d2dddef7SLois Curfman McInnes       /* Use the Brown/Saad method to compute h */
110d2dddef7SLois Curfman McInnes     } else {
111295f7fbbSLois Curfman McInnes       /* Compute error if desired */
112295f7fbbSLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
113145595abSSatish Balay       if ((ctx->need_err) || ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter-1)%ctx->compute_err_freq)))) {
114d2dddef7SLois Curfman McInnes         /* Use Jorge's method to compute noise */
1158b5db460SBarry Smith         ierr = SNESDiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr);
1161aa26658SKarl Rupp 
117369cc0aeSBarry Smith         ctx->error_rel = PetscSqrtReal(noise);
1181aa26658SKarl Rupp 
11957622a8eSBarry 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);
1201aa26658SKarl Rupp 
121295f7fbbSLois Curfman McInnes         ctx->compute_err_iter = iter;
1223b2a6d2fSBarry Smith         ctx->need_err         = PETSC_FALSE;
123d2dddef7SLois Curfman McInnes       }
124e6ed2b1fSLois Curfman McInnes 
125691598edSBarry Smith       ierr = VecDotBegin(U,a,&dot);CHKERRQ(ierr);
126691598edSBarry Smith       ierr = VecNormBegin(a,NORM_1,&sum);CHKERRQ(ierr);
127691598edSBarry Smith       ierr = VecNormBegin(a,NORM_2,&norm);CHKERRQ(ierr);
128691598edSBarry Smith       ierr = VecDotEnd(U,a,&dot);CHKERRQ(ierr);
129691598edSBarry Smith       ierr = VecNormEnd(a,NORM_1,&sum);CHKERRQ(ierr);
130691598edSBarry Smith       ierr = VecNormEnd(a,NORM_2,&norm);CHKERRQ(ierr);
131e6ed2b1fSLois Curfman McInnes 
132d2dddef7SLois Curfman McInnes 
133d2dddef7SLois Curfman McInnes       /* Safeguard for step sizes too small */
1341aa26658SKarl Rupp       if (sum == 0.0) {
1351aa26658SKarl Rupp         dot  = 1.0;
1361aa26658SKarl Rupp         norm = 1.0;
1371aa26658SKarl Rupp       } else if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum;
138145595abSSatish Balay       else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum;
139145595abSSatish Balay       h = PetscRealPart(ctx->error_rel*dot/(norm*norm));
140d2dddef7SLois Curfman McInnes     }
1411aa26658SKarl Rupp   } else h = ctx->h;
1421aa26658SKarl Rupp 
14357622a8eSBarry Smith   if (!ctx->jorge || !ctx->need_h) {ierr = PetscInfo1(snes,"h = %g\n",(double)h);CHKERRQ(ierr);}
144d2dddef7SLois Curfman McInnes 
145d2dddef7SLois Curfman McInnes   /* Evaluate function at F(u + ha) */
146db0b4b35SLois Curfman McInnes   hs   = h;
1472dcb1b2aSMatthew Knepley   ierr = VecWAXPY(w,hs,a,U);CHKERRQ(ierr);
148d2dddef7SLois Curfman McInnes   ierr = eval_fct(snes,w,y);CHKERRQ(ierr);
14979f36460SBarry Smith   ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr);
15079f36460SBarry Smith   ierr = VecScale(y,1.0/hs);CHKERRQ(ierr);
15139601f4eSBarry Smith   if (mat->nullsp) {ierr = MatNullSpaceRemove(mat->nullsp,y);CHKERRQ(ierr);}
152d2dddef7SLois Curfman McInnes 
1533ec795f1SBarry Smith   ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
15407a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
155d2dddef7SLois Curfman McInnes }
156d2dddef7SLois Curfman McInnes 
157d2dddef7SLois Curfman McInnes /*@C
15863dd3a1aSKris Buschelman    SNESMatrixFreeCreate2 - Creates a matrix-free matrix
159d2dddef7SLois Curfman McInnes    context for use with a SNES solver.  This matrix can be used as
160d2dddef7SLois Curfman McInnes    the Jacobian argument for the routine SNESSetJacobian().
161d2dddef7SLois Curfman McInnes 
162d2dddef7SLois Curfman McInnes    Input Parameters:
163d2dddef7SLois Curfman McInnes .  snes - the SNES context
164d2dddef7SLois Curfman McInnes .  x - vector where SNES solution is to be stored.
165d2dddef7SLois Curfman McInnes 
166d2dddef7SLois Curfman McInnes    Output Parameter:
167d2dddef7SLois Curfman McInnes .  J - the matrix-free matrix
168d2dddef7SLois Curfman McInnes 
16990c26489SBarry Smith    Level: advanced
17090c26489SBarry Smith 
171d2dddef7SLois Curfman McInnes    Notes:
172d2dddef7SLois Curfman McInnes    The matrix-free matrix context merely contains the function pointers
173d2dddef7SLois Curfman McInnes    and work space for performing finite difference approximations of
174d2dddef7SLois Curfman McInnes    Jacobian-vector products, J(u)*a, via
175d2dddef7SLois Curfman McInnes 
176d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h,
177d2dddef7SLois Curfman McInnes $   where by default
178d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
179d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
180d2dddef7SLois Curfman McInnes $   where
181d2dddef7SLois Curfman McInnes $        error_rel = square root of relative error in
182d2dddef7SLois Curfman McInnes $                    function evaluation
183d2dddef7SLois Curfman McInnes $        umin = minimum iterate parameter
184d2dddef7SLois Curfman McInnes $   Alternatively, the differencing parameter, h, can be set using
185e6ed2b1fSLois Curfman McInnes $   Jorge's nifty new strategy if one specifies the option
186e6ed2b1fSLois Curfman McInnes $          -snes_mf_jorge
187d2dddef7SLois Curfman McInnes 
1883ec795f1SBarry Smith    The user can set these parameters via MatMFFDSetFunctionError().
189a7f22e61SSatish Balay    See Users-Manual: ch_snes for details
190d2dddef7SLois Curfman McInnes 
191d2dddef7SLois Curfman McInnes    The user should call MatDestroy() when finished with the matrix-free
192d2dddef7SLois Curfman McInnes    matrix context.
193d2dddef7SLois Curfman McInnes 
194d2dddef7SLois Curfman McInnes    Options Database Keys:
195d2dddef7SLois Curfman McInnes $  -snes_mf_err <error_rel>
196d2dddef7SLois Curfman McInnes $  -snes_mf_unim <umin>
197295f7fbbSLois Curfman McInnes $  -snes_mf_compute_err
198295f7fbbSLois Curfman McInnes $  -snes_mf_freq_err <freq>
199e6ed2b1fSLois Curfman McInnes $  -snes_mf_jorge
200d2dddef7SLois Curfman McInnes 
201d2dddef7SLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix
202d2dddef7SLois Curfman McInnes 
2033ec795f1SBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError()
204d2dddef7SLois Curfman McInnes @*/
2057087cfbeSBarry Smith PetscErrorCode  SNESDefaultMatrixFreeCreate2(SNES snes,Vec x,Mat *J)
206d2dddef7SLois Curfman McInnes {
207d2dddef7SLois Curfman McInnes   MPI_Comm       comm;
208d2dddef7SLois Curfman McInnes   MFCtx_Private  *mfctx;
2096849ba73SBarry Smith   PetscErrorCode ierr;
210a7cc72afSBarry Smith   PetscInt       n,nloc;
211ace3abfcSBarry Smith   PetscBool      flg;
212d2dddef7SLois Curfman McInnes   char           p[64];
213d2dddef7SLois Curfman McInnes 
21407a3eed9SLois Curfman McInnes   PetscFunctionBegin;
215b00a9115SJed Brown   ierr                    = PetscNewLog(snes,&mfctx);CHKERRQ(ierr);
216d2dddef7SLois Curfman McInnes   mfctx->sp               = 0;
217d2dddef7SLois Curfman McInnes   mfctx->snes             = snes;
21877d8c4bbSBarry Smith   mfctx->error_rel        = PETSC_SQRT_MACHINE_EPSILON;
219d2dddef7SLois Curfman McInnes   mfctx->umin             = 1.e-6;
220d2dddef7SLois Curfman McInnes   mfctx->h                = 0.0;
221a7cc72afSBarry Smith   mfctx->need_h           = PETSC_TRUE;
222a7cc72afSBarry Smith   mfctx->need_err         = PETSC_FALSE;
223a7cc72afSBarry Smith   mfctx->compute_err      = PETSC_FALSE;
224295f7fbbSLois Curfman McInnes   mfctx->compute_err_freq = 0;
225295f7fbbSLois Curfman McInnes   mfctx->compute_err_iter = -1;
22690d69ab7SBarry Smith   mfctx->compute_err      = PETSC_FALSE;
22790d69ab7SBarry Smith   mfctx->jorge            = PETSC_FALSE;
2281aa26658SKarl Rupp 
229c5929fdfSBarry Smith   ierr = PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_err",&mfctx->error_rel,NULL);CHKERRQ(ierr);
230c5929fdfSBarry Smith   ierr = PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_umin",&mfctx->umin,NULL);CHKERRQ(ierr);
231c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_jorge",&mfctx->jorge,NULL);CHKERRQ(ierr);
232c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_compute_err",&mfctx->compute_err,NULL);CHKERRQ(ierr);
233c5929fdfSBarry Smith   ierr = PetscOptionsGetInt(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_freq_err",&mfctx->compute_err_freq,&flg);CHKERRQ(ierr);
234295f7fbbSLois Curfman McInnes   if (flg) {
235295f7fbbSLois Curfman McInnes     if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0;
2363b2a6d2fSBarry Smith     mfctx->compute_err = PETSC_TRUE;
237295f7fbbSLois Curfman McInnes   }
2383b2a6d2fSBarry Smith   if (mfctx->compute_err) mfctx->need_err = PETSC_TRUE;
239295f7fbbSLois Curfman McInnes   if (mfctx->jorge || mfctx->compute_err) {
2408b5db460SBarry Smith     ierr = SNESDiffParameterCreate_More(snes,x,&mfctx->data);CHKERRQ(ierr);
241d2dddef7SLois Curfman McInnes   } else mfctx->data = 0;
242d2dddef7SLois Curfman McInnes 
243c5929fdfSBarry Smith   ierr = PetscOptionsHasName(((PetscObject)snes)->options,NULL,"-help",&flg);CHKERRQ(ierr);
244549d3d68SSatish Balay   ierr = PetscStrcpy(p,"-");CHKERRQ(ierr);
2457adad957SLisandro Dalcin   if (((PetscObject)snes)->prefix) PetscStrcat(p,((PetscObject)snes)->prefix);
246d2dddef7SLois Curfman McInnes   if (flg) {
247ce94432eSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)snes)," Matrix-free Options (via SNES):\n");CHKERRQ(ierr);
24857622a8eSBarry 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);
24957622a8eSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"   %ssnes_mf_umin <umin>: see users manual (default %g)\n",p,(double)mfctx->umin);CHKERRQ(ierr);
250ce94432eSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"   %ssnes_mf_jorge: use Jorge More's method\n",p);CHKERRQ(ierr);
251ce94432eSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"   %ssnes_mf_compute_err: compute sqrt or relative error in function\n",p);CHKERRQ(ierr);
252ce94432eSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"   %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)\n",p);CHKERRQ(ierr);
253ce94432eSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"   %ssnes_mf_noise_file <file>: set file for printing noise info\n",p);CHKERRQ(ierr);
254d2dddef7SLois Curfman McInnes   }
255d2dddef7SLois Curfman McInnes   ierr = VecDuplicate(x,&mfctx->w);CHKERRQ(ierr);
256d2dddef7SLois Curfman McInnes   ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr);
257d2dddef7SLois Curfman McInnes   ierr = VecGetSize(x,&n);CHKERRQ(ierr);
258d2dddef7SLois Curfman McInnes   ierr = VecGetLocalSize(x,&nloc);CHKERRQ(ierr);
259f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,J);CHKERRQ(ierr);
260f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*J,nloc,n,n,n);CHKERRQ(ierr);
261f204ca49SKris Buschelman   ierr = MatSetType(*J,MATSHELL);CHKERRQ(ierr);
262f204ca49SKris Buschelman   ierr = MatShellSetContext(*J,mfctx);CHKERRQ(ierr);
263c134de8dSSatish Balay   ierr = MatShellSetOperation(*J,MATOP_MULT,(void (*)(void))SNESMatrixFreeMult2_Private);CHKERRQ(ierr);
264c134de8dSSatish Balay   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void (*)(void))SNESMatrixFreeDestroy2_Private);CHKERRQ(ierr);
265c134de8dSSatish Balay   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void (*)(void))SNESMatrixFreeView2_Private);CHKERRQ(ierr);
266*7831f9fcSJed Brown   ierr = MatSetUp(*J);CHKERRQ(ierr);
267d2dddef7SLois Curfman McInnes 
2683bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)*J,(PetscObject)mfctx->w);CHKERRQ(ierr);
2693bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)snes,(PetscObject)*J);CHKERRQ(ierr);
27007a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
271d2dddef7SLois Curfman McInnes }
272d2dddef7SLois Curfman McInnes 
27390c26489SBarry Smith /*@C
274758f395bSBarry Smith    SNESDefaultMatrixFreeSetParameters2 - Sets the parameters for the approximation of
275d2dddef7SLois Curfman McInnes    matrix-vector products using finite differences.
276d2dddef7SLois Curfman McInnes 
277d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h where
278d2dddef7SLois Curfman McInnes 
279d2dddef7SLois Curfman McInnes    either the user sets h directly here, or this parameter is computed via
280d2dddef7SLois Curfman McInnes 
281d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
282d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
283d2dddef7SLois Curfman McInnes $
284d2dddef7SLois Curfman McInnes 
285d2dddef7SLois Curfman McInnes    Input Parameters:
286758f395bSBarry Smith +  mat - the matrix
287d2dddef7SLois Curfman McInnes .  error_rel - relative error (should be set to the square root of
288d2dddef7SLois Curfman McInnes      the relative error in the function evaluations)
289d2dddef7SLois Curfman McInnes .  umin - minimum allowable u-value
290758f395bSBarry Smith -  h - differencing parameter
291d2dddef7SLois Curfman McInnes 
29290c26489SBarry Smith    Level: advanced
29390c26489SBarry Smith 
294d2dddef7SLois Curfman McInnes    Notes:
295d2dddef7SLois Curfman McInnes    If the user sets the parameter h directly, then this value will be used
296d2dddef7SLois Curfman McInnes    instead of the default computation indicated above.
297d2dddef7SLois Curfman McInnes 
298d2dddef7SLois Curfman McInnes .keywords: SNES, matrix-free, parameters
299d2dddef7SLois Curfman McInnes 
3005a655dc6SBarry Smith .seealso: MatCreateSNESMF()
301d2dddef7SLois Curfman McInnes @*/
3027087cfbeSBarry Smith PetscErrorCode  SNESDefaultMatrixFreeSetParameters2(Mat mat,PetscReal error,PetscReal umin,PetscReal h)
303d2dddef7SLois Curfman McInnes {
304d2dddef7SLois Curfman McInnes   MFCtx_Private  *ctx;
305dfbe8321SBarry Smith   PetscErrorCode ierr;
306d2dddef7SLois Curfman McInnes 
30707a3eed9SLois Curfman McInnes   PetscFunctionBegin;
308d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void**)&ctx);CHKERRQ(ierr);
309d2dddef7SLois Curfman McInnes   if (ctx) {
310d2dddef7SLois Curfman McInnes     if (error != PETSC_DEFAULT) ctx->error_rel = error;
311d2dddef7SLois Curfman McInnes     if (umin  != PETSC_DEFAULT) ctx->umin = umin;
312d2dddef7SLois Curfman McInnes     if (h     != PETSC_DEFAULT) {
313d2dddef7SLois Curfman McInnes       ctx->h      = h;
314a7cc72afSBarry Smith       ctx->need_h = PETSC_FALSE;
315d2dddef7SLois Curfman McInnes     }
316d2dddef7SLois Curfman McInnes   }
31707a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
318d2dddef7SLois Curfman McInnes }
319d2dddef7SLois Curfman McInnes 
3207087cfbeSBarry Smith PetscErrorCode  SNESUnSetMatrixFreeParameter(SNES snes)
321d2dddef7SLois Curfman McInnes {
322d2dddef7SLois Curfman McInnes   MFCtx_Private  *ctx;
323dfbe8321SBarry Smith   PetscErrorCode ierr;
324d2dddef7SLois Curfman McInnes   Mat            mat;
325d2dddef7SLois Curfman McInnes 
32607a3eed9SLois Curfman McInnes   PetscFunctionBegin;
3270298fd71SBarry Smith   ierr = SNESGetJacobian(snes,&mat,NULL,NULL,NULL);CHKERRQ(ierr);
328d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void**)&ctx);CHKERRQ(ierr);
329a7cc72afSBarry Smith   if (ctx) ctx->need_h = PETSC_TRUE;
33007a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
331d2dddef7SLois Curfman McInnes }
332d2dddef7SLois Curfman McInnes 
333