1d834d2e4SLois Curfman McInnes #ifdef PETSC_RCS_HEADER 2*f1d4f8a3SBarry Smith static char vcid[] = "$Id: snesmfj2.c,v 1.6 1997/10/21 03:15:31 curfman Exp bsmith $"; 3d2dddef7SLois Curfman McInnes #endif 4d2dddef7SLois Curfman McInnes 5d2dddef7SLois Curfman McInnes #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 6d2dddef7SLois Curfman McInnes #include "pinclude/pviewer.h" 7d2dddef7SLois Curfman McInnes #include <math.h> 8d2dddef7SLois Curfman McInnes 9d2dddef7SLois Curfman McInnes extern int DiffParameterCreate_More(SNES,Vec,void**); 10d2dddef7SLois Curfman McInnes extern int DiffParameterCompute_More(SNES,void*,Vec,Vec,double*,double*); 11d2dddef7SLois Curfman McInnes extern int DiffParameterDestroy_More(void*); 12d2dddef7SLois Curfman McInnes 13d2dddef7SLois Curfman McInnes typedef struct { /* default context for matrix-free SNES */ 14d2dddef7SLois Curfman McInnes SNES snes; /* SNES context */ 15d2dddef7SLois Curfman McInnes Vec w; /* work vector */ 16d2dddef7SLois Curfman McInnes PCNullSpace sp; /* null space context */ 17d2dddef7SLois Curfman McInnes double error_rel; /* square root of relative error in computing function */ 18d2dddef7SLois Curfman McInnes double umin; /* minimum allowable u'a value relative to |u|_1 */ 19d2dddef7SLois Curfman McInnes int jorge; /* flag indicating use of Jorge's method for determining 20d2dddef7SLois Curfman McInnes the differencing parameter */ 21db0b4b35SLois Curfman McInnes double h; /* differencing parameter */ 22d2dddef7SLois Curfman McInnes int need_h; /* flag indicating whether we must compute h */ 23295f7fbbSLois Curfman McInnes int need_err; /* flag indicating whether we must currently compute error_rel */ 24295f7fbbSLois Curfman McInnes int compute_err; /* flag indicating whether we must ever compute error_rel */ 25295f7fbbSLois Curfman McInnes int compute_err_iter; /* last iter where we've computer error_rel */ 26295f7fbbSLois Curfman McInnes int compute_err_freq; /* frequency of computing error_rel */ 27d2dddef7SLois Curfman McInnes void *data; /* implementation-specific data */ 28d2dddef7SLois Curfman McInnes } MFCtx_Private; 29d2dddef7SLois Curfman McInnes 30d2dddef7SLois Curfman McInnes #undef __FUNC__ 31d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeDestroy2_Private" /* ADIC Ignore */ 32*f1d4f8a3SBarry Smith int SNESMatrixFreeDestroy2_Private(Mat mat) 33d2dddef7SLois Curfman McInnes { 34d2dddef7SLois Curfman McInnes int ierr; 35d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 36d2dddef7SLois Curfman McInnes 37d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx); 38d2dddef7SLois Curfman McInnes ierr = VecDestroy(ctx->w); CHKERRQ(ierr); 39d2dddef7SLois Curfman McInnes if (ctx->sp) {ierr = PCNullSpaceDestroy(ctx->sp); CHKERRQ(ierr);} 40295f7fbbSLois Curfman McInnes if (ctx->jorge || ctx->compute_err) {ierr = DiffParameterDestroy_More(ctx->data); CHKERRQ(ierr);} 41d2dddef7SLois Curfman McInnes PetscFree(ctx); 42d2dddef7SLois Curfman McInnes return 0; 43d2dddef7SLois Curfman McInnes } 44d2dddef7SLois Curfman McInnes 45d2dddef7SLois Curfman McInnes #undef __FUNC__ 46d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeView2_Private" /* ADIC Ignore */ 47d2dddef7SLois Curfman McInnes /* 48d2dddef7SLois Curfman McInnes SNESMatrixFreeView2_Private - Views matrix-free parameters. 49d2dddef7SLois Curfman McInnes */ 50d2dddef7SLois Curfman McInnes int SNESMatrixFreeView2_Private(Mat J,Viewer viewer) 51d2dddef7SLois Curfman McInnes { 52d2dddef7SLois Curfman McInnes int ierr; 53d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 54d2dddef7SLois Curfman McInnes MPI_Comm comm; 55d2dddef7SLois Curfman McInnes FILE *fd; 56d2dddef7SLois Curfman McInnes ViewerType vtype; 57d2dddef7SLois Curfman McInnes 58d2dddef7SLois Curfman McInnes PetscObjectGetComm((PetscObject)J,&comm); 59d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 60d2dddef7SLois Curfman McInnes ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 61d2dddef7SLois Curfman McInnes ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 62d2dddef7SLois Curfman McInnes if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 63d2dddef7SLois Curfman McInnes PetscFPrintf(comm,fd," SNES matrix-free approximation:\n"); 64d2dddef7SLois Curfman McInnes if (ctx->jorge) 65d2dddef7SLois Curfman McInnes PetscFPrintf(comm,fd," using Jorge's method of determining differencing parameter\n"); 66d2dddef7SLois Curfman McInnes PetscFPrintf(comm,fd," err=%g (relative error in function evaluation)\n",ctx->error_rel); 67d2dddef7SLois Curfman McInnes PetscFPrintf(comm,fd," umin=%g (minimum iterate parameter)\n",ctx->umin); 68295f7fbbSLois Curfman McInnes if (ctx->compute_err) 69295f7fbbSLois Curfman McInnes PetscFPrintf(comm,fd," freq_err=%d (frequency for computing err)\n",ctx->compute_err_freq); 70d2dddef7SLois Curfman McInnes } 71d2dddef7SLois Curfman McInnes return 0; 72d2dddef7SLois Curfman McInnes } 73d2dddef7SLois Curfman McInnes 74e6ed2b1fSLois Curfman McInnes extern int VecDot_Seq(Vec,Vec,Scalar *); 75e6ed2b1fSLois Curfman McInnes extern int VecNorm_Seq(Vec,NormType,double *); 76e6ed2b1fSLois Curfman McInnes 77d2dddef7SLois Curfman McInnes #undef __FUNC__ 78d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESMatrixFreeMult2_Private" 79d2dddef7SLois Curfman McInnes /* 80d2dddef7SLois Curfman McInnes SNESMatrixFreeMult2_Private - Default matrix-free form for Jacobian-vector 81d2dddef7SLois Curfman McInnes product, y = F'(u)*a: 82d2dddef7SLois Curfman McInnes y = ( F(u + ha) - F(u) ) /h, 83d2dddef7SLois Curfman McInnes where F = nonlinear function, as set by SNESSetFunction() 84d2dddef7SLois Curfman McInnes u = current iterate 85d2dddef7SLois Curfman McInnes h = difference interval 86d2dddef7SLois Curfman McInnes */ 87d2dddef7SLois Curfman McInnes int SNESMatrixFreeMult2_Private(Mat mat,Vec a,Vec y) 88d2dddef7SLois Curfman McInnes { 89d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 90d2dddef7SLois Curfman McInnes SNES snes; 91db0b4b35SLois Curfman McInnes double h, norm, sum, umin, noise, ovalues[3],values[3]; 92db0b4b35SLois Curfman McInnes Scalar hs, dot, mone = -1.0; 93d2dddef7SLois Curfman McInnes Vec w,U,F; 94295f7fbbSLois Curfman McInnes int ierr, iter, (*eval_fct)(SNES,Vec,Vec); 95e6ed2b1fSLois Curfman McInnes MPI_Comm comm; 96d2dddef7SLois Curfman McInnes 97d2dddef7SLois Curfman McInnes /* We log matrix-free matrix-vector products separately, so that we can 98d2dddef7SLois Curfman McInnes separate the performance monitoring from the cases that use conventional 99d2dddef7SLois Curfman McInnes storage. We may eventually modify event logging to associate events 100d2dddef7SLois Curfman McInnes with particular objects, hence alleviating the more general problem. */ 101d2dddef7SLois Curfman McInnes PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0); 102d2dddef7SLois Curfman McInnes 103e6ed2b1fSLois Curfman McInnes PetscObjectGetComm((PetscObject)mat,&comm); 104e6ed2b1fSLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 105e6ed2b1fSLois Curfman McInnes snes = ctx->snes; 106e6ed2b1fSLois Curfman McInnes w = ctx->w; 107e6ed2b1fSLois Curfman McInnes umin = ctx->umin; 108e6ed2b1fSLois Curfman McInnes 109d2dddef7SLois Curfman McInnes ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr); 110d2dddef7SLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 111d2dddef7SLois Curfman McInnes eval_fct = SNESComputeFunction; 112d2dddef7SLois Curfman McInnes ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr); 113d2dddef7SLois Curfman McInnes } 114d2dddef7SLois Curfman McInnes else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 115d2dddef7SLois Curfman McInnes eval_fct = SNESComputeGradient; 116d2dddef7SLois Curfman McInnes ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr); 117d2dddef7SLois Curfman McInnes } 118d2dddef7SLois Curfman McInnes else SETERRQ(1,0,"Invalid method class"); 119d2dddef7SLois Curfman McInnes 120295f7fbbSLois Curfman McInnes 121d2dddef7SLois Curfman McInnes /* Determine a "good" step size, h */ 122d2dddef7SLois Curfman McInnes if (ctx->need_h) { 123d2dddef7SLois Curfman McInnes 124d2dddef7SLois Curfman McInnes /* Use Jorge's method to compute h */ 125d2dddef7SLois Curfman McInnes if (ctx->jorge) { 126d2dddef7SLois Curfman McInnes ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h); CHKERRQ(ierr); 127d2dddef7SLois Curfman McInnes 128d2dddef7SLois Curfman McInnes /* Use the Brown/Saad method to compute h */ 129d2dddef7SLois Curfman McInnes } else { 130295f7fbbSLois Curfman McInnes /* Compute error if desired */ 131295f7fbbSLois Curfman McInnes ierr = SNESGetIterationNumber(snes,&iter); CHKERRQ(ierr); 132295f7fbbSLois Curfman McInnes if ((ctx->need_err) || 133295f7fbbSLois Curfman McInnes ((ctx->compute_err_freq) && (ctx->compute_err_iter != iter) && (!((iter-1)%ctx->compute_err_freq)))) { 134d2dddef7SLois Curfman McInnes /* Use Jorge's method to compute noise */ 135d2dddef7SLois Curfman McInnes ierr = DiffParameterCompute_More(snes,ctx->data,U,a,&noise,&h); CHKERRQ(ierr); 136d2dddef7SLois Curfman McInnes ctx->error_rel = sqrt(noise); 137d2dddef7SLois Curfman McInnes PLogInfo(snes,"SNESMatrixFreeMult2_Private: Using Jorge's noise: noise=%g, sqrt(noise)=%g, h_more=%g\n", 138d2dddef7SLois Curfman McInnes noise,ctx->error_rel,h); 139295f7fbbSLois Curfman McInnes ctx->compute_err_iter = iter; 140d2dddef7SLois Curfman McInnes ctx->need_err = 0; 141d2dddef7SLois Curfman McInnes } 142e6ed2b1fSLois Curfman McInnes 143e6ed2b1fSLois Curfman McInnes /* 144d2dddef7SLois Curfman McInnes ierr = VecDot(U,a,&dot); CHKERRQ(ierr); 145d2dddef7SLois Curfman McInnes ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr); 146d2dddef7SLois Curfman McInnes ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr); 147e6ed2b1fSLois Curfman McInnes */ 148e6ed2b1fSLois Curfman McInnes 149e6ed2b1fSLois Curfman McInnes /* 150e6ed2b1fSLois Curfman McInnes Call the Seq Vector routines and then do a single reduction 151e6ed2b1fSLois Curfman McInnes to reduce the number of communications required 152e6ed2b1fSLois Curfman McInnes */ 153e6ed2b1fSLois Curfman McInnes 154db0b4b35SLois Curfman McInnes #if !defined(USE_PETSC_COMPLEX) 155e6ed2b1fSLois Curfman McInnes PLogEventBegin(VEC_Dot,U,a,0,0); 156e6ed2b1fSLois Curfman McInnes ierr = VecDot_Seq(U,a,ovalues); CHKERRQ(ierr); 157e6ed2b1fSLois Curfman McInnes PLogEventEnd(VEC_Dot,U,a,0,0); 158e6ed2b1fSLois Curfman McInnes PLogEventBegin(VEC_Norm,a,0,0,0); 159e6ed2b1fSLois Curfman McInnes ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr); 160e6ed2b1fSLois Curfman McInnes ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr); 161e6ed2b1fSLois Curfman McInnes ovalues[2] = ovalues[2]*ovalues[2]; 162e6ed2b1fSLois Curfman McInnes MPI_Allreduce(ovalues,values,3,MPI_DOUBLE,MPI_SUM,comm ); 163e6ed2b1fSLois Curfman McInnes dot = values[0]; sum = values[1]; norm = sqrt(values[2]); 164e6ed2b1fSLois Curfman McInnes PLogEventEnd(VEC_Norm,a,0,0,0); 165e6ed2b1fSLois Curfman McInnes #else 166e6ed2b1fSLois Curfman McInnes { 167e6ed2b1fSLois Curfman McInnes Scalar cvalues[3],covalues[3]; 168e6ed2b1fSLois Curfman McInnes 169e6ed2b1fSLois Curfman McInnes PLogEventBegin(VEC_Dot,U,a,0,0); 170e6ed2b1fSLois Curfman McInnes ierr = VecDot_Seq(U,a,covalues); CHKERRQ(ierr); 171e6ed2b1fSLois Curfman McInnes PLogEventEnd(VEC_Dot,U,a,0,0); 172e6ed2b1fSLois Curfman McInnes PLogEventBegin(VEC_Norm,a,0,0,0); 173e6ed2b1fSLois Curfman McInnes ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr); 174e6ed2b1fSLois Curfman McInnes ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr); 175e6ed2b1fSLois Curfman McInnes covalues[1] = ovalues[1]; 176e6ed2b1fSLois Curfman McInnes covalues[2] = ovalues[2]*ovalues[2]; 177e6ed2b1fSLois Curfman McInnes MPI_Allreduce(covalues,cvalues,6,MPI_DOUBLE,MPI_SUM,comm ); 178e6ed2b1fSLois Curfman McInnes dot = cvalues[0]; sum = PetscReal(cvalues[1]); norm = sqrt(PetscReal(cvalues[2])); 179e6ed2b1fSLois Curfman McInnes PLogEventEnd(VEC_Norm,a,0,0,0); 180e6ed2b1fSLois Curfman McInnes } 181e6ed2b1fSLois Curfman McInnes #endif 182d2dddef7SLois Curfman McInnes 183d2dddef7SLois Curfman McInnes /* Safeguard for step sizes too small */ 184d2dddef7SLois Curfman McInnes if (sum == 0.0) {dot = 1.0; norm = 1.0;} 185db0b4b35SLois Curfman McInnes #if defined(USE_PETSC_COMPLEX) 186d2dddef7SLois Curfman McInnes else if (abs(dot) < umin*sum && real(dot) >= 0.0) dot = umin*sum; 187d2dddef7SLois Curfman McInnes else if (abs(dot) < 0.0 && real(dot) > -umin*sum) dot = -umin*sum; 188d2dddef7SLois Curfman McInnes #else 189d2dddef7SLois Curfman McInnes else if (dot < umin*sum && dot >= 0.0) dot = umin*sum; 190d2dddef7SLois Curfman McInnes else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum; 191d2dddef7SLois Curfman McInnes #endif 192db0b4b35SLois Curfman McInnes h = PetscReal(ctx->error_rel*dot/(norm*norm)); 193d2dddef7SLois Curfman McInnes } 194d2dddef7SLois Curfman McInnes } else { 195d2dddef7SLois Curfman McInnes h = ctx->h; 196d2dddef7SLois Curfman McInnes } 197d2dddef7SLois Curfman McInnes if (!ctx->jorge || !ctx->need_h) PLogInfo(snes,"SNESMatrixFreeMult2_Private: h = %g\n",h); 198d2dddef7SLois Curfman McInnes 199d2dddef7SLois Curfman McInnes /* Evaluate function at F(u + ha) */ 200db0b4b35SLois Curfman McInnes hs = h; 201db0b4b35SLois Curfman McInnes ierr = VecWAXPY(&hs,a,U,w); CHKERRQ(ierr); 202d2dddef7SLois Curfman McInnes ierr = eval_fct(snes,w,y); CHKERRQ(ierr); 203d2dddef7SLois Curfman McInnes ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr); 204db0b4b35SLois Curfman McInnes hs = 1.0/hs; 205db0b4b35SLois Curfman McInnes ierr = VecScale(&hs,y); CHKERRQ(ierr); 206d2dddef7SLois Curfman McInnes if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);} 207d2dddef7SLois Curfman McInnes 208d2dddef7SLois Curfman McInnes PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0); 209d2dddef7SLois Curfman McInnes return 0; 210d2dddef7SLois Curfman McInnes } 211d2dddef7SLois Curfman McInnes 212d2dddef7SLois Curfman McInnes #undef __FUNC__ 213d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESDefaultMatrixFreeMatCreate2" 214d2dddef7SLois Curfman McInnes /*@C 215d2dddef7SLois Curfman McInnes SNESMatrixFreeMatCreate2 - Creates a matrix-free matrix 216d2dddef7SLois Curfman McInnes context for use with a SNES solver. This matrix can be used as 217d2dddef7SLois Curfman McInnes the Jacobian argument for the routine SNESSetJacobian(). 218d2dddef7SLois Curfman McInnes 219d2dddef7SLois Curfman McInnes Input Parameters: 220d2dddef7SLois Curfman McInnes . snes - the SNES context 221d2dddef7SLois Curfman McInnes . x - vector where SNES solution is to be stored. 222d2dddef7SLois Curfman McInnes 223d2dddef7SLois Curfman McInnes Output Parameter: 224d2dddef7SLois Curfman McInnes . J - the matrix-free matrix 225d2dddef7SLois Curfman McInnes 226d2dddef7SLois Curfman McInnes Notes: 227d2dddef7SLois Curfman McInnes The matrix-free matrix context merely contains the function pointers 228d2dddef7SLois Curfman McInnes and work space for performing finite difference approximations of 229d2dddef7SLois Curfman McInnes Jacobian-vector products, J(u)*a, via 230d2dddef7SLois Curfman McInnes 231d2dddef7SLois Curfman McInnes $ J(u)*a = [J(u+h*a) - J(u)]/h, 232d2dddef7SLois Curfman McInnes $ where by default 233d2dddef7SLois Curfman McInnes $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 234d2dddef7SLois Curfman McInnes $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 235d2dddef7SLois Curfman McInnes $ where 236d2dddef7SLois Curfman McInnes $ error_rel = square root of relative error in 237d2dddef7SLois Curfman McInnes $ function evaluation 238d2dddef7SLois Curfman McInnes $ umin = minimum iterate parameter 239d2dddef7SLois Curfman McInnes $ Alternatively, the differencing parameter, h, can be set using 240e6ed2b1fSLois Curfman McInnes $ Jorge's nifty new strategy if one specifies the option 241e6ed2b1fSLois Curfman McInnes $ -snes_mf_jorge 242d2dddef7SLois Curfman McInnes 243d2dddef7SLois Curfman McInnes The user can set these parameters via SNESSetMatrixFreeParameters(). 244d2dddef7SLois Curfman McInnes See the nonlinear solvers chapter of the users manual for details. 245d2dddef7SLois Curfman McInnes 246d2dddef7SLois Curfman McInnes The user should call MatDestroy() when finished with the matrix-free 247d2dddef7SLois Curfman McInnes matrix context. 248d2dddef7SLois Curfman McInnes 249d2dddef7SLois Curfman McInnes Options Database Keys: 250d2dddef7SLois Curfman McInnes $ -snes_mf_err <error_rel> 251d2dddef7SLois Curfman McInnes $ -snes_mf_unim <umin> 252295f7fbbSLois Curfman McInnes $ -snes_mf_compute_err 253295f7fbbSLois Curfman McInnes $ -snes_mf_freq_err <freq> 254e6ed2b1fSLois Curfman McInnes $ -snes_mf_jorge 255d2dddef7SLois Curfman McInnes 256d2dddef7SLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix 257d2dddef7SLois Curfman McInnes 258d2dddef7SLois Curfman McInnes .seealso: MatDestroy(), SNESSetMatrixFreeParameters() 259d2dddef7SLois Curfman McInnes @*/ 260d2dddef7SLois Curfman McInnes int SNESMatrixFreeMatCreate2(SNES snes,Vec x, Mat *J) 261d2dddef7SLois Curfman McInnes { 262d2dddef7SLois Curfman McInnes MPI_Comm comm; 263d2dddef7SLois Curfman McInnes MFCtx_Private *mfctx; 264d2dddef7SLois Curfman McInnes int n, nloc, ierr, flg; 265d2dddef7SLois Curfman McInnes char p[64]; 266d2dddef7SLois Curfman McInnes 267d2dddef7SLois Curfman McInnes mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx); 268d2dddef7SLois Curfman McInnes PLogObjectMemory(snes,sizeof(MFCtx_Private)); 269d2dddef7SLois Curfman McInnes mfctx->sp = 0; 270d2dddef7SLois Curfman McInnes mfctx->snes = snes; 271d2dddef7SLois Curfman McInnes mfctx->error_rel = 1.e-8; /* assumes double precision */ 272d2dddef7SLois Curfman McInnes mfctx->umin = 1.e-6; 273d2dddef7SLois Curfman McInnes mfctx->h = 0.0; 274d2dddef7SLois Curfman McInnes mfctx->need_h = 1; 275d2dddef7SLois Curfman McInnes mfctx->need_err = 0; 276295f7fbbSLois Curfman McInnes mfctx->compute_err = 0; 277295f7fbbSLois Curfman McInnes mfctx->compute_err_freq = 0; 278295f7fbbSLois Curfman McInnes mfctx->compute_err_iter = -1; 279d2dddef7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr); 280d2dddef7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr); 281e6ed2b1fSLois Curfman McInnes ierr = OptionsHasName(snes->prefix,"-snes_mf_jorge",&mfctx->jorge); CHKERRQ(ierr); 282295f7fbbSLois Curfman McInnes ierr = OptionsHasName(snes->prefix,"-snes_mf_compute_err",&mfctx->compute_err); CHKERRQ(ierr); 283295f7fbbSLois Curfman McInnes ierr = OptionsGetInt(snes->prefix,"-snes_mf_freq_err",&mfctx->compute_err_freq,&flg); CHKERRQ(ierr); 284295f7fbbSLois Curfman McInnes if (flg) { 285295f7fbbSLois Curfman McInnes if (mfctx->compute_err_freq < 0) mfctx->compute_err_freq = 0; 286295f7fbbSLois Curfman McInnes mfctx->compute_err = 1; 287295f7fbbSLois Curfman McInnes } 288295f7fbbSLois Curfman McInnes if (mfctx->compute_err == 1) mfctx->need_err = 1; 289295f7fbbSLois Curfman McInnes if (mfctx->jorge || mfctx->compute_err) { 290d2dddef7SLois Curfman McInnes ierr = DiffParameterCreate_More(snes,x,&mfctx->data); CHKERRQ(ierr); 291d2dddef7SLois Curfman McInnes } else mfctx->data = 0; 292d2dddef7SLois Curfman McInnes 293d2dddef7SLois Curfman McInnes ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 294d2dddef7SLois Curfman McInnes PetscStrcpy(p,"-"); 295d2dddef7SLois Curfman McInnes if (snes->prefix) PetscStrcat(p,snes->prefix); 296d2dddef7SLois Curfman McInnes if (flg) { 297295f7fbbSLois Curfman McInnes PetscPrintf(snes->comm," Matrix-free Options (via SNES):\n"); 298bd55573bSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_err <err>: set sqrt of relative error in function (default %g)\n",p,mfctx->error_rel); 299d2dddef7SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_umin <umin>: see users manual (default %g)\n",p,mfctx->umin); 300295f7fbbSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_jorge: use Jorge More's method\n",p); 301295f7fbbSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_compute_err: compute sqrt or relative error in function\n",p); 302295f7fbbSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_freq_err <freq>: frequency to recompute this (default only once)\n",p); 303295f7fbbSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_noise_file <file>: set file for printing noise info\n",p); 304d2dddef7SLois Curfman McInnes } 305d2dddef7SLois Curfman McInnes ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr); 306d2dddef7SLois Curfman McInnes ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr); 307d2dddef7SLois Curfman McInnes ierr = VecGetSize(x,&n); CHKERRQ(ierr); 308d2dddef7SLois Curfman McInnes ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr); 309d2dddef7SLois Curfman McInnes ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr); 310d2dddef7SLois Curfman McInnes ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult2_Private);CHKERRQ(ierr); 311d2dddef7SLois Curfman McInnes ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy2_Private);CHKERRQ(ierr); 312d2dddef7SLois Curfman McInnes ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView2_Private); CHKERRQ(ierr); 313d2dddef7SLois Curfman McInnes 314d2dddef7SLois Curfman McInnes PLogObjectParent(*J,mfctx->w); 315d2dddef7SLois Curfman McInnes PLogObjectParent(snes,*J); 316d2dddef7SLois Curfman McInnes return 0; 317d2dddef7SLois Curfman McInnes } 318d2dddef7SLois Curfman McInnes 319d2dddef7SLois Curfman McInnes #undef __FUNC__ 320d2dddef7SLois Curfman McInnes #define __FUNC__ "SNESSetMatrixFreeParameters2" 321d2dddef7SLois Curfman McInnes /*@ 322d2dddef7SLois Curfman McInnes SNESSetMatrixFreeParameters2 - Sets the parameters for the approximation of 323d2dddef7SLois Curfman McInnes matrix-vector products using finite differences. 324d2dddef7SLois Curfman McInnes 325d2dddef7SLois Curfman McInnes $ J(u)*a = [J(u+h*a) - J(u)]/h where 326d2dddef7SLois Curfman McInnes 327d2dddef7SLois Curfman McInnes either the user sets h directly here, or this parameter is computed via 328d2dddef7SLois Curfman McInnes 329d2dddef7SLois Curfman McInnes $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 330d2dddef7SLois Curfman McInnes $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 331d2dddef7SLois Curfman McInnes $ 332d2dddef7SLois Curfman McInnes 333d2dddef7SLois Curfman McInnes Input Parameters: 334d2dddef7SLois Curfman McInnes . snes - the SNES context 335d2dddef7SLois Curfman McInnes . error_rel - relative error (should be set to the square root of 336d2dddef7SLois Curfman McInnes the relative error in the function evaluations) 337d2dddef7SLois Curfman McInnes . umin - minimum allowable u-value 338d2dddef7SLois Curfman McInnes . h - differencing parameter 339d2dddef7SLois Curfman McInnes 340d2dddef7SLois Curfman McInnes Notes: 341d2dddef7SLois Curfman McInnes If the user sets the parameter h directly, then this value will be used 342d2dddef7SLois Curfman McInnes instead of the default computation indicated above. 343d2dddef7SLois Curfman McInnes 344d2dddef7SLois Curfman McInnes .keywords: SNES, matrix-free, parameters 345d2dddef7SLois Curfman McInnes 346d2dddef7SLois Curfman McInnes .seealso: SNESDefaultMatrixFreeMatCreate() 347d2dddef7SLois Curfman McInnes @*/ 348d2dddef7SLois Curfman McInnes int SNESSetMatrixFreeParameters2(SNES snes,double error,double umin,double h) 349d2dddef7SLois Curfman McInnes { 350d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 351d2dddef7SLois Curfman McInnes int ierr; 352d2dddef7SLois Curfman McInnes Mat mat; 353d2dddef7SLois Curfman McInnes 354d2dddef7SLois Curfman McInnes ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 355d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 356d2dddef7SLois Curfman McInnes if (ctx) { 357d2dddef7SLois Curfman McInnes if (error != PETSC_DEFAULT) ctx->error_rel = error; 358d2dddef7SLois Curfman McInnes if (umin != PETSC_DEFAULT) ctx->umin = umin; 359d2dddef7SLois Curfman McInnes if (h != PETSC_DEFAULT) { 360d2dddef7SLois Curfman McInnes ctx->h = h; 361d2dddef7SLois Curfman McInnes ctx->need_h = 0; 362d2dddef7SLois Curfman McInnes } 363d2dddef7SLois Curfman McInnes } 364d2dddef7SLois Curfman McInnes return 0; 365d2dddef7SLois Curfman McInnes } 366d2dddef7SLois Curfman McInnes 367d2dddef7SLois Curfman McInnes int SNESUnSetMatrixFreeParameter(SNES snes) 368d2dddef7SLois Curfman McInnes { 369d2dddef7SLois Curfman McInnes MFCtx_Private *ctx; 370d2dddef7SLois Curfman McInnes int ierr; 371d2dddef7SLois Curfman McInnes Mat mat; 372d2dddef7SLois Curfman McInnes 373d2dddef7SLois Curfman McInnes ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 374d2dddef7SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 375d2dddef7SLois Curfman McInnes if (ctx) ctx->need_h = 1; 376d2dddef7SLois Curfman McInnes return 0; 377d2dddef7SLois Curfman McInnes } 378d2dddef7SLois Curfman McInnes 379