1d2dddef7SLois Curfman McInnes 20a835dfdSSatish Balay #include "src/snes/snesimpl.h" /*I "petscsnes.h" I*/ 3d2dddef7SLois Curfman McInnes 4ca44d042SBarry Smith EXTERN int DiffParameterCreate_More(SNES,Vec,void**); 5ca44d042SBarry Smith EXTERN int DiffParameterCompute_More(SNES,void*,Vec,Vec,PetscReal*,PetscReal*); 6ca44d042SBarry Smith EXTERN int DiffParameterDestroy_More(void*); 7d2dddef7SLois Curfman McInnes 8d2dddef7SLois Curfman McInnes typedef struct { /* default context for matrix-free SNES */ 9d2dddef7SLois Curfman McInnes SNES snes; /* SNES context */ 10d2dddef7SLois Curfman McInnes Vec w; /* work vector */ 1174637425SBarry Smith MatNullSpace sp; /* null space context */ 12145595abSSatish Balay PetscReal error_rel; /* square root of relative error in computing function */ 13145595abSSatish Balay PetscReal umin; /* minimum allowable u'a value relative to |u|_1 */ 14d2dddef7SLois Curfman McInnes int jorge; /* flag indicating use of Jorge's method for determining 15d2dddef7SLois Curfman McInnes the differencing parameter */ 16145595abSSatish Balay PetscReal h; /* differencing parameter */ 17d2dddef7SLois Curfman McInnes int need_h; /* flag indicating whether we must compute h */ 18295f7fbbSLois Curfman McInnes int need_err; /* flag indicating whether we must currently compute error_rel */ 19295f7fbbSLois Curfman McInnes int compute_err; /* flag indicating whether we must ever compute error_rel */ 20295f7fbbSLois Curfman McInnes int compute_err_iter; /* last iter where we've computer error_rel */ 21295f7fbbSLois Curfman McInnes int compute_err_freq; /* frequency of computing error_rel */ 22d2dddef7SLois Curfman McInnes void *data; /* implementation-specific data */ 23d2dddef7SLois Curfman McInnes } MFCtx_Private; 24d2dddef7SLois Curfman McInnes 254a2ae208SSatish Balay #undef __FUNCT__ 264a2ae208SSatish Balay #define __FUNCT__ "SNESMatrixFreeDestroy2_Private" /* ADIC Ignore */ 27f1d4f8a3SBarry Smith int SNESMatrixFreeDestroy2_Private(Mat mat) 28d2dddef7SLois Curfman McInnes { 29d2dddef7SLois Curfman McInnes int ierr; 30d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 31d2dddef7SLois Curfman McInnes 3207a3eed9SLois Curfman McInnes PetscFunctionBegin; 33d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx); 34d2dddef7SLois Curfman McInnes ierr = VecDestroy(ctx->w);CHKERRQ(ierr); 3574637425SBarry Smith if (ctx->sp) {ierr = MatNullSpaceDestroy(ctx->sp);CHKERRQ(ierr);} 36295f7fbbSLois Curfman McInnes if (ctx->jorge || ctx->compute_err) {ierr = DiffParameterDestroy_More(ctx->data);CHKERRQ(ierr);} 37606d414cSSatish Balay ierr = PetscFree(ctx);CHKERRQ(ierr); 3807a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 39d2dddef7SLois Curfman McInnes } 40d2dddef7SLois Curfman McInnes 414a2ae208SSatish Balay #undef __FUNCT__ 424a2ae208SSatish Balay #define __FUNCT__ "SNESMatrixFreeView2_Private" /* ADIC Ignore */ 43d2dddef7SLois Curfman McInnes /* 44d2dddef7SLois Curfman McInnes SNESMatrixFreeView2_Private - Views matrix-free parameters. 45d2dddef7SLois Curfman McInnes */ 46b0a32e0cSBarry Smith int SNESMatrixFreeView2_Private(Mat J,PetscViewer viewer) 47d2dddef7SLois Curfman McInnes { 48d2dddef7SLois Curfman McInnes int ierr; 49d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 50*32077d6dSBarry Smith PetscTruth iascii; 51d2dddef7SLois Curfman McInnes 5207a3eed9SLois Curfman McInnes PetscFunctionBegin; 53d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(J,(void **)&ctx);CHKERRQ(ierr); 54*32077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 55*32077d6dSBarry Smith if (iascii) { 56b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," SNES matrix-free approximation:\n");CHKERRQ(ierr); 57758f395bSBarry Smith if (ctx->jorge) { 58b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," using Jorge's method of determining differencing parameter\n");CHKERRQ(ierr); 59758f395bSBarry Smith } 60b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," err=%g (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr); 61b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," umin=%g (minimum iterate parameter)\n",ctx->umin);CHKERRQ(ierr); 62758f395bSBarry Smith if (ctx->compute_err) { 63b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," freq_err=%d (frequency for computing err)\n",ctx->compute_err_freq);CHKERRQ(ierr); 64d2dddef7SLois Curfman McInnes } 65758f395bSBarry Smith } else { 6629bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by SNES matrix free Jorge",((PetscObject)viewer)->type_name); 67758f395bSBarry Smith } 6807a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 69d2dddef7SLois Curfman McInnes } 70d2dddef7SLois Curfman McInnes 714a2ae208SSatish Balay #undef __FUNCT__ 724a2ae208SSatish Balay #define __FUNCT__ "SNESMatrixFreeMult2_Private" 73d2dddef7SLois Curfman McInnes /* 74d2dddef7SLois Curfman McInnes SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector 75d2dddef7SLois Curfman McInnes product, y = F'(u)*a: 76d2dddef7SLois Curfman McInnes y = (F(u + ha) - F(u)) /h, 77d2dddef7SLois Curfman McInnes where F = nonlinear function, as set by SNESSetFunction() 78d2dddef7SLois Curfman McInnes u = current iterate 79d2dddef7SLois Curfman McInnes h = difference interval 80d2dddef7SLois Curfman McInnes */ 81d2dddef7SLois Curfman McInnes int SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y) 82d2dddef7SLois Curfman McInnes { 83d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 84d2dddef7SLois Curfman McInnes SNES snes; 85145595abSSatish Balay PetscReal h,norm,sum,umin,noise; 86ea709b57SSatish Balay PetscScalar hs,dot,mone = -1.0; 87d2dddef7SLois Curfman McInnes Vec w,U,F; 88295f7fbbSLois Curfman McInnes int ierr,iter,(*eval_fct)(SNES,Vec,Vec); 89e6ed2b1fSLois Curfman McInnes MPI_Comm comm; 90d2dddef7SLois Curfman McInnes 9107a3eed9SLois Curfman McInnes PetscFunctionBegin; 9207a3eed9SLois Curfman McInnes 93d2dddef7SLois Curfman McInnes /* We log matrix-free matrix-vector products separately, so that we can 94d2dddef7SLois Curfman McInnes separate the performance monitoring from the cases that use conventional 95d2dddef7SLois Curfman McInnes storage. We may eventually modify event logging to associate events 96d2dddef7SLois Curfman McInnes with particular objects, hence alleviating the more general problem. */ 9752ed857cSBarry Smith ierr = PetscLogEventBegin(MAT_MultMatrixFree,a,y,0,0);CHKERRQ(ierr); 98d2dddef7SLois Curfman McInnes 992d0c0e3bSBarry Smith ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 100e6ed2b1fSLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 101e6ed2b1fSLois Curfman McInnes snes = ctx->snes; 102e6ed2b1fSLois Curfman McInnes w = ctx->w; 103e6ed2b1fSLois Curfman McInnes umin = ctx->umin; 104e6ed2b1fSLois Curfman McInnes 105d2dddef7SLois Curfman McInnes ierr = SNESGetSolution(snes,&U);CHKERRQ(ierr); 106d2dddef7SLois Curfman McInnes eval_fct = SNESComputeFunction; 10774637425SBarry Smith ierr = SNESGetFunction(snes,&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 108295f7fbbSLois Curfman McInnes 109d2dddef7SLois Curfman McInnes /* Determine a "good" step size, h */ 110d2dddef7SLois Curfman McInnes if (ctx->need_h) { 111d2dddef7SLois Curfman McInnes 112d2dddef7SLois Curfman McInnes /* Use Jorge's method to compute h */ 113d2dddef7SLois Curfman McInnes if (ctx->jorge) { 114d2dddef7SLois Curfman McInnes ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr); 115d2dddef7SLois Curfman McInnes 116d2dddef7SLois Curfman McInnes /* Use the Brown/Saad method to compute h */ 117d2dddef7SLois Curfman McInnes } else { 118295f7fbbSLois Curfman McInnes /* Compute error if desired */ 119295f7fbbSLois Curfman McInnes ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr); 120145595abSSatish Balay if ((ctx->need_err) || ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter-1)%ctx->compute_err_freq)))) { 121d2dddef7SLois Curfman McInnes /* Use Jorge's method to compute noise */ 122d2dddef7SLois Curfman McInnes ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr); 123d2dddef7SLois Curfman McInnes ctx->error_rel = sqrt(noise); 124b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESMatrixFreeMult2_Private: Using Jorge's noise: noise=%g, sqrt(noise)=%g, h_more=%g\n",noise,ctx->error_rel,h); 125295f7fbbSLois Curfman McInnes ctx->compute_err_iter = iter; 126d2dddef7SLois Curfman McInnes ctx->need_err = 0; 127d2dddef7SLois Curfman McInnes } 128e6ed2b1fSLois Curfman McInnes 129691598edSBarry Smith ierr = VecDotBegin(U,a,&dot);CHKERRQ(ierr); 130691598edSBarry Smith ierr = VecNormBegin(a,NORM_1,&sum);CHKERRQ(ierr); 131691598edSBarry Smith ierr = VecNormBegin(a,NORM_2,&norm);CHKERRQ(ierr); 132691598edSBarry Smith ierr = VecDotEnd(U,a,&dot);CHKERRQ(ierr); 133691598edSBarry Smith ierr = VecNormEnd(a,NORM_1,&sum);CHKERRQ(ierr); 134691598edSBarry Smith ierr = VecNormEnd(a,NORM_2,&norm);CHKERRQ(ierr); 135e6ed2b1fSLois Curfman McInnes 136d2dddef7SLois Curfman McInnes 137d2dddef7SLois Curfman McInnes /* Safeguard for step sizes too small */ 138d2dddef7SLois Curfman McInnes if (sum == 0.0) {dot = 1.0; norm = 1.0;} 139aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 140145595abSSatish Balay else if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum; 141145595abSSatish Balay else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum; 142d2dddef7SLois Curfman McInnes #else 143d2dddef7SLois Curfman McInnes else if (dot < umin*sum && dot >= 0.0) dot = umin*sum; 144d2dddef7SLois Curfman McInnes else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum; 145d2dddef7SLois Curfman McInnes #endif 146145595abSSatish Balay h = PetscRealPart(ctx->error_rel*dot/(norm*norm)); 147d2dddef7SLois Curfman McInnes } 148d2dddef7SLois Curfman McInnes } else { 149d2dddef7SLois Curfman McInnes h = ctx->h; 150d2dddef7SLois Curfman McInnes } 151b0a32e0cSBarry Smith if (!ctx->jorge || !ctx->need_h) PetscLogInfo(snes,"SNESMatrixFreeMult2_Private: h = %g\n",h); 152d2dddef7SLois Curfman McInnes 153d2dddef7SLois Curfman McInnes /* Evaluate function at F(u + ha) */ 154db0b4b35SLois Curfman McInnes hs = h; 155db0b4b35SLois Curfman McInnes ierr = VecWAXPY(&hs,a,U,w);CHKERRQ(ierr); 156d2dddef7SLois Curfman McInnes ierr = eval_fct(snes,w,y);CHKERRQ(ierr); 157d2dddef7SLois Curfman McInnes ierr = VecAXPY(&mone,F,y);CHKERRQ(ierr); 158db0b4b35SLois Curfman McInnes hs = 1.0/hs; 159db0b4b35SLois Curfman McInnes ierr = VecScale(&hs,y);CHKERRQ(ierr); 16074637425SBarry Smith if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);CHKERRQ(ierr);} 161d2dddef7SLois Curfman McInnes 16252ed857cSBarry Smith ierr = PetscLogEventEnd(MAT_MultMatrixFree,a,y,0,0);CHKERRQ(ierr); 16307a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 164d2dddef7SLois Curfman McInnes } 165d2dddef7SLois Curfman McInnes 1664a2ae208SSatish Balay #undef __FUNCT__ 1674a2ae208SSatish Balay #define __FUNCT__ "SNESMatrixFreeMatCreate2" 168d2dddef7SLois Curfman McInnes /*@C 169d2dddef7SLois Curfman McInnes SNESMatrixFreeMatCreate2 - Creates a matrix-free matrix 170d2dddef7SLois Curfman McInnes context for use with a SNES solver. This matrix can be used as 171d2dddef7SLois Curfman McInnes the Jacobian argument for the routine SNESSetJacobian(). 172d2dddef7SLois Curfman McInnes 173d2dddef7SLois Curfman McInnes Input Parameters: 174d2dddef7SLois Curfman McInnes . snes - the SNES context 175d2dddef7SLois Curfman McInnes . x - vector where SNES solution is to be stored. 176d2dddef7SLois Curfman McInnes 177d2dddef7SLois Curfman McInnes Output Parameter: 178d2dddef7SLois Curfman McInnes . J - the matrix-free matrix 179d2dddef7SLois Curfman McInnes 18090c26489SBarry Smith Level: advanced 18190c26489SBarry Smith 182d2dddef7SLois Curfman McInnes Notes: 183d2dddef7SLois Curfman McInnes The matrix-free matrix context merely contains the function pointers 184d2dddef7SLois Curfman McInnes and work space for performing finite difference approximations of 185d2dddef7SLois Curfman McInnes Jacobian-vector products, J(u)*a, via 186d2dddef7SLois Curfman McInnes 187d2dddef7SLois Curfman McInnes $ J(u)*a = [J(u+h*a) - J(u)]/h, 188d2dddef7SLois Curfman McInnes $ where by default 189d2dddef7SLois Curfman McInnes $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 190d2dddef7SLois Curfman McInnes $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 191d2dddef7SLois Curfman McInnes $ where 192d2dddef7SLois Curfman McInnes $ error_rel = square root of relative error in 193d2dddef7SLois Curfman McInnes $ function evaluation 194d2dddef7SLois Curfman McInnes $ umin = minimum iterate parameter 195d2dddef7SLois Curfman McInnes $ Alternatively, the differencing parameter, h, can be set using 196e6ed2b1fSLois Curfman McInnes $ Jorge's nifty new strategy if one specifies the option 197e6ed2b1fSLois Curfman McInnes $ -snes_mf_jorge 198d2dddef7SLois Curfman McInnes 1995a655dc6SBarry Smith The user can set these parameters via MatSNESMFSetFunctionError(). 200d2dddef7SLois Curfman McInnes See the nonlinear solvers chapter of the users manual for details. 201d2dddef7SLois Curfman McInnes 202d2dddef7SLois Curfman McInnes The user should call MatDestroy() when finished with the matrix-free 203d2dddef7SLois Curfman McInnes matrix context. 204d2dddef7SLois Curfman McInnes 205d2dddef7SLois Curfman McInnes Options Database Keys: 206d2dddef7SLois Curfman McInnes $ -snes_mf_err <error_rel> 207d2dddef7SLois Curfman McInnes $ -snes_mf_unim <umin> 208295f7fbbSLois Curfman McInnes $ -snes_mf_compute_err 209295f7fbbSLois Curfman McInnes $ -snes_mf_freq_err <freq> 210e6ed2b1fSLois Curfman McInnes $ -snes_mf_jorge 211d2dddef7SLois Curfman McInnes 212d2dddef7SLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix 213d2dddef7SLois Curfman McInnes 2145a655dc6SBarry Smith .seealso: MatDestroy(), MatSNESMFSetFunctionError() 215d2dddef7SLois Curfman McInnes @*/ 216758f395bSBarry Smith int SNESDefaultMatrixFreeCreate2(SNES snes,Vec x,Mat *J) 217d2dddef7SLois Curfman McInnes { 218d2dddef7SLois Curfman McInnes MPI_Comm comm; 219d2dddef7SLois Curfman McInnes MFCtx_Private *mfctx; 220145595abSSatish Balay int n,nloc,ierr; 221145595abSSatish Balay PetscTruth flg; 222d2dddef7SLois Curfman McInnes char p[64]; 223d2dddef7SLois Curfman McInnes 22407a3eed9SLois Curfman McInnes PetscFunctionBegin; 22582502324SSatish Balay ierr = PetscNew(MFCtx_Private,&mfctx);CHKERRQ(ierr); 226549d3d68SSatish Balay ierr = PetscMemzero(mfctx,sizeof(MFCtx_Private));CHKERRQ(ierr); 227b0a32e0cSBarry Smith PetscLogObjectMemory(snes,sizeof(MFCtx_Private)); 228d2dddef7SLois Curfman McInnes mfctx->sp = 0; 229d2dddef7SLois Curfman McInnes mfctx->snes = snes; 23077d8c4bbSBarry Smith mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 231d2dddef7SLois Curfman McInnes mfctx->umin = 1.e-6; 232d2dddef7SLois Curfman McInnes mfctx->h = 0.0; 233d2dddef7SLois Curfman McInnes mfctx->need_h = 1; 234d2dddef7SLois Curfman McInnes mfctx->need_err = 0; 235295f7fbbSLois Curfman McInnes mfctx->compute_err = 0; 236295f7fbbSLois Curfman McInnes mfctx->compute_err_freq = 0; 237295f7fbbSLois Curfman McInnes mfctx->compute_err_iter = -1; 23887828ca2SBarry Smith ierr = PetscOptionsGetReal(snes->prefix,"-snes_mf_err",&mfctx->error_rel,PETSC_NULL);CHKERRQ(ierr); 23987828ca2SBarry Smith ierr = PetscOptionsGetReal(snes->prefix,"-snes_mf_umin",&mfctx->umin,PETSC_NULL);CHKERRQ(ierr); 24052ed857cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_jorge",(PetscTruth*)&mfctx->jorge);CHKERRQ(ierr); 24152ed857cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_compute_err",(PetscTruth*)&mfctx->compute_err);CHKERRQ(ierr); 242b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(snes->prefix,"-snes_mf_freq_err",&mfctx->compute_err_freq,&flg);CHKERRQ(ierr); 243295f7fbbSLois Curfman McInnes if (flg) { 244295f7fbbSLois Curfman McInnes if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0; 245295f7fbbSLois Curfman McInnes mfctx->compute_err = 1; 246295f7fbbSLois Curfman McInnes } 247295f7fbbSLois Curfman McInnes if (mfctx->compute_err == 1) mfctx->need_err = 1; 248295f7fbbSLois Curfman McInnes if (mfctx->jorge || mfctx->compute_err) { 249d2dddef7SLois Curfman McInnes ierr = DiffParameterCreate_More(snes,x,&mfctx->data);CHKERRQ(ierr); 250d2dddef7SLois Curfman McInnes } else mfctx->data = 0; 251d2dddef7SLois Curfman McInnes 252b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr); 253549d3d68SSatish Balay ierr = PetscStrcpy(p,"-");CHKERRQ(ierr); 254d2dddef7SLois Curfman McInnes if (snes->prefix) PetscStrcat(p,snes->prefix); 255d2dddef7SLois Curfman McInnes if (flg) { 256d132466eSBarry Smith ierr = PetscPrintf(snes->comm," Matrix-free Options (via SNES):\n");CHKERRQ(ierr); 257d132466eSBarry 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); 258d132466eSBarry Smith ierr = PetscPrintf(snes->comm," %ssnes_mf_umin <umin>: see users manual (default %g)\n",p,mfctx->umin);CHKERRQ(ierr); 259d132466eSBarry Smith ierr = PetscPrintf(snes->comm," %ssnes_mf_jorge: use Jorge More's method\n",p);CHKERRQ(ierr); 260d132466eSBarry Smith ierr = PetscPrintf(snes->comm," %ssnes_mf_compute_err: compute sqrt or relative error in function\n",p);CHKERRQ(ierr); 261d132466eSBarry Smith ierr = PetscPrintf(snes->comm," %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)\n",p);CHKERRQ(ierr); 262d132466eSBarry Smith ierr = PetscPrintf(snes->comm," %ssnes_mf_noise_file <file>: set file for printing noise info\n",p);CHKERRQ(ierr); 263d2dddef7SLois Curfman McInnes } 264d2dddef7SLois Curfman McInnes ierr = VecDuplicate(x,&mfctx->w);CHKERRQ(ierr); 265d2dddef7SLois Curfman McInnes ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr); 266d2dddef7SLois Curfman McInnes ierr = VecGetSize(x,&n);CHKERRQ(ierr); 267d2dddef7SLois Curfman McInnes ierr = VecGetLocalSize(x,&nloc);CHKERRQ(ierr); 268f204ca49SKris Buschelman ierr = MatCreate(comm,nloc,n,n,n,J);CHKERRQ(ierr); 269f204ca49SKris Buschelman ierr = MatSetType(*J,MATSHELL);CHKERRQ(ierr); 270f204ca49SKris Buschelman ierr = MatShellSetContext(*J,mfctx);CHKERRQ(ierr); 271c134de8dSSatish Balay ierr = MatShellSetOperation(*J,MATOP_MULT,(void(*)(void))SNESMatrixFreeMult2_Private);CHKERRQ(ierr); 272c134de8dSSatish Balay ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void(*)(void))SNESMatrixFreeDestroy2_Private);CHKERRQ(ierr); 273c134de8dSSatish Balay ierr = MatShellSetOperation(*J,MATOP_VIEW,(void(*)(void))SNESMatrixFreeView2_Private);CHKERRQ(ierr); 274d2dddef7SLois Curfman McInnes 275b0a32e0cSBarry Smith PetscLogObjectParent(*J,mfctx->w); 276b0a32e0cSBarry Smith PetscLogObjectParent(snes,*J); 27707a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 278d2dddef7SLois Curfman McInnes } 279d2dddef7SLois Curfman McInnes 2804a2ae208SSatish Balay #undef __FUNCT__ 2814a2ae208SSatish Balay #define __FUNCT__ "SNESDefaultMatrixFreeSetParameters2" 28290c26489SBarry Smith /*@C 283758f395bSBarry Smith SNESDefaultMatrixFreeSetParameters2 - Sets the parameters for the approximation of 284d2dddef7SLois Curfman McInnes matrix-vector products using finite differences. 285d2dddef7SLois Curfman McInnes 286d2dddef7SLois Curfman McInnes $ J(u)*a = [J(u+h*a) - J(u)]/h where 287d2dddef7SLois Curfman McInnes 288d2dddef7SLois Curfman McInnes either the user sets h directly here, or this parameter is computed via 289d2dddef7SLois Curfman McInnes 290d2dddef7SLois Curfman McInnes $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 291d2dddef7SLois Curfman McInnes $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 292d2dddef7SLois Curfman McInnes $ 293d2dddef7SLois Curfman McInnes 294d2dddef7SLois Curfman McInnes Input Parameters: 295758f395bSBarry Smith + mat - the matrix 296d2dddef7SLois Curfman McInnes . error_rel - relative error (should be set to the square root of 297d2dddef7SLois Curfman McInnes the relative error in the function evaluations) 298d2dddef7SLois Curfman McInnes . umin - minimum allowable u-value 299758f395bSBarry Smith - h - differencing parameter 300d2dddef7SLois Curfman McInnes 30190c26489SBarry Smith Level: advanced 30290c26489SBarry Smith 303d2dddef7SLois Curfman McInnes Notes: 304d2dddef7SLois Curfman McInnes If the user sets the parameter h directly, then this value will be used 305d2dddef7SLois Curfman McInnes instead of the default computation indicated above. 306d2dddef7SLois Curfman McInnes 307d2dddef7SLois Curfman McInnes .keywords: SNES, matrix-free, parameters 308d2dddef7SLois Curfman McInnes 3095a655dc6SBarry Smith .seealso: MatCreateSNESMF() 310d2dddef7SLois Curfman McInnes @*/ 311145595abSSatish Balay int SNESDefaultMatrixFreeSetParameters2(Mat mat,PetscReal error,PetscReal umin,PetscReal h) 312d2dddef7SLois Curfman McInnes { 313d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 314d2dddef7SLois Curfman McInnes int ierr; 315d2dddef7SLois Curfman McInnes 31607a3eed9SLois Curfman McInnes PetscFunctionBegin; 317d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 318d2dddef7SLois Curfman McInnes if (ctx) { 319d2dddef7SLois Curfman McInnes if (error != PETSC_DEFAULT) ctx->error_rel = error; 320d2dddef7SLois Curfman McInnes if (umin != PETSC_DEFAULT) ctx->umin = umin; 321d2dddef7SLois Curfman McInnes if (h != PETSC_DEFAULT) { 322d2dddef7SLois Curfman McInnes ctx->h = h; 323d2dddef7SLois Curfman McInnes ctx->need_h = 0; 324d2dddef7SLois Curfman McInnes } 325d2dddef7SLois Curfman McInnes } 32607a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 327d2dddef7SLois Curfman McInnes } 328d2dddef7SLois Curfman McInnes 329d2dddef7SLois Curfman McInnes int SNESUnSetMatrixFreeParameter(SNES snes) 330d2dddef7SLois Curfman McInnes { 331d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 332d2dddef7SLois Curfman McInnes int ierr; 333d2dddef7SLois Curfman McInnes Mat mat; 334d2dddef7SLois Curfman McInnes 33507a3eed9SLois Curfman McInnes PetscFunctionBegin; 33674637425SBarry Smith ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 337d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 338d2dddef7SLois Curfman McInnes if (ctx) ctx->need_h = 1; 33907a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 340d2dddef7SLois Curfman McInnes } 341d2dddef7SLois Curfman McInnes 342