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