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