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