1 2 #ifndef lint 3 static char vcid[] = "$Id: snesmfj.c,v 1.50 1997/06/12 22:52:40 bsmith Exp bsmith $"; 4 #endif 5 6 #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 7 #include "pinclude/pviewer.h" 8 #include <math.h> 9 10 typedef struct { /* default context for matrix-free SNES */ 11 SNES snes; /* SNES context */ 12 Vec w; /* work vector */ 13 PCNullSpace sp; /* null space context */ 14 double error_rel; /* square root of relative error in computing function */ 15 double umin; /* minimum allowable u'a value relative to |u|_1 */ 16 } MFCtx_Private; 17 18 #undef __FUNC__ 19 #define __FUNC__ "SNESMatrixFreeDestroy_Private" /* ADIC Ignore */ 20 int SNESMatrixFreeDestroy_Private(PetscObject obj) 21 { 22 int ierr; 23 Mat mat = (Mat) obj; 24 MFCtx_Private *ctx; 25 26 ierr = MatShellGetContext(mat,(void **)&ctx); 27 ierr = VecDestroy(ctx->w); CHKERRQ(ierr); 28 if (ctx->sp) {ierr = PCNullSpaceDestroy(ctx->sp); CHKERRQ(ierr);} 29 PetscFree(ctx); 30 return 0; 31 } 32 33 #undef __FUNC__ 34 #define __FUNC__ "SNESMatrixFreeView_Private" /* ADIC Ignore */ 35 /* 36 SNESMatrixFreeView_Private - Views matrix-free parameters. 37 */ 38 int SNESMatrixFreeView_Private(Mat J,Viewer viewer) 39 { 40 int ierr; 41 MFCtx_Private *ctx; 42 MPI_Comm comm; 43 FILE *fd; 44 ViewerType vtype; 45 46 PetscObjectGetComm((PetscObject)J,&comm); 47 ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 48 ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 49 ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 50 if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 51 PetscFPrintf(comm,fd," SNES matrix-free approximation:\n"); 52 PetscFPrintf(comm,fd," err=%g (relative error in function evaluation)\n",ctx->error_rel); 53 PetscFPrintf(comm,fd," umin=%g (minimum iterate parameter)\n",ctx->umin); 54 } 55 return 0; 56 } 57 58 extern int VecDot_Seq(Vec,Vec,Scalar *); 59 extern int VecNorm_Seq(Vec,NormType,double *); 60 61 #undef __FUNC__ 62 #define __FUNC__ "SNESMatrixFreeMult_Private" 63 /* 64 SNESMatrixFreeMult_Private - Default matrix-free form for Jacobian-vector 65 product, y = F'(u)*a: 66 y = ( F(u + ha) - F(u) ) /h, 67 where F = nonlinear function, as set by SNESSetFunction() 68 u = current iterate 69 h = difference interval 70 */ 71 int SNESMatrixFreeMult_Private(Mat mat,Vec a,Vec y) 72 { 73 MFCtx_Private *ctx; 74 SNES snes; 75 double ovalues[3],values[3],norm, sum, umin; 76 Scalar h, dot, mone = -1.0; 77 Vec w,U,F; 78 int ierr, (*eval_fct)(SNES,Vec,Vec); 79 MPI_Comm comm; 80 81 PetscObjectGetComm((PetscObject)mat,&comm); 82 ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 83 snes = ctx->snes; 84 w = ctx->w; 85 umin = ctx->umin; 86 87 /* We log matrix-free matrix-vector products separately, so that we can 88 separate the performance monitoring from the cases that use conventional 89 storage. We may eventually modify event logging to associate events 90 with particular objects, hence alleviating the more general problem. */ 91 92 ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr); 93 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 94 eval_fct = SNESComputeFunction; 95 ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr); 96 } 97 else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 98 eval_fct = SNESComputeGradient; 99 ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr); 100 } 101 else SETERRQ(1,0,"Invalid method class"); 102 103 /* Determine a "good" step size, h */ 104 105 /* 106 ierr = VecDot(U,a,&dot); CHKERRQ(ierr); 107 ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr); 108 ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr); 109 */ 110 111 /* 112 Call the Seq Vector routines and then do 113 a single reduction to reduce the number of communications 114 required 115 */ 116 117 PLogEventBegin(VEC_Dot,U,a,0,0); 118 ierr = VecDot_Seq(U,a,ovalues); CHKERRQ(ierr); 119 PLogEventEnd(VEC_Dot,U,a,0,0); 120 PLogEventBegin(VEC_Norm,a,0,0,0); 121 ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr); 122 ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr); 123 ovalues[2] = ovalues[2]*ovalues[2]; 124 MPI_Allreduce(ovalues,values,3,MPI_DOUBLE,MPI_SUM,comm ); 125 dot = values[0]; sum = values[1]; norm = sqrt(values[2]); 126 PLogEventEnd(VEC_Norm,a,0,0,0); 127 128 /* 129 the plogeventbegin() below should really be above, 130 but they cannot be nested so it excludes the time 131 to compute h 132 */ 133 PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0); 134 135 /* Safeguard for step sizes too small */ 136 if (sum == 0.0) {dot = 1.0; norm = 1.0;} 137 #if defined(PETSC_COMPLEX) 138 else if (abs(dot) < umin*sum && real(dot) >= 0.0) dot = umin*sum; 139 else if (abs(dot) < 0.0 && real(dot) > -umin*sum) dot = -umin*sum; 140 #else 141 else if (dot < umin*sum && dot >= 0.0) dot = umin*sum; 142 else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum; 143 #endif 144 h = ctx->error_rel*dot/(norm*norm); 145 146 /* Evaluate function at F(u + ha) */ 147 ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr); 148 ierr = eval_fct(snes,w,y); CHKERRQ(ierr); 149 ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr); 150 h = 1.0/h; 151 ierr = VecScale(&h,y); CHKERRQ(ierr); 152 if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);} 153 154 PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0); 155 return 0; 156 } 157 158 #undef __FUNC__ 159 #define __FUNC__ "SNESDefaultMatrixFreeMatCreate" 160 /*@C 161 SNESDefaultMatrixFreeMatCreate - Creates a matrix-free matrix 162 context for use with a SNES solver. This matrix can be used as 163 the Jacobian argument for the routine SNESSetJacobian(). 164 165 Input Parameters: 166 . snes - the SNES context 167 . x - vector where SNES solution is to be stored. 168 169 Output Parameter: 170 . J - the matrix-free matrix 171 172 Notes: 173 The matrix-free matrix context merely contains the function pointers 174 and work space for performing finite difference approximations of 175 Jacobian-vector products, J(u)*a, via 176 177 $ J(u)*a = [J(u+h*a) - J(u)]/h where 178 $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 179 $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 180 $ where 181 $ error_rel = square root of relative error in 182 $ function evaluation 183 $ umin = minimum iterate parameter 184 185 The user can set these parameters via SNESSetMatrixFreeParameters(). 186 See the nonlinear solvers chapter of the users manual for details. 187 188 The user should call MatDestroy() when finished with the matrix-free 189 matrix context. 190 191 Options Database Keys: 192 $ -snes_mf_err <error_rel> 193 $ -snes_mf_unim <umin> 194 195 .keywords: SNES, default, matrix-free, create, matrix 196 197 .seealso: MatDestroy(), SNESSetMatrixFreeParameters() 198 @*/ 199 int SNESDefaultMatrixFreeMatCreate(SNES snes,Vec x, Mat *J) 200 { 201 MPI_Comm comm; 202 MFCtx_Private *mfctx; 203 int n, nloc, ierr, flg; 204 char p[64]; 205 206 mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx); 207 PLogObjectMemory(snes,sizeof(MFCtx_Private)); 208 mfctx->sp = 0; 209 mfctx->snes = snes; 210 mfctx->error_rel = 1.e-8; /* assumes double precision */ 211 mfctx->umin = 1.e-6; 212 ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr); 213 ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr); 214 ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 215 PetscStrcpy(p,"-"); 216 if (snes->prefix) PetscStrcat(p,snes->prefix); 217 if (flg) { 218 PetscPrintf(snes->comm," %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel); 219 PetscPrintf(snes->comm," %ssnes_mf_umin <umin> see users manual (default %g)\n",p,mfctx->umin); 220 } 221 ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr); 222 ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr); 223 ierr = VecGetSize(x,&n); CHKERRQ(ierr); 224 ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr); 225 ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr); 226 ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult_Private);CHKERRQ(ierr); 227 ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy_Private);CHKERRQ(ierr); 228 ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView_Private); CHKERRQ(ierr); 229 PLogObjectParent(*J,mfctx->w); 230 PLogObjectParent(snes,*J); 231 return 0; 232 } 233 234 #undef __FUNC__ 235 #define __FUNC__ "SNESSetMatrixFreeParameters" 236 /*@ 237 SNESSetMatrixFreeParameters - Sets the parameters for the approximation of 238 matrix-vector products using finite differences. 239 240 $ J(u)*a = [J(u+h*a) - J(u)]/h where 241 $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 242 $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 243 $ 244 Input Parameters: 245 . snes - the SNES context 246 . error_rel - relative error (should be set to the square root of 247 the relative error in the function evaluations) 248 . umin - minimum allowable u-value 249 250 .keywords: SNES, matrix-free, parameters 251 252 .seealso: SNESDefaultMatrixFreeMatCreate() 253 @*/ 254 int SNESSetMatrixFreeParameters(SNES snes,double error,double umin) 255 { 256 MFCtx_Private *ctx; 257 int ierr; 258 Mat mat; 259 260 ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 261 ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 262 if (ctx) { 263 if (error != PETSC_DEFAULT) ctx->error_rel = error; 264 if (umin != PETSC_DEFAULT) ctx->umin = umin; 265 } 266 return 0; 267 } 268 269 #undef __FUNC__ 270 #define __FUNC__ "SNESDefaultMatrixFreeMatAddNullSpace" 271 /*@ 272 SNESDefaultMatrixFreeMatAddNullSpace - Provides a null space that 273 an operator is supposed to have. Since roundoff will create a 274 small component in the null space, if you know the null space 275 you may have it automatically removed. 276 277 Input Parameters: 278 . J - the matrix-free matrix context 279 . has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants 280 . n - number of vectors (excluding constant vector) in null space 281 . vecs - the vectors that span the null space (excluding the constant vector); 282 . these vectors must be orthonormal 283 284 .keywords: SNES, matrix-free, null space 285 @*/ 286 int SNESDefaultMatrixFreeMatAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs) 287 { 288 int ierr; 289 MFCtx_Private *ctx; 290 MPI_Comm comm; 291 292 PetscObjectGetComm((PetscObject)J,&comm); 293 294 ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 295 /* no context indicates that it is not the "matrix free" matrix type */ 296 if (!ctx) return 0; 297 ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr); 298 return 0; 299 } 300 301 302 303 304