1d2dddef7SLois Curfman McInnes 2af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 36370053bSBarry Smith /* matimpl.h is needed only for logging of matrix operation */ 4af0996ceSBarry Smith #include <petsc/private/matimpl.h> 5d2dddef7SLois Curfman McInnes 625acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESUnSetMatrixFreeParameter(SNES); 725acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 825acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDefaultMatrixFreeSetParameters2(Mat,PetscReal,PetscReal,PetscReal); 925acbd8eSLisandro Dalcin 1025acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterCreate_More(SNES,Vec,void**); 1125acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterCompute_More(SNES,void*,Vec,Vec,PetscReal*,PetscReal*); 1225acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterDestroy_More(void*); 13d2dddef7SLois Curfman McInnes 14d2dddef7SLois Curfman McInnes typedef struct { /* default context for matrix-free SNES */ 15d2dddef7SLois Curfman McInnes SNES snes; /* SNES context */ 16d2dddef7SLois Curfman McInnes Vec w; /* work vector */ 1774637425SBarry Smith MatNullSpace sp; /* null space context */ 18145595abSSatish Balay PetscReal error_rel; /* square root of relative error in computing function */ 19145595abSSatish Balay PetscReal umin; /* minimum allowable u'a value relative to |u|_1 */ 20f5af7f23SKarl Rupp PetscBool jorge; /* flag indicating use of Jorge's method for determining the differencing parameter */ 21145595abSSatish Balay PetscReal h; /* differencing parameter */ 22ace3abfcSBarry Smith PetscBool need_h; /* flag indicating whether we must compute h */ 23ace3abfcSBarry Smith PetscBool need_err; /* flag indicating whether we must currently compute error_rel */ 24ace3abfcSBarry Smith PetscBool compute_err; /* flag indicating whether we must ever compute error_rel */ 25a7cc72afSBarry Smith PetscInt compute_err_iter; /* last iter where we've computer error_rel */ 26a7cc72afSBarry Smith PetscInt compute_err_freq; /* frequency of computing error_rel */ 27d2dddef7SLois Curfman McInnes void *data; /* implementation-specific data */ 28d2dddef7SLois Curfman McInnes } MFCtx_Private; 29d2dddef7SLois Curfman McInnes 30dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeDestroy2_Private(Mat mat) 31d2dddef7SLois Curfman McInnes { 32d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 33d2dddef7SLois Curfman McInnes 3407a3eed9SLois Curfman McInnes PetscFunctionBegin; 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(mat,&ctx)); 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ctx->w)); 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceDestroy(&ctx->sp)); 38*5f80ce2aSJacob Faibussowitsch if (ctx->jorge || ctx->compute_err) CHKERRQ(SNESDiffParameterDestroy_More(ctx->data)); 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ctx)); 4007a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 41d2dddef7SLois Curfman McInnes } 42d2dddef7SLois Curfman McInnes 43d2dddef7SLois Curfman McInnes /* 44d2dddef7SLois Curfman McInnes SNESMatrixFreeView2_Private - Views matrix-free parameters. 45d2dddef7SLois Curfman McInnes */ 46dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeView2_Private(Mat J,PetscViewer viewer) 47d2dddef7SLois Curfman McInnes { 48d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 49ace3abfcSBarry Smith PetscBool iascii; 50d2dddef7SLois Curfman McInnes 5107a3eed9SLois Curfman McInnes PetscFunctionBegin; 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(J,&ctx)); 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 5432077d6dSBarry Smith if (iascii) { 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," SNES matrix-free approximation:\n")); 56758f395bSBarry Smith if (ctx->jorge) { 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," using Jorge's method of determining differencing parameter\n")); 58758f395bSBarry Smith } 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," err=%g (relative error in function evaluation)\n",(double)ctx->error_rel)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," umin=%g (minimum iterate parameter)\n",(double)ctx->umin)); 61758f395bSBarry Smith if (ctx->compute_err) { 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," freq_err=%D (frequency for computing err)\n",ctx->compute_err_freq)); 63d2dddef7SLois Curfman McInnes } 64758f395bSBarry Smith } 6507a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 66d2dddef7SLois Curfman McInnes } 67d2dddef7SLois Curfman McInnes 68d2dddef7SLois Curfman McInnes /* 69d2dddef7SLois Curfman McInnes SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector 70d2dddef7SLois Curfman McInnes product, y = F'(u)*a: 71d2dddef7SLois Curfman McInnes y = (F(u + ha) - F(u)) /h, 72d2dddef7SLois Curfman McInnes where F = nonlinear function, as set by SNESSetFunction() 73d2dddef7SLois Curfman McInnes u = current iterate 74d2dddef7SLois Curfman McInnes h = difference interval 75d2dddef7SLois Curfman McInnes */ 76dfbe8321SBarry Smith PetscErrorCode SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y) 77d2dddef7SLois Curfman McInnes { 78d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 79d2dddef7SLois Curfman McInnes SNES snes; 80145595abSSatish Balay PetscReal h,norm,sum,umin,noise; 8179f36460SBarry Smith PetscScalar hs,dot; 82d2dddef7SLois Curfman McInnes Vec w,U,F; 83e6ed2b1fSLois Curfman McInnes MPI_Comm comm; 8477431f27SBarry Smith PetscInt iter; 85*5f80ce2aSJacob Faibussowitsch PetscErrorCode (*eval_fct)(SNES,Vec,Vec); 86d2dddef7SLois Curfman McInnes 8707a3eed9SLois Curfman McInnes PetscFunctionBegin; 88d2dddef7SLois Curfman McInnes /* We log matrix-free matrix-vector products separately, so that we can 89d2dddef7SLois Curfman McInnes separate the performance monitoring from the cases that use conventional 90d2dddef7SLois Curfman McInnes storage. We may eventually modify event logging to associate events 91d2dddef7SLois Curfman McInnes with particular objects, hence alleviating the more general problem. */ 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventBegin(MATMFFD_Mult,a,y,0,0)); 93d2dddef7SLois Curfman McInnes 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)mat,&comm)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(mat,&ctx)); 96e6ed2b1fSLois Curfman McInnes snes = ctx->snes; 97e6ed2b1fSLois Curfman McInnes w = ctx->w; 98e6ed2b1fSLois Curfman McInnes umin = ctx->umin; 99e6ed2b1fSLois Curfman McInnes 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetSolution(snes,&U)); 101d2dddef7SLois Curfman McInnes eval_fct = SNESComputeFunction; 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetFunction(snes,&F,NULL,NULL)); 103295f7fbbSLois Curfman McInnes 104d2dddef7SLois Curfman McInnes /* Determine a "good" step size, h */ 105d2dddef7SLois Curfman McInnes if (ctx->need_h) { 106d2dddef7SLois Curfman McInnes 107d2dddef7SLois Curfman McInnes /* Use Jorge's method to compute h */ 108d2dddef7SLois Curfman McInnes if (ctx->jorge) { 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h)); 110d2dddef7SLois Curfman McInnes 111d2dddef7SLois Curfman McInnes /* Use the Brown/Saad method to compute h */ 112d2dddef7SLois Curfman McInnes } else { 113295f7fbbSLois Curfman McInnes /* Compute error if desired */ 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetIterationNumber(snes,&iter)); 115145595abSSatish Balay if ((ctx->need_err) || ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter-1)%ctx->compute_err_freq)))) { 116d2dddef7SLois Curfman McInnes /* Use Jorge's method to compute noise */ 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h)); 1181aa26658SKarl Rupp 119369cc0aeSBarry Smith ctx->error_rel = PetscSqrtReal(noise); 1201aa26658SKarl Rupp 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(snes,"Using Jorge's noise: noise=%g, sqrt(noise)=%g, h_more=%g\n",(double)noise,(double)ctx->error_rel,(double)h)); 1221aa26658SKarl Rupp 123295f7fbbSLois Curfman McInnes ctx->compute_err_iter = iter; 1243b2a6d2fSBarry Smith ctx->need_err = PETSC_FALSE; 125d2dddef7SLois Curfman McInnes } 126e6ed2b1fSLois Curfman McInnes 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDotBegin(U,a,&dot)); 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNormBegin(a,NORM_1,&sum)); 129*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNormBegin(a,NORM_2,&norm)); 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDotEnd(U,a,&dot)); 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNormEnd(a,NORM_1,&sum)); 132*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNormEnd(a,NORM_2,&norm)); 133e6ed2b1fSLois Curfman McInnes 134d2dddef7SLois Curfman McInnes /* Safeguard for step sizes too small */ 1351aa26658SKarl Rupp if (sum == 0.0) { 1361aa26658SKarl Rupp dot = 1.0; 1371aa26658SKarl Rupp norm = 1.0; 1381aa26658SKarl Rupp } else if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum; 139145595abSSatish Balay else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum; 140145595abSSatish Balay h = PetscRealPart(ctx->error_rel*dot/(norm*norm)); 141d2dddef7SLois Curfman McInnes } 1421aa26658SKarl Rupp } else h = ctx->h; 1431aa26658SKarl Rupp 144*5f80ce2aSJacob Faibussowitsch if (!ctx->jorge || !ctx->need_h) CHKERRQ(PetscInfo(snes,"h = %g\n",(double)h)); 145d2dddef7SLois Curfman McInnes 146d2dddef7SLois Curfman McInnes /* Evaluate function at F(u + ha) */ 147db0b4b35SLois Curfman McInnes hs = h; 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(w,hs,a,U)); 149*5f80ce2aSJacob Faibussowitsch CHKERRQ(eval_fct(snes,w,y)); 150*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,-1.0,F)); 151*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(y,1.0/hs)); 152*5f80ce2aSJacob Faibussowitsch if (mat->nullsp) CHKERRQ(MatNullSpaceRemove(mat->nullsp,y)); 153d2dddef7SLois Curfman McInnes 154*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventEnd(MATMFFD_Mult,a,y,0,0)); 15507a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 156d2dddef7SLois Curfman McInnes } 157d2dddef7SLois Curfman McInnes 158d2dddef7SLois Curfman McInnes /*@C 15963dd3a1aSKris Buschelman SNESMatrixFreeCreate2 - Creates a matrix-free matrix 160d2dddef7SLois Curfman McInnes context for use with a SNES solver. This matrix can be used as 161d2dddef7SLois Curfman McInnes the Jacobian argument for the routine SNESSetJacobian(). 162d2dddef7SLois Curfman McInnes 163d2dddef7SLois Curfman McInnes Input Parameters: 164a2b725a8SWilliam Gropp + snes - the SNES context 165a2b725a8SWilliam Gropp - x - vector where SNES solution is to be stored. 166d2dddef7SLois Curfman McInnes 167d2dddef7SLois Curfman McInnes Output Parameter: 168d2dddef7SLois Curfman McInnes . J - the matrix-free matrix 169d2dddef7SLois Curfman McInnes 17090c26489SBarry Smith Level: advanced 17190c26489SBarry Smith 172d2dddef7SLois Curfman McInnes Notes: 173d2dddef7SLois Curfman McInnes The matrix-free matrix context merely contains the function pointers 174d2dddef7SLois Curfman McInnes and work space for performing finite difference approximations of 175d2dddef7SLois Curfman McInnes Jacobian-vector products, J(u)*a, via 176d2dddef7SLois Curfman McInnes 177d2dddef7SLois Curfman McInnes $ J(u)*a = [J(u+h*a) - J(u)]/h, 178d2dddef7SLois Curfman McInnes $ where by default 179d2dddef7SLois Curfman McInnes $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 180d2dddef7SLois Curfman McInnes $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 181d2dddef7SLois Curfman McInnes $ where 182d2dddef7SLois Curfman McInnes $ error_rel = square root of relative error in 183d2dddef7SLois Curfman McInnes $ function evaluation 184d2dddef7SLois Curfman McInnes $ umin = minimum iterate parameter 185d2dddef7SLois Curfman McInnes $ Alternatively, the differencing parameter, h, can be set using 186e6ed2b1fSLois Curfman McInnes $ Jorge's nifty new strategy if one specifies the option 187e6ed2b1fSLois Curfman McInnes $ -snes_mf_jorge 188d2dddef7SLois Curfman McInnes 1893ec795f1SBarry Smith The user can set these parameters via MatMFFDSetFunctionError(). 190a7f22e61SSatish Balay See Users-Manual: ch_snes for details 191d2dddef7SLois Curfman McInnes 192d2dddef7SLois Curfman McInnes The user should call MatDestroy() when finished with the matrix-free 193d2dddef7SLois Curfman McInnes matrix context. 194d2dddef7SLois Curfman McInnes 195d2dddef7SLois Curfman McInnes Options Database Keys: 196d2dddef7SLois Curfman McInnes $ -snes_mf_err <error_rel> 197d2dddef7SLois Curfman McInnes $ -snes_mf_unim <umin> 198295f7fbbSLois Curfman McInnes $ -snes_mf_compute_err 199295f7fbbSLois Curfman McInnes $ -snes_mf_freq_err <freq> 200e6ed2b1fSLois Curfman McInnes $ -snes_mf_jorge 201d2dddef7SLois Curfman McInnes 2023ec795f1SBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError() 203d2dddef7SLois Curfman McInnes @*/ 2047087cfbeSBarry Smith PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES snes,Vec x,Mat *J) 205d2dddef7SLois Curfman McInnes { 206d2dddef7SLois Curfman McInnes MPI_Comm comm; 207d2dddef7SLois Curfman McInnes MFCtx_Private *mfctx; 208a7cc72afSBarry Smith PetscInt n,nloc; 209ace3abfcSBarry Smith PetscBool flg; 210d2dddef7SLois Curfman McInnes char p[64]; 211d2dddef7SLois Curfman McInnes 21207a3eed9SLois Curfman McInnes PetscFunctionBegin; 213*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(snes,&mfctx)); 2149e5d0892SLisandro Dalcin mfctx->sp = NULL; 215d2dddef7SLois Curfman McInnes mfctx->snes = snes; 21677d8c4bbSBarry Smith mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 217d2dddef7SLois Curfman McInnes mfctx->umin = 1.e-6; 218d2dddef7SLois Curfman McInnes mfctx->h = 0.0; 219a7cc72afSBarry Smith mfctx->need_h = PETSC_TRUE; 220a7cc72afSBarry Smith mfctx->need_err = PETSC_FALSE; 221a7cc72afSBarry Smith mfctx->compute_err = PETSC_FALSE; 222295f7fbbSLois Curfman McInnes mfctx->compute_err_freq = 0; 223295f7fbbSLois Curfman McInnes mfctx->compute_err_iter = -1; 22490d69ab7SBarry Smith mfctx->compute_err = PETSC_FALSE; 22590d69ab7SBarry Smith mfctx->jorge = PETSC_FALSE; 2261aa26658SKarl Rupp 227*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_err",&mfctx->error_rel,NULL)); 228*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_umin",&mfctx->umin,NULL)); 229*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_jorge",&mfctx->jorge,NULL)); 230*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_compute_err",&mfctx->compute_err,NULL)); 231*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_freq_err",&mfctx->compute_err_freq,&flg)); 232295f7fbbSLois Curfman McInnes if (flg) { 233295f7fbbSLois Curfman McInnes if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0; 2343b2a6d2fSBarry Smith mfctx->compute_err = PETSC_TRUE; 235295f7fbbSLois Curfman McInnes } 2363b2a6d2fSBarry Smith if (mfctx->compute_err) mfctx->need_err = PETSC_TRUE; 237295f7fbbSLois Curfman McInnes if (mfctx->jorge || mfctx->compute_err) { 238*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDiffParameterCreate_More(snes,x,&mfctx->data)); 2399e5d0892SLisandro Dalcin } else mfctx->data = NULL; 240d2dddef7SLois Curfman McInnes 241*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasHelp(((PetscObject)snes)->options,&flg)); 242*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncpy(p,"-",sizeof(p))); 243*5f80ce2aSJacob Faibussowitsch if (((PetscObject)snes)->prefix) CHKERRQ(PetscStrlcat(p,((PetscObject)snes)->prefix,sizeof(p))); 244d2dddef7SLois Curfman McInnes if (flg) { 245*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)snes)," Matrix-free Options (via SNES):\n")); 246*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)snes)," %ssnes_mf_err <err>: set sqrt of relative error in function (default %g)\n",p,(double)mfctx->error_rel)); 247*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)snes)," %ssnes_mf_umin <umin>: see users manual (default %g)\n",p,(double)mfctx->umin)); 248*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)snes)," %ssnes_mf_jorge: use Jorge More's method\n",p)); 249*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)snes)," %ssnes_mf_compute_err: compute sqrt or relative error in function\n",p)); 250*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)snes)," %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)\n",p)); 251*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)snes)," %ssnes_mf_noise_file <file>: set file for printing noise info\n",p)); 252d2dddef7SLois Curfman McInnes } 253*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&mfctx->w)); 254*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)x,&comm)); 255*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(x,&n)); 256*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(x,&nloc)); 257*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(comm,J)); 258*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(*J,nloc,n,n,n)); 259*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(*J,MATSHELL)); 260*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetContext(*J,mfctx)); 261*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(*J,MATOP_MULT,(void (*)(void))SNESMatrixFreeMult2_Private)); 262*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(*J,MATOP_DESTROY,(void (*)(void))SNESMatrixFreeDestroy2_Private)); 263*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(*J,MATOP_VIEW,(void (*)(void))SNESMatrixFreeView2_Private)); 264*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(*J)); 265d2dddef7SLois Curfman McInnes 266*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectParent((PetscObject)*J,(PetscObject)mfctx->w)); 267*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectParent((PetscObject)snes,(PetscObject)*J)); 26807a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 269d2dddef7SLois Curfman McInnes } 270d2dddef7SLois Curfman McInnes 27190c26489SBarry Smith /*@C 272758f395bSBarry Smith SNESDefaultMatrixFreeSetParameters2 - Sets the parameters for the approximation of 273d2dddef7SLois Curfman McInnes matrix-vector products using finite differences. 274d2dddef7SLois Curfman McInnes 275d2dddef7SLois Curfman McInnes $ J(u)*a = [J(u+h*a) - J(u)]/h where 276d2dddef7SLois Curfman McInnes 277d2dddef7SLois Curfman McInnes either the user sets h directly here, or this parameter is computed via 278d2dddef7SLois Curfman McInnes 279d2dddef7SLois Curfman McInnes $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 280d2dddef7SLois Curfman McInnes $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 281d2dddef7SLois Curfman McInnes $ 282d2dddef7SLois Curfman McInnes 283d2dddef7SLois Curfman McInnes Input Parameters: 284758f395bSBarry Smith + mat - the matrix 285d2dddef7SLois Curfman McInnes . error_rel - relative error (should be set to the square root of 286d2dddef7SLois Curfman McInnes the relative error in the function evaluations) 287d2dddef7SLois Curfman McInnes . umin - minimum allowable u-value 288758f395bSBarry Smith - h - differencing parameter 289d2dddef7SLois Curfman McInnes 29090c26489SBarry Smith Level: advanced 29190c26489SBarry Smith 292d2dddef7SLois Curfman McInnes Notes: 293d2dddef7SLois Curfman McInnes If the user sets the parameter h directly, then this value will be used 294d2dddef7SLois Curfman McInnes instead of the default computation indicated above. 295d2dddef7SLois Curfman McInnes 2965a655dc6SBarry Smith .seealso: MatCreateSNESMF() 297d2dddef7SLois Curfman McInnes @*/ 2987087cfbeSBarry Smith PetscErrorCode SNESDefaultMatrixFreeSetParameters2(Mat mat,PetscReal error,PetscReal umin,PetscReal h) 299d2dddef7SLois Curfman McInnes { 300d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 301d2dddef7SLois Curfman McInnes 30207a3eed9SLois Curfman McInnes PetscFunctionBegin; 303*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(mat,&ctx)); 304d2dddef7SLois Curfman McInnes if (ctx) { 305d2dddef7SLois Curfman McInnes if (error != PETSC_DEFAULT) ctx->error_rel = error; 306d2dddef7SLois Curfman McInnes if (umin != PETSC_DEFAULT) ctx->umin = umin; 307d2dddef7SLois Curfman McInnes if (h != PETSC_DEFAULT) { 308d2dddef7SLois Curfman McInnes ctx->h = h; 309a7cc72afSBarry Smith ctx->need_h = PETSC_FALSE; 310d2dddef7SLois Curfman McInnes } 311d2dddef7SLois Curfman McInnes } 31207a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 313d2dddef7SLois Curfman McInnes } 314d2dddef7SLois Curfman McInnes 3157087cfbeSBarry Smith PetscErrorCode SNESUnSetMatrixFreeParameter(SNES snes) 316d2dddef7SLois Curfman McInnes { 317d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 318d2dddef7SLois Curfman McInnes Mat mat; 319d2dddef7SLois Curfman McInnes 32007a3eed9SLois Curfman McInnes PetscFunctionBegin; 321*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetJacobian(snes,&mat,NULL,NULL,NULL)); 322*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(mat,&ctx)); 323a7cc72afSBarry Smith if (ctx) ctx->need_h = PETSC_TRUE; 32407a3eed9SLois Curfman McInnes PetscFunctionReturn(0); 325d2dddef7SLois Curfman McInnes } 326