1d2dddef7SLois Curfman McInnes 2b45d2f2cSJed Brown #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 36370053bSBarry Smith /* matimpl.h is needed only for logging of matrix operation */ 4b45d2f2cSJed Brown #include <petsc-private/matimpl.h> 5d2dddef7SLois Curfman McInnes 68b5db460SBarry Smith extern PetscErrorCode SNESDiffParameterCreate_More(SNES,Vec,void**); 78b5db460SBarry Smith extern PetscErrorCode SNESDiffParameterCompute_More(SNES,void*,Vec,Vec,PetscReal*,PetscReal*); 88b5db460SBarry Smith extern PetscErrorCode SNESDiffParameterDestroy_More(void*); 9d2dddef7SLois Curfman McInnes 10d2dddef7SLois Curfman McInnes typedef struct { /* default context for matrix-free SNES */ 11d2dddef7SLois Curfman McInnes SNES snes; /* SNES context */ 12d2dddef7SLois Curfman McInnes Vec w; /* work vector */ 1374637425SBarry Smith MatNullSpace sp; /* null space context */ 14145595abSSatish Balay PetscReal error_rel; /* square root of relative error in computing function */ 15145595abSSatish Balay PetscReal umin; /* minimum allowable u'a value relative to |u|_1 */ 16f5af7f23SKarl Rupp PetscBool jorge; /* flag indicating use of Jorge's method for determining the differencing parameter */ 17145595abSSatish Balay PetscReal h; /* differencing parameter */ 18ace3abfcSBarry Smith PetscBool need_h; /* flag indicating whether we must compute h */ 19ace3abfcSBarry Smith PetscBool need_err; /* flag indicating whether we must currently compute error_rel */ 20ace3abfcSBarry Smith PetscBool compute_err; /* flag indicating whether we must ever compute error_rel */ 21a7cc72afSBarry Smith PetscInt compute_err_iter; /* last iter where we've computer error_rel */ 22a7cc72afSBarry Smith PetscInt compute_err_freq; /* frequency of computing error_rel */ 23d2dddef7SLois Curfman McInnes void *data; /* implementation-specific data */ 24d2dddef7SLois Curfman McInnes } MFCtx_Private; 25d2dddef7SLois Curfman McInnes 264a2ae208SSatish Balay #undef __FUNCT__ 27bcaeba4dSBarry Smith #define __FUNCT__ "SNESMatrixFreeDestroy2_Private" 28dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeDestroy2_Private(Mat mat) 29d2dddef7SLois Curfman McInnes { 30dfbe8321SBarry Smith PetscErrorCode ierr; 31d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 32d2dddef7SLois Curfman McInnes 3307a3eed9SLois Curfman McInnes PetscFunctionBegin; 3422d28d08SBarry Smith ierr = MatShellGetContext(mat,(void**)&ctx);CHKERRQ(ierr); 356bf464f9SBarry Smith ierr = VecDestroy(&ctx->w);CHKERRQ(ierr); 366bf464f9SBarry Smith ierr = MatNullSpaceDestroy(&ctx->sp);CHKERRQ(ierr); 378b5db460SBarry Smith if (ctx->jorge || ctx->compute_err) {ierr = SNESDiffParameterDestroy_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__ 43bcaeba4dSBarry Smith #define __FUNCT__ "SNESMatrixFreeView2_Private" 44d2dddef7SLois Curfman McInnes /* 45d2dddef7SLois Curfman McInnes SNESMatrixFreeView2_Private - Views matrix-free parameters. 46d2dddef7SLois Curfman McInnes */ 47dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeView2_Private(Mat J,PetscViewer viewer) 48d2dddef7SLois Curfman McInnes { 49dfbe8321SBarry Smith PetscErrorCode ierr; 50d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 51ace3abfcSBarry Smith PetscBool iascii; 52d2dddef7SLois Curfman McInnes 5307a3eed9SLois Curfman McInnes PetscFunctionBegin; 54d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(J,(void**)&ctx);CHKERRQ(ierr); 55251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 5632077d6dSBarry Smith if (iascii) { 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 } 61a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," err=%G (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr); 62a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," umin=%G (minimum iterate parameter)\n",ctx->umin);CHKERRQ(ierr); 63758f395bSBarry Smith if (ctx->compute_err) { 6477431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," freq_err=%D (frequency for computing err)\n",ctx->compute_err_freq);CHKERRQ(ierr); 65d2dddef7SLois Curfman McInnes } 66758f395bSBarry Smith } 6707a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 68d2dddef7SLois Curfman McInnes } 69d2dddef7SLois Curfman McInnes 704a2ae208SSatish Balay #undef __FUNCT__ 714a2ae208SSatish Balay #define __FUNCT__ "SNESMatrixFreeMult2_Private" 72d2dddef7SLois Curfman McInnes /* 73d2dddef7SLois Curfman McInnes SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector 74d2dddef7SLois Curfman McInnes product, y = F'(u)*a: 75d2dddef7SLois Curfman McInnes y = (F(u + ha) - F(u)) /h, 76d2dddef7SLois Curfman McInnes where F = nonlinear function, as set by SNESSetFunction() 77d2dddef7SLois Curfman McInnes u = current iterate 78d2dddef7SLois Curfman McInnes h = difference interval 79d2dddef7SLois Curfman McInnes */ 80dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y) 81d2dddef7SLois Curfman McInnes { 82d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 83d2dddef7SLois Curfman McInnes SNES snes; 84145595abSSatish Balay PetscReal h,norm,sum,umin,noise; 8579f36460SBarry Smith PetscScalar hs,dot; 86d2dddef7SLois Curfman McInnes Vec w,U,F; 87dfbe8321SBarry Smith PetscErrorCode ierr,(*eval_fct)(SNES,Vec,Vec); 88e6ed2b1fSLois Curfman McInnes MPI_Comm comm; 8977431f27SBarry Smith PetscInt iter; 90d2dddef7SLois Curfman McInnes 9107a3eed9SLois Curfman McInnes PetscFunctionBegin; 92d2dddef7SLois Curfman McInnes /* We log matrix-free matrix-vector products separately, so that we can 93d2dddef7SLois Curfman McInnes separate the performance monitoring from the cases that use conventional 94d2dddef7SLois Curfman McInnes storage. We may eventually modify event logging to associate events 95d2dddef7SLois Curfman McInnes with particular objects, hence alleviating the more general problem. */ 963ec795f1SBarry Smith ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr); 97d2dddef7SLois Curfman McInnes 982d0c0e3bSBarry Smith ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 99e6ed2b1fSLois Curfman McInnes ierr = MatShellGetContext(mat,(void**)&ctx);CHKERRQ(ierr); 100e6ed2b1fSLois Curfman McInnes snes = ctx->snes; 101e6ed2b1fSLois Curfman McInnes w = ctx->w; 102e6ed2b1fSLois Curfman McInnes umin = ctx->umin; 103e6ed2b1fSLois Curfman McInnes 104d2dddef7SLois Curfman McInnes ierr = SNESGetSolution(snes,&U);CHKERRQ(ierr); 105d2dddef7SLois Curfman McInnes eval_fct = SNESComputeFunction; 1060298fd71SBarry Smith ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr); 107295f7fbbSLois Curfman McInnes 108d2dddef7SLois Curfman McInnes /* Determine a "good" step size, h */ 109d2dddef7SLois Curfman McInnes if (ctx->need_h) { 110d2dddef7SLois Curfman McInnes 111d2dddef7SLois Curfman McInnes /* Use Jorge's method to compute h */ 112d2dddef7SLois Curfman McInnes if (ctx->jorge) { 1138b5db460SBarry Smith ierr = SNESDiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr); 114d2dddef7SLois Curfman McInnes 115d2dddef7SLois Curfman McInnes /* Use the Brown/Saad method to compute h */ 116d2dddef7SLois Curfman McInnes } else { 117295f7fbbSLois Curfman McInnes /* Compute error if desired */ 118295f7fbbSLois Curfman McInnes ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr); 119145595abSSatish Balay if ((ctx->need_err) || ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter-1)%ctx->compute_err_freq)))) { 120d2dddef7SLois Curfman McInnes /* Use Jorge's method to compute noise */ 1218b5db460SBarry Smith ierr = SNESDiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr); 1221aa26658SKarl Rupp 123d2dddef7SLois Curfman McInnes ctx->error_rel = sqrt(noise); 1241aa26658SKarl Rupp 125ae15b995SBarry Smith ierr = PetscInfo3(snes,"Using Jorge's noise: noise=%G, sqrt(noise)=%G, h_more=%G\n",noise,ctx->error_rel,h);CHKERRQ(ierr); 1261aa26658SKarl Rupp 127295f7fbbSLois Curfman McInnes ctx->compute_err_iter = iter; 1283b2a6d2fSBarry Smith ctx->need_err = PETSC_FALSE; 129d2dddef7SLois Curfman McInnes } 130e6ed2b1fSLois Curfman McInnes 131691598edSBarry Smith ierr = VecDotBegin(U,a,&dot);CHKERRQ(ierr); 132691598edSBarry Smith ierr = VecNormBegin(a,NORM_1,&sum);CHKERRQ(ierr); 133691598edSBarry Smith ierr = VecNormBegin(a,NORM_2,&norm);CHKERRQ(ierr); 134691598edSBarry Smith ierr = VecDotEnd(U,a,&dot);CHKERRQ(ierr); 135691598edSBarry Smith ierr = VecNormEnd(a,NORM_1,&sum);CHKERRQ(ierr); 136691598edSBarry Smith ierr = VecNormEnd(a,NORM_2,&norm);CHKERRQ(ierr); 137e6ed2b1fSLois Curfman McInnes 138d2dddef7SLois Curfman McInnes 139d2dddef7SLois Curfman McInnes /* Safeguard for step sizes too small */ 1401aa26658SKarl Rupp if (sum == 0.0) { 1411aa26658SKarl Rupp dot = 1.0; 1421aa26658SKarl Rupp norm = 1.0; 1431aa26658SKarl Rupp } else if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum; 144145595abSSatish Balay else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum; 145145595abSSatish Balay h = PetscRealPart(ctx->error_rel*dot/(norm*norm)); 146d2dddef7SLois Curfman McInnes } 1471aa26658SKarl Rupp } else h = ctx->h; 1481aa26658SKarl Rupp 149ae15b995SBarry Smith if (!ctx->jorge || !ctx->need_h) {ierr = PetscInfo1(snes,"h = %G\n",h);CHKERRQ(ierr);} 150d2dddef7SLois Curfman McInnes 151d2dddef7SLois Curfman McInnes /* Evaluate function at F(u + ha) */ 152db0b4b35SLois Curfman McInnes hs = h; 1532dcb1b2aSMatthew Knepley ierr = VecWAXPY(w,hs,a,U);CHKERRQ(ierr); 154d2dddef7SLois Curfman McInnes ierr = eval_fct(snes,w,y);CHKERRQ(ierr); 15579f36460SBarry Smith ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr); 15679f36460SBarry Smith ierr = VecScale(y,1.0/hs);CHKERRQ(ierr); 1570298fd71SBarry Smith if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,NULL);CHKERRQ(ierr);} 158d2dddef7SLois Curfman McInnes 1593ec795f1SBarry Smith ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr); 16007a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 161d2dddef7SLois Curfman McInnes } 162d2dddef7SLois Curfman McInnes 1634a2ae208SSatish Balay #undef __FUNCT__ 164d8f46077SPeter Brune #define __FUNCT__ "SNESDefaultMatrixFreeCreate2" 165d2dddef7SLois Curfman McInnes /*@C 16663dd3a1aSKris Buschelman SNESMatrixFreeCreate2 - Creates a matrix-free matrix 167d2dddef7SLois Curfman McInnes context for use with a SNES solver. This matrix can be used as 168d2dddef7SLois Curfman McInnes the Jacobian argument for the routine SNESSetJacobian(). 169d2dddef7SLois Curfman McInnes 170d2dddef7SLois Curfman McInnes Input Parameters: 171d2dddef7SLois Curfman McInnes . snes - the SNES context 172d2dddef7SLois Curfman McInnes . x - vector where SNES solution is to be stored. 173d2dddef7SLois Curfman McInnes 174d2dddef7SLois Curfman McInnes Output Parameter: 175d2dddef7SLois Curfman McInnes . J - the matrix-free matrix 176d2dddef7SLois Curfman McInnes 17790c26489SBarry Smith Level: advanced 17890c26489SBarry Smith 179d2dddef7SLois Curfman McInnes Notes: 180d2dddef7SLois Curfman McInnes The matrix-free matrix context merely contains the function pointers 181d2dddef7SLois Curfman McInnes and work space for performing finite difference approximations of 182d2dddef7SLois Curfman McInnes Jacobian-vector products, J(u)*a, via 183d2dddef7SLois Curfman McInnes 184d2dddef7SLois Curfman McInnes $ J(u)*a = [J(u+h*a) - J(u)]/h, 185d2dddef7SLois Curfman McInnes $ where by default 186d2dddef7SLois Curfman McInnes $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 187d2dddef7SLois Curfman McInnes $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 188d2dddef7SLois Curfman McInnes $ where 189d2dddef7SLois Curfman McInnes $ error_rel = square root of relative error in 190d2dddef7SLois Curfman McInnes $ function evaluation 191d2dddef7SLois Curfman McInnes $ umin = minimum iterate parameter 192d2dddef7SLois Curfman McInnes $ Alternatively, the differencing parameter, h, can be set using 193e6ed2b1fSLois Curfman McInnes $ Jorge's nifty new strategy if one specifies the option 194e6ed2b1fSLois Curfman McInnes $ -snes_mf_jorge 195d2dddef7SLois Curfman McInnes 1963ec795f1SBarry Smith The user can set these parameters via MatMFFDSetFunctionError(). 1970598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_snes">SNES chapter of the users manual</A>. 198d2dddef7SLois Curfman McInnes 199d2dddef7SLois Curfman McInnes The user should call MatDestroy() when finished with the matrix-free 200d2dddef7SLois Curfman McInnes matrix context. 201d2dddef7SLois Curfman McInnes 202d2dddef7SLois Curfman McInnes Options Database Keys: 203d2dddef7SLois Curfman McInnes $ -snes_mf_err <error_rel> 204d2dddef7SLois Curfman McInnes $ -snes_mf_unim <umin> 205295f7fbbSLois Curfman McInnes $ -snes_mf_compute_err 206295f7fbbSLois Curfman McInnes $ -snes_mf_freq_err <freq> 207e6ed2b1fSLois Curfman McInnes $ -snes_mf_jorge 208d2dddef7SLois Curfman McInnes 209d2dddef7SLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix 210d2dddef7SLois Curfman McInnes 2113ec795f1SBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError() 212d2dddef7SLois Curfman McInnes @*/ 2137087cfbeSBarry Smith PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES snes,Vec x,Mat *J) 214d2dddef7SLois Curfman McInnes { 215d2dddef7SLois Curfman McInnes MPI_Comm comm; 216d2dddef7SLois Curfman McInnes MFCtx_Private *mfctx; 2176849ba73SBarry Smith PetscErrorCode ierr; 218a7cc72afSBarry Smith PetscInt n,nloc; 219ace3abfcSBarry Smith PetscBool flg; 220d2dddef7SLois Curfman McInnes char p[64]; 221d2dddef7SLois Curfman McInnes 22207a3eed9SLois Curfman McInnes PetscFunctionBegin; 22338f2d2fdSLisandro Dalcin ierr = PetscNewLog(snes,MFCtx_Private,&mfctx);CHKERRQ(ierr); 224d2dddef7SLois Curfman McInnes mfctx->sp = 0; 225d2dddef7SLois Curfman McInnes mfctx->snes = snes; 22677d8c4bbSBarry Smith mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 227d2dddef7SLois Curfman McInnes mfctx->umin = 1.e-6; 228d2dddef7SLois Curfman McInnes mfctx->h = 0.0; 229a7cc72afSBarry Smith mfctx->need_h = PETSC_TRUE; 230a7cc72afSBarry Smith mfctx->need_err = PETSC_FALSE; 231a7cc72afSBarry Smith mfctx->compute_err = PETSC_FALSE; 232295f7fbbSLois Curfman McInnes mfctx->compute_err_freq = 0; 233295f7fbbSLois Curfman McInnes mfctx->compute_err_iter = -1; 23490d69ab7SBarry Smith mfctx->compute_err = PETSC_FALSE; 23590d69ab7SBarry Smith mfctx->jorge = PETSC_FALSE; 2361aa26658SKarl Rupp 2370298fd71SBarry Smith ierr = PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_mf_err",&mfctx->error_rel,NULL);CHKERRQ(ierr); 2380298fd71SBarry Smith ierr = PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_mf_umin",&mfctx->umin,NULL);CHKERRQ(ierr); 2390298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_mf_jorge",&mfctx->jorge,NULL);CHKERRQ(ierr); 2400298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_mf_compute_err",&mfctx->compute_err,NULL);CHKERRQ(ierr); 2417adad957SLisandro Dalcin ierr = PetscOptionsGetInt(((PetscObject)snes)->prefix,"-snes_mf_freq_err",&mfctx->compute_err_freq,&flg);CHKERRQ(ierr); 242295f7fbbSLois Curfman McInnes if (flg) { 243295f7fbbSLois Curfman McInnes if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0; 2443b2a6d2fSBarry Smith mfctx->compute_err = PETSC_TRUE; 245295f7fbbSLois Curfman McInnes } 2463b2a6d2fSBarry Smith if (mfctx->compute_err) mfctx->need_err = PETSC_TRUE; 247295f7fbbSLois Curfman McInnes if (mfctx->jorge || mfctx->compute_err) { 2488b5db460SBarry Smith ierr = SNESDiffParameterCreate_More(snes,x,&mfctx->data);CHKERRQ(ierr); 249d2dddef7SLois Curfman McInnes } else mfctx->data = 0; 250d2dddef7SLois Curfman McInnes 2510298fd71SBarry Smith ierr = PetscOptionsHasName(NULL,"-help",&flg);CHKERRQ(ierr); 252549d3d68SSatish Balay ierr = PetscStrcpy(p,"-");CHKERRQ(ierr); 2537adad957SLisandro Dalcin if (((PetscObject)snes)->prefix) PetscStrcat(p,((PetscObject)snes)->prefix); 254d2dddef7SLois Curfman McInnes if (flg) { 255*ce94432eSBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)snes)," Matrix-free Options (via SNES):\n");CHKERRQ(ierr); 256*ce94432eSBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)snes)," %ssnes_mf_err <err>: set sqrt of relative error in function (default %G)\n",p,mfctx->error_rel);CHKERRQ(ierr); 257*ce94432eSBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)snes)," %ssnes_mf_umin <umin>: see users manual (default %G)\n",p,mfctx->umin);CHKERRQ(ierr); 258*ce94432eSBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)snes)," %ssnes_mf_jorge: use Jorge More's method\n",p);CHKERRQ(ierr); 259*ce94432eSBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)snes)," %ssnes_mf_compute_err: compute sqrt or relative error in function\n",p);CHKERRQ(ierr); 260*ce94432eSBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)snes)," %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)\n",p);CHKERRQ(ierr); 261*ce94432eSBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)snes)," %ssnes_mf_noise_file <file>: set file for printing noise info\n",p);CHKERRQ(ierr); 262d2dddef7SLois Curfman McInnes } 263d2dddef7SLois Curfman McInnes ierr = VecDuplicate(x,&mfctx->w);CHKERRQ(ierr); 264d2dddef7SLois Curfman McInnes ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr); 265d2dddef7SLois Curfman McInnes ierr = VecGetSize(x,&n);CHKERRQ(ierr); 266d2dddef7SLois Curfman McInnes ierr = VecGetLocalSize(x,&nloc);CHKERRQ(ierr); 267f69a0ea3SMatthew Knepley ierr = MatCreate(comm,J);CHKERRQ(ierr); 268f69a0ea3SMatthew Knepley ierr = MatSetSizes(*J,nloc,n,n,n);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 27552e6d16bSBarry Smith ierr = PetscLogObjectParent(*J,mfctx->w);CHKERRQ(ierr); 27652e6d16bSBarry Smith ierr = PetscLogObjectParent(snes,*J);CHKERRQ(ierr); 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 @*/ 3117087cfbeSBarry Smith PetscErrorCode SNESDefaultMatrixFreeSetParameters2(Mat mat,PetscReal error,PetscReal umin,PetscReal h) 312d2dddef7SLois Curfman McInnes { 313d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 314dfbe8321SBarry Smith PetscErrorCode 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; 323a7cc72afSBarry Smith ctx->need_h = PETSC_FALSE; 324d2dddef7SLois Curfman McInnes } 325d2dddef7SLois Curfman McInnes } 32607a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 327d2dddef7SLois Curfman McInnes } 328d2dddef7SLois Curfman McInnes 3297087cfbeSBarry Smith PetscErrorCode SNESUnSetMatrixFreeParameter(SNES snes) 330d2dddef7SLois Curfman McInnes { 331d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 332dfbe8321SBarry Smith PetscErrorCode ierr; 333d2dddef7SLois Curfman McInnes Mat mat; 334d2dddef7SLois Curfman McInnes 33507a3eed9SLois Curfman McInnes PetscFunctionBegin; 3360298fd71SBarry Smith ierr = SNESGetJacobian(snes,&mat,NULL,NULL,NULL);CHKERRQ(ierr); 337d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(mat,(void**)&ctx);CHKERRQ(ierr); 338a7cc72afSBarry Smith if (ctx) ctx->need_h = PETSC_TRUE; 33907a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 340d2dddef7SLois Curfman McInnes } 341d2dddef7SLois Curfman McInnes 342