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 */ 16ace3abfcSBarry Smith PetscBool jorge; /* flag indicating use of Jorge's method for determining 17d2dddef7SLois Curfman McInnes the differencing parameter */ 18145595abSSatish Balay PetscReal h; /* differencing parameter */ 19ace3abfcSBarry Smith PetscBool need_h; /* flag indicating whether we must compute h */ 20ace3abfcSBarry Smith PetscBool need_err; /* flag indicating whether we must currently compute error_rel */ 21ace3abfcSBarry Smith PetscBool compute_err; /* flag indicating whether we must ever compute error_rel */ 22a7cc72afSBarry Smith PetscInt compute_err_iter; /* last iter where we've computer error_rel */ 23a7cc72afSBarry Smith PetscInt compute_err_freq; /* frequency of computing error_rel */ 24d2dddef7SLois Curfman McInnes void *data; /* implementation-specific data */ 25d2dddef7SLois Curfman McInnes } MFCtx_Private; 26d2dddef7SLois Curfman McInnes 274a2ae208SSatish Balay #undef __FUNCT__ 28bcaeba4dSBarry Smith #define __FUNCT__ "SNESMatrixFreeDestroy2_Private" 29dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeDestroy2_Private(Mat mat) 30d2dddef7SLois Curfman McInnes { 31dfbe8321SBarry Smith PetscErrorCode ierr; 32d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 33d2dddef7SLois Curfman McInnes 3407a3eed9SLois Curfman McInnes PetscFunctionBegin; 35*22d28d08SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 366bf464f9SBarry Smith ierr = VecDestroy(&ctx->w);CHKERRQ(ierr); 376bf464f9SBarry Smith ierr = MatNullSpaceDestroy(&ctx->sp);CHKERRQ(ierr); 388b5db460SBarry Smith if (ctx->jorge || ctx->compute_err) {ierr = SNESDiffParameterDestroy_More(ctx->data);CHKERRQ(ierr);} 39606d414cSSatish Balay ierr = PetscFree(ctx);CHKERRQ(ierr); 4007a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 41d2dddef7SLois Curfman McInnes } 42d2dddef7SLois Curfman McInnes 434a2ae208SSatish Balay #undef __FUNCT__ 44bcaeba4dSBarry Smith #define __FUNCT__ "SNESMatrixFreeView2_Private" 45d2dddef7SLois Curfman McInnes /* 46d2dddef7SLois Curfman McInnes SNESMatrixFreeView2_Private - Views matrix-free parameters. 47d2dddef7SLois Curfman McInnes */ 48dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeView2_Private(Mat J,PetscViewer viewer) 49d2dddef7SLois Curfman McInnes { 50dfbe8321SBarry Smith PetscErrorCode ierr; 51d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 52ace3abfcSBarry Smith PetscBool iascii; 53d2dddef7SLois Curfman McInnes 5407a3eed9SLois Curfman McInnes PetscFunctionBegin; 55d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(J,(void **)&ctx);CHKERRQ(ierr); 56251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 5732077d6dSBarry Smith if (iascii) { 58b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," SNES matrix-free approximation:\n");CHKERRQ(ierr); 59758f395bSBarry Smith if (ctx->jorge) { 60b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," using Jorge's method of determining differencing parameter\n");CHKERRQ(ierr); 61758f395bSBarry Smith } 62a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," err=%G (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr); 63a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," umin=%G (minimum iterate parameter)\n",ctx->umin);CHKERRQ(ierr); 64758f395bSBarry Smith if (ctx->compute_err) { 6577431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," freq_err=%D (frequency for computing err)\n",ctx->compute_err_freq);CHKERRQ(ierr); 66d2dddef7SLois Curfman McInnes } 67758f395bSBarry Smith } else { 68e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by SNES matrix free Jorge",((PetscObject)viewer)->type_name); 69758f395bSBarry Smith } 7007a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 71d2dddef7SLois Curfman McInnes } 72d2dddef7SLois Curfman McInnes 734a2ae208SSatish Balay #undef __FUNCT__ 744a2ae208SSatish Balay #define __FUNCT__ "SNESMatrixFreeMult2_Private" 75d2dddef7SLois Curfman McInnes /* 76d2dddef7SLois Curfman McInnes SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector 77d2dddef7SLois Curfman McInnes product, y = F'(u)*a: 78d2dddef7SLois Curfman McInnes y = (F(u + ha) - F(u)) /h, 79d2dddef7SLois Curfman McInnes where F = nonlinear function, as set by SNESSetFunction() 80d2dddef7SLois Curfman McInnes u = current iterate 81d2dddef7SLois Curfman McInnes h = difference interval 82d2dddef7SLois Curfman McInnes */ 83dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y) 84d2dddef7SLois Curfman McInnes { 85d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 86d2dddef7SLois Curfman McInnes SNES snes; 87145595abSSatish Balay PetscReal h,norm,sum,umin,noise; 8879f36460SBarry Smith PetscScalar hs,dot; 89d2dddef7SLois Curfman McInnes Vec w,U,F; 90dfbe8321SBarry Smith PetscErrorCode ierr,(*eval_fct)(SNES,Vec,Vec); 91e6ed2b1fSLois Curfman McInnes MPI_Comm comm; 9277431f27SBarry Smith PetscInt iter; 93d2dddef7SLois Curfman McInnes 9407a3eed9SLois Curfman McInnes PetscFunctionBegin; 9507a3eed9SLois Curfman McInnes 96d2dddef7SLois Curfman McInnes /* We log matrix-free matrix-vector products separately, so that we can 97d2dddef7SLois Curfman McInnes separate the performance monitoring from the cases that use conventional 98d2dddef7SLois Curfman McInnes storage. We may eventually modify event logging to associate events 99d2dddef7SLois Curfman McInnes with particular objects, hence alleviating the more general problem. */ 1003ec795f1SBarry Smith ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr); 101d2dddef7SLois Curfman McInnes 1022d0c0e3bSBarry Smith ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 103e6ed2b1fSLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 104e6ed2b1fSLois Curfman McInnes snes = ctx->snes; 105e6ed2b1fSLois Curfman McInnes w = ctx->w; 106e6ed2b1fSLois Curfman McInnes umin = ctx->umin; 107e6ed2b1fSLois Curfman McInnes 108d2dddef7SLois Curfman McInnes ierr = SNESGetSolution(snes,&U);CHKERRQ(ierr); 109d2dddef7SLois Curfman McInnes eval_fct = SNESComputeFunction; 11074637425SBarry Smith ierr = SNESGetFunction(snes,&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 111295f7fbbSLois Curfman McInnes 112d2dddef7SLois Curfman McInnes /* Determine a "good" step size, h */ 113d2dddef7SLois Curfman McInnes if (ctx->need_h) { 114d2dddef7SLois Curfman McInnes 115d2dddef7SLois Curfman McInnes /* Use Jorge's method to compute h */ 116d2dddef7SLois Curfman McInnes if (ctx->jorge) { 1178b5db460SBarry Smith ierr = SNESDiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr); 118d2dddef7SLois Curfman McInnes 119d2dddef7SLois Curfman McInnes /* Use the Brown/Saad method to compute h */ 120d2dddef7SLois Curfman McInnes } else { 121295f7fbbSLois Curfman McInnes /* Compute error if desired */ 122295f7fbbSLois Curfman McInnes ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr); 123145595abSSatish Balay if ((ctx->need_err) || ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter-1)%ctx->compute_err_freq)))) { 124d2dddef7SLois Curfman McInnes /* Use Jorge's method to compute noise */ 1258b5db460SBarry Smith ierr = SNESDiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h);CHKERRQ(ierr); 126d2dddef7SLois Curfman McInnes ctx->error_rel = sqrt(noise); 127ae15b995SBarry Smith ierr = PetscInfo3(snes,"Using Jorge's noise: noise=%G, sqrt(noise)=%G, h_more=%G\n",noise,ctx->error_rel,h);CHKERRQ(ierr); 128295f7fbbSLois Curfman McInnes ctx->compute_err_iter = iter; 1293b2a6d2fSBarry Smith ctx->need_err = PETSC_FALSE; 130d2dddef7SLois Curfman McInnes } 131e6ed2b1fSLois Curfman McInnes 132691598edSBarry Smith ierr = VecDotBegin(U,a,&dot);CHKERRQ(ierr); 133691598edSBarry Smith ierr = VecNormBegin(a,NORM_1,&sum);CHKERRQ(ierr); 134691598edSBarry Smith ierr = VecNormBegin(a,NORM_2,&norm);CHKERRQ(ierr); 135691598edSBarry Smith ierr = VecDotEnd(U,a,&dot);CHKERRQ(ierr); 136691598edSBarry Smith ierr = VecNormEnd(a,NORM_1,&sum);CHKERRQ(ierr); 137691598edSBarry Smith ierr = VecNormEnd(a,NORM_2,&norm);CHKERRQ(ierr); 138e6ed2b1fSLois Curfman McInnes 139d2dddef7SLois Curfman McInnes 140d2dddef7SLois Curfman McInnes /* Safeguard for step sizes too small */ 141d2dddef7SLois Curfman McInnes if (sum == 0.0) {dot = 1.0; norm = 1.0;} 142145595abSSatish Balay else if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum; 143145595abSSatish Balay else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum; 144145595abSSatish Balay h = PetscRealPart(ctx->error_rel*dot/(norm*norm)); 145d2dddef7SLois Curfman McInnes } 146d2dddef7SLois Curfman McInnes } else { 147d2dddef7SLois Curfman McInnes h = ctx->h; 148d2dddef7SLois Curfman McInnes } 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); 15774637425SBarry Smith if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,PETSC_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; 2367adad957SLisandro Dalcin ierr = PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_mf_err",&mfctx->error_rel,PETSC_NULL);CHKERRQ(ierr); 2377adad957SLisandro Dalcin ierr = PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_mf_umin",&mfctx->umin,PETSC_NULL);CHKERRQ(ierr); 238acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_mf_jorge",&mfctx->jorge,PETSC_NULL);CHKERRQ(ierr); 239acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_mf_compute_err",&mfctx->compute_err,PETSC_NULL);CHKERRQ(ierr); 2407adad957SLisandro Dalcin ierr = PetscOptionsGetInt(((PetscObject)snes)->prefix,"-snes_mf_freq_err",&mfctx->compute_err_freq,&flg);CHKERRQ(ierr); 241295f7fbbSLois Curfman McInnes if (flg) { 242295f7fbbSLois Curfman McInnes if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0; 2433b2a6d2fSBarry Smith mfctx->compute_err = PETSC_TRUE; 244295f7fbbSLois Curfman McInnes } 2453b2a6d2fSBarry Smith if (mfctx->compute_err) mfctx->need_err = PETSC_TRUE; 246295f7fbbSLois Curfman McInnes if (mfctx->jorge || mfctx->compute_err) { 2478b5db460SBarry Smith ierr = SNESDiffParameterCreate_More(snes,x,&mfctx->data);CHKERRQ(ierr); 248d2dddef7SLois Curfman McInnes } else mfctx->data = 0; 249d2dddef7SLois Curfman McInnes 250b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr); 251549d3d68SSatish Balay ierr = PetscStrcpy(p,"-");CHKERRQ(ierr); 2527adad957SLisandro Dalcin if (((PetscObject)snes)->prefix) PetscStrcat(p,((PetscObject)snes)->prefix); 253d2dddef7SLois Curfman McInnes if (flg) { 2547adad957SLisandro Dalcin ierr = PetscPrintf(((PetscObject)snes)->comm," Matrix-free Options (via SNES):\n");CHKERRQ(ierr); 2557adad957SLisandro Dalcin ierr = PetscPrintf(((PetscObject)snes)->comm," %ssnes_mf_err <err>: set sqrt of relative error in function (default %G)\n",p,mfctx->error_rel);CHKERRQ(ierr); 2567adad957SLisandro Dalcin ierr = PetscPrintf(((PetscObject)snes)->comm," %ssnes_mf_umin <umin>: see users manual (default %G)\n",p,mfctx->umin);CHKERRQ(ierr); 2577adad957SLisandro Dalcin ierr = PetscPrintf(((PetscObject)snes)->comm," %ssnes_mf_jorge: use Jorge More's method\n",p);CHKERRQ(ierr); 2587adad957SLisandro Dalcin ierr = PetscPrintf(((PetscObject)snes)->comm," %ssnes_mf_compute_err: compute sqrt or relative error in function\n",p);CHKERRQ(ierr); 2597adad957SLisandro Dalcin ierr = PetscPrintf(((PetscObject)snes)->comm," %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)\n",p);CHKERRQ(ierr); 2607adad957SLisandro Dalcin ierr = PetscPrintf(((PetscObject)snes)->comm," %ssnes_mf_noise_file <file>: set file for printing noise info\n",p);CHKERRQ(ierr); 261d2dddef7SLois Curfman McInnes } 262d2dddef7SLois Curfman McInnes ierr = VecDuplicate(x,&mfctx->w);CHKERRQ(ierr); 263d2dddef7SLois Curfman McInnes ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr); 264d2dddef7SLois Curfman McInnes ierr = VecGetSize(x,&n);CHKERRQ(ierr); 265d2dddef7SLois Curfman McInnes ierr = VecGetLocalSize(x,&nloc);CHKERRQ(ierr); 266f69a0ea3SMatthew Knepley ierr = MatCreate(comm,J);CHKERRQ(ierr); 267f69a0ea3SMatthew Knepley ierr = MatSetSizes(*J,nloc,n,n,n);CHKERRQ(ierr); 268f204ca49SKris Buschelman ierr = MatSetType(*J,MATSHELL);CHKERRQ(ierr); 269f204ca49SKris Buschelman ierr = MatShellSetContext(*J,mfctx);CHKERRQ(ierr); 270c134de8dSSatish Balay ierr = MatShellSetOperation(*J,MATOP_MULT,(void(*)(void))SNESMatrixFreeMult2_Private);CHKERRQ(ierr); 271c134de8dSSatish Balay ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void(*)(void))SNESMatrixFreeDestroy2_Private);CHKERRQ(ierr); 272c134de8dSSatish Balay ierr = MatShellSetOperation(*J,MATOP_VIEW,(void(*)(void))SNESMatrixFreeView2_Private);CHKERRQ(ierr); 273d2dddef7SLois Curfman McInnes 27452e6d16bSBarry Smith ierr = PetscLogObjectParent(*J,mfctx->w);CHKERRQ(ierr); 27552e6d16bSBarry Smith ierr = PetscLogObjectParent(snes,*J);CHKERRQ(ierr); 27607a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 277d2dddef7SLois Curfman McInnes } 278d2dddef7SLois Curfman McInnes 2794a2ae208SSatish Balay #undef __FUNCT__ 2804a2ae208SSatish Balay #define __FUNCT__ "SNESDefaultMatrixFreeSetParameters2" 28190c26489SBarry Smith /*@C 282758f395bSBarry Smith SNESDefaultMatrixFreeSetParameters2 - Sets the parameters for the approximation of 283d2dddef7SLois Curfman McInnes matrix-vector products using finite differences. 284d2dddef7SLois Curfman McInnes 285d2dddef7SLois Curfman McInnes $ J(u)*a = [J(u+h*a) - J(u)]/h where 286d2dddef7SLois Curfman McInnes 287d2dddef7SLois Curfman McInnes either the user sets h directly here, or this parameter is computed via 288d2dddef7SLois Curfman McInnes 289d2dddef7SLois Curfman McInnes $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 290d2dddef7SLois Curfman McInnes $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 291d2dddef7SLois Curfman McInnes $ 292d2dddef7SLois Curfman McInnes 293d2dddef7SLois Curfman McInnes Input Parameters: 294758f395bSBarry Smith + mat - the matrix 295d2dddef7SLois Curfman McInnes . error_rel - relative error (should be set to the square root of 296d2dddef7SLois Curfman McInnes the relative error in the function evaluations) 297d2dddef7SLois Curfman McInnes . umin - minimum allowable u-value 298758f395bSBarry Smith - h - differencing parameter 299d2dddef7SLois Curfman McInnes 30090c26489SBarry Smith Level: advanced 30190c26489SBarry Smith 302d2dddef7SLois Curfman McInnes Notes: 303d2dddef7SLois Curfman McInnes If the user sets the parameter h directly, then this value will be used 304d2dddef7SLois Curfman McInnes instead of the default computation indicated above. 305d2dddef7SLois Curfman McInnes 306d2dddef7SLois Curfman McInnes .keywords: SNES, matrix-free, parameters 307d2dddef7SLois Curfman McInnes 3085a655dc6SBarry Smith .seealso: MatCreateSNESMF() 309d2dddef7SLois Curfman McInnes @*/ 3107087cfbeSBarry Smith PetscErrorCode SNESDefaultMatrixFreeSetParameters2(Mat mat,PetscReal error,PetscReal umin,PetscReal h) 311d2dddef7SLois Curfman McInnes { 312d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 313dfbe8321SBarry Smith PetscErrorCode ierr; 314d2dddef7SLois Curfman McInnes 31507a3eed9SLois Curfman McInnes PetscFunctionBegin; 316d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 317d2dddef7SLois Curfman McInnes if (ctx) { 318d2dddef7SLois Curfman McInnes if (error != PETSC_DEFAULT) ctx->error_rel = error; 319d2dddef7SLois Curfman McInnes if (umin != PETSC_DEFAULT) ctx->umin = umin; 320d2dddef7SLois Curfman McInnes if (h != PETSC_DEFAULT) { 321d2dddef7SLois Curfman McInnes ctx->h = h; 322a7cc72afSBarry Smith ctx->need_h = PETSC_FALSE; 323d2dddef7SLois Curfman McInnes } 324d2dddef7SLois Curfman McInnes } 32507a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 326d2dddef7SLois Curfman McInnes } 327d2dddef7SLois Curfman McInnes 3287087cfbeSBarry Smith PetscErrorCode SNESUnSetMatrixFreeParameter(SNES snes) 329d2dddef7SLois Curfman McInnes { 330d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 331dfbe8321SBarry Smith PetscErrorCode ierr; 332d2dddef7SLois Curfman McInnes Mat mat; 333d2dddef7SLois Curfman McInnes 33407a3eed9SLois Curfman McInnes PetscFunctionBegin; 33574637425SBarry Smith ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 336d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 337a7cc72afSBarry Smith if (ctx) ctx->need_h = PETSC_TRUE; 33807a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 339d2dddef7SLois Curfman McInnes } 340d2dddef7SLois Curfman McInnes 341