xref: /petsc/src/snes/interface/noise/snesmfj2.c (revision 758f395b27ac8d6dbc276772ed0d06fa12e89e10)
1d834d2e4SLois Curfman McInnes #ifdef PETSC_RCS_HEADER
2*758f395bSBarry Smith static char vcid[] = "$Id: snesmfj2.c,v 1.10 1998/09/14 17:26:28 curfman Exp bsmith $";
3d2dddef7SLois Curfman McInnes #endif
4d2dddef7SLois Curfman McInnes 
5d2dddef7SLois Curfman McInnes #include "src/snes/snesimpl.h"   /*I  "snes.h"   I*/
6d2dddef7SLois Curfman McInnes 
7d2dddef7SLois Curfman McInnes extern int DiffParameterCreate_More(SNES,Vec,void**);
8d2dddef7SLois Curfman McInnes extern int DiffParameterCompute_More(SNES,void*,Vec,Vec,double*,double*);
9d2dddef7SLois Curfman McInnes extern int DiffParameterDestroy_More(void*);
10d2dddef7SLois Curfman McInnes 
11d2dddef7SLois Curfman McInnes typedef struct {  /* default context for matrix-free SNES */
12d2dddef7SLois Curfman McInnes   SNES        snes;             /* SNES context */
13d2dddef7SLois Curfman McInnes   Vec         w;                /* work vector */
14d2dddef7SLois Curfman McInnes   PCNullSpace sp;               /* null space context */
15d2dddef7SLois Curfman McInnes   double      error_rel;        /* square root of relative error in computing function */
16d2dddef7SLois Curfman McInnes   double      umin;             /* minimum allowable u'a value relative to |u|_1 */
17d2dddef7SLois Curfman McInnes   int         jorge;            /* flag indicating use of Jorge's method for determining
18d2dddef7SLois Curfman McInnes                                    the differencing parameter */
19db0b4b35SLois Curfman McInnes   double      h;                /* differencing parameter */
20d2dddef7SLois Curfman McInnes   int         need_h;           /* flag indicating whether we must compute h */
21295f7fbbSLois Curfman McInnes   int         need_err;         /* flag indicating whether we must currently compute error_rel */
22295f7fbbSLois Curfman McInnes   int         compute_err;      /* flag indicating whether we must ever compute error_rel */
23295f7fbbSLois Curfman McInnes   int         compute_err_iter; /* last iter where we've computer error_rel */
24295f7fbbSLois Curfman McInnes   int         compute_err_freq; /* frequency of computing error_rel */
25d2dddef7SLois Curfman McInnes   void        *data;            /* implementation-specific data */
26d2dddef7SLois Curfman McInnes } MFCtx_Private;
27d2dddef7SLois Curfman McInnes 
28d2dddef7SLois Curfman McInnes #undef __FUNC__
29d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeDestroy2_Private" /* ADIC Ignore */
30f1d4f8a3SBarry Smith int SNESMatrixFreeDestroy2_Private(Mat mat)
31d2dddef7SLois Curfman McInnes {
32d2dddef7SLois Curfman McInnes   int           ierr;
33d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
34d2dddef7SLois Curfman McInnes 
3507a3eed9SLois Curfman McInnes   PetscFunctionBegin;
36d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx);
37d2dddef7SLois Curfman McInnes   ierr = VecDestroy(ctx->w); CHKERRQ(ierr);
38d2dddef7SLois Curfman McInnes   if (ctx->sp) {ierr = PCNullSpaceDestroy(ctx->sp); CHKERRQ(ierr);}
39295f7fbbSLois Curfman McInnes   if (ctx->jorge || ctx->compute_err) {ierr = DiffParameterDestroy_More(ctx->data); CHKERRQ(ierr);}
40d2dddef7SLois Curfman McInnes   PetscFree(ctx);
4107a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
42d2dddef7SLois Curfman McInnes }
43d2dddef7SLois Curfman McInnes 
44d2dddef7SLois Curfman McInnes #undef __FUNC__
45d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeView2_Private" /* ADIC Ignore */
46d2dddef7SLois Curfman McInnes /*
47d2dddef7SLois Curfman McInnes    SNESMatrixFreeView2_Private - Views matrix-free parameters.
48d2dddef7SLois Curfman McInnes  */
49d2dddef7SLois Curfman McInnes int SNESMatrixFreeView2_Private(Mat J,Viewer viewer)
50d2dddef7SLois Curfman McInnes {
51d2dddef7SLois Curfman McInnes   int           ierr;
52d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
53d2dddef7SLois Curfman McInnes   MPI_Comm      comm;
54d2dddef7SLois Curfman McInnes   FILE          *fd;
55d2dddef7SLois Curfman McInnes   ViewerType    vtype;
56d2dddef7SLois Curfman McInnes 
5707a3eed9SLois Curfman McInnes   PetscFunctionBegin;
58d2dddef7SLois Curfman McInnes   PetscObjectGetComm((PetscObject)J,&comm);
59d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
60d2dddef7SLois Curfman McInnes   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
61d2dddef7SLois Curfman McInnes   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
62*758f395bSBarry Smith   if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
63d2dddef7SLois Curfman McInnes      PetscFPrintf(comm,fd,"  SNES matrix-free approximation:\n");
64*758f395bSBarry Smith      if (ctx->jorge) {
65d2dddef7SLois Curfman McInnes        PetscFPrintf(comm,fd,"    using Jorge's method of determining differencing parameter\n");
66*758f395bSBarry Smith      }
67d2dddef7SLois Curfman McInnes      PetscFPrintf(comm,fd,"    err=%g (relative error in function evaluation)\n",ctx->error_rel);
68d2dddef7SLois Curfman McInnes      PetscFPrintf(comm,fd,"    umin=%g (minimum iterate parameter)\n",ctx->umin);
69*758f395bSBarry Smith      if (ctx->compute_err) {
70295f7fbbSLois Curfman McInnes        PetscFPrintf(comm,fd,"    freq_err=%d (frequency for computing err)\n",ctx->compute_err_freq);
71d2dddef7SLois Curfman McInnes      }
72*758f395bSBarry Smith   } else {
73*758f395bSBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
74*758f395bSBarry Smith   }
7507a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
76d2dddef7SLois Curfman McInnes }
77d2dddef7SLois Curfman McInnes 
78e6ed2b1fSLois Curfman McInnes extern int VecDot_Seq(Vec,Vec,Scalar *);
79e6ed2b1fSLois Curfman McInnes extern int VecNorm_Seq(Vec,NormType,double *);
80e6ed2b1fSLois Curfman McInnes 
81d2dddef7SLois Curfman McInnes #undef __FUNC__
82d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeMult2_Private"
83d2dddef7SLois Curfman McInnes /*
84d2dddef7SLois Curfman McInnes   SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector
85d2dddef7SLois Curfman McInnes   product, y = F'(u)*a:
86d2dddef7SLois Curfman McInnes         y = ( F(u + ha) - F(u) ) /h,
87d2dddef7SLois Curfman McInnes   where F = nonlinear function, as set by SNESSetFunction()
88d2dddef7SLois Curfman McInnes         u = current iterate
89d2dddef7SLois Curfman McInnes         h = difference interval
90d2dddef7SLois Curfman McInnes */
91d2dddef7SLois Curfman McInnes int SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y)
92d2dddef7SLois Curfman McInnes {
93d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
94d2dddef7SLois Curfman McInnes   SNES          snes;
95db0b4b35SLois Curfman McInnes   double        h, norm, sum, umin, noise, ovalues[3],values[3];
96db0b4b35SLois Curfman McInnes   Scalar        hs, dot, mone = -1.0;
97d2dddef7SLois Curfman McInnes   Vec           w,U,F;
98295f7fbbSLois Curfman McInnes   int           ierr, iter, (*eval_fct)(SNES,Vec,Vec);
99e6ed2b1fSLois Curfman McInnes   MPI_Comm      comm;
100d2dddef7SLois Curfman McInnes 
10107a3eed9SLois Curfman McInnes   PetscFunctionBegin;
10207a3eed9SLois Curfman McInnes 
103d2dddef7SLois Curfman McInnes   /* We log matrix-free matrix-vector products separately, so that we can
104d2dddef7SLois Curfman McInnes      separate the performance monitoring from the cases that use conventional
105d2dddef7SLois Curfman McInnes      storage.  We may eventually modify event logging to associate events
106d2dddef7SLois Curfman McInnes      with particular objects, hence alleviating the more general problem. */
107d2dddef7SLois Curfman McInnes   PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0);
108d2dddef7SLois Curfman McInnes 
109e6ed2b1fSLois Curfman McInnes   PetscObjectGetComm((PetscObject)mat,&comm);
110e6ed2b1fSLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
111e6ed2b1fSLois Curfman McInnes   snes = ctx->snes;
112e6ed2b1fSLois Curfman McInnes   w    = ctx->w;
113e6ed2b1fSLois Curfman McInnes   umin = ctx->umin;
114e6ed2b1fSLois Curfman McInnes 
115d2dddef7SLois Curfman McInnes   ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr);
116d2dddef7SLois Curfman McInnes   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
117d2dddef7SLois Curfman McInnes     eval_fct = SNESComputeFunction;
118d2dddef7SLois Curfman McInnes     ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr);
119d2dddef7SLois Curfman McInnes   }
120d2dddef7SLois Curfman McInnes   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
121d2dddef7SLois Curfman McInnes     eval_fct = SNESComputeGradient;
122d2dddef7SLois Curfman McInnes     ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr);
123d2dddef7SLois Curfman McInnes   }
124d2dddef7SLois Curfman McInnes   else SETERRQ(1,0,"Invalid method class");
125d2dddef7SLois Curfman McInnes 
126295f7fbbSLois Curfman McInnes 
127d2dddef7SLois Curfman McInnes   /* Determine a "good" step size, h */
128d2dddef7SLois Curfman McInnes   if (ctx->need_h) {
129d2dddef7SLois Curfman McInnes 
130d2dddef7SLois Curfman McInnes     /* Use Jorge's method to compute h */
131d2dddef7SLois Curfman McInnes     if (ctx->jorge) {
132d2dddef7SLois Curfman McInnes       ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h); CHKERRQ(ierr);
133d2dddef7SLois Curfman McInnes 
134d2dddef7SLois Curfman McInnes     /* Use the Brown/Saad method to compute h */
135d2dddef7SLois Curfman McInnes     } else {
136295f7fbbSLois Curfman McInnes       /* Compute error if desired */
137295f7fbbSLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&iter); CHKERRQ(ierr);
138295f7fbbSLois Curfman McInnes       if ((ctx->need_err) ||
139295f7fbbSLois Curfman McInnes           ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter-1)%ctx->compute_err_freq)))) {
140d2dddef7SLois Curfman McInnes         /* Use Jorge's method to compute noise */
141d2dddef7SLois Curfman McInnes         ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h); CHKERRQ(ierr);
142d2dddef7SLois Curfman McInnes         ctx->error_rel = sqrt(noise);
143d2dddef7SLois Curfman McInnes         PLogInfo(snes,"SNESMatrixFreeMult2_Private: Using Jorge's noise: noise=%g, sqrt(noise)=%g, h_more=%g\n",
144d2dddef7SLois Curfman McInnes             noise,ctx->error_rel,h);
145295f7fbbSLois Curfman McInnes         ctx->compute_err_iter = iter;
146d2dddef7SLois Curfman McInnes         ctx->need_err = 0;
147d2dddef7SLois Curfman McInnes       }
148e6ed2b1fSLois Curfman McInnes 
149e6ed2b1fSLois Curfman McInnes       /*
150d2dddef7SLois Curfman McInnes       ierr = VecDot(U,a,&dot); CHKERRQ(ierr);
151d2dddef7SLois Curfman McInnes       ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr);
152d2dddef7SLois Curfman McInnes       ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr);
153e6ed2b1fSLois Curfman McInnes       */
154e6ed2b1fSLois Curfman McInnes 
155e6ed2b1fSLois Curfman McInnes      /*
156e6ed2b1fSLois Curfman McInnes         Call the Seq Vector routines and then do a single reduction
157e6ed2b1fSLois Curfman McInnes         to reduce the number of communications required
158e6ed2b1fSLois Curfman McInnes      */
159e6ed2b1fSLois Curfman McInnes 
160db0b4b35SLois Curfman McInnes #if !defined(USE_PETSC_COMPLEX)
161e6ed2b1fSLois Curfman McInnes      PLogEventBegin(VEC_Dot,U,a,0,0);
162e6ed2b1fSLois Curfman McInnes      ierr = VecDot_Seq(U,a,ovalues); CHKERRQ(ierr);
163e6ed2b1fSLois Curfman McInnes      PLogEventEnd(VEC_Dot,U,a,0,0);
164e6ed2b1fSLois Curfman McInnes      PLogEventBegin(VEC_Norm,a,0,0,0);
165e6ed2b1fSLois Curfman McInnes      ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr);
166e6ed2b1fSLois Curfman McInnes      ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr);
167e6ed2b1fSLois Curfman McInnes      ovalues[2] = ovalues[2]*ovalues[2];
168e6ed2b1fSLois Curfman McInnes      MPI_Allreduce(ovalues,values,3,MPI_DOUBLE,MPI_SUM,comm );
169e6ed2b1fSLois Curfman McInnes      dot = values[0]; sum = values[1]; norm = sqrt(values[2]);
170e6ed2b1fSLois Curfman McInnes      PLogEventEnd(VEC_Norm,a,0,0,0);
171e6ed2b1fSLois Curfman McInnes #else
172e6ed2b1fSLois Curfman McInnes      {
173e6ed2b1fSLois Curfman McInnes        Scalar cvalues[3],covalues[3];
174e6ed2b1fSLois Curfman McInnes 
175e6ed2b1fSLois Curfman McInnes        PLogEventBegin(VEC_Dot,U,a,0,0);
176e6ed2b1fSLois Curfman McInnes        ierr = VecDot_Seq(U,a,covalues); CHKERRQ(ierr);
177e6ed2b1fSLois Curfman McInnes        PLogEventEnd(VEC_Dot,U,a,0,0);
178e6ed2b1fSLois Curfman McInnes        PLogEventBegin(VEC_Norm,a,0,0,0);
179e6ed2b1fSLois Curfman McInnes        ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr);
180e6ed2b1fSLois Curfman McInnes        ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr);
181e6ed2b1fSLois Curfman McInnes        covalues[1] = ovalues[1];
182e6ed2b1fSLois Curfman McInnes        covalues[2] = ovalues[2]*ovalues[2];
183e6ed2b1fSLois Curfman McInnes        MPI_Allreduce(covalues,cvalues,6,MPI_DOUBLE,MPI_SUM,comm );
184e6ed2b1fSLois Curfman McInnes        dot = cvalues[0]; sum = PetscReal(cvalues[1]); norm = sqrt(PetscReal(cvalues[2]));
185e6ed2b1fSLois Curfman McInnes        PLogEventEnd(VEC_Norm,a,0,0,0);
186e6ed2b1fSLois Curfman McInnes      }
187e6ed2b1fSLois Curfman McInnes #endif
188d2dddef7SLois Curfman McInnes 
189d2dddef7SLois Curfman McInnes       /* Safeguard for step sizes too small */
190d2dddef7SLois Curfman McInnes       if (sum == 0.0) {dot = 1.0; norm = 1.0;}
191db0b4b35SLois Curfman McInnes #if defined(USE_PETSC_COMPLEX)
192e20fef11SSatish Balay       else if (PetscAbsScalar(dot) < umin*sum && PetscReal(dot) >= 0.0) dot = umin*sum;
193e20fef11SSatish Balay       else if (PetscAbsScalar(dot) < 0.0 && PetscReal(dot) > -umin*sum) dot = -umin*sum;
194d2dddef7SLois Curfman McInnes #else
195d2dddef7SLois Curfman McInnes       else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
196d2dddef7SLois Curfman McInnes       else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
197d2dddef7SLois Curfman McInnes #endif
198db0b4b35SLois Curfman McInnes       h = PetscReal(ctx->error_rel*dot/(norm*norm));
199d2dddef7SLois Curfman McInnes     }
200d2dddef7SLois Curfman McInnes   } else {
201d2dddef7SLois Curfman McInnes     h = ctx->h;
202d2dddef7SLois Curfman McInnes   }
203d2dddef7SLois Curfman McInnes   if (!ctx->jorge || !ctx->need_h) PLogInfo(snes,"SNESMatrixFreeMult2_Private: h = %g\n",h);
204d2dddef7SLois Curfman McInnes 
205d2dddef7SLois Curfman McInnes   /* Evaluate function at F(u + ha) */
206db0b4b35SLois Curfman McInnes   hs = h;
207db0b4b35SLois Curfman McInnes   ierr = VecWAXPY(&hs,a,U,w); CHKERRQ(ierr);
208d2dddef7SLois Curfman McInnes   ierr = eval_fct(snes,w,y); CHKERRQ(ierr);
209d2dddef7SLois Curfman McInnes   ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr);
210db0b4b35SLois Curfman McInnes   hs = 1.0/hs;
211db0b4b35SLois Curfman McInnes   ierr = VecScale(&hs,y); CHKERRQ(ierr);
212d2dddef7SLois Curfman McInnes   if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);}
213d2dddef7SLois Curfman McInnes 
214d2dddef7SLois Curfman McInnes   PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0);
21507a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
216d2dddef7SLois Curfman McInnes }
217d2dddef7SLois Curfman McInnes 
218d2dddef7SLois Curfman McInnes #undef __FUNC__
2191a7c28d5SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeMatCreate2"
220d2dddef7SLois Curfman McInnes /*@C
221d2dddef7SLois Curfman McInnes    SNESMatrixFreeMatCreate2 - Creates a matrix-free matrix
222d2dddef7SLois Curfman McInnes    context for use with a SNES solver.  This matrix can be used as
223d2dddef7SLois Curfman McInnes    the Jacobian argument for the routine SNESSetJacobian().
224d2dddef7SLois Curfman McInnes 
225d2dddef7SLois Curfman McInnes    Input Parameters:
226d2dddef7SLois Curfman McInnes .  snes - the SNES context
227d2dddef7SLois Curfman McInnes .  x - vector where SNES solution is to be stored.
228d2dddef7SLois Curfman McInnes 
229d2dddef7SLois Curfman McInnes    Output Parameter:
230d2dddef7SLois Curfman McInnes .  J - the matrix-free matrix
231d2dddef7SLois Curfman McInnes 
232d2dddef7SLois Curfman McInnes    Notes:
233d2dddef7SLois Curfman McInnes    The matrix-free matrix context merely contains the function pointers
234d2dddef7SLois Curfman McInnes    and work space for performing finite difference approximations of
235d2dddef7SLois Curfman McInnes    Jacobian-vector products, J(u)*a, via
236d2dddef7SLois Curfman McInnes 
237d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h,
238d2dddef7SLois Curfman McInnes $   where by default
239d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
240d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
241d2dddef7SLois Curfman McInnes $   where
242d2dddef7SLois Curfman McInnes $        error_rel = square root of relative error in
243d2dddef7SLois Curfman McInnes $                    function evaluation
244d2dddef7SLois Curfman McInnes $        umin = minimum iterate parameter
245d2dddef7SLois Curfman McInnes $   Alternatively, the differencing parameter, h, can be set using
246e6ed2b1fSLois Curfman McInnes $   Jorge's nifty new strategy if one specifies the option
247e6ed2b1fSLois Curfman McInnes $          -snes_mf_jorge
248d2dddef7SLois Curfman McInnes 
249*758f395bSBarry Smith    The user can set these parameters via MatSNESFDMFSetFunctionError().
250d2dddef7SLois Curfman McInnes    See the nonlinear solvers chapter of the users manual for details.
251d2dddef7SLois Curfman McInnes 
252d2dddef7SLois Curfman McInnes    The user should call MatDestroy() when finished with the matrix-free
253d2dddef7SLois Curfman McInnes    matrix context.
254d2dddef7SLois Curfman McInnes 
255d2dddef7SLois Curfman McInnes    Options Database Keys:
256d2dddef7SLois Curfman McInnes $  -snes_mf_err <error_rel>
257d2dddef7SLois Curfman McInnes $  -snes_mf_unim <umin>
258295f7fbbSLois Curfman McInnes $  -snes_mf_compute_err
259295f7fbbSLois Curfman McInnes $  -snes_mf_freq_err <freq>
260e6ed2b1fSLois Curfman McInnes $  -snes_mf_jorge
261d2dddef7SLois Curfman McInnes 
262d2dddef7SLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix
263d2dddef7SLois Curfman McInnes 
264*758f395bSBarry Smith .seealso: MatDestroy(), MatSNESFDMFSetFunctionError()
265d2dddef7SLois Curfman McInnes @*/
266*758f395bSBarry Smith int SNESDefaultMatrixFreeCreate2(SNES snes,Vec x, Mat *J)
267d2dddef7SLois Curfman McInnes {
268d2dddef7SLois Curfman McInnes   MPI_Comm      comm;
269d2dddef7SLois Curfman McInnes   MFCtx_Private *mfctx;
270d2dddef7SLois Curfman McInnes   int           n, nloc, ierr, flg;
271d2dddef7SLois Curfman McInnes   char          p[64];
272d2dddef7SLois Curfman McInnes 
27307a3eed9SLois Curfman McInnes   PetscFunctionBegin;
274d2dddef7SLois Curfman McInnes   mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx);
27507a3eed9SLois Curfman McInnes   PetscMemzero(mfctx,sizeof(MFCtx_Private));
276d2dddef7SLois Curfman McInnes   PLogObjectMemory(snes,sizeof(MFCtx_Private));
277d2dddef7SLois Curfman McInnes   mfctx->sp   = 0;
278d2dddef7SLois Curfman McInnes   mfctx->snes = snes;
279d2dddef7SLois Curfman McInnes   mfctx->error_rel        = 1.e-8; /* assumes double precision */
280d2dddef7SLois Curfman McInnes   mfctx->umin             = 1.e-6;
281d2dddef7SLois Curfman McInnes   mfctx->h                = 0.0;
282d2dddef7SLois Curfman McInnes   mfctx->need_h           = 1;
283d2dddef7SLois Curfman McInnes   mfctx->need_err         = 0;
284295f7fbbSLois Curfman McInnes   mfctx->compute_err      = 0;
285295f7fbbSLois Curfman McInnes   mfctx->compute_err_freq = 0;
286295f7fbbSLois Curfman McInnes   mfctx->compute_err_iter = -1;
287d2dddef7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr);
288d2dddef7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr);
289e6ed2b1fSLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf_jorge",&mfctx->jorge); CHKERRQ(ierr);
290295f7fbbSLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf_compute_err",&mfctx->compute_err); CHKERRQ(ierr);
291295f7fbbSLois Curfman McInnes   ierr = OptionsGetInt(snes->prefix,"-snes_mf_freq_err",&mfctx->compute_err_freq,&flg); CHKERRQ(ierr);
292295f7fbbSLois Curfman McInnes   if (flg) {
293295f7fbbSLois Curfman McInnes     if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0;
294295f7fbbSLois Curfman McInnes     mfctx->compute_err = 1;
295295f7fbbSLois Curfman McInnes   }
296295f7fbbSLois Curfman McInnes   if (mfctx->compute_err == 1) mfctx->need_err = 1;
297295f7fbbSLois Curfman McInnes   if (mfctx->jorge || mfctx->compute_err) {
298d2dddef7SLois Curfman McInnes     ierr = DiffParameterCreate_More(snes,x,&mfctx->data); CHKERRQ(ierr);
299d2dddef7SLois Curfman McInnes   } else mfctx->data = 0;
300d2dddef7SLois Curfman McInnes 
301d2dddef7SLois Curfman McInnes   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
302d2dddef7SLois Curfman McInnes   PetscStrcpy(p,"-");
303d2dddef7SLois Curfman McInnes   if (snes->prefix) PetscStrcat(p,snes->prefix);
304d2dddef7SLois Curfman McInnes   if (flg) {
305295f7fbbSLois Curfman McInnes     PetscPrintf(snes->comm," Matrix-free Options (via SNES):\n");
306bd55573bSLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_err <err>: set sqrt of relative error in function (default %g)\n",p,mfctx->error_rel);
307d2dddef7SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_umin <umin>: see users manual (default %g)\n",p,mfctx->umin);
308295f7fbbSLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_jorge: use Jorge More's method\n",p);
309295f7fbbSLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_compute_err: compute sqrt or relative error in function\n",p);
310295f7fbbSLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)\n",p);
311295f7fbbSLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_noise_file <file>: set file for printing noise info\n",p);
312d2dddef7SLois Curfman McInnes   }
313d2dddef7SLois Curfman McInnes   ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr);
314d2dddef7SLois Curfman McInnes   ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr);
315d2dddef7SLois Curfman McInnes   ierr = VecGetSize(x,&n); CHKERRQ(ierr);
316d2dddef7SLois Curfman McInnes   ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr);
317d2dddef7SLois Curfman McInnes   ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr);
318d2dddef7SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult2_Private);CHKERRQ(ierr);
319d2dddef7SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy2_Private);CHKERRQ(ierr);
320d2dddef7SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView2_Private); CHKERRQ(ierr);
321d2dddef7SLois Curfman McInnes 
322d2dddef7SLois Curfman McInnes   PLogObjectParent(*J,mfctx->w);
323d2dddef7SLois Curfman McInnes   PLogObjectParent(snes,*J);
32407a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
325d2dddef7SLois Curfman McInnes }
326d2dddef7SLois Curfman McInnes 
327d2dddef7SLois Curfman McInnes #undef __FUNC__
328*758f395bSBarry Smith #define __FUNC__ "SNESDefaultMatrixFreeSetParameters2"
329d2dddef7SLois Curfman McInnes /*@
330*758f395bSBarry Smith    SNESDefaultMatrixFreeSetParameters2 - Sets the parameters for the approximation of
331d2dddef7SLois Curfman McInnes    matrix-vector products using finite differences.
332d2dddef7SLois Curfman McInnes 
333d2dddef7SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h where
334d2dddef7SLois Curfman McInnes 
335d2dddef7SLois Curfman McInnes    either the user sets h directly here, or this parameter is computed via
336d2dddef7SLois Curfman McInnes 
337d2dddef7SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
338d2dddef7SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
339d2dddef7SLois Curfman McInnes $
340d2dddef7SLois Curfman McInnes 
341d2dddef7SLois Curfman McInnes    Input Parameters:
342*758f395bSBarry Smith +  mat - the matrix
343d2dddef7SLois Curfman McInnes .  error_rel - relative error (should be set to the square root of
344d2dddef7SLois Curfman McInnes      the relative error in the function evaluations)
345d2dddef7SLois Curfman McInnes .  umin - minimum allowable u-value
346*758f395bSBarry Smith -  h - differencing parameter
347d2dddef7SLois Curfman McInnes 
348d2dddef7SLois Curfman McInnes    Notes:
349d2dddef7SLois Curfman McInnes    If the user sets the parameter h directly, then this value will be used
350d2dddef7SLois Curfman McInnes    instead of the default computation indicated above.
351d2dddef7SLois Curfman McInnes 
352d2dddef7SLois Curfman McInnes .keywords: SNES, matrix-free, parameters
353d2dddef7SLois Curfman McInnes 
354*758f395bSBarry Smith .seealso: MatCreateSNESFDMF()
355d2dddef7SLois Curfman McInnes @*/
356*758f395bSBarry Smith int SNESDefaultMatrixFreeSetParameters2(Mat mat,double error,double umin,double h)
357d2dddef7SLois Curfman McInnes {
358d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
359d2dddef7SLois Curfman McInnes   int           ierr;
360d2dddef7SLois Curfman McInnes 
36107a3eed9SLois Curfman McInnes   PetscFunctionBegin;
362d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
363d2dddef7SLois Curfman McInnes   if (ctx) {
364d2dddef7SLois Curfman McInnes     if (error != PETSC_DEFAULT) ctx->error_rel = error;
365d2dddef7SLois Curfman McInnes     if (umin  != PETSC_DEFAULT) ctx->umin = umin;
366d2dddef7SLois Curfman McInnes     if (h     != PETSC_DEFAULT) {
367d2dddef7SLois Curfman McInnes       ctx->h = h;
368d2dddef7SLois Curfman McInnes       ctx->need_h = 0;
369d2dddef7SLois Curfman McInnes     }
370d2dddef7SLois Curfman McInnes   }
37107a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
372d2dddef7SLois Curfman McInnes }
373d2dddef7SLois Curfman McInnes 
374d2dddef7SLois Curfman McInnes int SNESUnSetMatrixFreeParameter(SNES snes)
375d2dddef7SLois Curfman McInnes {
376d2dddef7SLois Curfman McInnes   MFCtx_Private *ctx;
377d2dddef7SLois Curfman McInnes   int           ierr;
378d2dddef7SLois Curfman McInnes   Mat           mat;
379d2dddef7SLois Curfman McInnes 
38007a3eed9SLois Curfman McInnes   PetscFunctionBegin;
381d2dddef7SLois Curfman McInnes   ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr);
382d2dddef7SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
383d2dddef7SLois Curfman McInnes   if (ctx) ctx->need_h = 1;
38407a3eed9SLois Curfman McInnes   PetscFunctionReturn(0);
385d2dddef7SLois Curfman McInnes }
386d2dddef7SLois Curfman McInnes 
387