1 #ifdef PETSC_RCS_HEADER 2 static char vcid[] = "$Id: snesut.c,v 1.43 1999/01/13 21:46:21 bsmith Exp curfman $"; 3 #endif 4 5 #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 6 7 #undef __FUNC__ 8 #define __FUNC__ "SNESVecViewMonitor" 9 /*@C 10 SNESVecViewMonitor - Monitors progress of the SNES solvers by calling 11 VecView() for the approximate solution at each iteration. 12 13 Collective on SNES 14 15 Input Parameters: 16 + snes - the SNES context 17 . its - iteration number 18 . fgnorm - 2-norm of residual (or gradient) 19 - dummy - either a viewer or PETSC_NULL 20 21 Level: intermediate 22 23 .keywords: SNES, nonlinear, vector, monitor, view 24 25 .seealso: SNESSetMonitor(), SNESDefaultMonitor(), VecView() 26 @*/ 27 int SNESVecViewMonitor(SNES snes,int its,double fgnorm,void *dummy) 28 { 29 int ierr; 30 Vec x; 31 Viewer viewer = (Viewer) dummy; 32 33 PetscFunctionBegin; 34 ierr = SNESGetSolution(snes,&x);CHKERRQ(ierr); 35 if (!viewer) { 36 MPI_Comm comm; 37 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 38 viewer = VIEWER_DRAW_(comm); 39 } 40 ierr = VecView(x,viewer);CHKERRQ(ierr); 41 42 PetscFunctionReturn(0); 43 } 44 45 #undef __FUNC__ 46 #define __FUNC__ "SNESDefaultMonitor" 47 /*@C 48 SNESDefaultMonitor - Monitoring progress of the SNES solvers (default). 49 50 Collective on SNES 51 52 Input Parameters: 53 + snes - the SNES context 54 . its - iteration number 55 . fgnorm - 2-norm of residual (or gradient) 56 - dummy - unused context 57 58 Notes: 59 For SNES_NONLINEAR_EQUATIONS methods the routine prints the 60 residual norm at each iteration. 61 62 For SNES_UNCONSTRAINED_MINIMIZATION methods the routine prints the 63 function value and gradient norm at each iteration. 64 65 Level: intermediate 66 67 .keywords: SNES, nonlinear, default, monitor, norm 68 69 .seealso: SNESSetMonitor(), SNESVecViewMonitor() 70 @*/ 71 int SNESDefaultMonitor(SNES snes,int its,double fgnorm,void *dummy) 72 { 73 PetscFunctionBegin; 74 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 75 PetscPrintf(snes->comm, "iter = %d, SNES Function norm %g \n",its,fgnorm); 76 } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 77 PetscPrintf(snes->comm,"iter = %d, SNES Function value %g, Gradient norm %g \n",its,snes->fc,fgnorm); 78 } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class"); 79 PetscFunctionReturn(0); 80 } 81 82 /* ---------------------------------------------------------------- */ 83 #undef __FUNC__ 84 #define __FUNC__ "SNESDefaultSMonitor" 85 /* 86 Default (short) SNES Monitor, same as SNESDefaultMonitor() except 87 it prints fewer digits of the residual as the residual gets smaller. 88 This is because the later digits are meaningless and are often 89 different on different machines; by using this routine different 90 machines will usually generate the same output. 91 */ 92 int SNESDefaultSMonitor(SNES snes,int its, double fgnorm,void *dummy) 93 { 94 PetscFunctionBegin; 95 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 96 if (fgnorm > 1.e-9) { 97 PetscPrintf(snes->comm, "iter = %d, SNES Function norm %g \n",its,fgnorm); 98 } else if (fgnorm > 1.e-11){ 99 PetscPrintf(snes->comm, "iter = %d, SNES Function norm %5.3e \n",its,fgnorm); 100 } else { 101 PetscPrintf(snes->comm, "iter = %d, SNES Function norm < 1.e-11\n",its); 102 } 103 } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 104 if (fgnorm > 1.e-9) { 105 PetscPrintf(snes->comm, 106 "iter = %d, SNES Function value %g, Gradient norm %g \n",its,snes->fc,fgnorm); 107 } else if (fgnorm > 1.e-11) { 108 PetscPrintf(snes->comm, 109 "iter = %d, SNES Function value %g, Gradient norm %5.3e \n",its,snes->fc,fgnorm); 110 } else { 111 PetscPrintf(snes->comm, 112 "iter = %d, SNES Function value %g, Gradient norm < 1.e-11\n",its,snes->fc); 113 } 114 } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class"); 115 PetscFunctionReturn(0); 116 } 117 /* ---------------------------------------------------------------- */ 118 #undef __FUNC__ 119 #define __FUNC__ "SNESConverged_EQ_LS" 120 /*@C 121 SNESConverged_EQ_LS - Monitors the convergence of the solvers for 122 systems of nonlinear equations (default). 123 124 Collective on SNES 125 126 Input Parameters: 127 + snes - the SNES context 128 . xnorm - 2-norm of current iterate 129 . pnorm - 2-norm of current step 130 . fnorm - 2-norm of function 131 - dummy - unused context 132 133 Returns: 134 + 2 - if ( fnorm < atol ), 135 . 3 - if ( pnorm < xtol*xnorm ), 136 . 4 - if ( fnorm < rtol*fnorm0 ), 137 . -2 - if ( nfct > maxf ), 138 - 0 - otherwise, 139 140 where 141 + maxf - maximum number of function evaluations, 142 set with SNESSetTolerances() 143 . nfct - number of function evaluations, 144 . atol - absolute function norm tolerance, 145 set with SNESSetTolerances() 146 - rtol - relative function norm tolerance, set with SNESSetTolerances() 147 148 Level: intermediate 149 150 .keywords: SNES, nonlinear, default, converged, convergence 151 152 .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged() 153 @*/ 154 int SNESConverged_EQ_LS(SNES snes,double xnorm,double pnorm,double fnorm,void *dummy) 155 { 156 PetscFunctionBegin; 157 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 158 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 159 } 160 /* Note: Reserve return code 1, -1 for compatibility with SNESConverged_EQ_TR */ 161 if (fnorm != fnorm) { 162 PLogInfo(snes,"SNESConverged_EQ_LS:Failed to converged, function norm is NaN\n"); 163 PetscFunctionReturn(-3); 164 } 165 if (fnorm <= snes->ttol) { 166 PLogInfo(snes, 167 "SNESConverged_EQ_LS:Converged due to function norm %g < %g (relative tolerance)\n",fnorm,snes->ttol); 168 PetscFunctionReturn(4); 169 } 170 171 if (fnorm < snes->atol) { 172 PLogInfo(snes, 173 "SNESConverged_EQ_LS: Converged due to function norm %g < %g\n",fnorm,snes->atol); 174 PetscFunctionReturn(2); 175 } 176 if (pnorm < snes->xtol*(xnorm)) { 177 PLogInfo(snes, 178 "SNESConverged_EQ_LS: Converged due to small update length: %g < %g * %g\n", 179 pnorm,snes->xtol,xnorm); 180 PetscFunctionReturn(3); 181 } 182 if (snes->nfuncs > snes->max_funcs) { 183 PLogInfo(snes,"SNESConverged_EQ_LS: Exceeded maximum number of function evaluations: %d > %d\n", 184 snes->nfuncs, snes->max_funcs ); 185 PetscFunctionReturn(-2); 186 } 187 PetscFunctionReturn(0); 188 } 189 /* ------------------------------------------------------------ */ 190 #undef __FUNC__ 191 #define __FUNC__ "SNES_KSP_SetConvergenceTestEW" 192 /*@ 193 SNES_KSP_SetConvergenceTestEW - Sets alternative convergence test 194 for the linear solvers within an inexact Newton method. 195 196 Collective on SNES 197 198 Input Parameter: 199 . snes - SNES context 200 201 Notes: 202 Currently, the default is to use a constant relative tolerance for 203 the inner linear solvers. Alternatively, one can use the 204 Eisenstat-Walker method, where the relative convergence tolerance 205 is reset at each Newton iteration according progress of the nonlinear 206 solver. 207 208 Level: advanced 209 210 Reference: 211 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 212 inexact Newton method", SISC 17 (1), pp.16-32, 1996. 213 214 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton 215 @*/ 216 int SNES_KSP_SetConvergenceTestEW(SNES snes) 217 { 218 PetscFunctionBegin; 219 snes->ksp_ewconv = 1; 220 PetscFunctionReturn(0); 221 } 222 223 #undef __FUNC__ 224 #define __FUNC__ "SNES_KSP_SetParametersEW" 225 /*@ 226 SNES_KSP_SetParametersEW - Sets parameters for Eisenstat-Walker 227 convergence criteria for the linear solvers within an inexact 228 Newton method. 229 230 Collective on SNES 231 232 Input Parameters: 233 + snes - SNES context 234 . version - version 1 or 2 (default is 2) 235 . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1) 236 . rtol_max - maximum relative tolerance (0 <= rtol_max < 1) 237 . alpha - power for version 2 rtol computation (1 < alpha <= 2) 238 . alpha2 - power for safeguard 239 . gamma2 - multiplicative factor for version 2 rtol computation 240 (0 <= gamma2 <= 1) 241 - threshold - threshold for imposing safeguard (0 < threshold < 1) 242 243 Note: 244 Use PETSC_DEFAULT to retain the default for any of the parameters. 245 246 Level: advanced 247 248 Reference: 249 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 250 inexact Newton method", Utah State University Math. Stat. Dept. Res. 251 Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput. 252 253 .keywords: SNES, KSP, Eisenstat, Walker, set, parameters 254 255 .seealso: SNES_KSP_SetConvergenceTestEW() 256 @*/ 257 int SNES_KSP_SetParametersEW(SNES snes,int version,double rtol_0, 258 double rtol_max,double gamma2,double alpha, 259 double alpha2,double threshold) 260 { 261 SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx; 262 263 PetscFunctionBegin; 264 if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"No Eisenstat-Walker context existing"); 265 if (version != PETSC_DEFAULT) kctx->version = version; 266 if (rtol_0 != PETSC_DEFAULT) kctx->rtol_0 = rtol_0; 267 if (rtol_max != PETSC_DEFAULT) kctx->rtol_max = rtol_max; 268 if (gamma2 != PETSC_DEFAULT) kctx->gamma = gamma2; 269 if (alpha != PETSC_DEFAULT) kctx->alpha = alpha; 270 if (alpha2 != PETSC_DEFAULT) kctx->alpha2 = alpha2; 271 if (threshold != PETSC_DEFAULT) kctx->threshold = threshold; 272 if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) { 273 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"0.0 <= rtol_0 < 1.0: %g",kctx->rtol_0); 274 } 275 if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) { 276 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"0.0 <= rtol_max < 1.0\n",kctx->rtol_max); 277 } 278 if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) { 279 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"0.0 < threshold < 1.0\n",kctx->threshold); 280 } 281 if (kctx->gamma < 0.0 || kctx->gamma > 1.0) { 282 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"0.0 <= alpha <= 1.0\n",kctx->gamma); 283 } 284 if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) { 285 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"1.0 < alpha <= 2.0\n",kctx->alpha); 286 } 287 if (kctx->version != 1 && kctx->version !=2) { 288 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"Only versions 1 and 2 are supported: %d",kctx->version); 289 } 290 PetscFunctionReturn(0); 291 } 292 293 #undef __FUNC__ 294 #define __FUNC__ "SNES_KSP_EW_ComputeRelativeTolerance_Private" 295 int SNES_KSP_EW_ComputeRelativeTolerance_Private(SNES snes,KSP ksp) 296 { 297 SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx; 298 double rtol = 0.0, stol; 299 int ierr; 300 301 PetscFunctionBegin; 302 if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"No Eisenstat-Walker context exists"); 303 if (snes->iter == 1) { 304 rtol = kctx->rtol_0; 305 } else { 306 if (kctx->version == 1) { 307 rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last; 308 if (rtol < 0.0) rtol = -rtol; 309 stol = pow(kctx->rtol_last,kctx->alpha2); 310 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 311 } else if (kctx->version == 2) { 312 rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha); 313 stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha); 314 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 315 } else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"Only versions 1 or 2 are supported: %d",kctx->version); 316 } 317 rtol = PetscMin(rtol,kctx->rtol_max); 318 kctx->rtol_last = rtol; 319 PLogInfo(snes,"SNESConverged_EQ_LS: iter %d, Eisenstat-Walker (version %d) KSP rtol = %g\n",snes->iter,kctx->version,rtol); 320 ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); CHKERRQ(ierr); 321 kctx->norm_last = snes->norm; 322 PetscFunctionReturn(0); 323 } 324 325 #undef __FUNC__ 326 #define __FUNC__ "SNES_KSP_EW_Converged_Private" 327 int SNES_KSP_EW_Converged_Private(KSP ksp,int n,double rnorm,void *ctx) 328 { 329 SNES snes = (SNES)ctx; 330 SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx; 331 int convinfo; 332 333 PetscFunctionBegin; 334 if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"No Eisenstat-Walker context set"); 335 if (n == 0) SNES_KSP_EW_ComputeRelativeTolerance_Private(snes,ksp); 336 convinfo = KSPDefaultConverged(ksp,n,rnorm,ctx); 337 kctx->lresid_last = rnorm; 338 if (convinfo) { 339 PLogInfo(snes,"SNES_KSP_EW_Converged_Private: KSP iterations=%d, rnorm=%g\n",n,rnorm); 340 } 341 PetscFunctionReturn(convinfo); 342 } 343 344 345