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