1*d2dddef7SLois Curfman McInnes 2*d2dddef7SLois Curfman McInnes #ifndef lint 3*d2dddef7SLois Curfman McInnes static char vcid[] = "$Id: utils.c,v 1.8 1997/05/07 16:25:49 curfman Exp $"; 4*d2dddef7SLois Curfman McInnes #endif 5*d2dddef7SLois Curfman McInnes 6*d2dddef7SLois Curfman McInnes #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 7*d2dddef7SLois Curfman McInnes #include "pinclude/pviewer.h" 8*d2dddef7SLois Curfman McInnes #include <math.h> 9*d2dddef7SLois Curfman McInnes 10*d2dddef7SLois Curfman McInnes extern int DiffParameterCreate_More(SNES,Vec,void**); 11*d2dddef7SLois Curfman McInnes extern int DiffParameterCompute_More(SNES,void*,Vec,Vec,double*,double*); 12*d2dddef7SLois Curfman McInnes extern int DiffParameterDestroy_More(void*); 13*d2dddef7SLois Curfman McInnes 14*d2dddef7SLois Curfman McInnes typedef struct { /* default context for matrix-free SNES */ 15*d2dddef7SLois Curfman McInnes SNES snes; /* SNES context */ 16*d2dddef7SLois Curfman McInnes Vec w; /* work vector */ 17*d2dddef7SLois Curfman McInnes PCNullSpace sp; /* null space context */ 18*d2dddef7SLois Curfman McInnes double error_rel; /* square root of relative error in computing function */ 19*d2dddef7SLois Curfman McInnes double umin; /* minimum allowable u'a value relative to |u|_1 */ 20*d2dddef7SLois Curfman McInnes int jorge; /* flag indicating use of Jorge's method for determining 21*d2dddef7SLois Curfman McInnes the differencing parameter */ 22*d2dddef7SLois Curfman McInnes Scalar h; /* differencing parameter */ 23*d2dddef7SLois Curfman McInnes int need_h; /* flag indicating whether we must compute h */ 24*d2dddef7SLois Curfman McInnes int need_err; /* flag indicating whether we must compute error_rel */ 25*d2dddef7SLois Curfman McInnes void *data; /* implementation-specific data */ 26*d2dddef7SLois Curfman McInnes } MFCtx_Private; 27*d2dddef7SLois Curfman McInnes 28*d2dddef7SLois Curfman McInnes #undef __FUNC__ 29*d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeDestroy2_Private" /* ADIC Ignore */ 30*d2dddef7SLois Curfman McInnes int SNESMatrixFreeDestroy2_Private(PetscObject obj) 31*d2dddef7SLois Curfman McInnes { 32*d2dddef7SLois Curfman McInnes int ierr; 33*d2dddef7SLois Curfman McInnes Mat mat = (Mat) obj; 34*d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 35*d2dddef7SLois Curfman McInnes 36*d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx); 37*d2dddef7SLois Curfman McInnes ierr = VecDestroy(ctx->w); CHKERRQ(ierr); 38*d2dddef7SLois Curfman McInnes if (ctx->sp) {ierr = PCNullSpaceDestroy(ctx->sp); CHKERRQ(ierr);} 39*d2dddef7SLois Curfman McInnes if (ctx->jorge || ctx->need_err) {ierr = DiffParameterDestroy_More(ctx->data); CHKERRQ(ierr);} 40*d2dddef7SLois Curfman McInnes PetscFree(ctx); 41*d2dddef7SLois Curfman McInnes return 0; 42*d2dddef7SLois Curfman McInnes } 43*d2dddef7SLois Curfman McInnes 44*d2dddef7SLois Curfman McInnes #undef __FUNC__ 45*d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeView2_Private" /* ADIC Ignore */ 46*d2dddef7SLois Curfman McInnes /* 47*d2dddef7SLois Curfman McInnes SNESMatrixFreeView2_Private - Views matrix-free parameters. 48*d2dddef7SLois Curfman McInnes */ 49*d2dddef7SLois Curfman McInnes int SNESMatrixFreeView2_Private(Mat J,Viewer viewer) 50*d2dddef7SLois Curfman McInnes { 51*d2dddef7SLois Curfman McInnes int ierr; 52*d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 53*d2dddef7SLois Curfman McInnes MPI_Comm comm; 54*d2dddef7SLois Curfman McInnes FILE *fd; 55*d2dddef7SLois Curfman McInnes ViewerType vtype; 56*d2dddef7SLois Curfman McInnes 57*d2dddef7SLois Curfman McInnes PetscObjectGetComm((PetscObject)J,&comm); 58*d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 59*d2dddef7SLois Curfman McInnes ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 60*d2dddef7SLois Curfman McInnes ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 61*d2dddef7SLois Curfman McInnes if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 62*d2dddef7SLois Curfman McInnes PetscFPrintf(comm,fd," SNES matrix-free approximation:\n"); 63*d2dddef7SLois Curfman McInnes if (ctx->jorge) 64*d2dddef7SLois Curfman McInnes PetscFPrintf(comm,fd," using Jorge's method of determining differencing parameter\n"); 65*d2dddef7SLois Curfman McInnes PetscFPrintf(comm,fd," err=%g (relative error in function evaluation)\n",ctx->error_rel); 66*d2dddef7SLois Curfman McInnes PetscFPrintf(comm,fd," umin=%g (minimum iterate parameter)\n",ctx->umin); 67*d2dddef7SLois Curfman McInnes } 68*d2dddef7SLois Curfman McInnes return 0; 69*d2dddef7SLois Curfman McInnes } 70*d2dddef7SLois Curfman McInnes 71*d2dddef7SLois Curfman McInnes #undef __FUNC__ 72*d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeMult2_Private" 73*d2dddef7SLois Curfman McInnes /* 74*d2dddef7SLois Curfman McInnes SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector 75*d2dddef7SLois Curfman McInnes product, y = F'(u)*a: 76*d2dddef7SLois Curfman McInnes y = ( F(u + ha) - F(u) ) /h, 77*d2dddef7SLois Curfman McInnes where F = nonlinear function, as set by SNESSetFunction() 78*d2dddef7SLois Curfman McInnes u = current iterate 79*d2dddef7SLois Curfman McInnes h = difference interval 80*d2dddef7SLois Curfman McInnes */ 81*d2dddef7SLois Curfman McInnes int SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y) 82*d2dddef7SLois Curfman McInnes { 83*d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 84*d2dddef7SLois Curfman McInnes SNES snes; 85*d2dddef7SLois Curfman McInnes double norm, sum, umin, noise; 86*d2dddef7SLois Curfman McInnes Scalar h, dot, mone = -1.0; 87*d2dddef7SLois Curfman McInnes Vec w,U,F; 88*d2dddef7SLois Curfman McInnes int ierr, (*eval_fct)(SNES,Vec,Vec); 89*d2dddef7SLois Curfman McInnes 90*d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 91*d2dddef7SLois Curfman McInnes snes = ctx->snes; 92*d2dddef7SLois Curfman McInnes w = ctx->w; 93*d2dddef7SLois Curfman McInnes umin = ctx->umin; 94*d2dddef7SLois Curfman McInnes 95*d2dddef7SLois Curfman McInnes /* We log matrix-free matrix-vector products separately, so that we can 96*d2dddef7SLois Curfman McInnes separate the performance monitoring from the cases that use conventional 97*d2dddef7SLois Curfman McInnes storage. We may eventually modify event logging to associate events 98*d2dddef7SLois Curfman McInnes with particular objects, hence alleviating the more general problem. */ 99*d2dddef7SLois Curfman McInnes PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0); 100*d2dddef7SLois Curfman McInnes 101*d2dddef7SLois Curfman McInnes ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr); 102*d2dddef7SLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 103*d2dddef7SLois Curfman McInnes eval_fct = SNESComputeFunction; 104*d2dddef7SLois Curfman McInnes ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr); 105*d2dddef7SLois Curfman McInnes } 106*d2dddef7SLois Curfman McInnes else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 107*d2dddef7SLois Curfman McInnes eval_fct = SNESComputeGradient; 108*d2dddef7SLois Curfman McInnes ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr); 109*d2dddef7SLois Curfman McInnes } 110*d2dddef7SLois Curfman McInnes else SETERRQ(1,0,"Invalid method class"); 111*d2dddef7SLois Curfman McInnes 112*d2dddef7SLois Curfman McInnes /* Determine a "good" step size, h */ 113*d2dddef7SLois Curfman McInnes if (ctx->need_h) { 114*d2dddef7SLois Curfman McInnes 115*d2dddef7SLois Curfman McInnes /* Use Jorge's method to compute h */ 116*d2dddef7SLois Curfman McInnes if (ctx->jorge) { 117*d2dddef7SLois Curfman McInnes ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h); CHKERRQ(ierr); 118*d2dddef7SLois Curfman McInnes 119*d2dddef7SLois Curfman McInnes /* Use the Brown/Saad method to compute h */ 120*d2dddef7SLois Curfman McInnes } else { 121*d2dddef7SLois Curfman McInnes if (ctx->need_err) { 122*d2dddef7SLois Curfman McInnes /* Use Jorge's method to compute noise */ 123*d2dddef7SLois Curfman McInnes ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h); CHKERRQ(ierr); 124*d2dddef7SLois Curfman McInnes ctx->error_rel = sqrt(noise); 125*d2dddef7SLois Curfman McInnes PLogInfo(snes,"SNESMatrixFreeMult2_Private: Using Jorge's noise: noise=%g, sqrt(noise)=%g, h_more=%g\n", 126*d2dddef7SLois Curfman McInnes noise,ctx->error_rel,h); 127*d2dddef7SLois Curfman McInnes 128*d2dddef7SLois Curfman McInnes ctx->need_err = 0; 129*d2dddef7SLois Curfman McInnes } 130*d2dddef7SLois Curfman McInnes ierr = VecDot(U,a,&dot); CHKERRQ(ierr); 131*d2dddef7SLois Curfman McInnes ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr); 132*d2dddef7SLois Curfman McInnes ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr); 133*d2dddef7SLois Curfman McInnes 134*d2dddef7SLois Curfman McInnes /* Safeguard for step sizes too small */ 135*d2dddef7SLois Curfman McInnes if (sum == 0.0) {dot = 1.0; norm = 1.0;} 136*d2dddef7SLois Curfman McInnes #if defined(PETSC_COMPLEX) 137*d2dddef7SLois Curfman McInnes else if (abs(dot) < umin*sum && real(dot) >= 0.0) dot = umin*sum; 138*d2dddef7SLois Curfman McInnes else if (abs(dot) < 0.0 && real(dot) > -umin*sum) dot = -umin*sum; 139*d2dddef7SLois Curfman McInnes #else 140*d2dddef7SLois Curfman McInnes else if (dot < umin*sum && dot >= 0.0) dot = umin*sum; 141*d2dddef7SLois Curfman McInnes else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum; 142*d2dddef7SLois Curfman McInnes #endif 143*d2dddef7SLois Curfman McInnes h = ctx->error_rel*dot/(norm*norm); 144*d2dddef7SLois Curfman McInnes } 145*d2dddef7SLois Curfman McInnes } else { 146*d2dddef7SLois Curfman McInnes h = ctx->h; 147*d2dddef7SLois Curfman McInnes } 148*d2dddef7SLois Curfman McInnes if (!ctx->jorge || !ctx->need_h) PLogInfo(snes,"SNESMatrixFreeMult2_Private: h = %g\n",h); 149*d2dddef7SLois Curfman McInnes 150*d2dddef7SLois Curfman McInnes /* Evaluate function at F(u + ha) */ 151*d2dddef7SLois Curfman McInnes ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr); 152*d2dddef7SLois Curfman McInnes ierr = eval_fct(snes,w,y); CHKERRQ(ierr); 153*d2dddef7SLois Curfman McInnes ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr); 154*d2dddef7SLois Curfman McInnes h = 1.0/h; 155*d2dddef7SLois Curfman McInnes ierr = VecScale(&h,y); CHKERRQ(ierr); 156*d2dddef7SLois Curfman McInnes if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);} 157*d2dddef7SLois Curfman McInnes 158*d2dddef7SLois Curfman McInnes PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0); 159*d2dddef7SLois Curfman McInnes return 0; 160*d2dddef7SLois Curfman McInnes } 161*d2dddef7SLois Curfman McInnes 162*d2dddef7SLois Curfman McInnes #undef __FUNC__ 163*d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESDefaultMatrixFreeMatCreate2" 164*d2dddef7SLois Curfman McInnes /*@C 165*d2dddef7SLois Curfman McInnes SNESMatrixFreeMatCreate2 - Creates a matrix-free matrix 166*d2dddef7SLois Curfman McInnes context for use with a SNES solver. This matrix can be used as 167*d2dddef7SLois Curfman McInnes the Jacobian argument for the routine SNESSetJacobian(). 168*d2dddef7SLois Curfman McInnes 169*d2dddef7SLois Curfman McInnes Input Parameters: 170*d2dddef7SLois Curfman McInnes . snes - the SNES context 171*d2dddef7SLois Curfman McInnes . x - vector where SNES solution is to be stored. 172*d2dddef7SLois Curfman McInnes 173*d2dddef7SLois Curfman McInnes Output Parameter: 174*d2dddef7SLois Curfman McInnes . J - the matrix-free matrix 175*d2dddef7SLois Curfman McInnes 176*d2dddef7SLois Curfman McInnes Notes: 177*d2dddef7SLois Curfman McInnes The matrix-free matrix context merely contains the function pointers 178*d2dddef7SLois Curfman McInnes and work space for performing finite difference approximations of 179*d2dddef7SLois Curfman McInnes Jacobian-vector products, J(u)*a, via 180*d2dddef7SLois Curfman McInnes 181*d2dddef7SLois Curfman McInnes $ J(u)*a = [J(u+h*a) - J(u)]/h, 182*d2dddef7SLois Curfman McInnes $ where by default 183*d2dddef7SLois Curfman McInnes $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 184*d2dddef7SLois Curfman McInnes $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 185*d2dddef7SLois Curfman McInnes $ where 186*d2dddef7SLois Curfman McInnes $ error_rel = square root of relative error in 187*d2dddef7SLois Curfman McInnes $ function evaluation 188*d2dddef7SLois Curfman McInnes $ umin = minimum iterate parameter 189*d2dddef7SLois Curfman McInnes $ Alternatively, the differencing parameter, h, can be set using 190*d2dddef7SLois Curfman McInnes $ Jorge's nifty new strategy if the option 191*d2dddef7SLois Curfman McInnes $ -matrix_free_jorge 192*d2dddef7SLois Curfman McInnes # is used. 193*d2dddef7SLois Curfman McInnes 194*d2dddef7SLois Curfman McInnes The user can set these parameters via SNESSetMatrixFreeParameters(). 195*d2dddef7SLois Curfman McInnes See the nonlinear solvers chapter of the users manual for details. 196*d2dddef7SLois Curfman McInnes 197*d2dddef7SLois Curfman McInnes The user should call MatDestroy() when finished with the matrix-free 198*d2dddef7SLois Curfman McInnes matrix context. 199*d2dddef7SLois Curfman McInnes 200*d2dddef7SLois Curfman McInnes Options Database Keys: 201*d2dddef7SLois Curfman McInnes $ -snes_mf_err <error_rel> 202*d2dddef7SLois Curfman McInnes $ -snes_mf_unim <umin> 203*d2dddef7SLois Curfman McInnes $ -matrix_free_jorge 204*d2dddef7SLois Curfman McInnes 205*d2dddef7SLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix 206*d2dddef7SLois Curfman McInnes 207*d2dddef7SLois Curfman McInnes .seealso: MatDestroy(), SNESSetMatrixFreeParameters() 208*d2dddef7SLois Curfman McInnes @*/ 209*d2dddef7SLois Curfman McInnes int SNESMatrixFreeMatCreate2(SNES snes,Vec x, Mat *J) 210*d2dddef7SLois Curfman McInnes { 211*d2dddef7SLois Curfman McInnes MPI_Comm comm; 212*d2dddef7SLois Curfman McInnes MFCtx_Private *mfctx; 213*d2dddef7SLois Curfman McInnes int n, nloc, ierr, flg; 214*d2dddef7SLois Curfman McInnes char p[64]; 215*d2dddef7SLois Curfman McInnes 216*d2dddef7SLois Curfman McInnes mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx); 217*d2dddef7SLois Curfman McInnes PLogObjectMemory(snes,sizeof(MFCtx_Private)); 218*d2dddef7SLois Curfman McInnes mfctx->sp = 0; 219*d2dddef7SLois Curfman McInnes mfctx->snes = snes; 220*d2dddef7SLois Curfman McInnes mfctx->error_rel = 1.e-8; /* assumes double precision */ 221*d2dddef7SLois Curfman McInnes mfctx->umin = 1.e-6; 222*d2dddef7SLois Curfman McInnes mfctx->h = 0.0; 223*d2dddef7SLois Curfman McInnes mfctx->need_h = 1; 224*d2dddef7SLois Curfman McInnes mfctx->need_err = 0; 225*d2dddef7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr); 226*d2dddef7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr); 227*d2dddef7SLois Curfman McInnes ierr = OptionsHasName(snes->prefix,"-matrix_free_jorge",&mfctx->jorge); CHKERRQ(ierr); 228*d2dddef7SLois Curfman McInnes ierr = OptionsHasName(snes->prefix,"-compute_err",&mfctx->need_err); CHKERRQ(ierr); 229*d2dddef7SLois Curfman McInnes if (mfctx->jorge || mfctx->need_err) { 230*d2dddef7SLois Curfman McInnes ierr = DiffParameterCreate_More(snes,x,&mfctx->data); CHKERRQ(ierr); 231*d2dddef7SLois Curfman McInnes } else mfctx->data = 0; 232*d2dddef7SLois Curfman McInnes 233*d2dddef7SLois Curfman McInnes ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 234*d2dddef7SLois Curfman McInnes PetscStrcpy(p,"-"); 235*d2dddef7SLois Curfman McInnes if (snes->prefix) PetscStrcat(p,snes->prefix); 236*d2dddef7SLois Curfman McInnes if (flg) { 237*d2dddef7SLois Curfman McInnes PetscPrintf(snes->comm," -matrix_free_jorge: use Jorge More's method\n"); 238*d2dddef7SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel); 239*d2dddef7SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_umin <umin>: see users manual (default %g)\n",p,mfctx->umin); 240*d2dddef7SLois Curfman McInnes } 241*d2dddef7SLois Curfman McInnes ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr); 242*d2dddef7SLois Curfman McInnes ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr); 243*d2dddef7SLois Curfman McInnes ierr = VecGetSize(x,&n); CHKERRQ(ierr); 244*d2dddef7SLois Curfman McInnes ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr); 245*d2dddef7SLois Curfman McInnes ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr); 246*d2dddef7SLois Curfman McInnes ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult2_Private);CHKERRQ(ierr); 247*d2dddef7SLois Curfman McInnes ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy2_Private);CHKERRQ(ierr); 248*d2dddef7SLois Curfman McInnes ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView2_Private); CHKERRQ(ierr); 249*d2dddef7SLois Curfman McInnes 250*d2dddef7SLois Curfman McInnes PLogObjectParent(*J,mfctx->w); 251*d2dddef7SLois Curfman McInnes PLogObjectParent(snes,*J); 252*d2dddef7SLois Curfman McInnes return 0; 253*d2dddef7SLois Curfman McInnes } 254*d2dddef7SLois Curfman McInnes 255*d2dddef7SLois Curfman McInnes #undef __FUNC__ 256*d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESSetMatrixFreeParameters2" 257*d2dddef7SLois Curfman McInnes /*@ 258*d2dddef7SLois Curfman McInnes SNESSetMatrixFreeParameters2 - Sets the parameters for the approximation of 259*d2dddef7SLois Curfman McInnes matrix-vector products using finite differences. 260*d2dddef7SLois Curfman McInnes 261*d2dddef7SLois Curfman McInnes $ J(u)*a = [J(u+h*a) - J(u)]/h where 262*d2dddef7SLois Curfman McInnes 263*d2dddef7SLois Curfman McInnes either the user sets h directly here, or this parameter is computed via 264*d2dddef7SLois Curfman McInnes 265*d2dddef7SLois Curfman McInnes $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 266*d2dddef7SLois Curfman McInnes $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 267*d2dddef7SLois Curfman McInnes $ 268*d2dddef7SLois Curfman McInnes 269*d2dddef7SLois Curfman McInnes Input Parameters: 270*d2dddef7SLois Curfman McInnes . snes - the SNES context 271*d2dddef7SLois Curfman McInnes . error_rel - relative error (should be set to the square root of 272*d2dddef7SLois Curfman McInnes the relative error in the function evaluations) 273*d2dddef7SLois Curfman McInnes . umin - minimum allowable u-value 274*d2dddef7SLois Curfman McInnes . h - differencing parameter 275*d2dddef7SLois Curfman McInnes 276*d2dddef7SLois Curfman McInnes Notes: 277*d2dddef7SLois Curfman McInnes If the user sets the parameter h directly, then this value will be used 278*d2dddef7SLois Curfman McInnes instead of the default computation indicated above. 279*d2dddef7SLois Curfman McInnes 280*d2dddef7SLois Curfman McInnes .keywords: SNES, matrix-free, parameters 281*d2dddef7SLois Curfman McInnes 282*d2dddef7SLois Curfman McInnes .seealso: SNESDefaultMatrixFreeMatCreate() 283*d2dddef7SLois Curfman McInnes @*/ 284*d2dddef7SLois Curfman McInnes int SNESSetMatrixFreeParameters2(SNES snes,double error,double umin,double h) 285*d2dddef7SLois Curfman McInnes { 286*d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 287*d2dddef7SLois Curfman McInnes int ierr; 288*d2dddef7SLois Curfman McInnes Mat mat; 289*d2dddef7SLois Curfman McInnes 290*d2dddef7SLois Curfman McInnes ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 291*d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 292*d2dddef7SLois Curfman McInnes if (ctx) { 293*d2dddef7SLois Curfman McInnes if (error != PETSC_DEFAULT) ctx->error_rel = error; 294*d2dddef7SLois Curfman McInnes if (umin != PETSC_DEFAULT) ctx->umin = umin; 295*d2dddef7SLois Curfman McInnes if (h != PETSC_DEFAULT) { 296*d2dddef7SLois Curfman McInnes ctx->h = h; 297*d2dddef7SLois Curfman McInnes ctx->need_h = 0; 298*d2dddef7SLois Curfman McInnes } 299*d2dddef7SLois Curfman McInnes } 300*d2dddef7SLois Curfman McInnes return 0; 301*d2dddef7SLois Curfman McInnes } 302*d2dddef7SLois Curfman McInnes 303*d2dddef7SLois Curfman McInnes int SNESUnSetMatrixFreeParameter(SNES snes) 304*d2dddef7SLois Curfman McInnes { 305*d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 306*d2dddef7SLois Curfman McInnes int ierr; 307*d2dddef7SLois Curfman McInnes Mat mat; 308*d2dddef7SLois Curfman McInnes 309*d2dddef7SLois Curfman McInnes ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 310*d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 311*d2dddef7SLois Curfman McInnes if (ctx) ctx->need_h = 1; 312*d2dddef7SLois Curfman McInnes return 0; 313*d2dddef7SLois Curfman McInnes } 314*d2dddef7SLois Curfman McInnes 315