1 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 2 3 static SNESMSType SNESMSDefault = SNESMSM62; 4 static PetscBool SNESMSRegisterAllCalled; 5 static PetscBool SNESMSPackageInitialized; 6 7 typedef struct _SNESMSTableau *SNESMSTableau; 8 struct _SNESMSTableau { 9 char *name; 10 PetscInt nstages; /* Number of stages */ 11 PetscInt nregisters; /* Number of registers */ 12 PetscReal stability; /* Scaled stability region */ 13 PetscReal *gamma; /* Coefficients of 3S* method */ 14 PetscReal *delta; /* Coefficients of 3S* method */ 15 PetscReal *betasub; /* Subdiagonal of beta in Shu-Osher form */ 16 }; 17 18 typedef struct _SNESMSTableauLink *SNESMSTableauLink; 19 struct _SNESMSTableauLink { 20 struct _SNESMSTableau tab; 21 SNESMSTableauLink next; 22 }; 23 static SNESMSTableauLink SNESMSTableauList; 24 25 typedef struct { 26 SNESMSTableau tableau; /* Tableau in low-storage form */ 27 PetscReal damping; /* Damping parameter, like length of (pseudo) time step */ 28 PetscBool norms; /* Compute norms, usually only for monitoring purposes */ 29 } SNES_MS; 30 31 /*@C 32 SNESMSRegisterAll - Registers all of the multi-stage methods in `SNESMS` 33 34 Logically Collective 35 36 Level: developer 37 38 .seealso: `SNES`, `SNESMS`, `SNESMSRegisterDestroy()` 39 @*/ 40 PetscErrorCode SNESMSRegisterAll(void) 41 { 42 PetscFunctionBegin; 43 if (SNESMSRegisterAllCalled) PetscFunctionReturn(0); 44 SNESMSRegisterAllCalled = PETSC_TRUE; 45 46 { 47 PetscReal alpha[1] = {1.0}; 48 PetscCall(SNESMSRegister(SNESMSEULER, 1, 1, 1.0, NULL, NULL, alpha)); 49 } 50 51 { 52 PetscReal gamma[3][6] = { 53 {0.0000000000000000E+00, -7.0304722367110606E-01, -1.9836719667506464E-01, -1.6023843981863788E+00, 9.4483822882855284E-02, -1.4204296130641869E-01}, 54 {1.0000000000000000E+00, 1.1111025767083920E+00, 5.6150921583923230E-01, 7.4151723494934041E-01, 3.1714538168600587E-01, 4.6479276238548706E-01 }, 55 {0.0000000000000000E+00, 0.0000000000000000E+00, 0.0000000000000000E+00, 6.7968174970583317E-01, -4.1755042846051737E-03, -1.9115668129923846E-01} 56 }; 57 PetscReal delta[6] = {1.0000000000000000E+00, 5.3275427433201750E-01, 6.0143544663985238E-01, 4.5874077053842177E-01, 2.7544386906104651E-01, 0.0000000000000000E+00}; 58 PetscReal betasub[6] = {8.4753115429481929E-01, 7.4018896368655618E-01, 6.5963574086583309E-03, 4.6747795645517759E-01, 1.3314545813643919E-01, 5.3260800028018784E-01}; 59 PetscCall(SNESMSRegister(SNESMSM62, 6, 3, 1.0, &gamma[0][0], delta, betasub)); 60 } 61 62 { /* Jameson (1983) */ 63 PetscReal alpha[4] = {0.25, 0.5, 0.55, 1.0}; 64 PetscCall(SNESMSRegister(SNESMSJAMESON83, 4, 1, 1.0, NULL, NULL, alpha)); 65 } 66 67 { /* Van Leer, Tai, and Powell (1989) 1 stage, order 1 */ 68 PetscReal alpha[1] = {1.0}; 69 PetscCall(SNESMSRegister(SNESMSVLTP11, 1, 1, 0.5, NULL, NULL, alpha)); 70 } 71 { /* Van Leer, Tai, and Powell (1989) 2 stage, order 1 */ 72 PetscReal alpha[2] = {0.3333, 1.0}; 73 PetscCall(SNESMSRegister(SNESMSVLTP21, 2, 1, 1.0, NULL, NULL, alpha)); 74 } 75 { /* Van Leer, Tai, and Powell (1989) 3 stage, order 1 */ 76 PetscReal alpha[3] = {0.1481, 0.4000, 1.0}; 77 PetscCall(SNESMSRegister(SNESMSVLTP31, 3, 1, 1.5, NULL, NULL, alpha)); 78 } 79 { /* Van Leer, Tai, and Powell (1989) 4 stage, order 1 */ 80 PetscReal alpha[4] = {0.0833, 0.2069, 0.4265, 1.0}; 81 PetscCall(SNESMSRegister(SNESMSVLTP41, 4, 1, 2.0, NULL, NULL, alpha)); 82 } 83 { /* Van Leer, Tai, and Powell (1989) 5 stage, order 1 */ 84 PetscReal alpha[5] = {0.0533, 0.1263, 0.2375, 0.4414, 1.0}; 85 PetscCall(SNESMSRegister(SNESMSVLTP51, 5, 1, 2.5, NULL, NULL, alpha)); 86 } 87 { /* Van Leer, Tai, and Powell (1989) 6 stage, order 1 */ 88 PetscReal alpha[6] = {0.0370, 0.0851, 0.1521, 0.2562, 0.4512, 1.0}; 89 PetscCall(SNESMSRegister(SNESMSVLTP61, 6, 1, 3.0, NULL, NULL, alpha)); 90 } 91 PetscFunctionReturn(0); 92 } 93 94 /*@C 95 SNESMSRegisterDestroy - Frees the list of schemes that were registered by `SNESMSRegister()`. 96 97 Not Collective 98 99 Level: developer 100 101 .seealso: `SNESMSRegister()`, `SNESMSRegisterAll()` 102 @*/ 103 PetscErrorCode SNESMSRegisterDestroy(void) 104 { 105 SNESMSTableauLink link; 106 107 PetscFunctionBegin; 108 while ((link = SNESMSTableauList)) { 109 SNESMSTableau t = &link->tab; 110 SNESMSTableauList = link->next; 111 112 PetscCall(PetscFree(t->name)); 113 PetscCall(PetscFree(t->gamma)); 114 PetscCall(PetscFree(t->delta)); 115 PetscCall(PetscFree(t->betasub)); 116 PetscCall(PetscFree(link)); 117 } 118 SNESMSRegisterAllCalled = PETSC_FALSE; 119 PetscFunctionReturn(0); 120 } 121 122 /*@C 123 SNESMSInitializePackage - This function initializes everything in the `SNESMS` package. It is called 124 from `SNESInitializePackage()`. 125 126 Level: developer 127 128 .seealso: `PetscInitialize()` 129 @*/ 130 PetscErrorCode SNESMSInitializePackage(void) 131 { 132 PetscFunctionBegin; 133 if (SNESMSPackageInitialized) PetscFunctionReturn(0); 134 SNESMSPackageInitialized = PETSC_TRUE; 135 136 PetscCall(SNESMSRegisterAll()); 137 PetscCall(PetscRegisterFinalize(SNESMSFinalizePackage)); 138 PetscFunctionReturn(0); 139 } 140 141 /*@C 142 SNESMSFinalizePackage - This function destroys everything in the `SNESMS` package. It is 143 called from `PetscFinalize()`. 144 145 Level: developer 146 147 .seealso: `PetscFinalize()` 148 @*/ 149 PetscErrorCode SNESMSFinalizePackage(void) 150 { 151 PetscFunctionBegin; 152 SNESMSPackageInitialized = PETSC_FALSE; 153 154 PetscCall(SNESMSRegisterDestroy()); 155 PetscFunctionReturn(0); 156 } 157 158 /*@C 159 SNESMSRegister - register a multistage scheme for `SNESMS` 160 161 Logically Collective 162 163 Input Parameters: 164 + name - identifier for method 165 . nstages - number of stages 166 . nregisters - number of registers used by low-storage implementation 167 . stability - scaled stability region 168 . gamma - coefficients, see Ketcheson's paper 169 . delta - coefficients, see Ketcheson's paper 170 - betasub - subdiagonal of Shu-Osher form 171 172 Notes: 173 The notation is described in Ketcheson (2010) Runge-Kutta methods with minimum storage implementations. 174 175 Many multistage schemes are of the form 176 $ X_0 = X^{(n)} 177 $ X_k = X_0 + \alpha_k * F(X_{k-1}), k = 1,\ldots,s 178 $ X^{(n+1)} = X_s 179 These methods can be registered with 180 .vb 181 SNESMSRegister("name",s,1,stability,NULL,NULL,alpha); 182 .ve 183 184 Level: advanced 185 186 .seealso: `SNESMS` 187 @*/ 188 PetscErrorCode SNESMSRegister(SNESMSType name, PetscInt nstages, PetscInt nregisters, PetscReal stability, const PetscReal gamma[], const PetscReal delta[], const PetscReal betasub[]) 189 { 190 SNESMSTableauLink link; 191 SNESMSTableau t; 192 193 PetscFunctionBegin; 194 PetscValidCharPointer(name, 1); 195 PetscCheck(nstages >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have at least one stage"); 196 if (gamma || delta) { 197 PetscCheck(nregisters == 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support for methods written in 3-register form"); 198 PetscValidRealPointer(gamma, 5); 199 PetscValidRealPointer(delta, 6); 200 } else { 201 PetscCheck(nregisters == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support for methods written in 1-register form"); 202 } 203 PetscValidRealPointer(betasub, 7); 204 205 PetscCall(SNESMSInitializePackage()); 206 PetscCall(PetscNew(&link)); 207 t = &link->tab; 208 PetscCall(PetscStrallocpy(name, &t->name)); 209 t->nstages = nstages; 210 t->nregisters = nregisters; 211 t->stability = stability; 212 213 if (gamma && delta) { 214 PetscCall(PetscMalloc1(nstages * nregisters, &t->gamma)); 215 PetscCall(PetscMalloc1(nstages, &t->delta)); 216 PetscCall(PetscArraycpy(t->gamma, gamma, nstages * nregisters)); 217 PetscCall(PetscArraycpy(t->delta, delta, nstages)); 218 } 219 PetscCall(PetscMalloc1(nstages, &t->betasub)); 220 PetscCall(PetscArraycpy(t->betasub, betasub, nstages)); 221 222 link->next = SNESMSTableauList; 223 SNESMSTableauList = link; 224 PetscFunctionReturn(0); 225 } 226 227 /* 228 X - initial state, updated in-place. 229 F - residual, computed at the initial X on input 230 */ 231 static PetscErrorCode SNESMSStep_3Sstar(SNES snes, Vec X, Vec F) 232 { 233 SNES_MS *ms = (SNES_MS *)snes->data; 234 SNESMSTableau t = ms->tableau; 235 const PetscReal *gamma = t->gamma, *delta = t->delta, *betasub = t->betasub; 236 Vec S1, S2, S3, Y; 237 PetscInt i, nstages = t->nstages; 238 239 PetscFunctionBegin; 240 Y = snes->work[0]; 241 S1 = X; 242 S2 = snes->work[1]; 243 S3 = snes->work[2]; 244 PetscCall(VecZeroEntries(S2)); 245 PetscCall(VecCopy(X, S3)); 246 for (i = 0; i < nstages; i++) { 247 Vec Ss[4]; 248 PetscScalar scoeff[4]; 249 250 Ss[0] = S1; 251 Ss[1] = S2; 252 Ss[2] = S3; 253 Ss[3] = Y; 254 255 scoeff[0] = gamma[0 * nstages + i] - 1; 256 scoeff[1] = gamma[1 * nstages + i]; 257 scoeff[2] = gamma[2 * nstages + i]; 258 scoeff[3] = -betasub[i] * ms->damping; 259 260 PetscCall(VecAXPY(S2, delta[i], S1)); 261 if (i > 0) PetscCall(SNESComputeFunction(snes, S1, F)); 262 PetscCall(KSPSolve(snes->ksp, F, Y)); 263 PetscCall(VecMAXPY(S1, 4, scoeff, Ss)); 264 } 265 PetscFunctionReturn(0); 266 } 267 268 /* 269 X - initial state, updated in-place. 270 F - residual, computed at the initial X on input 271 */ 272 static PetscErrorCode SNESMSStep_Basic(SNES snes, Vec X, Vec F) 273 { 274 SNES_MS *ms = (SNES_MS *)snes->data; 275 SNESMSTableau tab = ms->tableau; 276 const PetscReal *alpha = tab->betasub, h = ms->damping; 277 PetscInt i, nstages = tab->nstages; 278 Vec X0 = snes->work[0]; 279 280 PetscFunctionBegin; 281 PetscCall(VecCopy(X, X0)); 282 for (i = 0; i < nstages; i++) { 283 if (i > 0) PetscCall(SNESComputeFunction(snes, X, F)); 284 PetscCall(KSPSolve(snes->ksp, F, X)); 285 PetscCall(VecAYPX(X, -alpha[i] * h, X0)); 286 } 287 PetscFunctionReturn(0); 288 } 289 290 static PetscErrorCode SNESMSStep_Step(SNES snes, Vec X, Vec F) 291 { 292 SNES_MS *ms = (SNES_MS *)snes->data; 293 SNESMSTableau tab = ms->tableau; 294 295 PetscFunctionBegin; 296 if (tab->gamma && tab->delta) { 297 PetscCall(SNESMSStep_3Sstar(snes, X, F)); 298 } else { 299 PetscCall(SNESMSStep_Basic(snes, X, F)); 300 } 301 PetscFunctionReturn(0); 302 } 303 304 static PetscErrorCode SNESMSStep_Norms(SNES snes, PetscInt iter, Vec F) 305 { 306 SNES_MS *ms = (SNES_MS *)snes->data; 307 PetscReal fnorm; 308 309 PetscFunctionBegin; 310 if (ms->norms) { 311 PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 312 SNESCheckFunctionNorm(snes, fnorm); 313 /* Monitor convergence */ 314 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 315 snes->iter = iter; 316 snes->norm = fnorm; 317 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 318 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 319 PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 320 /* Test for convergence */ 321 PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 322 } else if (iter > 0) { 323 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 324 snes->iter = iter; 325 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 326 } 327 PetscFunctionReturn(0); 328 } 329 330 static PetscErrorCode SNESSolve_MS(SNES snes) 331 { 332 SNES_MS *ms = (SNES_MS *)snes->data; 333 Vec X = snes->vec_sol, F = snes->vec_func; 334 PetscInt i; 335 336 PetscFunctionBegin; 337 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); 338 PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 339 340 snes->reason = SNES_CONVERGED_ITERATING; 341 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 342 snes->iter = 0; 343 snes->norm = 0; 344 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 345 346 if (!snes->vec_func_init_set) { 347 PetscCall(SNESComputeFunction(snes, X, F)); 348 } else snes->vec_func_init_set = PETSC_FALSE; 349 350 PetscCall(SNESMSStep_Norms(snes, 0, F)); 351 if (snes->reason) PetscFunctionReturn(0); 352 353 for (i = 0; i < snes->max_its; i++) { 354 /* Call general purpose update function */ 355 PetscTryTypeMethod(snes, update, snes->iter); 356 357 if (i == 0 && snes->jacobian) { 358 /* This method does not require a Jacobian, but it is usually preconditioned by PBJacobi */ 359 PetscCall(SNESComputeJacobian(snes, snes->vec_sol, snes->jacobian, snes->jacobian_pre)); 360 SNESCheckJacobianDomainerror(snes); 361 PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 362 } 363 364 PetscCall(SNESMSStep_Step(snes, X, F)); 365 366 if (i < snes->max_its - 1 || ms->norms) PetscCall(SNESComputeFunction(snes, X, F)); 367 368 PetscCall(SNESMSStep_Norms(snes, i + 1, F)); 369 if (snes->reason) PetscFunctionReturn(0); 370 } 371 if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; 372 PetscFunctionReturn(0); 373 } 374 375 static PetscErrorCode SNESSetUp_MS(SNES snes) 376 { 377 SNES_MS *ms = (SNES_MS *)snes->data; 378 SNESMSTableau tab = ms->tableau; 379 PetscInt nwork = tab->nregisters; 380 381 PetscFunctionBegin; 382 PetscCall(SNESSetWorkVecs(snes, nwork)); 383 PetscCall(SNESSetUpMatrices(snes)); 384 PetscFunctionReturn(0); 385 } 386 387 static PetscErrorCode SNESReset_MS(SNES snes) 388 { 389 PetscFunctionBegin; 390 PetscFunctionReturn(0); 391 } 392 393 static PetscErrorCode SNESDestroy_MS(SNES snes) 394 { 395 PetscFunctionBegin; 396 PetscCall(SNESReset_MS(snes)); 397 PetscCall(PetscFree(snes->data)); 398 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetType_C", NULL)); 399 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetType_C", NULL)); 400 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetDamping_C", NULL)); 401 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetDamping_C", NULL)); 402 PetscFunctionReturn(0); 403 } 404 405 static PetscErrorCode SNESView_MS(SNES snes, PetscViewer viewer) 406 { 407 PetscBool iascii; 408 SNES_MS *ms = (SNES_MS *)snes->data; 409 SNESMSTableau tab = ms->tableau; 410 411 PetscFunctionBegin; 412 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 413 if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, " multi-stage method type: %s\n", tab->name)); 414 PetscFunctionReturn(0); 415 } 416 417 static PetscErrorCode SNESSetFromOptions_MS(SNES snes, PetscOptionItems *PetscOptionsObject) 418 { 419 SNES_MS *ms = (SNES_MS *)snes->data; 420 421 PetscFunctionBegin; 422 PetscOptionsHeadBegin(PetscOptionsObject, "SNES MS options"); 423 { 424 SNESMSTableauLink link; 425 PetscInt count, choice; 426 PetscBool flg; 427 const char **namelist; 428 SNESMSType mstype; 429 PetscReal damping; 430 431 PetscCall(SNESMSGetType(snes, &mstype)); 432 for (link = SNESMSTableauList, count = 0; link; link = link->next, count++) 433 ; 434 PetscCall(PetscMalloc1(count, (char ***)&namelist)); 435 for (link = SNESMSTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name; 436 PetscCall(PetscOptionsEList("-snes_ms_type", "Multistage smoother type", "SNESMSSetType", (const char *const *)namelist, count, mstype, &choice, &flg)); 437 if (flg) PetscCall(SNESMSSetType(snes, namelist[choice])); 438 PetscCall(PetscFree(namelist)); 439 PetscCall(SNESMSGetDamping(snes, &damping)); 440 PetscCall(PetscOptionsReal("-snes_ms_damping", "Damping for multistage method", "SNESMSSetDamping", damping, &damping, &flg)); 441 if (flg) PetscCall(SNESMSSetDamping(snes, damping)); 442 PetscCall(PetscOptionsBool("-snes_ms_norms", "Compute norms for monitoring", "none", ms->norms, &ms->norms, NULL)); 443 } 444 PetscOptionsHeadEnd(); 445 PetscFunctionReturn(0); 446 } 447 448 static PetscErrorCode SNESMSGetType_MS(SNES snes, SNESMSType *mstype) 449 { 450 SNES_MS *ms = (SNES_MS *)snes->data; 451 SNESMSTableau tab = ms->tableau; 452 453 PetscFunctionBegin; 454 *mstype = tab->name; 455 PetscFunctionReturn(0); 456 } 457 458 static PetscErrorCode SNESMSSetType_MS(SNES snes, SNESMSType mstype) 459 { 460 SNES_MS *ms = (SNES_MS *)snes->data; 461 SNESMSTableauLink link; 462 PetscBool match; 463 464 PetscFunctionBegin; 465 if (ms->tableau) { 466 PetscCall(PetscStrcmp(ms->tableau->name, mstype, &match)); 467 if (match) PetscFunctionReturn(0); 468 } 469 for (link = SNESMSTableauList; link; link = link->next) { 470 PetscCall(PetscStrcmp(link->tab.name, mstype, &match)); 471 if (match) { 472 if (snes->setupcalled) PetscCall(SNESReset_MS(snes)); 473 ms->tableau = &link->tab; 474 if (snes->setupcalled) PetscCall(SNESSetUp_MS(snes)); 475 PetscFunctionReturn(0); 476 } 477 } 478 SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", mstype); 479 } 480 481 /*@C 482 SNESMSGetType - Get the type of multistage smoother `SNESMS` 483 484 Not collective 485 486 Input Parameter: 487 . snes - nonlinear solver context 488 489 Output Parameter: 490 . mstype - type of multistage method 491 492 Level: advanced 493 494 .seealso: `SNESMS`, `SNESMSSetType()`, `SNESMSType`, `SNESMS` 495 @*/ 496 PetscErrorCode SNESMSGetType(SNES snes, SNESMSType *mstype) 497 { 498 PetscFunctionBegin; 499 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 500 PetscValidPointer(mstype, 2); 501 PetscUseMethod(snes, "SNESMSGetType_C", (SNES, SNESMSType *), (snes, mstype)); 502 PetscFunctionReturn(0); 503 } 504 505 /*@C 506 SNESMSSetType - Set the type of multistage smoother `SNESMS` 507 508 Logically collective 509 510 Input Parameters: 511 + snes - nonlinear solver context 512 - mstype - type of multistage method 513 514 Level: advanced 515 516 .seealso: `SNESMS`, `SNESMSGetType()`, `SNESMSType`, `SNESMS` 517 @*/ 518 PetscErrorCode SNESMSSetType(SNES snes, SNESMSType mstype) 519 { 520 PetscFunctionBegin; 521 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 522 PetscValidCharPointer(mstype, 2); 523 PetscTryMethod(snes, "SNESMSSetType_C", (SNES, SNESMSType), (snes, mstype)); 524 PetscFunctionReturn(0); 525 } 526 527 static PetscErrorCode SNESMSGetDamping_MS(SNES snes, PetscReal *damping) 528 { 529 SNES_MS *ms = (SNES_MS *)snes->data; 530 531 PetscFunctionBegin; 532 *damping = ms->damping; 533 PetscFunctionReturn(0); 534 } 535 536 static PetscErrorCode SNESMSSetDamping_MS(SNES snes, PetscReal damping) 537 { 538 SNES_MS *ms = (SNES_MS *)snes->data; 539 540 PetscFunctionBegin; 541 ms->damping = damping; 542 PetscFunctionReturn(0); 543 } 544 545 /*@C 546 SNESMSGetDamping - Get the damping parameter of `SNESMS` multistage scheme 547 548 Not collective 549 550 Input Parameter: 551 . snes - nonlinear solver context 552 553 Output Parameter: 554 . damping - damping parameter 555 556 Level: advanced 557 558 .seealso: `SNESMSSetDamping()`, `SNESMS` 559 @*/ 560 PetscErrorCode SNESMSGetDamping(SNES snes, PetscReal *damping) 561 { 562 PetscFunctionBegin; 563 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 564 PetscValidRealPointer(damping, 2); 565 PetscUseMethod(snes, "SNESMSGetDamping_C", (SNES, PetscReal *), (snes, damping)); 566 PetscFunctionReturn(0); 567 } 568 569 /*@C 570 SNESMSSetDamping - Set the damping parameter for a `SNESMS` multistage scheme 571 572 Logically collective 573 574 Input Parameters: 575 + snes - nonlinear solver context 576 - damping - damping parameter 577 578 Level: advanced 579 580 .seealso: `SNESMSGetDamping()`, `SNESMS` 581 @*/ 582 PetscErrorCode SNESMSSetDamping(SNES snes, PetscReal damping) 583 { 584 PetscFunctionBegin; 585 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 586 PetscValidLogicalCollectiveReal(snes, damping, 2); 587 PetscTryMethod(snes, "SNESMSSetDamping_C", (SNES, PetscReal), (snes, damping)); 588 PetscFunctionReturn(0); 589 } 590 591 /*MC 592 SNESMS - multi-stage smoothers 593 594 Options Database Keys: 595 + -snes_ms_type - type of multi-stage smoother 596 - -snes_ms_damping - damping for multi-stage method 597 598 Notes: 599 These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for 600 FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebyshev). 601 602 Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds. 603 604 The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with `SNESMSRegister()`. 605 606 References: 607 + * - Ketcheson (2010) Runge Kutta methods with minimum storage implementations (https://doi.org/10.1016/j.jcp.2009.11.006). 608 . * - Jameson (1983) Solution of the Euler equations for two dimensional transonic flow by a multigrid method (https://doi.org/10.1016/0096-3003(83)90019-X). 609 . * - Pierce and Giles (1997) Preconditioned multigrid methods for compressible flow calculations on stretched meshes (https://doi.org/10.1006/jcph.1997.5772). 610 - * - Van Leer, Tai, and Powell (1989) Design of optimally smoothing multi-stage schemes for the Euler equations (https://doi.org/10.2514/6.1989-1933). 611 612 Level: advanced 613 614 .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESMS`, `SNESFAS`, `KSPCHEBYSHEV` 615 616 M*/ 617 PETSC_EXTERN PetscErrorCode SNESCreate_MS(SNES snes) 618 { 619 SNES_MS *ms; 620 621 PetscFunctionBegin; 622 PetscCall(SNESMSInitializePackage()); 623 624 snes->ops->setup = SNESSetUp_MS; 625 snes->ops->solve = SNESSolve_MS; 626 snes->ops->destroy = SNESDestroy_MS; 627 snes->ops->setfromoptions = SNESSetFromOptions_MS; 628 snes->ops->view = SNESView_MS; 629 snes->ops->reset = SNESReset_MS; 630 631 snes->usesnpc = PETSC_FALSE; 632 snes->usesksp = PETSC_TRUE; 633 634 snes->alwayscomputesfinalresidual = PETSC_FALSE; 635 636 PetscCall(PetscNew(&ms)); 637 snes->data = (void *)ms; 638 ms->damping = 0.9; 639 ms->norms = PETSC_FALSE; 640 641 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetType_C", SNESMSGetType_MS)); 642 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetType_C", SNESMSSetType_MS)); 643 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetDamping_C", SNESMSGetDamping_MS)); 644 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetDamping_C", SNESMSSetDamping_MS)); 645 646 PetscCall(SNESMSSetType(snes, SNESMSDefault)); 647 PetscFunctionReturn(0); 648 } 649