1 #ifdef PETSC_RCS_HEADER 2 static char vcid[] = "$Id: snesmfj.c,v 1.69 1998/10/26 03:26:40 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 Scalar currenth; /* last differencing parameter used */ 16 Scalar *historyh; /* history of h */ 17 int ncurrenth,maxcurrenth; 18 } MFCtx_Private; 19 20 #undef __FUNC__ 21 #define __FUNC__ "SNESMatrixFreeDestroy_Private" 22 int SNESMatrixFreeDestroy_Private(Mat mat) 23 { 24 int ierr; 25 MFCtx_Private *ctx; 26 27 PetscFunctionBegin; 28 ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 29 ierr = VecDestroy(ctx->w); CHKERRQ(ierr); 30 if (ctx->sp) {ierr = PCNullSpaceDestroy(ctx->sp);CHKERRQ(ierr);} 31 PetscFree(ctx); 32 PetscFunctionReturn(0); 33 } 34 35 #undef __FUNC__ 36 #define __FUNC__ "SNESMatrixFreeView_Private" 37 /* 38 SNESMatrixFreeView_Private - Views matrix-free parameters. 39 40 */ 41 int SNESMatrixFreeView_Private(Mat J,Viewer viewer) 42 { 43 int ierr; 44 MFCtx_Private *ctx; 45 MPI_Comm comm; 46 FILE *fd; 47 ViewerType vtype; 48 49 PetscFunctionBegin; 50 PetscObjectGetComm((PetscObject)J,&comm); 51 ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 52 ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 53 ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 54 if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 55 PetscFPrintf(comm,fd," SNES matrix-free approximation:\n"); 56 PetscFPrintf(comm,fd," err=%g (relative error in function evaluation)\n",ctx->error_rel); 57 PetscFPrintf(comm,fd," umin=%g (minimum iterate parameter)\n",ctx->umin); 58 } else { 59 SETERRQ(1,1,"Viewer type not supported for this object"); 60 } 61 PetscFunctionReturn(0); 62 } 63 64 extern int VecDot_Seq(Vec,Vec,Scalar *); 65 extern int VecNorm_Seq(Vec,NormType,double *); 66 67 #undef __FUNC__ 68 #define __FUNC__ "SNESMatrixFreeMult_Private" 69 /* 70 SNESMatrixFreeMult_Private - Default matrix-free form for Jacobian-vector 71 product, y = F'(u)*a: 72 y = ( F(u + ha) - F(u) ) /h, 73 where F = nonlinear function, as set by SNESSetFunction() 74 u = current iterate 75 h = difference interval 76 */ 77 int SNESMatrixFreeMult_Private(Mat mat,Vec a,Vec y) 78 { 79 MFCtx_Private *ctx; 80 SNES snes; 81 double ovalues[3],norm, sum, umin; 82 Scalar h, dot, mone = -1.0; 83 Vec w,U,F; 84 int ierr, (*eval_fct)(SNES,Vec,Vec); 85 MPI_Comm comm; 86 #if !defined(USE_PETSC_COMPLEX) 87 double values[3]; 88 #endif 89 90 PetscFunctionBegin; 91 PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0); 92 93 PetscObjectGetComm((PetscObject)mat,&comm); 94 ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 95 snes = ctx->snes; 96 w = ctx->w; 97 umin = ctx->umin; 98 99 /* We log matrix-free matrix-vector products separately, so that we can 100 separate the performance monitoring from the cases that use conventional 101 storage. We may eventually modify event logging to associate events 102 with particular objects, hence alleviating the more general problem. */ 103 104 ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr); 105 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 106 eval_fct = SNESComputeFunction; 107 ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr); 108 } 109 else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 110 eval_fct = SNESComputeGradient; 111 ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr); 112 } 113 else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid method class"); 114 115 /* Determine a "good" step size, h */ 116 117 /* 118 ierr = VecDot(U,a,&dot); CHKERRQ(ierr); 119 ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr); 120 ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr); 121 */ 122 123 /* 124 Call the Seq Vector routines and then do a single reduction 125 to reduce the number of communications required 126 */ 127 128 #if !defined(USE_PETSC_COMPLEX) 129 PLogEventBegin(VEC_Dot,U,a,0,0); 130 ierr = VecDot_Seq(U,a,ovalues); CHKERRQ(ierr); 131 PLogEventEnd(VEC_Dot,U,a,0,0); 132 PLogEventBegin(VEC_Norm,a,0,0,0); 133 ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr); 134 ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr); 135 ovalues[2] = ovalues[2]*ovalues[2]; 136 PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm); 137 ierr = MPI_Allreduce(ovalues,values,3,MPI_DOUBLE,MPI_SUM,comm );CHKERRQ(ierr); 138 PLogEventBarrierEnd(VEC_NormBarrier,0,0,0,0,comm); 139 dot = values[0]; sum = values[1]; norm = sqrt(values[2]); 140 PLogEventEnd(VEC_Norm,a,0,0,0); 141 #else 142 { 143 Scalar cvalues[3],covalues[3]; 144 145 PLogEventBegin(VEC_Dot,U,a,0,0); 146 ierr = VecDot_Seq(U,a,covalues); CHKERRQ(ierr); 147 PLogEventEnd(VEC_Dot,U,a,0,0); 148 PLogEventBegin(VEC_Norm,a,0,0,0); 149 ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr); 150 ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr); 151 covalues[1] = ovalues[1]; 152 covalues[2] = ovalues[2]*ovalues[2]; 153 PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm); 154 ierr = MPI_Allreduce(covalues,cvalues,6,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 155 PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm); 156 dot = cvalues[0]; sum = PetscReal(cvalues[1]); norm = sqrt(PetscReal(cvalues[2])); 157 PLogEventEnd(VEC_Norm,a,0,0,0); 158 } 159 #endif 160 161 162 /* Safeguard for step sizes too small */ 163 if (sum == 0.0) {dot = 1.0; norm = 1.0;} 164 #if defined(USE_PETSC_COMPLEX) 165 else if (PetscAbsScalar(dot) < umin*sum && PetscReal(dot) >= 0.0) dot = umin*sum; 166 else if (PetscAbsScalar(dot) < 0.0 && PetscReal(dot) > -umin*sum) dot = -umin*sum; 167 #else 168 else if (dot < umin*sum && dot >= 0.0) dot = umin*sum; 169 else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum; 170 #endif 171 h = ctx->error_rel*dot/(norm*norm); 172 173 /* keep a record of the current differencing parameter h */ 174 ctx->currenth = h; 175 #if defined(USE_PETSC_COMPLEX) 176 PLogInfo(mat,"Current differencing parameter: %g + %g i\n",PetscReal(h),PetscImaginary(h)); 177 #else 178 PLogInfo(mat,"Current differencing parameter: %g\n",h); 179 #endif 180 if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 181 ctx->historyh[ctx->ncurrenth++] = h; 182 } 183 184 /* Evaluate function at F(u + ha) */ 185 ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr); 186 ierr = eval_fct(snes,w,y); CHKERRQ(ierr); 187 188 ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr); 189 h = 1.0/h; 190 ierr = VecScale(&h,y); CHKERRQ(ierr); 191 if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);} 192 193 PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0); 194 PetscFunctionReturn(0); 195 } 196 197 #undef __FUNC__ 198 #define __FUNC__ "SNESDefaultMatrixFreeCreate" 199 /*@C 200 SNESDefaultMatrixFreeCreate - Creates a matrix-free matrix 201 context for use with a SNES solver. This matrix can be used as 202 the Jacobian argument for the routine SNESSetJacobian(). 203 204 Collective on SNES and Vec 205 206 Input Parameters: 207 + snes - the SNES context 208 - x - vector where SNES solution is to be stored. 209 210 Output Parameter: 211 . J - the matrix-free matrix 212 213 Notes: 214 The matrix-free matrix context merely contains the function pointers 215 and work space for performing finite difference approximations of 216 Jacobian-vector products, J(u)*a, via 217 218 .vb 219 J(u)*a = [J(u+h*a) - J(u)]/h where 220 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 221 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 222 where 223 error_rel = square root of relative error in function evaluation 224 umin = minimum iterate parameter 225 .ve 226 227 The user can set these parameters via SNESDefaultMatrixFreeSetParameters(). 228 See the nonlinear solvers chapter of the users manual for details. 229 230 The user should call MatDestroy() when finished with the matrix-free 231 matrix context. 232 233 Options Database Keys: 234 + -snes_mf_err <error_rel> - Sets error_rel 235 - -snes_mf_unim <umin> - Sets umin 236 237 .keywords: SNES, default, matrix-free, create, matrix 238 239 .seealso: MatDestroy(), SNESDefaultMatrixFreeSetParameters() 240 @*/ 241 int SNESDefaultMatrixFreeCreate(SNES snes,Vec x, Mat *J) 242 { 243 MPI_Comm comm; 244 MFCtx_Private *mfctx; 245 int n, nloc, ierr, flg; 246 char p[64]; 247 248 PetscFunctionBegin; 249 mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx); 250 PLogObjectMemory(snes,sizeof(MFCtx_Private)); 251 mfctx->sp = 0; 252 mfctx->snes = snes; 253 mfctx->error_rel = 1.e-8; /* assumes double precision */ 254 mfctx->umin = 1.e-6; 255 mfctx->currenth = 0.0; 256 mfctx->historyh = PETSC_NULL; 257 mfctx->ncurrenth = 0; 258 mfctx->maxcurrenth = 0; 259 260 ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr); 261 ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr); 262 ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 263 PetscStrcpy(p,"-"); 264 if (snes->prefix) PetscStrcat(p,snes->prefix); 265 if (flg) { 266 (*PetscHelpPrintf)(snes->comm," %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel); 267 (*PetscHelpPrintf)(snes->comm," %ssnes_mf_umin <umin> see users manual (default %g)\n",p,mfctx->umin); 268 } 269 ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr); 270 ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr); 271 ierr = VecGetSize(x,&n); CHKERRQ(ierr); 272 ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr); 273 ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr); 274 ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult_Private);CHKERRQ(ierr); 275 ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy_Private);CHKERRQ(ierr); 276 ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView_Private); CHKERRQ(ierr); 277 PLogObjectParent(*J,mfctx->w); 278 PLogObjectParent(snes,*J); 279 PetscFunctionReturn(0); 280 } 281 282 #undef __FUNC__ 283 #define __FUNC__ "SNESDefaultMatrixFreeGetH" 284 /*@ 285 SNESDefaultMatrixFreeGetH - Gets the last h that was used as the differencing 286 parameter. 287 288 Not Collective 289 290 Input Parameters: 291 . mat - the matrix obtained with SNESDefaultMatrixFreeCreate() 292 293 Output Paramter: 294 . h - the differencing step size 295 296 .keywords: SNES, matrix-free, parameters 297 298 .seealso: SNESDefaultMatrixFreeCreate() 299 @*/ 300 int SNESDefaultMatrixFreeGetH(Mat mat,Scalar *h) 301 { 302 MFCtx_Private *ctx; 303 int ierr; 304 305 PetscFunctionBegin; 306 ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 307 if (ctx) { 308 *h = ctx->currenth; 309 } 310 PetscFunctionReturn(0); 311 } 312 313 #undef __FUNC__ 314 #define __FUNC__ "SNESDefaultMatrixFreeKSPMonitor" 315 /* 316 SNESDefaultMatrixFreeKSPMonitor - A KSP monitor for use with the default PETSc 317 SNES matrix free routines. Prints the h differencing parameter used at each 318 timestep. 319 320 */ 321 int SNESDefaultMatrixFreeKSPMonitor(KSP ksp,int n,double rnorm,void *dummy) 322 { 323 PC pc; 324 MFCtx_Private *ctx; 325 int ierr; 326 Mat mat; 327 MPI_Comm comm; 328 PetscTruth nonzeroinitialguess; 329 330 PetscFunctionBegin; 331 ierr = PetscObjectGetComm((PetscObject)ksp,&comm);CHKERRQ(ierr); 332 ierr = KSPGetPC(ksp,&pc); CHKERRQ(ierr); 333 ierr = KSPGetInitialGuessNonzero(ksp,&nonzeroinitialguess);CHKERRQ(ierr); 334 ierr = PCGetOperators(pc,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 335 ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 336 if (n > 0 || nonzeroinitialguess) { 337 #if defined(USE_PETSC_COMPLEX) 338 PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g + %g i\n",n,rnorm, 339 PetscReal(ctx->currenth),PetscImaginary(ctx->currenth)); 340 #else 341 PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g \n",n,rnorm,ctx->currenth); 342 #endif 343 } else { 344 PetscPrintf(comm,"%d KSP Residual norm %14.12e\n",n,rnorm); 345 } 346 PetscFunctionReturn(0); 347 } 348 349 #undef __FUNC__ 350 #define __FUNC__ "SNESDefaultMatrixFreeSetParameters" 351 /*@ 352 SNESDefaultMatrixFreeSetParameters - Sets the parameters for the approximation of 353 matrix-vector products using finite differences. 354 355 Collective on Mat 356 357 Input Parameters: 358 + mat - the matrix free matrix created via SNESDefaultMatrixFreeCreate() 359 . error_rel - relative error (should be set to the square root of 360 the relative error in the function evaluations) 361 - umin - minimum allowable u-value 362 363 Notes: 364 The default matrix-free matrix-vector product routine computes 365 .vb 366 J(u)*a = [J(u+h*a) - J(u)]/h where 367 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 368 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 369 .ve 370 371 Options Database Keys: 372 + -snes_mf_err <error_rel> - Sets error_rel 373 - -snes_mf_unim <umin> - Sets umin 374 375 .keywords: SNES, matrix-free, parameters 376 377 .seealso: SNESDefaultMatrixFreeCreate() 378 @*/ 379 int SNESDefaultMatrixFreeSetParameters(Mat mat,double error,double umin) 380 { 381 MFCtx_Private *ctx; 382 int ierr; 383 384 PetscFunctionBegin; 385 ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 386 if (ctx) { 387 if (error != PETSC_DEFAULT) ctx->error_rel = error; 388 if (umin != PETSC_DEFAULT) ctx->umin = umin; 389 } 390 PetscFunctionReturn(0); 391 } 392 393 #undef __FUNC__ 394 #define __FUNC__ "SNESDefaultMatrixFreeAddNullSpace" 395 /*@ 396 SNESDefaultMatrixFreeAddNullSpace - Provides a null space that 397 an operator is supposed to have. Since roundoff will create a 398 small component in the null space, if you know the null space 399 you may have it automatically removed. 400 401 Collective on Mat 402 403 Input Parameters: 404 + J - the matrix-free matrix context 405 . has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants 406 . n - number of vectors (excluding constant vector) in null space 407 - vecs - the vectors that span the null space (excluding the constant vector); 408 these vectors must be orthonormal 409 410 .keywords: SNES, matrix-free, null space 411 @*/ 412 int SNESDefaultMatrixFreeAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs) 413 { 414 int ierr; 415 MFCtx_Private *ctx; 416 MPI_Comm comm; 417 418 PetscFunctionBegin; 419 PetscObjectGetComm((PetscObject)J,&comm); 420 421 ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 422 /* no context indicates that it is not the "matrix free" matrix type */ 423 if (!ctx) PetscFunctionReturn(0); 424 ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr); 425 PetscFunctionReturn(0); 426 } 427 428 429 430 431