1*ea709b57SSatish Balay /*$Id: snesmfj2.c,v 1.33 2001/08/06 21:17:10 bsmith Exp balay $*/ 2d2dddef7SLois Curfman McInnes 30a835dfdSSatish Balay #include "src/snes/snesimpl.h" /*I "petscsnes.h" I*/ 4d2dddef7SLois Curfman McInnes 5ca44d042SBarry Smith EXTERN int DiffParameterCreate_More(SNES,Vec,void**); 6ca44d042SBarry Smith EXTERN int DiffParameterCompute_More(SNES,void*,Vec,Vec,PetscReal*,PetscReal*); 7ca44d042SBarry Smith EXTERN int DiffParameterDestroy_More(void*); 8d2dddef7SLois Curfman McInnes 9d2dddef7SLois Curfman McInnes typedef struct { /* default context for matrix-free SNES */ 10d2dddef7SLois Curfman McInnes SNES snes; /* SNES context */ 11d2dddef7SLois Curfman McInnes Vec w; /* work vector */ 1274637425SBarry Smith MatNullSpace sp; /* null space context */ 13145595abSSatish Balay PetscReal error_rel; /* square root of relative error in computing function */ 14145595abSSatish Balay PetscReal umin; /* minimum allowable u'a value relative to |u|_1 */ 15d2dddef7SLois Curfman McInnes int jorge; /* flag indicating use of Jorge's method for determining 16d2dddef7SLois Curfman McInnes the differencing parameter */ 17145595abSSatish Balay PetscReal h; /* differencing parameter */ 18d2dddef7SLois Curfman McInnes int need_h; /* flag indicating whether we must compute h */ 19295f7fbbSLois Curfman McInnes int need_err; /* flag indicating whether we must currently compute error_rel */ 20295f7fbbSLois Curfman McInnes int compute_err; /* flag indicating whether we must ever compute error_rel */ 21295f7fbbSLois Curfman McInnes int compute_err_iter; /* last iter where we've computer error_rel */ 22295f7fbbSLois Curfman McInnes int compute_err_freq; /* frequency of computing error_rel */ 23d2dddef7SLois Curfman McInnes void *data; /* implementation-specific data */ 24d2dddef7SLois Curfman McInnes } MFCtx_Private; 25d2dddef7SLois Curfman McInnes 264a2ae208SSatish Balay #undef __FUNCT__ 274a2ae208SSatish Balay #define __FUNCT__ "SNESMatrixFreeDestroy2_Private" /* ADIC Ignore */ 28f1d4f8a3SBarry Smith int SNESMatrixFreeDestroy2_Private(Mat mat) 29d2dddef7SLois Curfman McInnes { 30d2dddef7SLois Curfman McInnes int ierr; 31d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 32d2dddef7SLois Curfman McInnes 3307a3eed9SLois Curfman McInnes PetscFunctionBegin; 34d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx); 35d2dddef7SLois Curfman McInnes ierr = VecDestroy(ctx->w);CHKERRQ(ierr); 3674637425SBarry Smith if (ctx->sp) {ierr = MatNullSpaceDestroy(ctx->sp);CHKERRQ(ierr);} 37295f7fbbSLois Curfman McInnes if (ctx->jorge || ctx->compute_err) {ierr = DiffParameterDestroy_More(ctx->data);CHKERRQ(ierr);} 38606d414cSSatish Balay ierr = PetscFree(ctx);CHKERRQ(ierr); 3907a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 40d2dddef7SLois Curfman McInnes } 41d2dddef7SLois Curfman McInnes 424a2ae208SSatish Balay #undef __FUNCT__ 434a2ae208SSatish Balay #define __FUNCT__ "SNESMatrixFreeView2_Private" /* ADIC Ignore */ 44d2dddef7SLois Curfman McInnes /* 45d2dddef7SLois Curfman McInnes SNESMatrixFreeView2_Private - Views matrix-free parameters. 46d2dddef7SLois Curfman McInnes */ 47b0a32e0cSBarry Smith int SNESMatrixFreeView2_Private(Mat J,PetscViewer viewer) 48d2dddef7SLois Curfman McInnes { 49d2dddef7SLois Curfman McInnes int ierr; 50d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 516831982aSBarry Smith PetscTruth isascii; 52d2dddef7SLois Curfman McInnes 5307a3eed9SLois Curfman McInnes PetscFunctionBegin; 54d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(J,(void **)&ctx);CHKERRQ(ierr); 55b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 560f5bd95cSBarry Smith if (isascii) { 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 } 61b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," err=%g (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr); 62b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," umin=%g (minimum iterate parameter)\n",ctx->umin);CHKERRQ(ierr); 63758f395bSBarry Smith if (ctx->compute_err) { 64b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," freq_err=%d (frequency for computing err)\n",ctx->compute_err_freq);CHKERRQ(ierr); 65d2dddef7SLois Curfman McInnes } 66758f395bSBarry Smith } else { 6729bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by SNES matrix free Jorge",((PetscObject)viewer)->type_name); 68758f395bSBarry Smith } 6907a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 70d2dddef7SLois Curfman McInnes } 71d2dddef7SLois Curfman McInnes 724a2ae208SSatish Balay #undef __FUNCT__ 734a2ae208SSatish Balay #define __FUNCT__ "SNESMatrixFreeMult2_Private" 74d2dddef7SLois Curfman McInnes /* 75d2dddef7SLois Curfman McInnes SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector 76d2dddef7SLois Curfman McInnes product, y = F'(u)*a: 77d2dddef7SLois Curfman McInnes y = (F(u + ha) - F(u)) /h, 78d2dddef7SLois Curfman McInnes where F = nonlinear function, as set by SNESSetFunction() 79d2dddef7SLois Curfman McInnes u = current iterate 80d2dddef7SLois Curfman McInnes h = difference interval 81d2dddef7SLois Curfman McInnes */ 82d2dddef7SLois Curfman McInnes int SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y) 83d2dddef7SLois Curfman McInnes { 84d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 85d2dddef7SLois Curfman McInnes SNES snes; 86145595abSSatish Balay PetscReal h,norm,sum,umin,noise; 87*ea709b57SSatish Balay PetscScalar hs,dot,mone = -1.0; 88d2dddef7SLois Curfman McInnes Vec w,U,F; 89295f7fbbSLois Curfman McInnes int ierr,iter,(*eval_fct)(SNES,Vec,Vec); 90e6ed2b1fSLois Curfman McInnes MPI_Comm comm; 91d2dddef7SLois Curfman McInnes 9207a3eed9SLois Curfman McInnes PetscFunctionBegin; 9307a3eed9SLois Curfman McInnes 94d2dddef7SLois Curfman McInnes /* We log matrix-free matrix-vector products separately, so that we can 95d2dddef7SLois Curfman McInnes separate the performance monitoring from the cases that use conventional 96d2dddef7SLois Curfman McInnes storage. We may eventually modify event logging to associate events 97d2dddef7SLois Curfman McInnes with particular objects, hence alleviating the more general problem. */ 98b0a32e0cSBarry Smith ierr = PetscLogEventBegin(MAT_MatrixFreeMult,a,y,0,0);CHKERRQ(ierr); 99d2dddef7SLois Curfman McInnes 1002d0c0e3bSBarry Smith ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 101e6ed2b1fSLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 102e6ed2b1fSLois Curfman McInnes snes = ctx->snes; 103e6ed2b1fSLois Curfman McInnes w = ctx->w; 104e6ed2b1fSLois Curfman McInnes umin = ctx->umin; 105e6ed2b1fSLois Curfman McInnes 106d2dddef7SLois Curfman McInnes ierr = SNESGetSolution(snes,&U);CHKERRQ(ierr); 107d2dddef7SLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 108d2dddef7SLois Curfman McInnes eval_fct = SNESComputeFunction; 10974637425SBarry Smith ierr = SNESGetFunction(snes,&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 110d2dddef7SLois Curfman McInnes } 111d2dddef7SLois Curfman McInnes else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 112d2dddef7SLois Curfman McInnes eval_fct = SNESComputeGradient; 113184914b5SBarry Smith ierr = SNESGetGradient(snes,&F,PETSC_NULL);CHKERRQ(ierr); 114d2dddef7SLois Curfman McInnes } 11529bbc08cSBarry Smith else SETERRQ(1,"Invalid method class"); 116d2dddef7SLois Curfman McInnes 117295f7fbbSLois Curfman McInnes 118d2dddef7SLois Curfman McInnes /* Determine a "good" step size, h */ 119d2dddef7SLois Curfman McInnes if (ctx->need_h) { 120d2dddef7SLois Curfman McInnes 121d2dddef7SLois Curfman McInnes /* Use Jorge's method to compute h */ 122d2dddef7SLois Curfman McInnes if (ctx->jorge) { 123d2dddef7SLois Curfman McInnes ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr); 124d2dddef7SLois Curfman McInnes 125d2dddef7SLois Curfman McInnes /* Use the Brown/Saad method to compute h */ 126d2dddef7SLois Curfman McInnes } else { 127295f7fbbSLois Curfman McInnes /* Compute error if desired */ 128295f7fbbSLois Curfman McInnes ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr); 129145595abSSatish Balay if ((ctx->need_err) || ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter-1)%ctx->compute_err_freq)))) { 130d2dddef7SLois Curfman McInnes /* Use Jorge's method to compute noise */ 131d2dddef7SLois Curfman McInnes ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr); 132d2dddef7SLois Curfman McInnes ctx->error_rel = sqrt(noise); 133b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESMatrixFreeMult2_Private: Using Jorge's noise: noise=%g, sqrt(noise)=%g, h_more=%g\n",noise,ctx->error_rel,h); 134295f7fbbSLois Curfman McInnes ctx->compute_err_iter = iter; 135d2dddef7SLois Curfman McInnes ctx->need_err = 0; 136d2dddef7SLois Curfman McInnes } 137e6ed2b1fSLois Curfman McInnes 138691598edSBarry Smith ierr = VecDotBegin(U,a,&dot);CHKERRQ(ierr); 139691598edSBarry Smith ierr = VecNormBegin(a,NORM_1,&sum);CHKERRQ(ierr); 140691598edSBarry Smith ierr = VecNormBegin(a,NORM_2,&norm);CHKERRQ(ierr); 141691598edSBarry Smith ierr = VecDotEnd(U,a,&dot);CHKERRQ(ierr); 142691598edSBarry Smith ierr = VecNormEnd(a,NORM_1,&sum);CHKERRQ(ierr); 143691598edSBarry Smith ierr = VecNormEnd(a,NORM_2,&norm);CHKERRQ(ierr); 144e6ed2b1fSLois Curfman McInnes 145d2dddef7SLois Curfman McInnes 146d2dddef7SLois Curfman McInnes /* Safeguard for step sizes too small */ 147d2dddef7SLois Curfman McInnes if (sum == 0.0) {dot = 1.0; norm = 1.0;} 148aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 149145595abSSatish Balay else if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum; 150145595abSSatish Balay else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum; 151d2dddef7SLois Curfman McInnes #else 152d2dddef7SLois Curfman McInnes else if (dot < umin*sum && dot >= 0.0) dot = umin*sum; 153d2dddef7SLois Curfman McInnes else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum; 154d2dddef7SLois Curfman McInnes #endif 155145595abSSatish Balay h = PetscRealPart(ctx->error_rel*dot/(norm*norm)); 156d2dddef7SLois Curfman McInnes } 157d2dddef7SLois Curfman McInnes } else { 158d2dddef7SLois Curfman McInnes h = ctx->h; 159d2dddef7SLois Curfman McInnes } 160b0a32e0cSBarry Smith if (!ctx->jorge || !ctx->need_h) PetscLogInfo(snes,"SNESMatrixFreeMult2_Private: h = %g\n",h); 161d2dddef7SLois Curfman McInnes 162d2dddef7SLois Curfman McInnes /* Evaluate function at F(u + ha) */ 163db0b4b35SLois Curfman McInnes hs = h; 164db0b4b35SLois Curfman McInnes ierr = VecWAXPY(&hs,a,U,w);CHKERRQ(ierr); 165d2dddef7SLois Curfman McInnes ierr = eval_fct(snes,w,y);CHKERRQ(ierr); 166d2dddef7SLois Curfman McInnes ierr = VecAXPY(&mone,F,y);CHKERRQ(ierr); 167db0b4b35SLois Curfman McInnes hs = 1.0/hs; 168db0b4b35SLois Curfman McInnes ierr = VecScale(&hs,y);CHKERRQ(ierr); 16974637425SBarry Smith if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);CHKERRQ(ierr);} 170d2dddef7SLois Curfman McInnes 171b0a32e0cSBarry Smith ierr = PetscLogEventEnd(MAT_MatrixFreeMult,a,y,0,0);CHKERRQ(ierr); 17207a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 173d2dddef7SLois Curfman McInnes } 174d2dddef7SLois Curfman McInnes 1754a2ae208SSatish Balay #undef __FUNCT__ 1764a2ae208SSatish Balay #define __FUNCT__ "SNESMatrixFreeMatCreate2" 177d2dddef7SLois Curfman McInnes /*@C 178d2dddef7SLois Curfman McInnes SNESMatrixFreeMatCreate2 - Creates a matrix-free matrix 179d2dddef7SLois Curfman McInnes context for use with a SNES solver. This matrix can be used as 180d2dddef7SLois Curfman McInnes the Jacobian argument for the routine SNESSetJacobian(). 181d2dddef7SLois Curfman McInnes 182d2dddef7SLois Curfman McInnes Input Parameters: 183d2dddef7SLois Curfman McInnes . snes - the SNES context 184d2dddef7SLois Curfman McInnes . x - vector where SNES solution is to be stored. 185d2dddef7SLois Curfman McInnes 186d2dddef7SLois Curfman McInnes Output Parameter: 187d2dddef7SLois Curfman McInnes . J - the matrix-free matrix 188d2dddef7SLois Curfman McInnes 189d2dddef7SLois Curfman McInnes Notes: 190d2dddef7SLois Curfman McInnes The matrix-free matrix context merely contains the function pointers 191d2dddef7SLois Curfman McInnes and work space for performing finite difference approximations of 192d2dddef7SLois Curfman McInnes Jacobian-vector products, J(u)*a, via 193d2dddef7SLois Curfman McInnes 194d2dddef7SLois Curfman McInnes $ J(u)*a = [J(u+h*a) - J(u)]/h, 195d2dddef7SLois Curfman McInnes $ where by default 196d2dddef7SLois Curfman McInnes $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 197d2dddef7SLois Curfman McInnes $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 198d2dddef7SLois Curfman McInnes $ where 199d2dddef7SLois Curfman McInnes $ error_rel = square root of relative error in 200d2dddef7SLois Curfman McInnes $ function evaluation 201d2dddef7SLois Curfman McInnes $ umin = minimum iterate parameter 202d2dddef7SLois Curfman McInnes $ Alternatively, the differencing parameter, h, can be set using 203e6ed2b1fSLois Curfman McInnes $ Jorge's nifty new strategy if one specifies the option 204e6ed2b1fSLois Curfman McInnes $ -snes_mf_jorge 205d2dddef7SLois Curfman McInnes 2065a655dc6SBarry Smith The user can set these parameters via MatSNESMFSetFunctionError(). 207d2dddef7SLois Curfman McInnes See the nonlinear solvers chapter of the users manual for details. 208d2dddef7SLois Curfman McInnes 209d2dddef7SLois Curfman McInnes The user should call MatDestroy() when finished with the matrix-free 210d2dddef7SLois Curfman McInnes matrix context. 211d2dddef7SLois Curfman McInnes 212d2dddef7SLois Curfman McInnes Options Database Keys: 213d2dddef7SLois Curfman McInnes $ -snes_mf_err <error_rel> 214d2dddef7SLois Curfman McInnes $ -snes_mf_unim <umin> 215295f7fbbSLois Curfman McInnes $ -snes_mf_compute_err 216295f7fbbSLois Curfman McInnes $ -snes_mf_freq_err <freq> 217e6ed2b1fSLois Curfman McInnes $ -snes_mf_jorge 218d2dddef7SLois Curfman McInnes 219d2dddef7SLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix 220d2dddef7SLois Curfman McInnes 2215a655dc6SBarry Smith .seealso: MatDestroy(), MatSNESMFSetFunctionError() 222d2dddef7SLois Curfman McInnes @*/ 223758f395bSBarry Smith int SNESDefaultMatrixFreeCreate2(SNES snes,Vec x,Mat *J) 224d2dddef7SLois Curfman McInnes { 225d2dddef7SLois Curfman McInnes MPI_Comm comm; 226d2dddef7SLois Curfman McInnes MFCtx_Private *mfctx; 227145595abSSatish Balay int n,nloc,ierr; 228145595abSSatish Balay PetscTruth flg; 229d2dddef7SLois Curfman McInnes char p[64]; 230d2dddef7SLois Curfman McInnes 23107a3eed9SLois Curfman McInnes PetscFunctionBegin; 23282502324SSatish Balay ierr = PetscNew(MFCtx_Private,&mfctx);CHKERRQ(ierr); 233549d3d68SSatish Balay ierr = PetscMemzero(mfctx,sizeof(MFCtx_Private));CHKERRQ(ierr); 234b0a32e0cSBarry Smith PetscLogObjectMemory(snes,sizeof(MFCtx_Private)); 235d2dddef7SLois Curfman McInnes mfctx->sp = 0; 236d2dddef7SLois Curfman McInnes mfctx->snes = snes; 237145595abSSatish Balay mfctx->error_rel = 1.e-8; /* assumes PetscReal precision */ 238d2dddef7SLois Curfman McInnes mfctx->umin = 1.e-6; 239d2dddef7SLois Curfman McInnes mfctx->h = 0.0; 240d2dddef7SLois Curfman McInnes mfctx->need_h = 1; 241d2dddef7SLois Curfman McInnes mfctx->need_err = 0; 242295f7fbbSLois Curfman McInnes mfctx->compute_err = 0; 243295f7fbbSLois Curfman McInnes mfctx->compute_err_freq = 0; 244295f7fbbSLois Curfman McInnes mfctx->compute_err_iter = -1; 24587828ca2SBarry Smith ierr = PetscOptionsGetReal(snes->prefix,"-snes_mf_err",&mfctx->error_rel,PETSC_NULL);CHKERRQ(ierr); 24687828ca2SBarry Smith ierr = PetscOptionsGetReal(snes->prefix,"-snes_mf_umin",&mfctx->umin,PETSC_NULL);CHKERRQ(ierr); 247b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_jorge",&mfctx->jorge);CHKERRQ(ierr); 248b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_compute_err",&mfctx->compute_err);CHKERRQ(ierr); 249b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(snes->prefix,"-snes_mf_freq_err",&mfctx->compute_err_freq,&flg);CHKERRQ(ierr); 250295f7fbbSLois Curfman McInnes if (flg) { 251295f7fbbSLois Curfman McInnes if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0; 252295f7fbbSLois Curfman McInnes mfctx->compute_err = 1; 253295f7fbbSLois Curfman McInnes } 254295f7fbbSLois Curfman McInnes if (mfctx->compute_err == 1) mfctx->need_err = 1; 255295f7fbbSLois Curfman McInnes if (mfctx->jorge || mfctx->compute_err) { 256d2dddef7SLois Curfman McInnes ierr = DiffParameterCreate_More(snes,x,&mfctx->data);CHKERRQ(ierr); 257d2dddef7SLois Curfman McInnes } else mfctx->data = 0; 258d2dddef7SLois Curfman McInnes 259b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr); 260549d3d68SSatish Balay ierr = PetscStrcpy(p,"-");CHKERRQ(ierr); 261d2dddef7SLois Curfman McInnes if (snes->prefix) PetscStrcat(p,snes->prefix); 262d2dddef7SLois Curfman McInnes if (flg) { 263d132466eSBarry Smith ierr = PetscPrintf(snes->comm," Matrix-free Options (via SNES):\n");CHKERRQ(ierr); 264d132466eSBarry Smith ierr = PetscPrintf(snes->comm," %ssnes_mf_err <err>: set sqrt of relative error in function (default %g)\n",p,mfctx->error_rel);CHKERRQ(ierr); 265d132466eSBarry Smith ierr = PetscPrintf(snes->comm," %ssnes_mf_umin <umin>: see users manual (default %g)\n",p,mfctx->umin);CHKERRQ(ierr); 266d132466eSBarry Smith ierr = PetscPrintf(snes->comm," %ssnes_mf_jorge: use Jorge More's method\n",p);CHKERRQ(ierr); 267d132466eSBarry Smith ierr = PetscPrintf(snes->comm," %ssnes_mf_compute_err: compute sqrt or relative error in function\n",p);CHKERRQ(ierr); 268d132466eSBarry Smith ierr = PetscPrintf(snes->comm," %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)\n",p);CHKERRQ(ierr); 269d132466eSBarry Smith ierr = PetscPrintf(snes->comm," %ssnes_mf_noise_file <file>: set file for printing noise info\n",p);CHKERRQ(ierr); 270d2dddef7SLois Curfman McInnes } 271d2dddef7SLois Curfman McInnes ierr = VecDuplicate(x,&mfctx->w);CHKERRQ(ierr); 272d2dddef7SLois Curfman McInnes ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr); 273d2dddef7SLois Curfman McInnes ierr = VecGetSize(x,&n);CHKERRQ(ierr); 274d2dddef7SLois Curfman McInnes ierr = VecGetLocalSize(x,&nloc);CHKERRQ(ierr); 275d2dddef7SLois Curfman McInnes ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J);CHKERRQ(ierr); 27637bd1cefSSatish Balay ierr = MatShellSetOperation(*J,MATOP_MULT,(void(*)())SNESMatrixFreeMult2_Private);CHKERRQ(ierr); 27737bd1cefSSatish Balay ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void(*)())SNESMatrixFreeDestroy2_Private);CHKERRQ(ierr); 27837bd1cefSSatish Balay ierr = MatShellSetOperation(*J,MATOP_VIEW,(void(*)())SNESMatrixFreeView2_Private);CHKERRQ(ierr); 279d2dddef7SLois Curfman McInnes 280b0a32e0cSBarry Smith PetscLogObjectParent(*J,mfctx->w); 281b0a32e0cSBarry Smith PetscLogObjectParent(snes,*J); 28207a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 283d2dddef7SLois Curfman McInnes } 284d2dddef7SLois Curfman McInnes 2854a2ae208SSatish Balay #undef __FUNCT__ 2864a2ae208SSatish Balay #define __FUNCT__ "SNESDefaultMatrixFreeSetParameters2" 287d2dddef7SLois Curfman McInnes /*@ 288758f395bSBarry Smith SNESDefaultMatrixFreeSetParameters2 - Sets the parameters for the approximation of 289d2dddef7SLois Curfman McInnes matrix-vector products using finite differences. 290d2dddef7SLois Curfman McInnes 291d2dddef7SLois Curfman McInnes $ J(u)*a = [J(u+h*a) - J(u)]/h where 292d2dddef7SLois Curfman McInnes 293d2dddef7SLois Curfman McInnes either the user sets h directly here, or this parameter is computed via 294d2dddef7SLois Curfman McInnes 295d2dddef7SLois Curfman McInnes $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 296d2dddef7SLois Curfman McInnes $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 297d2dddef7SLois Curfman McInnes $ 298d2dddef7SLois Curfman McInnes 299d2dddef7SLois Curfman McInnes Input Parameters: 300758f395bSBarry Smith + mat - the matrix 301d2dddef7SLois Curfman McInnes . error_rel - relative error (should be set to the square root of 302d2dddef7SLois Curfman McInnes the relative error in the function evaluations) 303d2dddef7SLois Curfman McInnes . umin - minimum allowable u-value 304758f395bSBarry Smith - h - differencing parameter 305d2dddef7SLois Curfman McInnes 306d2dddef7SLois Curfman McInnes Notes: 307d2dddef7SLois Curfman McInnes If the user sets the parameter h directly, then this value will be used 308d2dddef7SLois Curfman McInnes instead of the default computation indicated above. 309d2dddef7SLois Curfman McInnes 310d2dddef7SLois Curfman McInnes .keywords: SNES, matrix-free, parameters 311d2dddef7SLois Curfman McInnes 3125a655dc6SBarry Smith .seealso: MatCreateSNESMF() 313d2dddef7SLois Curfman McInnes @*/ 314145595abSSatish Balay int SNESDefaultMatrixFreeSetParameters2(Mat mat,PetscReal error,PetscReal umin,PetscReal h) 315d2dddef7SLois Curfman McInnes { 316d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 317d2dddef7SLois Curfman McInnes int ierr; 318d2dddef7SLois Curfman McInnes 31907a3eed9SLois Curfman McInnes PetscFunctionBegin; 320d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 321d2dddef7SLois Curfman McInnes if (ctx) { 322d2dddef7SLois Curfman McInnes if (error != PETSC_DEFAULT) ctx->error_rel = error; 323d2dddef7SLois Curfman McInnes if (umin != PETSC_DEFAULT) ctx->umin = umin; 324d2dddef7SLois Curfman McInnes if (h != PETSC_DEFAULT) { 325d2dddef7SLois Curfman McInnes ctx->h = h; 326d2dddef7SLois Curfman McInnes ctx->need_h = 0; 327d2dddef7SLois Curfman McInnes } 328d2dddef7SLois Curfman McInnes } 32907a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 330d2dddef7SLois Curfman McInnes } 331d2dddef7SLois Curfman McInnes 332d2dddef7SLois Curfman McInnes int SNESUnSetMatrixFreeParameter(SNES snes) 333d2dddef7SLois Curfman McInnes { 334d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 335d2dddef7SLois Curfman McInnes int ierr; 336d2dddef7SLois Curfman McInnes Mat mat; 337d2dddef7SLois Curfman McInnes 33807a3eed9SLois Curfman McInnes PetscFunctionBegin; 33974637425SBarry Smith ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 340d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 341d2dddef7SLois Curfman McInnes if (ctx) ctx->need_h = 1; 34207a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 343d2dddef7SLois Curfman McInnes } 344d2dddef7SLois Curfman McInnes 345