1 #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 2 3 typedef struct { 4 PetscInt sweeps; /* number of sweeps through the local subdomain before neighbor communication */ 5 PetscInt max_its; /* maximum iterations of the inner pointblock solver */ 6 PetscReal rtol; /* relative tolerance of the inner pointblock solver */ 7 PetscReal abstol; /* absolute tolerance of the inner pointblock solver */ 8 PetscReal stol; /* step tolerance of the inner pointblock solver */ 9 } SNES_GS; 10 11 12 #undef __FUNCT__ 13 #define __FUNCT__ "SNESGSSetTolerances" 14 /*@ 15 SNESGSSetTolerances - Sets various parameters used in convergence tests. 16 17 Logically Collective on SNES 18 19 Input Parameters: 20 + snes - the SNES context 21 . abstol - absolute convergence tolerance 22 . rtol - relative convergence tolerance 23 . stol - convergence tolerance in terms of the norm of the change in the solution between steps, || delta x || < stol*|| x || 24 - maxit - maximum number of iterations 25 26 27 Options Database Keys: 28 + -snes_gs_atol <abstol> - Sets abstol 29 . -snes_gs_rtol <rtol> - Sets rtol 30 . -snes_gs_stol <stol> - Sets stol 31 - -snes_max_it <maxit> - Sets maxit 32 33 Level: intermediate 34 35 .keywords: SNES, nonlinear, gauss-seidel, set, convergence, tolerances 36 37 .seealso: SNESSetTrustRegionTolerance() 38 @*/ 39 PetscErrorCode SNESGSSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit) 40 { 41 SNES_GS *gs = (SNES_GS*)snes->data; 42 43 PetscFunctionBegin; 44 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 45 46 if (abstol != PETSC_DEFAULT) { 47 if (abstol < 0.0) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %G must be non-negative",abstol); 48 gs->abstol = abstol; 49 } 50 if (rtol != PETSC_DEFAULT) { 51 if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %G must be non-negative and less than 1.0",rtol); 52 gs->rtol = rtol; 53 } 54 if (stol != PETSC_DEFAULT) { 55 if (stol < 0.0) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %G must be non-negative",stol); 56 gs->stol = stol; 57 } 58 if (maxit != PETSC_DEFAULT) { 59 if (maxit < 0) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit); 60 gs->max_its = maxit; 61 } 62 PetscFunctionReturn(0); 63 } 64 65 #undef __FUNCT__ 66 #define __FUNCT__ "SNESGSGetTolerances" 67 /*@ 68 SNESGSGetTolerances - Gets various parameters used in convergence tests. 69 70 Not Collective 71 72 Input Parameters: 73 + snes - the SNES context 74 . atol - absolute convergence tolerance 75 . rtol - relative convergence tolerance 76 . stol - convergence tolerance in terms of the norm 77 of the change in the solution between steps 78 - maxit - maximum number of iterations 79 80 Notes: 81 The user can specify PETSC_NULL for any parameter that is not needed. 82 83 Level: intermediate 84 85 .keywords: SNES, nonlinear, get, convergence, tolerances 86 87 .seealso: SNESSetTolerances() 88 @*/ 89 PetscErrorCode SNESGSGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit) 90 { 91 SNES_GS *gs = (SNES_GS*)snes->data; 92 93 PetscFunctionBegin; 94 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 95 if (atol) *atol = gs->abstol; 96 if (rtol) *rtol = gs->rtol; 97 if (stol) *stol = gs->stol; 98 if (maxit) *maxit = gs->max_its; 99 PetscFunctionReturn(0); 100 } 101 102 #undef __FUNCT__ 103 #define __FUNCT__ "SNESGSSetSweeps" 104 /*@ 105 SNESGSSetSweeps - Sets the number of sweeps of GS to use. 106 107 Input Parameters: 108 + snes - the SNES context 109 - sweeps - the number of sweeps of GS to perform. 110 111 Level: intermediate 112 113 .keywords: SNES, nonlinear, set, Gauss-Siedel 114 115 .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESGSGetSweeps() 116 @*/ 117 118 PetscErrorCode SNESGSSetSweeps(SNES snes, PetscInt sweeps) 119 { 120 SNES_GS *gs = (SNES_GS*)snes->data; 121 122 PetscFunctionBegin; 123 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 124 gs->sweeps = sweeps; 125 PetscFunctionReturn(0); 126 } 127 128 #undef __FUNCT__ 129 #define __FUNCT__ "SNESGSGetSweeps" 130 /*@ 131 SNESGSGetSweeps - Gets the number of sweeps GS will use. 132 133 Input Parameters: 134 . snes - the SNES context 135 136 Output Parameters: 137 . sweeps - the number of sweeps of GS to perform. 138 139 Level: intermediate 140 141 .keywords: SNES, nonlinear, set, Gauss-Siedel 142 143 .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESGSSetSweeps() 144 @*/ 145 PetscErrorCode SNESGSGetSweeps(SNES snes, PetscInt * sweeps) 146 { 147 SNES_GS *gs = (SNES_GS*)snes->data; 148 149 PetscFunctionBegin; 150 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 151 *sweeps = gs->sweeps; 152 PetscFunctionReturn(0); 153 } 154 155 156 #undef __FUNCT__ 157 #define __FUNCT__ "SNESDefaultApplyGS" 158 PetscErrorCode SNESDefaultApplyGS(SNES snes,Vec X,Vec F,void *ctx) 159 { 160 PetscFunctionBegin; 161 /* see if there's a coloring on the DM */ 162 163 PetscFunctionReturn(0); 164 } 165 166 #undef __FUNCT__ 167 #define __FUNCT__ "SNESReset_GS" 168 PetscErrorCode SNESReset_GS(SNES snes) 169 { 170 PetscFunctionBegin; 171 PetscFunctionReturn(0); 172 } 173 174 #undef __FUNCT__ 175 #define __FUNCT__ "SNESDestroy_GS" 176 PetscErrorCode SNESDestroy_GS(SNES snes) 177 { 178 PetscErrorCode ierr; 179 180 PetscFunctionBegin; 181 ierr = SNESReset_GS(snes);CHKERRQ(ierr); 182 ierr = PetscFree(snes->data);CHKERRQ(ierr); 183 PetscFunctionReturn(0); 184 } 185 186 #undef __FUNCT__ 187 #define __FUNCT__ "SNESSetUp_GS" 188 PetscErrorCode SNESSetUp_GS(SNES snes) 189 { 190 PetscFunctionBegin; 191 PetscFunctionReturn(0); 192 } 193 194 #undef __FUNCT__ 195 #define __FUNCT__ "SNESSetFromOptions_GS" 196 PetscErrorCode SNESSetFromOptions_GS(SNES snes) 197 { 198 SNES_GS *gs = (SNES_GS*)snes->data; 199 PetscErrorCode ierr; 200 PetscInt sweeps,max_its=PETSC_DEFAULT; 201 PetscReal rtol=PETSC_DEFAULT,atol=PETSC_DEFAULT,stol=PETSC_DEFAULT; 202 PetscBool flg,flg1,flg2,flg3; 203 204 PetscFunctionBegin; 205 ierr = PetscOptionsHead("SNES GS options");CHKERRQ(ierr); 206 /* GS Options */ 207 ierr = PetscOptionsInt("-snes_gs_sweeps","Number of sweeps of GS to apply","SNESComputeGS",gs->sweeps,&sweeps,&flg);CHKERRQ(ierr); 208 if (flg) { 209 ierr = SNESGSSetSweeps(snes,sweeps);CHKERRQ(ierr); 210 } 211 ierr = PetscOptionsReal("-snes_gs_atol","Number of sweeps of GS to apply","SNESComputeGS",gs->abstol,&atol,&flg);CHKERRQ(ierr); 212 ierr = PetscOptionsReal("-snes_gs_rtol","Number of sweeps of GS to apply","SNESComputeGS",gs->rtol,&rtol,&flg1);CHKERRQ(ierr); 213 ierr = PetscOptionsReal("-snes_gs_stol","Number of sweeps of GS to apply","SNESComputeGS",gs->stol,&stol,&flg2);CHKERRQ(ierr); 214 ierr = PetscOptionsInt("-snes_gs_max_it","Number of sweeps of GS to apply","SNESComputeGS",gs->max_its,&max_its,&flg3);CHKERRQ(ierr); 215 if (flg || flg1 || flg2 || flg3) { 216 ierr = SNESGSSetTolerances(snes,atol,rtol,stol,max_its);CHKERRQ(ierr); 217 } 218 ierr = PetscOptionsTail();CHKERRQ(ierr); 219 PetscFunctionReturn(0); 220 } 221 222 #undef __FUNCT__ 223 #define __FUNCT__ "SNESView_GS" 224 PetscErrorCode SNESView_GS(SNES snes, PetscViewer viewer) 225 { 226 PetscFunctionBegin; 227 PetscFunctionReturn(0); 228 } 229 230 #undef __FUNCT__ 231 #define __FUNCT__ "SNESSolve_GS" 232 PetscErrorCode SNESSolve_GS(SNES snes) 233 { 234 Vec F; 235 Vec X; 236 Vec B; 237 PetscInt i; 238 PetscReal fnorm; 239 PetscErrorCode ierr; 240 SNESNormType normtype; 241 242 PetscFunctionBegin; 243 244 X = snes->vec_sol; 245 F = snes->vec_func; 246 B = snes->vec_rhs; 247 248 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 249 snes->iter = 0; 250 snes->norm = 0.; 251 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 252 snes->reason = SNES_CONVERGED_ITERATING; 253 254 ierr = SNESGetNormType(snes, &normtype);CHKERRQ(ierr); 255 if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 256 /* compute the initial function and preconditioned update delX */ 257 if (!snes->vec_func_init_set) { 258 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 259 if (snes->domainerror) { 260 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 261 PetscFunctionReturn(0); 262 } 263 } else { 264 snes->vec_func_init_set = PETSC_FALSE; 265 } 266 267 /* convergence test */ 268 if (!snes->norm_init_set) { 269 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 270 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 271 } else { 272 fnorm = snes->norm_init; 273 snes->norm_init_set = PETSC_FALSE; 274 } 275 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 276 snes->iter = 0; 277 snes->norm = fnorm; 278 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 279 SNESLogConvHistory(snes,snes->norm,0); 280 ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 281 282 /* set parameter for default relative tolerance convergence test */ 283 snes->ttol = fnorm*snes->rtol; 284 285 /* test convergence */ 286 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 287 if (snes->reason) PetscFunctionReturn(0); 288 } else { 289 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 290 SNESLogConvHistory(snes,snes->norm,0); 291 ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 292 } 293 294 295 /* Call general purpose update function */ 296 if (snes->ops->update) { 297 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 298 } 299 300 for (i = 0; i < snes->max_its; i++) { 301 ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 302 /* only compute norms if requested or about to exit due to maximum iterations */ 303 if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) { 304 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 305 if (snes->domainerror) { 306 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 307 PetscFunctionReturn(0); 308 } 309 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 310 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 311 } 312 /* Monitor convergence */ 313 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 314 snes->iter = i+1; 315 snes->norm = fnorm; 316 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 317 SNESLogConvHistory(snes,snes->norm,0); 318 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 319 /* Test for convergence */ 320 if (normtype == SNES_NORM_FUNCTION) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 321 if (snes->reason) PetscFunctionReturn(0); 322 /* Call general purpose update function */ 323 if (snes->ops->update) { 324 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 325 } 326 } 327 if (normtype == SNES_NORM_FUNCTION) { 328 if (i == snes->max_its) { 329 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 330 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 331 } 332 } else { 333 if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */ 334 } 335 PetscFunctionReturn(0); 336 } 337 338 /*MC 339 SNESGS - Just calls the user-provided solution routine provided with SNESSetGS() 340 341 Level: advanced 342 343 Notes: 344 the Gauss-Seidel smoother is inherited through composition. If a solver has been created with SNESGetPC(), it will have 345 its parent's Gauss-Seidel routine associated with it. 346 347 .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetGS(), SNESType (for list of available types) 348 M*/ 349 350 EXTERN_C_BEGIN 351 #undef __FUNCT__ 352 #define __FUNCT__ "SNESCreate_GS" 353 PetscErrorCode SNESCreate_GS(SNES snes) 354 { 355 SNES_GS *gs; 356 PetscErrorCode ierr; 357 358 PetscFunctionBegin; 359 snes->ops->destroy = SNESDestroy_GS; 360 snes->ops->setup = SNESSetUp_GS; 361 snes->ops->setfromoptions = SNESSetFromOptions_GS; 362 snes->ops->view = SNESView_GS; 363 snes->ops->solve = SNESSolve_GS; 364 snes->ops->reset = SNESReset_GS; 365 366 snes->usesksp = PETSC_FALSE; 367 snes->usespc = PETSC_FALSE; 368 369 if (!snes->tolerancesset) { 370 snes->max_its = 10000; 371 snes->max_funcs = 10000; 372 } 373 374 ierr = PetscNewLog(snes, SNES_GS, &gs);CHKERRQ(ierr); 375 376 gs->sweeps = 1; 377 gs->rtol = 1e-5; 378 gs->abstol = 1e-15; 379 gs->stol = 1e-12; 380 gs->max_its = 50; 381 382 snes->data = (void*) gs; 383 384 PetscFunctionReturn(0); 385 } 386 EXTERN_C_END 387