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