1 2 #ifndef lint 3 static char vcid[] = "$Id: snesmfj.c,v 1.51 1997/06/15 14:07:37 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 a single reduction 113 to reduce the number of communications required 114 */ 115 116 #if !defined(PETSC_COMPLEX) 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 #else 128 { 129 Scalar cvalues[3],covalues[3]; 130 131 PLogEventBegin(VEC_Dot,U,a,0,0); 132 ierr = VecDot_Seq(U,a,covalues); CHKERRQ(ierr); 133 PLogEventEnd(VEC_Dot,U,a,0,0); 134 PLogEventBegin(VEC_Norm,a,0,0,0); 135 ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr); 136 ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr); 137 covalues[1] = ovalues[1]; 138 covalues[2] = ovalues[2]*ovalues[2]; 139 MPI_Allreduce(covalues,cvalues,6,MPI_DOUBLE,MPI_SUM,comm ); 140 dot = cvalues[0]; sum = PetscReal(cvalues[1]); norm = sqrt(PetscReal(cvalues[2])); 141 PLogEventEnd(VEC_Norm,a,0,0,0); 142 } 143 #endif 144 145 /* 146 the plogeventbegin() below should really be above, 147 but they cannot be nested so it excludes the time 148 to compute h 149 */ 150 PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0); 151 152 /* Safeguard for step sizes too small */ 153 if (sum == 0.0) {dot = 1.0; norm = 1.0;} 154 #if defined(PETSC_COMPLEX) 155 else if (abs(dot) < umin*sum && real(dot) >= 0.0) dot = umin*sum; 156 else if (abs(dot) < 0.0 && real(dot) > -umin*sum) dot = -umin*sum; 157 #else 158 else if (dot < umin*sum && dot >= 0.0) dot = umin*sum; 159 else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum; 160 #endif 161 h = ctx->error_rel*dot/(norm*norm); 162 163 /* Evaluate function at F(u + ha) */ 164 ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr); 165 ierr = eval_fct(snes,w,y); CHKERRQ(ierr); 166 ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr); 167 h = 1.0/h; 168 ierr = VecScale(&h,y); CHKERRQ(ierr); 169 if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);} 170 171 PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0); 172 return 0; 173 } 174 175 #undef __FUNC__ 176 #define __FUNC__ "SNESDefaultMatrixFreeMatCreate" 177 /*@C 178 SNESDefaultMatrixFreeMatCreate - Creates a matrix-free matrix 179 context for use with a SNES solver. This matrix can be used as 180 the Jacobian argument for the routine SNESSetJacobian(). 181 182 Input Parameters: 183 . snes - the SNES context 184 . x - vector where SNES solution is to be stored. 185 186 Output Parameter: 187 . J - the matrix-free matrix 188 189 Notes: 190 The matrix-free matrix context merely contains the function pointers 191 and work space for performing finite difference approximations of 192 Jacobian-vector products, J(u)*a, via 193 194 $ J(u)*a = [J(u+h*a) - J(u)]/h where 195 $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 196 $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 197 $ where 198 $ error_rel = square root of relative error in 199 $ function evaluation 200 $ umin = minimum iterate parameter 201 202 The user can set these parameters via SNESSetMatrixFreeParameters(). 203 See the nonlinear solvers chapter of the users manual for details. 204 205 The user should call MatDestroy() when finished with the matrix-free 206 matrix context. 207 208 Options Database Keys: 209 $ -snes_mf_err <error_rel> 210 $ -snes_mf_unim <umin> 211 212 .keywords: SNES, default, matrix-free, create, matrix 213 214 .seealso: MatDestroy(), SNESSetMatrixFreeParameters() 215 @*/ 216 int SNESDefaultMatrixFreeMatCreate(SNES snes,Vec x, Mat *J) 217 { 218 MPI_Comm comm; 219 MFCtx_Private *mfctx; 220 int n, nloc, ierr, flg; 221 char p[64]; 222 223 mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx); 224 PLogObjectMemory(snes,sizeof(MFCtx_Private)); 225 mfctx->sp = 0; 226 mfctx->snes = snes; 227 mfctx->error_rel = 1.e-8; /* assumes double precision */ 228 mfctx->umin = 1.e-6; 229 ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr); 230 ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr); 231 ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 232 PetscStrcpy(p,"-"); 233 if (snes->prefix) PetscStrcat(p,snes->prefix); 234 if (flg) { 235 PetscPrintf(snes->comm," %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel); 236 PetscPrintf(snes->comm," %ssnes_mf_umin <umin> see users manual (default %g)\n",p,mfctx->umin); 237 } 238 ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr); 239 ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr); 240 ierr = VecGetSize(x,&n); CHKERRQ(ierr); 241 ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr); 242 ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr); 243 ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult_Private);CHKERRQ(ierr); 244 ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy_Private);CHKERRQ(ierr); 245 ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView_Private); CHKERRQ(ierr); 246 PLogObjectParent(*J,mfctx->w); 247 PLogObjectParent(snes,*J); 248 return 0; 249 } 250 251 #undef __FUNC__ 252 #define __FUNC__ "SNESSetMatrixFreeParameters" 253 /*@ 254 SNESSetMatrixFreeParameters - Sets the parameters for the approximation of 255 matrix-vector products using finite differences. 256 257 $ J(u)*a = [J(u+h*a) - J(u)]/h where 258 $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 259 $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 260 $ 261 Input Parameters: 262 . snes - the SNES context 263 . error_rel - relative error (should be set to the square root of 264 the relative error in the function evaluations) 265 . umin - minimum allowable u-value 266 267 .keywords: SNES, matrix-free, parameters 268 269 .seealso: SNESDefaultMatrixFreeMatCreate() 270 @*/ 271 int SNESSetMatrixFreeParameters(SNES snes,double error,double umin) 272 { 273 MFCtx_Private *ctx; 274 int ierr; 275 Mat mat; 276 277 ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 278 ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 279 if (ctx) { 280 if (error != PETSC_DEFAULT) ctx->error_rel = error; 281 if (umin != PETSC_DEFAULT) ctx->umin = umin; 282 } 283 return 0; 284 } 285 286 #undef __FUNC__ 287 #define __FUNC__ "SNESDefaultMatrixFreeMatAddNullSpace" 288 /*@ 289 SNESDefaultMatrixFreeMatAddNullSpace - Provides a null space that 290 an operator is supposed to have. Since roundoff will create a 291 small component in the null space, if you know the null space 292 you may have it automatically removed. 293 294 Input Parameters: 295 . J - the matrix-free matrix context 296 . has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants 297 . n - number of vectors (excluding constant vector) in null space 298 . vecs - the vectors that span the null space (excluding the constant vector); 299 . these vectors must be orthonormal 300 301 .keywords: SNES, matrix-free, null space 302 @*/ 303 int SNESDefaultMatrixFreeMatAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs) 304 { 305 int ierr; 306 MFCtx_Private *ctx; 307 MPI_Comm comm; 308 309 PetscObjectGetComm((PetscObject)J,&comm); 310 311 ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 312 /* no context indicates that it is not the "matrix free" matrix type */ 313 if (!ctx) return 0; 314 ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr); 315 return 0; 316 } 317 318 319 320 321