1 #include <../src/snes/impls/gs/gsimpl.h> /*I "petscsnes.h" I*/ 2 3 /*@ 4 SNESNGSSetTolerances - Sets various parameters used in convergence tests. 5 6 Logically Collective on SNES 7 8 Input Parameters: 9 + snes - the SNES context 10 . abstol - absolute convergence tolerance 11 . rtol - relative convergence tolerance 12 . stol - convergence tolerance in terms of the norm of the change in the solution between steps, || delta x || < stol*|| x || 13 - maxit - maximum number of iterations 14 15 Options Database Keys: 16 + -snes_ngs_atol <abstol> - Sets abstol 17 . -snes_ngs_rtol <rtol> - Sets rtol 18 . -snes_ngs_stol <stol> - Sets stol 19 - -snes_max_it <maxit> - Sets maxit 20 21 Level: intermediate 22 23 .seealso: `SNESSetTrustRegionTolerance()` 24 @*/ 25 PetscErrorCode SNESNGSSetTolerances(SNES snes, PetscReal abstol, PetscReal rtol, PetscReal stol, PetscInt maxit) { 26 SNES_NGS *gs = (SNES_NGS *)snes->data; 27 28 PetscFunctionBegin; 29 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 30 31 if (abstol != PETSC_DEFAULT) { 32 PetscCheck(abstol >= 0.0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %g must be non-negative", (double)abstol); 33 gs->abstol = abstol; 34 } 35 if (rtol != PETSC_DEFAULT) { 36 PetscCheck(rtol >= 0.0 && rtol < 1.0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Relative tolerance %g must be non-negative and less than 1.0", (double)rtol); 37 gs->rtol = rtol; 38 } 39 if (stol != PETSC_DEFAULT) { 40 PetscCheck(stol >= 0.0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Step tolerance %g must be non-negative", (double)stol); 41 gs->stol = stol; 42 } 43 if (maxit != PETSC_DEFAULT) { 44 PetscCheck(maxit >= 0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", maxit); 45 gs->max_its = maxit; 46 } 47 PetscFunctionReturn(0); 48 } 49 50 /*@ 51 SNESNGSGetTolerances - Gets various parameters used in convergence tests. 52 53 Not Collective 54 55 Input Parameters: 56 + snes - the SNES context 57 . atol - absolute convergence tolerance 58 . rtol - relative convergence tolerance 59 . stol - convergence tolerance in terms of the norm 60 of the change in the solution between steps 61 - maxit - maximum number of iterations 62 63 Notes: 64 The user can specify NULL for any parameter that is not needed. 65 66 Level: intermediate 67 68 .seealso: `SNESSetTolerances()` 69 @*/ 70 PetscErrorCode SNESNGSGetTolerances(SNES snes, PetscReal *atol, PetscReal *rtol, PetscReal *stol, PetscInt *maxit) { 71 SNES_NGS *gs = (SNES_NGS *)snes->data; 72 73 PetscFunctionBegin; 74 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 75 if (atol) *atol = gs->abstol; 76 if (rtol) *rtol = gs->rtol; 77 if (stol) *stol = gs->stol; 78 if (maxit) *maxit = gs->max_its; 79 PetscFunctionReturn(0); 80 } 81 82 /*@ 83 SNESNGSSetSweeps - Sets the number of sweeps of GS to use. 84 85 Input Parameters: 86 + snes - the SNES context 87 - sweeps - the number of sweeps of GS to perform. 88 89 Level: intermediate 90 91 .seealso: `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSGetSweeps()` 92 @*/ 93 94 PetscErrorCode SNESNGSSetSweeps(SNES snes, PetscInt sweeps) { 95 SNES_NGS *gs = (SNES_NGS *)snes->data; 96 97 PetscFunctionBegin; 98 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 99 gs->sweeps = sweeps; 100 PetscFunctionReturn(0); 101 } 102 103 /*@ 104 SNESNGSGetSweeps - Gets the number of sweeps GS will use. 105 106 Input Parameters: 107 . snes - the SNES context 108 109 Output Parameters: 110 . sweeps - the number of sweeps of GS to perform. 111 112 Level: intermediate 113 114 .seealso: `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSSetSweeps()` 115 @*/ 116 PetscErrorCode SNESNGSGetSweeps(SNES snes, PetscInt *sweeps) { 117 SNES_NGS *gs = (SNES_NGS *)snes->data; 118 119 PetscFunctionBegin; 120 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 121 *sweeps = gs->sweeps; 122 PetscFunctionReturn(0); 123 } 124 125 PetscErrorCode SNESReset_NGS(SNES snes) { 126 SNES_NGS *gs = (SNES_NGS *)snes->data; 127 128 PetscFunctionBegin; 129 PetscCall(ISColoringDestroy(&gs->coloring)); 130 PetscFunctionReturn(0); 131 } 132 133 PetscErrorCode SNESDestroy_NGS(SNES snes) { 134 PetscFunctionBegin; 135 PetscCall(SNESReset_NGS(snes)); 136 PetscCall(PetscFree(snes->data)); 137 PetscFunctionReturn(0); 138 } 139 140 PetscErrorCode SNESSetUp_NGS(SNES snes) { 141 PetscErrorCode (*f)(SNES, Vec, Vec, void *); 142 143 PetscFunctionBegin; 144 PetscCall(SNESGetNGS(snes, &f, NULL)); 145 if (!f) { PetscCall(SNESSetNGS(snes, SNESComputeNGSDefaultSecant, NULL)); } 146 PetscFunctionReturn(0); 147 } 148 149 PetscErrorCode SNESSetFromOptions_NGS(SNES snes, PetscOptionItems *PetscOptionsObject) { 150 SNES_NGS *gs = (SNES_NGS *)snes->data; 151 PetscInt sweeps, max_its = PETSC_DEFAULT; 152 PetscReal rtol = PETSC_DEFAULT, atol = PETSC_DEFAULT, stol = PETSC_DEFAULT; 153 PetscBool flg, flg1, flg2, flg3; 154 155 PetscFunctionBegin; 156 PetscOptionsHeadBegin(PetscOptionsObject, "SNES GS options"); 157 /* GS Options */ 158 PetscCall(PetscOptionsInt("-snes_ngs_sweeps", "Number of sweeps of GS to apply", "SNESComputeGS", gs->sweeps, &sweeps, &flg)); 159 if (flg) PetscCall(SNESNGSSetSweeps(snes, sweeps)); 160 PetscCall(PetscOptionsReal("-snes_ngs_atol", "Absolute residual tolerance for GS iteration", "SNESComputeGS", gs->abstol, &atol, &flg)); 161 PetscCall(PetscOptionsReal("-snes_ngs_rtol", "Relative residual tolerance for GS iteration", "SNESComputeGS", gs->rtol, &rtol, &flg1)); 162 PetscCall(PetscOptionsReal("-snes_ngs_stol", "Absolute update tolerance for GS iteration", "SNESComputeGS", gs->stol, &stol, &flg2)); 163 PetscCall(PetscOptionsInt("-snes_ngs_max_it", "Maximum number of sweeps of GS to apply", "SNESComputeGS", gs->max_its, &max_its, &flg3)); 164 if (flg || flg1 || flg2 || flg3) { PetscCall(SNESNGSSetTolerances(snes, atol, rtol, stol, max_its)); } 165 flg = PETSC_FALSE; 166 PetscCall(PetscOptionsBool("-snes_ngs_secant", "Use finite difference secant approximation with coloring", "", flg, &flg, NULL)); 167 if (flg) { 168 PetscCall(SNESSetNGS(snes, SNESComputeNGSDefaultSecant, NULL)); 169 PetscCall(PetscInfo(snes, "Setting default finite difference secant approximation with coloring\n")); 170 } 171 PetscCall(PetscOptionsReal("-snes_ngs_secant_h", "Differencing parameter for secant search", "", gs->h, &gs->h, NULL)); 172 PetscCall(PetscOptionsBool("-snes_ngs_secant_mat_coloring", "Use the graph coloring of the Jacobian for the secant GS", "", gs->secant_mat, &gs->secant_mat, &flg)); 173 174 PetscOptionsHeadEnd(); 175 PetscFunctionReturn(0); 176 } 177 178 PetscErrorCode SNESView_NGS(SNES snes, PetscViewer viewer) { 179 PetscErrorCode (*f)(SNES, Vec, Vec, void *); 180 SNES_NGS *gs = (SNES_NGS *)snes->data; 181 PetscBool iascii; 182 183 PetscFunctionBegin; 184 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 185 if (iascii) { 186 PetscCall(DMSNESGetNGS(snes->dm, &f, NULL)); 187 if (f == SNESComputeNGSDefaultSecant) { PetscCall(PetscViewerASCIIPrintf(viewer, " Use finite difference secant approximation with coloring with h = %g \n", (double)gs->h)); } 188 } 189 PetscFunctionReturn(0); 190 } 191 192 PetscErrorCode SNESSolve_NGS(SNES snes) { 193 Vec F; 194 Vec X; 195 Vec B; 196 PetscInt i; 197 PetscReal fnorm; 198 SNESNormSchedule normschedule; 199 200 PetscFunctionBegin; 201 202 PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 203 204 PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 205 X = snes->vec_sol; 206 F = snes->vec_func; 207 B = snes->vec_rhs; 208 209 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 210 snes->iter = 0; 211 snes->norm = 0.; 212 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 213 snes->reason = SNES_CONVERGED_ITERATING; 214 215 PetscCall(SNESGetNormSchedule(snes, &normschedule)); 216 if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { 217 /* compute the initial function and preconditioned update delX */ 218 if (!snes->vec_func_init_set) { 219 PetscCall(SNESComputeFunction(snes, X, F)); 220 } else snes->vec_func_init_set = PETSC_FALSE; 221 222 PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 223 SNESCheckFunctionNorm(snes, fnorm); 224 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 225 snes->iter = 0; 226 snes->norm = fnorm; 227 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 228 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 229 PetscCall(SNESMonitor(snes, 0, snes->norm)); 230 231 /* test convergence */ 232 PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 233 if (snes->reason) PetscFunctionReturn(0); 234 } else { 235 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 236 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 237 } 238 239 /* Call general purpose update function */ 240 PetscTryTypeMethod(snes, update, snes->iter); 241 242 for (i = 0; i < snes->max_its; i++) { 243 PetscCall(SNESComputeNGS(snes, B, X)); 244 /* only compute norms if requested or about to exit due to maximum iterations */ 245 if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 246 PetscCall(SNESComputeFunction(snes, X, F)); 247 PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 248 SNESCheckFunctionNorm(snes, fnorm); 249 /* Monitor convergence */ 250 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 251 snes->iter = i + 1; 252 snes->norm = fnorm; 253 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 254 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 255 PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 256 } 257 /* Test for convergence */ 258 if (normschedule == SNES_NORM_ALWAYS) PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 259 if (snes->reason) PetscFunctionReturn(0); 260 /* Call general purpose update function */ 261 PetscTryTypeMethod(snes, update, snes->iter); 262 } 263 if (normschedule == SNES_NORM_ALWAYS) { 264 if (i == snes->max_its) { 265 PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", snes->max_its)); 266 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 267 } 268 } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */ 269 PetscFunctionReturn(0); 270 } 271 272 /*MC 273 SNESNGS - Either calls the user-provided solution routine provided with SNESSetNGS() or does a finite difference secant approximation 274 using coloring. 275 276 Level: advanced 277 278 Options Database: 279 + -snes_ngs_sweeps <n> - Number of sweeps of GS to apply 280 . -snes_ngs_atol <atol> - Absolute residual tolerance for GS iteration 281 . -snes_ngs_rtol <rtol> - Relative residual tolerance for GS iteration 282 . -snes_ngs_stol <stol> - Absolute update tolerance for GS iteration 283 . -snes_ngs_max_it <maxit> - Maximum number of sweeps of GS to apply 284 . -snes_ngs_secant - Use pointwise secant local Jacobian approximation with coloring instead of user provided Gauss-Seidel routine, this is 285 used by default if no user provided Gauss-Seidel routine is available. Requires either that a DM that can compute a coloring 286 is available or a Jacobian sparse matrix is provided (from which to get the coloring). 287 . -snes_ngs_secant_h <h> - Differencing parameter for secant approximation 288 . -snes_ngs_secant_mat_coloring - Use the graph coloring of the Jacobian for the secant GS even if a DM is available. 289 - -snes_norm_schedule <none, always, initialonly, finalonly, initialfinalonly> - how often the residual norms are computed 290 291 Notes: 292 the Gauss-Seidel smoother is inherited through composition. If a solver has been created with SNESGetNPC(), it will have 293 its parent's Gauss-Seidel routine associated with it. 294 295 By default this routine computes the solution norm at each iteration, this can be time consuming, you can turn this off with SNESSetNormSchedule() or -snes_norm_schedule 296 References: 297 . * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 298 SIAM Review, 57(4), 2015 299 300 .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESSetNGS()`, `SNESType`, `SNESNGSSetSweeps()`, `SNESNGSSetTolerances()`, 301 `SNESSetNormSchedule()` 302 M*/ 303 304 PETSC_EXTERN PetscErrorCode SNESCreate_NGS(SNES snes) { 305 SNES_NGS *gs; 306 307 PetscFunctionBegin; 308 snes->ops->destroy = SNESDestroy_NGS; 309 snes->ops->setup = SNESSetUp_NGS; 310 snes->ops->setfromoptions = SNESSetFromOptions_NGS; 311 snes->ops->view = SNESView_NGS; 312 snes->ops->solve = SNESSolve_NGS; 313 snes->ops->reset = SNESReset_NGS; 314 315 snes->usesksp = PETSC_FALSE; 316 snes->usesnpc = PETSC_FALSE; 317 318 snes->alwayscomputesfinalresidual = PETSC_FALSE; 319 320 if (!snes->tolerancesset) { 321 snes->max_its = 10000; 322 snes->max_funcs = 10000; 323 } 324 325 PetscCall(PetscNewLog(snes, &gs)); 326 327 gs->sweeps = 1; 328 gs->rtol = 1e-5; 329 gs->abstol = PETSC_MACHINE_EPSILON; 330 gs->stol = 1000 * PETSC_MACHINE_EPSILON; 331 gs->max_its = 50; 332 gs->h = PETSC_SQRT_MACHINE_EPSILON; 333 334 snes->data = (void *)gs; 335 PetscFunctionReturn(0); 336 } 337