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 Not Collective, but should be called by all processes which will need the schemes to be registered 35 36 Level: advanced 37 38 .seealso: `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: advanced 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 160 161 Not Collective, but the same schemes should be registered on all processes on which they will be used 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; Ss[1] = S2; Ss[2] = S3; Ss[3] = Y; 251 252 scoeff[0] = gamma[0*nstages+i] - 1; 253 scoeff[1] = gamma[1*nstages+i]; 254 scoeff[2] = gamma[2*nstages+i]; 255 scoeff[3] = -betasub[i]*ms->damping; 256 257 PetscCall(VecAXPY(S2,delta[i],S1)); 258 if (i > 0) { 259 PetscCall(SNESComputeFunction(snes,S1,F)); 260 } 261 PetscCall(KSPSolve(snes->ksp,F,Y)); 262 PetscCall(VecMAXPY(S1,4,scoeff,Ss)); 263 } 264 PetscFunctionReturn(0); 265 } 266 267 /* 268 X - initial state, updated in-place. 269 F - residual, computed at the initial X on input 270 */ 271 static PetscErrorCode SNESMSStep_Basic(SNES snes,Vec X,Vec F) 272 { 273 SNES_MS *ms = (SNES_MS*)snes->data; 274 SNESMSTableau tab = ms->tableau; 275 const PetscReal *alpha = tab->betasub, h = ms->damping; 276 PetscInt i,nstages = tab->nstages; 277 Vec X0 = snes->work[0]; 278 279 PetscFunctionBegin; 280 PetscCall(VecCopy(X,X0)); 281 for (i = 0; i < nstages; i++) { 282 if (i > 0) { 283 PetscCall(SNESComputeFunction(snes,X,F)); 284 } 285 PetscCall(KSPSolve(snes->ksp,F,X)); 286 PetscCall(VecAYPX(X,-alpha[i]*h,X0)); 287 } 288 PetscFunctionReturn(0); 289 } 290 291 static PetscErrorCode SNESMSStep_Step(SNES snes,Vec X,Vec F) 292 { 293 SNES_MS *ms = (SNES_MS*)snes->data; 294 SNESMSTableau tab = ms->tableau; 295 296 PetscFunctionBegin; 297 if (tab->gamma && tab->delta) { 298 PetscCall(SNESMSStep_3Sstar(snes,X,F)); 299 } else { 300 PetscCall(SNESMSStep_Basic(snes,X,F)); 301 } 302 PetscFunctionReturn(0); 303 } 304 305 static PetscErrorCode SNESMSStep_Norms(SNES snes,PetscInt iter,Vec F) 306 { 307 SNES_MS *ms = (SNES_MS*)snes->data; 308 PetscReal fnorm; 309 310 PetscFunctionBegin; 311 if (ms->norms) { 312 PetscCall(VecNorm(F,NORM_2,&fnorm)); /* fnorm <- ||F|| */ 313 SNESCheckFunctionNorm(snes,fnorm); 314 /* Monitor convergence */ 315 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 316 snes->iter = iter; 317 snes->norm = fnorm; 318 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 319 PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0)); 320 PetscCall(SNESMonitor(snes,snes->iter,snes->norm)); 321 /* Test for convergence */ 322 PetscCall((*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP)); 323 } else if (iter > 0) { 324 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 325 snes->iter = iter; 326 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 327 } 328 PetscFunctionReturn(0); 329 } 330 331 static PetscErrorCode SNESSolve_MS(SNES snes) 332 { 333 SNES_MS *ms = (SNES_MS*)snes->data; 334 Vec X = snes->vec_sol,F = snes->vec_func; 335 PetscInt i; 336 337 PetscFunctionBegin; 338 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); 339 PetscCall(PetscCitationsRegister(SNESCitation,&SNEScite)); 340 341 snes->reason = SNES_CONVERGED_ITERATING; 342 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 343 snes->iter = 0; 344 snes->norm = 0; 345 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 346 347 if (!snes->vec_func_init_set) { 348 PetscCall(SNESComputeFunction(snes,X,F)); 349 } else snes->vec_func_init_set = PETSC_FALSE; 350 351 PetscCall(SNESMSStep_Norms(snes,0,F)); 352 if (snes->reason) PetscFunctionReturn(0); 353 354 for (i = 0; i < snes->max_its; i++) { 355 356 /* Call general purpose update function */ 357 if (snes->ops->update) PetscCall((*snes->ops->update)(snes,snes->iter)); 358 359 if (i == 0 && snes->jacobian) { 360 /* This method does not require a Jacobian, but it is usually preconditioned by PBJacobi */ 361 PetscCall(SNESComputeJacobian(snes,snes->vec_sol,snes->jacobian,snes->jacobian_pre)); 362 SNESCheckJacobianDomainerror(snes); 363 PetscCall(KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre)); 364 } 365 366 PetscCall(SNESMSStep_Step(snes,X,F)); 367 368 if (i < snes->max_its-1 || ms->norms) { 369 PetscCall(SNESComputeFunction(snes,X,F)); 370 } 371 372 PetscCall(SNESMSStep_Norms(snes,i+1,F)); 373 if (snes->reason) PetscFunctionReturn(0); 374 } 375 if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; 376 PetscFunctionReturn(0); 377 } 378 379 static PetscErrorCode SNESSetUp_MS(SNES snes) 380 { 381 SNES_MS *ms = (SNES_MS*)snes->data; 382 SNESMSTableau tab = ms->tableau; 383 PetscInt nwork = tab->nregisters; 384 385 PetscFunctionBegin; 386 PetscCall(SNESSetWorkVecs(snes,nwork)); 387 PetscCall(SNESSetUpMatrices(snes)); 388 PetscFunctionReturn(0); 389 } 390 391 static PetscErrorCode SNESReset_MS(SNES snes) 392 { 393 PetscFunctionBegin; 394 PetscFunctionReturn(0); 395 } 396 397 static PetscErrorCode SNESDestroy_MS(SNES snes) 398 { 399 PetscFunctionBegin; 400 PetscCall(SNESReset_MS(snes)); 401 PetscCall(PetscFree(snes->data)); 402 PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetType_C",NULL)); 403 PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetType_C",NULL)); 404 PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetDamping_C",NULL)); 405 PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetDamping_C",NULL)); 406 PetscFunctionReturn(0); 407 } 408 409 static PetscErrorCode SNESView_MS(SNES snes,PetscViewer viewer) 410 { 411 PetscBool iascii; 412 SNES_MS *ms = (SNES_MS*)snes->data; 413 SNESMSTableau tab = ms->tableau; 414 415 PetscFunctionBegin; 416 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 417 if (iascii) { 418 PetscCall(PetscViewerASCIIPrintf(viewer," multi-stage method type: %s\n",tab->name)); 419 } 420 PetscFunctionReturn(0); 421 } 422 423 static PetscErrorCode SNESSetFromOptions_MS(PetscOptionItems *PetscOptionsObject,SNES snes) 424 { 425 SNES_MS *ms = (SNES_MS*)snes->data; 426 427 PetscFunctionBegin; 428 PetscOptionsHeadBegin(PetscOptionsObject,"SNES MS options"); 429 { 430 SNESMSTableauLink link; 431 PetscInt count,choice; 432 PetscBool flg; 433 const char **namelist; 434 SNESMSType mstype; 435 PetscReal damping; 436 437 PetscCall(SNESMSGetType(snes,&mstype)); 438 for (link=SNESMSTableauList,count=0; link; link=link->next,count++) ; 439 PetscCall(PetscMalloc1(count,(char***)&namelist)); 440 for (link=SNESMSTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 441 PetscCall(PetscOptionsEList("-snes_ms_type","Multistage smoother type","SNESMSSetType",(const char*const*)namelist,count,mstype,&choice,&flg)); 442 if (flg) PetscCall(SNESMSSetType(snes,namelist[choice])); 443 PetscCall(PetscFree(namelist)); 444 PetscCall(SNESMSGetDamping(snes,&damping)); 445 PetscCall(PetscOptionsReal("-snes_ms_damping","Damping for multistage method","SNESMSSetDamping",damping,&damping,&flg)); 446 if (flg) PetscCall(SNESMSSetDamping(snes,damping)); 447 PetscCall(PetscOptionsBool("-snes_ms_norms","Compute norms for monitoring","none",ms->norms,&ms->norms,NULL)); 448 } 449 PetscOptionsHeadEnd(); 450 PetscFunctionReturn(0); 451 } 452 453 static PetscErrorCode SNESMSGetType_MS(SNES snes,SNESMSType *mstype) 454 { 455 SNES_MS *ms = (SNES_MS*)snes->data; 456 SNESMSTableau tab = ms->tableau; 457 458 PetscFunctionBegin; 459 *mstype = tab->name; 460 PetscFunctionReturn(0); 461 } 462 463 static PetscErrorCode SNESMSSetType_MS(SNES snes,SNESMSType mstype) 464 { 465 SNES_MS *ms = (SNES_MS*)snes->data; 466 SNESMSTableauLink link; 467 PetscBool match; 468 469 PetscFunctionBegin; 470 if (ms->tableau) { 471 PetscCall(PetscStrcmp(ms->tableau->name,mstype,&match)); 472 if (match) PetscFunctionReturn(0); 473 } 474 for (link = SNESMSTableauList; link; link=link->next) { 475 PetscCall(PetscStrcmp(link->tab.name,mstype,&match)); 476 if (match) { 477 if (snes->setupcalled) PetscCall(SNESReset_MS(snes)); 478 ms->tableau = &link->tab; 479 if (snes->setupcalled) PetscCall(SNESSetUp_MS(snes)); 480 PetscFunctionReturn(0); 481 } 482 } 483 SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mstype); 484 } 485 486 /*@C 487 SNESMSGetType - Get the type of multistage smoother 488 489 Not collective 490 491 Input Parameter: 492 . snes - nonlinear solver context 493 494 Output Parameter: 495 . mstype - type of multistage method 496 497 Level: beginner 498 499 .seealso: `SNESMSSetType()`, `SNESMSType`, `SNESMS` 500 @*/ 501 PetscErrorCode SNESMSGetType(SNES snes,SNESMSType *mstype) 502 { 503 PetscFunctionBegin; 504 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 505 PetscValidPointer(mstype,2); 506 PetscUseMethod(snes,"SNESMSGetType_C",(SNES,SNESMSType*),(snes,mstype)); 507 PetscFunctionReturn(0); 508 } 509 510 /*@C 511 SNESMSSetType - Set the type of multistage smoother 512 513 Logically collective 514 515 Input Parameters: 516 + snes - nonlinear solver context 517 - mstype - type of multistage method 518 519 Level: beginner 520 521 .seealso: `SNESMSGetType()`, `SNESMSType`, `SNESMS` 522 @*/ 523 PetscErrorCode SNESMSSetType(SNES snes,SNESMSType mstype) 524 { 525 PetscFunctionBegin; 526 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 527 PetscValidCharPointer(mstype,2); 528 PetscTryMethod(snes,"SNESMSSetType_C",(SNES,SNESMSType),(snes,mstype)); 529 PetscFunctionReturn(0); 530 } 531 532 static PetscErrorCode SNESMSGetDamping_MS(SNES snes,PetscReal *damping) 533 { 534 SNES_MS *ms = (SNES_MS*)snes->data; 535 536 PetscFunctionBegin; 537 *damping = ms->damping; 538 PetscFunctionReturn(0); 539 } 540 541 static PetscErrorCode SNESMSSetDamping_MS(SNES snes,PetscReal damping) 542 { 543 SNES_MS *ms = (SNES_MS*)snes->data; 544 545 PetscFunctionBegin; 546 ms->damping = damping; 547 PetscFunctionReturn(0); 548 } 549 550 /*@C 551 SNESMSGetDamping - Get the damping parameter 552 553 Not collective 554 555 Input Parameter: 556 . snes - nonlinear solver context 557 558 Output Parameter: 559 . damping - damping parameter 560 561 Level: advanced 562 563 .seealso: `SNESMSSetDamping()`, `SNESMS` 564 @*/ 565 PetscErrorCode SNESMSGetDamping(SNES snes,PetscReal *damping) 566 { 567 PetscFunctionBegin; 568 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 569 PetscValidRealPointer(damping,2); 570 PetscUseMethod(snes,"SNESMSGetDamping_C",(SNES,PetscReal*),(snes,damping)); 571 PetscFunctionReturn(0); 572 } 573 574 /*@C 575 SNESMSSetDamping - Set the damping parameter 576 577 Logically collective 578 579 Input Parameters: 580 + snes - nonlinear solver context 581 - damping - damping parameter 582 583 Level: advanced 584 585 .seealso: `SNESMSGetDamping()`, `SNESMS` 586 @*/ 587 PetscErrorCode SNESMSSetDamping(SNES snes,PetscReal damping) 588 { 589 PetscFunctionBegin; 590 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 591 PetscValidLogicalCollectiveReal(snes,damping,2); 592 PetscTryMethod(snes,"SNESMSSetDamping_C",(SNES,PetscReal),(snes,damping)); 593 PetscFunctionReturn(0); 594 } 595 596 /* -------------------------------------------------------------------------- */ 597 /*MC 598 SNESMS - multi-stage smoothers 599 600 Options Database: 601 602 + -snes_ms_type - type of multi-stage smoother 603 - -snes_ms_damping - damping for multi-stage method 604 605 Notes: 606 These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for 607 FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebyshev). 608 609 Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds. 610 611 The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with SNESMSRegister(). 612 613 References: 614 + * - Ketcheson (2010) Runge Kutta methods with minimum storage implementations (https://doi.org/10.1016/j.jcp.2009.11.006). 615 . * - 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). 616 . * - Pierce and Giles (1997) Preconditioned multigrid methods for compressible flow calculations on stretched meshes (https://doi.org/10.1006/jcph.1997.5772). 617 - * - 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). 618 619 Level: beginner 620 621 .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESMS`, `SNESFAS`, `KSPCHEBYSHEV` 622 623 M*/ 624 PETSC_EXTERN PetscErrorCode SNESCreate_MS(SNES snes) 625 { 626 SNES_MS *ms; 627 628 PetscFunctionBegin; 629 PetscCall(SNESMSInitializePackage()); 630 631 snes->ops->setup = SNESSetUp_MS; 632 snes->ops->solve = SNESSolve_MS; 633 snes->ops->destroy = SNESDestroy_MS; 634 snes->ops->setfromoptions = SNESSetFromOptions_MS; 635 snes->ops->view = SNESView_MS; 636 snes->ops->reset = SNESReset_MS; 637 638 snes->usesnpc = PETSC_FALSE; 639 snes->usesksp = PETSC_TRUE; 640 641 snes->alwayscomputesfinalresidual = PETSC_FALSE; 642 643 PetscCall(PetscNewLog(snes,&ms)); 644 snes->data = (void*)ms; 645 ms->damping = 0.9; 646 ms->norms = PETSC_FALSE; 647 648 PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetType_C",SNESMSGetType_MS)); 649 PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetType_C",SNESMSSetType_MS)); 650 PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetDamping_C",SNESMSGetDamping_MS)); 651 PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetDamping_C",SNESMSSetDamping_MS)); 652 653 PetscCall(SNESMSSetType(snes,SNESMSDefault)); 654 PetscFunctionReturn(0); 655 } 656