1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 237e1895aSJed Brown 319fd82e9SBarry Smith static SNESMSType SNESMSDefault = SNESMSM62; 437e1895aSJed Brown static PetscBool SNESMSRegisterAllCalled; 537e1895aSJed Brown static PetscBool SNESMSPackageInitialized; 637e1895aSJed Brown 737e1895aSJed Brown typedef struct _SNESMSTableau *SNESMSTableau; 837e1895aSJed Brown struct _SNESMSTableau { 937e1895aSJed Brown char *name; 1037e1895aSJed Brown PetscInt nstages; /* Number of stages */ 1137e1895aSJed Brown PetscInt nregisters; /* Number of registers */ 1237e1895aSJed Brown PetscReal stability; /* Scaled stability region */ 1337e1895aSJed Brown PetscReal *gamma; /* Coefficients of 3S* method */ 1437e1895aSJed Brown PetscReal *delta; /* Coefficients of 3S* method */ 1537e1895aSJed Brown PetscReal *betasub; /* Subdiagonal of beta in Shu-Osher form */ 1637e1895aSJed Brown }; 1737e1895aSJed Brown 1837e1895aSJed Brown typedef struct _SNESMSTableauLink *SNESMSTableauLink; 1937e1895aSJed Brown struct _SNESMSTableauLink { 2037e1895aSJed Brown struct _SNESMSTableau tab; 2137e1895aSJed Brown SNESMSTableauLink next; 2237e1895aSJed Brown }; 2337e1895aSJed Brown static SNESMSTableauLink SNESMSTableauList; 2437e1895aSJed Brown 2537e1895aSJed Brown typedef struct { 2637e1895aSJed Brown SNESMSTableau tableau; /* Tableau in low-storage form */ 2737e1895aSJed Brown PetscReal damping; /* Damping parameter, like length of (pseudo) time step */ 2837e1895aSJed Brown PetscBool norms; /* Compute norms, usually only for monitoring purposes */ 2937e1895aSJed Brown } SNES_MS; 3037e1895aSJed Brown 3137e1895aSJed Brown /*@C 32f6dfbefdSBarry Smith SNESMSRegisterAll - Registers all of the multi-stage methods in `SNESMS` 3337e1895aSJed Brown 34f6dfbefdSBarry Smith Logically Collective 3537e1895aSJed Brown 36f6dfbefdSBarry Smith Level: developer 3737e1895aSJed Brown 38f6dfbefdSBarry Smith .seealso: `SNES`, `SNESMS`, `SNESMSRegisterDestroy()` 3937e1895aSJed Brown @*/ 409371c9d4SSatish Balay PetscErrorCode SNESMSRegisterAll(void) { 4137e1895aSJed Brown PetscFunctionBegin; 4237e1895aSJed Brown if (SNESMSRegisterAllCalled) PetscFunctionReturn(0); 4337e1895aSJed Brown SNESMSRegisterAllCalled = PETSC_TRUE; 4437e1895aSJed Brown 4537e1895aSJed Brown { 469ce202e0SLisandro Dalcin PetscReal alpha[1] = {1.0}; 479566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSEULER, 1, 1, 1.0, NULL, NULL, alpha)); 4837e1895aSJed Brown } 4937e1895aSJed Brown 5037e1895aSJed Brown { 5137e1895aSJed Brown PetscReal gamma[3][6] = { 5237e1895aSJed Brown {0.0000000000000000E+00, -7.0304722367110606E-01, -1.9836719667506464E-01, -1.6023843981863788E+00, 9.4483822882855284E-02, -1.4204296130641869E-01}, 5337e1895aSJed Brown {1.0000000000000000E+00, 1.1111025767083920E+00, 5.6150921583923230E-01, 7.4151723494934041E-01, 3.1714538168600587E-01, 4.6479276238548706E-01 }, 5437e1895aSJed Brown {0.0000000000000000E+00, 0.0000000000000000E+00, 0.0000000000000000E+00, 6.7968174970583317E-01, -4.1755042846051737E-03, -1.9115668129923846E-01} 5537e1895aSJed Brown }; 5637e1895aSJed Brown PetscReal delta[6] = {1.0000000000000000E+00, 5.3275427433201750E-01, 6.0143544663985238E-01, 4.5874077053842177E-01, 2.7544386906104651E-01, 0.0000000000000000E+00}; 5737e1895aSJed Brown PetscReal betasub[6] = {8.4753115429481929E-01, 7.4018896368655618E-01, 6.5963574086583309E-03, 4.6747795645517759E-01, 1.3314545813643919E-01, 5.3260800028018784E-01}; 589566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSM62, 6, 3, 1.0, &gamma[0][0], delta, betasub)); 5937e1895aSJed Brown } 609ce202e0SLisandro Dalcin 619ce202e0SLisandro Dalcin { /* Jameson (1983) */ 629ce202e0SLisandro Dalcin PetscReal alpha[4] = {0.25, 0.5, 0.55, 1.0}; 639566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSJAMESON83, 4, 1, 1.0, NULL, NULL, alpha)); 64a97cb6bcSJed Brown } 653847c725SLisandro Dalcin 663847c725SLisandro Dalcin { /* Van Leer, Tai, and Powell (1989) 1 stage, order 1 */ 679ce202e0SLisandro Dalcin PetscReal alpha[1] = {1.0}; 689566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP11, 1, 1, 0.5, NULL, NULL, alpha)); 693847c725SLisandro Dalcin } 70a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 2 stage, order 1 */ 719ce202e0SLisandro Dalcin PetscReal alpha[2] = {0.3333, 1.0}; 729566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP21, 2, 1, 1.0, NULL, NULL, alpha)); 73a97cb6bcSJed Brown } 74a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 3 stage, order 1 */ 759ce202e0SLisandro Dalcin PetscReal alpha[3] = {0.1481, 0.4000, 1.0}; 769566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP31, 3, 1, 1.5, NULL, NULL, alpha)); 77a97cb6bcSJed Brown } 78a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 4 stage, order 1 */ 799ce202e0SLisandro Dalcin PetscReal alpha[4] = {0.0833, 0.2069, 0.4265, 1.0}; 809566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP41, 4, 1, 2.0, NULL, NULL, alpha)); 81a97cb6bcSJed Brown } 82a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 5 stage, order 1 */ 839ce202e0SLisandro Dalcin PetscReal alpha[5] = {0.0533, 0.1263, 0.2375, 0.4414, 1.0}; 849566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP51, 5, 1, 2.5, NULL, NULL, alpha)); 85a97cb6bcSJed Brown } 86a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 6 stage, order 1 */ 879ce202e0SLisandro Dalcin PetscReal alpha[6] = {0.0370, 0.0851, 0.1521, 0.2562, 0.4512, 1.0}; 889566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP61, 6, 1, 3.0, NULL, NULL, alpha)); 89a97cb6bcSJed Brown } 9037e1895aSJed Brown PetscFunctionReturn(0); 9137e1895aSJed Brown } 9237e1895aSJed Brown 9337e1895aSJed Brown /*@C 94f6dfbefdSBarry Smith SNESMSRegisterDestroy - Frees the list of schemes that were registered by `SNESMSRegister()`. 9537e1895aSJed Brown 9637e1895aSJed Brown Not Collective 9737e1895aSJed Brown 98f6dfbefdSBarry Smith Level: developer 9937e1895aSJed Brown 100db781477SPatrick Sanan .seealso: `SNESMSRegister()`, `SNESMSRegisterAll()` 10137e1895aSJed Brown @*/ 1029371c9d4SSatish Balay PetscErrorCode SNESMSRegisterDestroy(void) { 10337e1895aSJed Brown SNESMSTableauLink link; 10437e1895aSJed Brown 10537e1895aSJed Brown PetscFunctionBegin; 10637e1895aSJed Brown while ((link = SNESMSTableauList)) { 10737e1895aSJed Brown SNESMSTableau t = &link->tab; 10837e1895aSJed Brown SNESMSTableauList = link->next; 1091aa26658SKarl Rupp 1109566063dSJacob Faibussowitsch PetscCall(PetscFree(t->name)); 1119566063dSJacob Faibussowitsch PetscCall(PetscFree(t->gamma)); 1129566063dSJacob Faibussowitsch PetscCall(PetscFree(t->delta)); 1139566063dSJacob Faibussowitsch PetscCall(PetscFree(t->betasub)); 1149566063dSJacob Faibussowitsch PetscCall(PetscFree(link)); 11537e1895aSJed Brown } 11637e1895aSJed Brown SNESMSRegisterAllCalled = PETSC_FALSE; 11737e1895aSJed Brown PetscFunctionReturn(0); 11837e1895aSJed Brown } 11937e1895aSJed Brown 12037e1895aSJed Brown /*@C 121f6dfbefdSBarry Smith SNESMSInitializePackage - This function initializes everything in the `SNESMS` package. It is called 122f6dfbefdSBarry Smith from `SNESInitializePackage()`. 12337e1895aSJed Brown 12437e1895aSJed Brown Level: developer 12537e1895aSJed Brown 126db781477SPatrick Sanan .seealso: `PetscInitialize()` 12737e1895aSJed Brown @*/ 1289371c9d4SSatish Balay PetscErrorCode SNESMSInitializePackage(void) { 12937e1895aSJed Brown PetscFunctionBegin; 13037e1895aSJed Brown if (SNESMSPackageInitialized) PetscFunctionReturn(0); 13137e1895aSJed Brown SNESMSPackageInitialized = PETSC_TRUE; 1321aa26658SKarl Rupp 1339566063dSJacob Faibussowitsch PetscCall(SNESMSRegisterAll()); 1349566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(SNESMSFinalizePackage)); 13537e1895aSJed Brown PetscFunctionReturn(0); 13637e1895aSJed Brown } 13737e1895aSJed Brown 13837e1895aSJed Brown /*@C 139f6dfbefdSBarry Smith SNESMSFinalizePackage - This function destroys everything in the `SNESMS` package. It is 140f6dfbefdSBarry Smith called from `PetscFinalize()`. 14137e1895aSJed Brown 14237e1895aSJed Brown Level: developer 14337e1895aSJed Brown 144db781477SPatrick Sanan .seealso: `PetscFinalize()` 14537e1895aSJed Brown @*/ 1469371c9d4SSatish Balay PetscErrorCode SNESMSFinalizePackage(void) { 14737e1895aSJed Brown PetscFunctionBegin; 14837e1895aSJed Brown SNESMSPackageInitialized = PETSC_FALSE; 1491aa26658SKarl Rupp 1509566063dSJacob Faibussowitsch PetscCall(SNESMSRegisterDestroy()); 15137e1895aSJed Brown PetscFunctionReturn(0); 15237e1895aSJed Brown } 15337e1895aSJed Brown 15437e1895aSJed Brown /*@C 155f6dfbefdSBarry Smith SNESMSRegister - register a multistage scheme for `SNESMS` 15637e1895aSJed Brown 157f6dfbefdSBarry Smith Logically Collective 15837e1895aSJed Brown 15937e1895aSJed Brown Input Parameters: 16037e1895aSJed Brown + name - identifier for method 16137e1895aSJed Brown . nstages - number of stages 16237e1895aSJed Brown . nregisters - number of registers used by low-storage implementation 1639ce202e0SLisandro Dalcin . stability - scaled stability region 16437e1895aSJed Brown . gamma - coefficients, see Ketcheson's paper 16537e1895aSJed Brown . delta - coefficients, see Ketcheson's paper 16637e1895aSJed Brown - betasub - subdiagonal of Shu-Osher form 16737e1895aSJed Brown 16837e1895aSJed Brown Notes: 16937e1895aSJed Brown The notation is described in Ketcheson (2010) Runge-Kutta methods with minimum storage implementations. 17037e1895aSJed Brown 1719ce202e0SLisandro Dalcin Many multistage schemes are of the form 1729ce202e0SLisandro Dalcin $ X_0 = X^{(n)} 1739ce202e0SLisandro Dalcin $ X_k = X_0 + \alpha_k * F(X_{k-1}), k = 1,\ldots,s 1749ce202e0SLisandro Dalcin $ X^{(n+1)} = X_s 1759ce202e0SLisandro Dalcin These methods can be registered with 1769ce202e0SLisandro Dalcin .vb 1779ce202e0SLisandro Dalcin SNESMSRegister("name",s,1,stability,NULL,NULL,alpha); 1789ce202e0SLisandro Dalcin .ve 1799ce202e0SLisandro Dalcin 18037e1895aSJed Brown Level: advanced 18137e1895aSJed Brown 182db781477SPatrick Sanan .seealso: `SNESMS` 18337e1895aSJed Brown @*/ 1849371c9d4SSatish Balay PetscErrorCode SNESMSRegister(SNESMSType name, PetscInt nstages, PetscInt nregisters, PetscReal stability, const PetscReal gamma[], const PetscReal delta[], const PetscReal betasub[]) { 18537e1895aSJed Brown SNESMSTableauLink link; 18637e1895aSJed Brown SNESMSTableau t; 18737e1895aSJed Brown 18837e1895aSJed Brown PetscFunctionBegin; 18937e1895aSJed Brown PetscValidCharPointer(name, 1); 19008401ef6SPierre Jolivet PetscCheck(nstages >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have at least one stage"); 1919ce202e0SLisandro Dalcin if (gamma || delta) { 19208401ef6SPierre Jolivet PetscCheck(nregisters == 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support for methods written in 3-register form"); 193064a246eSJacob Faibussowitsch PetscValidRealPointer(gamma, 5); 194064a246eSJacob Faibussowitsch PetscValidRealPointer(delta, 6); 1959ce202e0SLisandro Dalcin } else { 19608401ef6SPierre Jolivet PetscCheck(nregisters == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support for methods written in 1-register form"); 1979ce202e0SLisandro Dalcin } 198064a246eSJacob Faibussowitsch PetscValidRealPointer(betasub, 7); 19937e1895aSJed Brown 2009566063dSJacob Faibussowitsch PetscCall(SNESMSInitializePackage()); 2019566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 20237e1895aSJed Brown t = &link->tab; 2039566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &t->name)); 20437e1895aSJed Brown t->nstages = nstages; 20537e1895aSJed Brown t->nregisters = nregisters; 20637e1895aSJed Brown t->stability = stability; 2071aa26658SKarl Rupp 2089ce202e0SLisandro Dalcin if (gamma && delta) { 2099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nstages * nregisters, &t->gamma)); 2109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nstages, &t->delta)); 2119566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->gamma, gamma, nstages * nregisters)); 2129566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->delta, delta, nstages)); 2139ce202e0SLisandro Dalcin } 2149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nstages, &t->betasub)); 2159566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->betasub, betasub, nstages)); 21637e1895aSJed Brown 21737e1895aSJed Brown link->next = SNESMSTableauList; 21837e1895aSJed Brown SNESMSTableauList = link; 21937e1895aSJed Brown PetscFunctionReturn(0); 22037e1895aSJed Brown } 22137e1895aSJed Brown 22237e1895aSJed Brown /* 22337e1895aSJed Brown X - initial state, updated in-place. 22437e1895aSJed Brown F - residual, computed at the initial X on input 22537e1895aSJed Brown */ 2269371c9d4SSatish Balay static PetscErrorCode SNESMSStep_3Sstar(SNES snes, Vec X, Vec F) { 22737e1895aSJed Brown SNES_MS *ms = (SNES_MS *)snes->data; 22837e1895aSJed Brown SNESMSTableau t = ms->tableau; 22937e1895aSJed Brown const PetscReal *gamma = t->gamma, *delta = t->delta, *betasub = t->betasub; 23037e1895aSJed Brown Vec S1, S2, S3, Y; 231bf80fd75SLisandro Dalcin PetscInt i, nstages = t->nstages; 23237e1895aSJed Brown 23337e1895aSJed Brown PetscFunctionBegin; 23437e1895aSJed Brown Y = snes->work[0]; 23537e1895aSJed Brown S1 = X; 23637e1895aSJed Brown S2 = snes->work[1]; 23737e1895aSJed Brown S3 = snes->work[2]; 2389566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(S2)); 2399566063dSJacob Faibussowitsch PetscCall(VecCopy(X, S3)); 24037e1895aSJed Brown for (i = 0; i < nstages; i++) { 241da80777bSKarl Rupp Vec Ss[4]; 242da80777bSKarl Rupp PetscScalar scoeff[4]; 243da80777bSKarl Rupp 2449371c9d4SSatish Balay Ss[0] = S1; 2459371c9d4SSatish Balay Ss[1] = S2; 2469371c9d4SSatish Balay Ss[2] = S3; 2479371c9d4SSatish Balay Ss[3] = Y; 2481aa26658SKarl Rupp 249bf80fd75SLisandro Dalcin scoeff[0] = gamma[0 * nstages + i] - 1; 250da80777bSKarl Rupp scoeff[1] = gamma[1 * nstages + i]; 251da80777bSKarl Rupp scoeff[2] = gamma[2 * nstages + i]; 252da80777bSKarl Rupp scoeff[3] = -betasub[i] * ms->damping; 2531aa26658SKarl Rupp 2549566063dSJacob Faibussowitsch PetscCall(VecAXPY(S2, delta[i], S1)); 25548a46eb9SPierre Jolivet if (i > 0) PetscCall(SNESComputeFunction(snes, S1, F)); 2569566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, Y)); 2579566063dSJacob Faibussowitsch PetscCall(VecMAXPY(S1, 4, scoeff, Ss)); 25837e1895aSJed Brown } 25937e1895aSJed Brown PetscFunctionReturn(0); 26037e1895aSJed Brown } 26137e1895aSJed Brown 2629ce202e0SLisandro Dalcin /* 2639ce202e0SLisandro Dalcin X - initial state, updated in-place. 2649ce202e0SLisandro Dalcin F - residual, computed at the initial X on input 2659ce202e0SLisandro Dalcin */ 2669371c9d4SSatish Balay static PetscErrorCode SNESMSStep_Basic(SNES snes, Vec X, Vec F) { 2679ce202e0SLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 2689ce202e0SLisandro Dalcin SNESMSTableau tab = ms->tableau; 2699ce202e0SLisandro Dalcin const PetscReal *alpha = tab->betasub, h = ms->damping; 2709ce202e0SLisandro Dalcin PetscInt i, nstages = tab->nstages; 2719ce202e0SLisandro Dalcin Vec X0 = snes->work[0]; 2729ce202e0SLisandro Dalcin 2739ce202e0SLisandro Dalcin PetscFunctionBegin; 2749566063dSJacob Faibussowitsch PetscCall(VecCopy(X, X0)); 2759ce202e0SLisandro Dalcin for (i = 0; i < nstages; i++) { 27648a46eb9SPierre Jolivet if (i > 0) PetscCall(SNESComputeFunction(snes, X, F)); 2779566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, X)); 2789566063dSJacob Faibussowitsch PetscCall(VecAYPX(X, -alpha[i] * h, X0)); 2799ce202e0SLisandro Dalcin } 2809ce202e0SLisandro Dalcin PetscFunctionReturn(0); 2819ce202e0SLisandro Dalcin } 2829ce202e0SLisandro Dalcin 2839371c9d4SSatish Balay static PetscErrorCode SNESMSStep_Step(SNES snes, Vec X, Vec F) { 2849ce202e0SLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 2859ce202e0SLisandro Dalcin SNESMSTableau tab = ms->tableau; 2869ce202e0SLisandro Dalcin 2879ce202e0SLisandro Dalcin PetscFunctionBegin; 2889ce202e0SLisandro Dalcin if (tab->gamma && tab->delta) { 2899566063dSJacob Faibussowitsch PetscCall(SNESMSStep_3Sstar(snes, X, F)); 2909ce202e0SLisandro Dalcin } else { 2919566063dSJacob Faibussowitsch PetscCall(SNESMSStep_Basic(snes, X, F)); 2929ce202e0SLisandro Dalcin } 2939ce202e0SLisandro Dalcin PetscFunctionReturn(0); 2949ce202e0SLisandro Dalcin } 2959ce202e0SLisandro Dalcin 2969371c9d4SSatish Balay static PetscErrorCode SNESMSStep_Norms(SNES snes, PetscInt iter, Vec F) { 297bf80fd75SLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 298bf80fd75SLisandro Dalcin PetscReal fnorm; 299bf80fd75SLisandro Dalcin 300bf80fd75SLisandro Dalcin PetscFunctionBegin; 301bf80fd75SLisandro Dalcin if (ms->norms) { 3029566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 303bf80fd75SLisandro Dalcin SNESCheckFunctionNorm(snes, fnorm); 304bf80fd75SLisandro Dalcin /* Monitor convergence */ 3059566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 306bf80fd75SLisandro Dalcin snes->iter = iter; 307bf80fd75SLisandro Dalcin snes->norm = fnorm; 3089566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3099566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 3109566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 311bf80fd75SLisandro Dalcin /* Test for convergence */ 312dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 313bf80fd75SLisandro Dalcin } else if (iter > 0) { 3149566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 315bf80fd75SLisandro Dalcin snes->iter = iter; 3169566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 317bf80fd75SLisandro Dalcin } 318bf80fd75SLisandro Dalcin PetscFunctionReturn(0); 319bf80fd75SLisandro Dalcin } 320bf80fd75SLisandro Dalcin 3219371c9d4SSatish Balay static PetscErrorCode SNESSolve_MS(SNES snes) { 32237e1895aSJed Brown SNES_MS *ms = (SNES_MS *)snes->data; 32337e1895aSJed Brown Vec X = snes->vec_sol, F = snes->vec_func; 32437e1895aSJed Brown PetscInt i; 32537e1895aSJed Brown 32637e1895aSJed Brown PetscFunctionBegin; 3270b121fc5SBarry Smith 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); 3289566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 329bf80fd75SLisandro Dalcin 33037e1895aSJed Brown snes->reason = SNES_CONVERGED_ITERATING; 3319566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 33237e1895aSJed Brown snes->iter = 0; 333bf80fd75SLisandro Dalcin snes->norm = 0; 3349566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 335bf80fd75SLisandro Dalcin 336e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 3379566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 3381aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 3391aa26658SKarl Rupp 3409566063dSJacob Faibussowitsch PetscCall(SNESMSStep_Norms(snes, 0, F)); 34137e1895aSJed Brown if (snes->reason) PetscFunctionReturn(0); 342bf80fd75SLisandro Dalcin 343bf80fd75SLisandro Dalcin for (i = 0; i < snes->max_its; i++) { 34437e1895aSJed Brown /* Call general purpose update function */ 345dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 346bf80fd75SLisandro Dalcin 347bf80fd75SLisandro Dalcin if (i == 0 && snes->jacobian) { 348bf80fd75SLisandro Dalcin /* This method does not require a Jacobian, but it is usually preconditioned by PBJacobi */ 3499566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, snes->vec_sol, snes->jacobian, snes->jacobian_pre)); 350bf80fd75SLisandro Dalcin SNESCheckJacobianDomainerror(snes); 3519566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 352bf80fd75SLisandro Dalcin } 353bf80fd75SLisandro Dalcin 3549566063dSJacob Faibussowitsch PetscCall(SNESMSStep_Step(snes, X, F)); 35537e1895aSJed Brown 35648a46eb9SPierre Jolivet if (i < snes->max_its - 1 || ms->norms) PetscCall(SNESComputeFunction(snes, X, F)); 35737e1895aSJed Brown 3589566063dSJacob Faibussowitsch PetscCall(SNESMSStep_Norms(snes, i + 1, F)); 35937e1895aSJed Brown if (snes->reason) PetscFunctionReturn(0); 36037e1895aSJed Brown } 36137e1895aSJed Brown if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; 36237e1895aSJed Brown PetscFunctionReturn(0); 36337e1895aSJed Brown } 36437e1895aSJed Brown 3659371c9d4SSatish Balay static PetscErrorCode SNESSetUp_MS(SNES snes) { 3669ce202e0SLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 3679ce202e0SLisandro Dalcin SNESMSTableau tab = ms->tableau; 3689ce202e0SLisandro Dalcin PetscInt nwork = tab->nregisters; 36937e1895aSJed Brown 37037e1895aSJed Brown PetscFunctionBegin; 3719566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, nwork)); 3729566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 37337e1895aSJed Brown PetscFunctionReturn(0); 37437e1895aSJed Brown } 37537e1895aSJed Brown 3769371c9d4SSatish Balay static PetscErrorCode SNESReset_MS(SNES snes) { 37737e1895aSJed Brown PetscFunctionBegin; 37837e1895aSJed Brown PetscFunctionReturn(0); 37937e1895aSJed Brown } 38037e1895aSJed Brown 3819371c9d4SSatish Balay static PetscErrorCode SNESDestroy_MS(SNES snes) { 38237e1895aSJed Brown PetscFunctionBegin; 3839566063dSJacob Faibussowitsch PetscCall(SNESReset_MS(snes)); 3849566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 3859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetType_C", NULL)); 3869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetType_C", NULL)); 3879566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetDamping_C", NULL)); 3889566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetDamping_C", NULL)); 38937e1895aSJed Brown PetscFunctionReturn(0); 39037e1895aSJed Brown } 39137e1895aSJed Brown 3929371c9d4SSatish Balay static PetscErrorCode SNESView_MS(SNES snes, PetscViewer viewer) { 39337e1895aSJed Brown PetscBool iascii; 39437e1895aSJed Brown SNES_MS *ms = (SNES_MS *)snes->data; 39557715debSLisandro Dalcin SNESMSTableau tab = ms->tableau; 39637e1895aSJed Brown 39737e1895aSJed Brown PetscFunctionBegin; 3989566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 39948a46eb9SPierre Jolivet if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, " multi-stage method type: %s\n", tab->name)); 40037e1895aSJed Brown PetscFunctionReturn(0); 40137e1895aSJed Brown } 40237e1895aSJed Brown 4039371c9d4SSatish Balay static PetscErrorCode SNESSetFromOptions_MS(SNES snes, PetscOptionItems *PetscOptionsObject) { 404725b86d8SJed Brown SNES_MS *ms = (SNES_MS *)snes->data; 40537e1895aSJed Brown 40637e1895aSJed Brown PetscFunctionBegin; 407d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES MS options"); 40837e1895aSJed Brown { 40937e1895aSJed Brown SNESMSTableauLink link; 41037e1895aSJed Brown PetscInt count, choice; 41137e1895aSJed Brown PetscBool flg; 41237e1895aSJed Brown const char **namelist; 41357715debSLisandro Dalcin SNESMSType mstype; 41457715debSLisandro Dalcin PetscReal damping; 41537e1895aSJed Brown 4169566063dSJacob Faibussowitsch PetscCall(SNESMSGetType(snes, &mstype)); 4179371c9d4SSatish Balay for (link = SNESMSTableauList, count = 0; link; link = link->next, count++) 4189371c9d4SSatish Balay ; 4199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, (char ***)&namelist)); 42037e1895aSJed Brown for (link = SNESMSTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name; 4219566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-snes_ms_type", "Multistage smoother type", "SNESMSSetType", (const char *const *)namelist, count, mstype, &choice, &flg)); 4229566063dSJacob Faibussowitsch if (flg) PetscCall(SNESMSSetType(snes, namelist[choice])); 4239566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 4249566063dSJacob Faibussowitsch PetscCall(SNESMSGetDamping(snes, &damping)); 4259566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ms_damping", "Damping for multistage method", "SNESMSSetDamping", damping, &damping, &flg)); 4269566063dSJacob Faibussowitsch if (flg) PetscCall(SNESMSSetDamping(snes, damping)); 4279566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ms_norms", "Compute norms for monitoring", "none", ms->norms, &ms->norms, NULL)); 42837e1895aSJed Brown } 429d0609cedSBarry Smith PetscOptionsHeadEnd(); 43037e1895aSJed Brown PetscFunctionReturn(0); 43137e1895aSJed Brown } 43237e1895aSJed Brown 4339371c9d4SSatish Balay static PetscErrorCode SNESMSGetType_MS(SNES snes, SNESMSType *mstype) { 43457715debSLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 43557715debSLisandro Dalcin SNESMSTableau tab = ms->tableau; 43657715debSLisandro Dalcin 43757715debSLisandro Dalcin PetscFunctionBegin; 43857715debSLisandro Dalcin *mstype = tab->name; 43957715debSLisandro Dalcin PetscFunctionReturn(0); 44057715debSLisandro Dalcin } 44157715debSLisandro Dalcin 4429371c9d4SSatish Balay static PetscErrorCode SNESMSSetType_MS(SNES snes, SNESMSType mstype) { 44337e1895aSJed Brown SNES_MS *ms = (SNES_MS *)snes->data; 44437e1895aSJed Brown SNESMSTableauLink link; 44537e1895aSJed Brown PetscBool match; 44637e1895aSJed Brown 44737e1895aSJed Brown PetscFunctionBegin; 44837e1895aSJed Brown if (ms->tableau) { 4499566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ms->tableau->name, mstype, &match)); 45037e1895aSJed Brown if (match) PetscFunctionReturn(0); 45137e1895aSJed Brown } 45237e1895aSJed Brown for (link = SNESMSTableauList; link; link = link->next) { 4539566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name, mstype, &match)); 45437e1895aSJed Brown if (match) { 4559566063dSJacob Faibussowitsch if (snes->setupcalled) PetscCall(SNESReset_MS(snes)); 45637e1895aSJed Brown ms->tableau = &link->tab; 4579566063dSJacob Faibussowitsch if (snes->setupcalled) PetscCall(SNESSetUp_MS(snes)); 45837e1895aSJed Brown PetscFunctionReturn(0); 45937e1895aSJed Brown } 46037e1895aSJed Brown } 46198921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", mstype); 46237e1895aSJed Brown } 46337e1895aSJed Brown 46437e1895aSJed Brown /*@C 465f6dfbefdSBarry Smith SNESMSGetType - Get the type of multistage smoother `SNESMS` 46657715debSLisandro Dalcin 46757715debSLisandro Dalcin Not collective 46857715debSLisandro Dalcin 46957715debSLisandro Dalcin Input Parameter: 47057715debSLisandro Dalcin . snes - nonlinear solver context 47157715debSLisandro Dalcin 47257715debSLisandro Dalcin Output Parameter: 47357715debSLisandro Dalcin . mstype - type of multistage method 47457715debSLisandro Dalcin 475f6dfbefdSBarry Smith Level: advanced 47657715debSLisandro Dalcin 477f6dfbefdSBarry Smith .seealso: `SNESMS`, `SNESMSSetType()`, `SNESMSType`, `SNESMS` 47857715debSLisandro Dalcin @*/ 4799371c9d4SSatish Balay PetscErrorCode SNESMSGetType(SNES snes, SNESMSType *mstype) { 48057715debSLisandro Dalcin PetscFunctionBegin; 48157715debSLisandro Dalcin PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 48257715debSLisandro Dalcin PetscValidPointer(mstype, 2); 483cac4c232SBarry Smith PetscUseMethod(snes, "SNESMSGetType_C", (SNES, SNESMSType *), (snes, mstype)); 48457715debSLisandro Dalcin PetscFunctionReturn(0); 48557715debSLisandro Dalcin } 48657715debSLisandro Dalcin 48757715debSLisandro Dalcin /*@C 488f6dfbefdSBarry Smith SNESMSSetType - Set the type of multistage smoother `SNESMS` 48937e1895aSJed Brown 49037e1895aSJed Brown Logically collective 49137e1895aSJed Brown 492d8d19677SJose E. Roman Input Parameters: 49337e1895aSJed Brown + snes - nonlinear solver context 49437e1895aSJed Brown - mstype - type of multistage method 49537e1895aSJed Brown 496f6dfbefdSBarry Smith Level: advanced 49737e1895aSJed Brown 498f6dfbefdSBarry Smith .seealso: `SNESMS`, `SNESMSGetType()`, `SNESMSType`, `SNESMS` 49937e1895aSJed Brown @*/ 5009371c9d4SSatish Balay PetscErrorCode SNESMSSetType(SNES snes, SNESMSType mstype) { 50137e1895aSJed Brown PetscFunctionBegin; 50237e1895aSJed Brown PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 503dadcf809SJacob Faibussowitsch PetscValidCharPointer(mstype, 2); 504cac4c232SBarry Smith PetscTryMethod(snes, "SNESMSSetType_C", (SNES, SNESMSType), (snes, mstype)); 50557715debSLisandro Dalcin PetscFunctionReturn(0); 50657715debSLisandro Dalcin } 50757715debSLisandro Dalcin 5089371c9d4SSatish Balay static PetscErrorCode SNESMSGetDamping_MS(SNES snes, PetscReal *damping) { 50957715debSLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 51057715debSLisandro Dalcin 51157715debSLisandro Dalcin PetscFunctionBegin; 51257715debSLisandro Dalcin *damping = ms->damping; 51357715debSLisandro Dalcin PetscFunctionReturn(0); 51457715debSLisandro Dalcin } 51557715debSLisandro Dalcin 5169371c9d4SSatish Balay static PetscErrorCode SNESMSSetDamping_MS(SNES snes, PetscReal damping) { 51757715debSLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 51857715debSLisandro Dalcin 51957715debSLisandro Dalcin PetscFunctionBegin; 52057715debSLisandro Dalcin ms->damping = damping; 52157715debSLisandro Dalcin PetscFunctionReturn(0); 52257715debSLisandro Dalcin } 52357715debSLisandro Dalcin 52457715debSLisandro Dalcin /*@C 525f6dfbefdSBarry Smith SNESMSGetDamping - Get the damping parameter of `SNESMS` multistage scheme 52657715debSLisandro Dalcin 52757715debSLisandro Dalcin Not collective 52857715debSLisandro Dalcin 52957715debSLisandro Dalcin Input Parameter: 53057715debSLisandro Dalcin . snes - nonlinear solver context 53157715debSLisandro Dalcin 53257715debSLisandro Dalcin Output Parameter: 53357715debSLisandro Dalcin . damping - damping parameter 53457715debSLisandro Dalcin 53557715debSLisandro Dalcin Level: advanced 53657715debSLisandro Dalcin 537db781477SPatrick Sanan .seealso: `SNESMSSetDamping()`, `SNESMS` 53857715debSLisandro Dalcin @*/ 5399371c9d4SSatish Balay PetscErrorCode SNESMSGetDamping(SNES snes, PetscReal *damping) { 54057715debSLisandro Dalcin PetscFunctionBegin; 54157715debSLisandro Dalcin PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 542dadcf809SJacob Faibussowitsch PetscValidRealPointer(damping, 2); 543cac4c232SBarry Smith PetscUseMethod(snes, "SNESMSGetDamping_C", (SNES, PetscReal *), (snes, damping)); 54457715debSLisandro Dalcin PetscFunctionReturn(0); 54557715debSLisandro Dalcin } 54657715debSLisandro Dalcin 54757715debSLisandro Dalcin /*@C 548f6dfbefdSBarry Smith SNESMSSetDamping - Set the damping parameter for a `SNESMS` multistage scheme 54957715debSLisandro Dalcin 55057715debSLisandro Dalcin Logically collective 55157715debSLisandro Dalcin 552d8d19677SJose E. Roman Input Parameters: 55357715debSLisandro Dalcin + snes - nonlinear solver context 55457715debSLisandro Dalcin - damping - damping parameter 55557715debSLisandro Dalcin 55657715debSLisandro Dalcin Level: advanced 55757715debSLisandro Dalcin 558db781477SPatrick Sanan .seealso: `SNESMSGetDamping()`, `SNESMS` 55957715debSLisandro Dalcin @*/ 5609371c9d4SSatish Balay PetscErrorCode SNESMSSetDamping(SNES snes, PetscReal damping) { 56157715debSLisandro Dalcin PetscFunctionBegin; 56257715debSLisandro Dalcin PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 56357715debSLisandro Dalcin PetscValidLogicalCollectiveReal(snes, damping, 2); 564cac4c232SBarry Smith PetscTryMethod(snes, "SNESMSSetDamping_C", (SNES, PetscReal), (snes, damping)); 56537e1895aSJed Brown PetscFunctionReturn(0); 56637e1895aSJed Brown } 56737e1895aSJed Brown 56837e1895aSJed Brown /*MC 56937e1895aSJed Brown SNESMS - multi-stage smoothers 57037e1895aSJed Brown 571f6dfbefdSBarry Smith Options Database Keys: 57237e1895aSJed Brown + -snes_ms_type - type of multi-stage smoother 57337e1895aSJed Brown - -snes_ms_damping - damping for multi-stage method 57437e1895aSJed Brown 57537e1895aSJed Brown Notes: 57637e1895aSJed Brown These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for 5776c9de887SHong Zhang FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebyshev). 57837e1895aSJed Brown 57937e1895aSJed Brown Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds. 58037e1895aSJed Brown 581f6dfbefdSBarry Smith The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with `SNESMSRegister()`. 58237e1895aSJed Brown 58337e1895aSJed Brown References: 584606c0280SSatish Balay + * - Ketcheson (2010) Runge Kutta methods with minimum storage implementations (https://doi.org/10.1016/j.jcp.2009.11.006). 585606c0280SSatish Balay . * - 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). 586606c0280SSatish Balay . * - Pierce and Giles (1997) Preconditioned multigrid methods for compressible flow calculations on stretched meshes (https://doi.org/10.1006/jcph.1997.5772). 587606c0280SSatish Balay - * - 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). 58837e1895aSJed Brown 589f6dfbefdSBarry Smith Level: advanced 59037e1895aSJed Brown 591db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESMS`, `SNESFAS`, `KSPCHEBYSHEV` 59237e1895aSJed Brown 59337e1895aSJed Brown M*/ 5949371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_MS(SNES snes) { 59537e1895aSJed Brown SNES_MS *ms; 59637e1895aSJed Brown 59737e1895aSJed Brown PetscFunctionBegin; 5989566063dSJacob Faibussowitsch PetscCall(SNESMSInitializePackage()); 59937e1895aSJed Brown 60037e1895aSJed Brown snes->ops->setup = SNESSetUp_MS; 60137e1895aSJed Brown snes->ops->solve = SNESSolve_MS; 60237e1895aSJed Brown snes->ops->destroy = SNESDestroy_MS; 60337e1895aSJed Brown snes->ops->setfromoptions = SNESSetFromOptions_MS; 60437e1895aSJed Brown snes->ops->view = SNESView_MS; 60537e1895aSJed Brown snes->ops->reset = SNESReset_MS; 60637e1895aSJed Brown 607efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 60837e1895aSJed Brown snes->usesksp = PETSC_TRUE; 60937e1895aSJed Brown 6104fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 6114fc747eaSLawrence Mitchell 612*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ms)); 61337e1895aSJed Brown snes->data = (void *)ms; 61437e1895aSJed Brown ms->damping = 0.9; 61537e1895aSJed Brown ms->norms = PETSC_FALSE; 61637e1895aSJed Brown 6179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetType_C", SNESMSGetType_MS)); 6189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetType_C", SNESMSSetType_MS)); 6199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetDamping_C", SNESMSGetDamping_MS)); 6209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetDamping_C", SNESMSSetDamping_MS)); 62157715debSLisandro Dalcin 6229566063dSJacob Faibussowitsch PetscCall(SNESMSSetType(snes, SNESMSDefault)); 62337e1895aSJed Brown PetscFunctionReturn(0); 62437e1895aSJed Brown } 625