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 @*/ 40*d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSRegisterAll(void) 41*d71ae5a4SJacob Faibussowitsch { 4237e1895aSJed Brown PetscFunctionBegin; 4337e1895aSJed Brown if (SNESMSRegisterAllCalled) PetscFunctionReturn(0); 4437e1895aSJed Brown SNESMSRegisterAllCalled = PETSC_TRUE; 4537e1895aSJed Brown 4637e1895aSJed Brown { 479ce202e0SLisandro Dalcin PetscReal alpha[1] = {1.0}; 489566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSEULER, 1, 1, 1.0, NULL, NULL, alpha)); 4937e1895aSJed Brown } 5037e1895aSJed Brown 5137e1895aSJed Brown { 5237e1895aSJed Brown PetscReal gamma[3][6] = { 5337e1895aSJed Brown {0.0000000000000000E+00, -7.0304722367110606E-01, -1.9836719667506464E-01, -1.6023843981863788E+00, 9.4483822882855284E-02, -1.4204296130641869E-01}, 5437e1895aSJed Brown {1.0000000000000000E+00, 1.1111025767083920E+00, 5.6150921583923230E-01, 7.4151723494934041E-01, 3.1714538168600587E-01, 4.6479276238548706E-01 }, 5537e1895aSJed Brown {0.0000000000000000E+00, 0.0000000000000000E+00, 0.0000000000000000E+00, 6.7968174970583317E-01, -4.1755042846051737E-03, -1.9115668129923846E-01} 5637e1895aSJed Brown }; 5737e1895aSJed Brown PetscReal delta[6] = {1.0000000000000000E+00, 5.3275427433201750E-01, 6.0143544663985238E-01, 4.5874077053842177E-01, 2.7544386906104651E-01, 0.0000000000000000E+00}; 5837e1895aSJed Brown PetscReal betasub[6] = {8.4753115429481929E-01, 7.4018896368655618E-01, 6.5963574086583309E-03, 4.6747795645517759E-01, 1.3314545813643919E-01, 5.3260800028018784E-01}; 599566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSM62, 6, 3, 1.0, &gamma[0][0], delta, betasub)); 6037e1895aSJed Brown } 619ce202e0SLisandro Dalcin 629ce202e0SLisandro Dalcin { /* Jameson (1983) */ 639ce202e0SLisandro Dalcin PetscReal alpha[4] = {0.25, 0.5, 0.55, 1.0}; 649566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSJAMESON83, 4, 1, 1.0, NULL, NULL, alpha)); 65a97cb6bcSJed Brown } 663847c725SLisandro Dalcin 673847c725SLisandro Dalcin { /* Van Leer, Tai, and Powell (1989) 1 stage, order 1 */ 689ce202e0SLisandro Dalcin PetscReal alpha[1] = {1.0}; 699566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP11, 1, 1, 0.5, NULL, NULL, alpha)); 703847c725SLisandro Dalcin } 71a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 2 stage, order 1 */ 729ce202e0SLisandro Dalcin PetscReal alpha[2] = {0.3333, 1.0}; 739566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP21, 2, 1, 1.0, NULL, NULL, alpha)); 74a97cb6bcSJed Brown } 75a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 3 stage, order 1 */ 769ce202e0SLisandro Dalcin PetscReal alpha[3] = {0.1481, 0.4000, 1.0}; 779566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP31, 3, 1, 1.5, NULL, NULL, alpha)); 78a97cb6bcSJed Brown } 79a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 4 stage, order 1 */ 809ce202e0SLisandro Dalcin PetscReal alpha[4] = {0.0833, 0.2069, 0.4265, 1.0}; 819566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP41, 4, 1, 2.0, NULL, NULL, alpha)); 82a97cb6bcSJed Brown } 83a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 5 stage, order 1 */ 849ce202e0SLisandro Dalcin PetscReal alpha[5] = {0.0533, 0.1263, 0.2375, 0.4414, 1.0}; 859566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP51, 5, 1, 2.5, NULL, NULL, alpha)); 86a97cb6bcSJed Brown } 87a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 6 stage, order 1 */ 889ce202e0SLisandro Dalcin PetscReal alpha[6] = {0.0370, 0.0851, 0.1521, 0.2562, 0.4512, 1.0}; 899566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP61, 6, 1, 3.0, NULL, NULL, alpha)); 90a97cb6bcSJed Brown } 9137e1895aSJed Brown PetscFunctionReturn(0); 9237e1895aSJed Brown } 9337e1895aSJed Brown 9437e1895aSJed Brown /*@C 95f6dfbefdSBarry Smith SNESMSRegisterDestroy - Frees the list of schemes that were registered by `SNESMSRegister()`. 9637e1895aSJed Brown 9737e1895aSJed Brown Not Collective 9837e1895aSJed Brown 99f6dfbefdSBarry Smith Level: developer 10037e1895aSJed Brown 101db781477SPatrick Sanan .seealso: `SNESMSRegister()`, `SNESMSRegisterAll()` 10237e1895aSJed Brown @*/ 103*d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSRegisterDestroy(void) 104*d71ae5a4SJacob Faibussowitsch { 10537e1895aSJed Brown SNESMSTableauLink link; 10637e1895aSJed Brown 10737e1895aSJed Brown PetscFunctionBegin; 10837e1895aSJed Brown while ((link = SNESMSTableauList)) { 10937e1895aSJed Brown SNESMSTableau t = &link->tab; 11037e1895aSJed Brown SNESMSTableauList = link->next; 1111aa26658SKarl Rupp 1129566063dSJacob Faibussowitsch PetscCall(PetscFree(t->name)); 1139566063dSJacob Faibussowitsch PetscCall(PetscFree(t->gamma)); 1149566063dSJacob Faibussowitsch PetscCall(PetscFree(t->delta)); 1159566063dSJacob Faibussowitsch PetscCall(PetscFree(t->betasub)); 1169566063dSJacob Faibussowitsch PetscCall(PetscFree(link)); 11737e1895aSJed Brown } 11837e1895aSJed Brown SNESMSRegisterAllCalled = PETSC_FALSE; 11937e1895aSJed Brown PetscFunctionReturn(0); 12037e1895aSJed Brown } 12137e1895aSJed Brown 12237e1895aSJed Brown /*@C 123f6dfbefdSBarry Smith SNESMSInitializePackage - This function initializes everything in the `SNESMS` package. It is called 124f6dfbefdSBarry Smith from `SNESInitializePackage()`. 12537e1895aSJed Brown 12637e1895aSJed Brown Level: developer 12737e1895aSJed Brown 128db781477SPatrick Sanan .seealso: `PetscInitialize()` 12937e1895aSJed Brown @*/ 130*d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSInitializePackage(void) 131*d71ae5a4SJacob Faibussowitsch { 13237e1895aSJed Brown PetscFunctionBegin; 13337e1895aSJed Brown if (SNESMSPackageInitialized) PetscFunctionReturn(0); 13437e1895aSJed Brown SNESMSPackageInitialized = PETSC_TRUE; 1351aa26658SKarl Rupp 1369566063dSJacob Faibussowitsch PetscCall(SNESMSRegisterAll()); 1379566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(SNESMSFinalizePackage)); 13837e1895aSJed Brown PetscFunctionReturn(0); 13937e1895aSJed Brown } 14037e1895aSJed Brown 14137e1895aSJed Brown /*@C 142f6dfbefdSBarry Smith SNESMSFinalizePackage - This function destroys everything in the `SNESMS` package. It is 143f6dfbefdSBarry Smith called from `PetscFinalize()`. 14437e1895aSJed Brown 14537e1895aSJed Brown Level: developer 14637e1895aSJed Brown 147db781477SPatrick Sanan .seealso: `PetscFinalize()` 14837e1895aSJed Brown @*/ 149*d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSFinalizePackage(void) 150*d71ae5a4SJacob Faibussowitsch { 15137e1895aSJed Brown PetscFunctionBegin; 15237e1895aSJed Brown SNESMSPackageInitialized = PETSC_FALSE; 1531aa26658SKarl Rupp 1549566063dSJacob Faibussowitsch PetscCall(SNESMSRegisterDestroy()); 15537e1895aSJed Brown PetscFunctionReturn(0); 15637e1895aSJed Brown } 15737e1895aSJed Brown 15837e1895aSJed Brown /*@C 159f6dfbefdSBarry Smith SNESMSRegister - register a multistage scheme for `SNESMS` 16037e1895aSJed Brown 161f6dfbefdSBarry Smith Logically Collective 16237e1895aSJed Brown 16337e1895aSJed Brown Input Parameters: 16437e1895aSJed Brown + name - identifier for method 16537e1895aSJed Brown . nstages - number of stages 16637e1895aSJed Brown . nregisters - number of registers used by low-storage implementation 1679ce202e0SLisandro Dalcin . stability - scaled stability region 16837e1895aSJed Brown . gamma - coefficients, see Ketcheson's paper 16937e1895aSJed Brown . delta - coefficients, see Ketcheson's paper 17037e1895aSJed Brown - betasub - subdiagonal of Shu-Osher form 17137e1895aSJed Brown 17237e1895aSJed Brown Notes: 17337e1895aSJed Brown The notation is described in Ketcheson (2010) Runge-Kutta methods with minimum storage implementations. 17437e1895aSJed Brown 1759ce202e0SLisandro Dalcin Many multistage schemes are of the form 1769ce202e0SLisandro Dalcin $ X_0 = X^{(n)} 1779ce202e0SLisandro Dalcin $ X_k = X_0 + \alpha_k * F(X_{k-1}), k = 1,\ldots,s 1789ce202e0SLisandro Dalcin $ X^{(n+1)} = X_s 1799ce202e0SLisandro Dalcin These methods can be registered with 1809ce202e0SLisandro Dalcin .vb 1819ce202e0SLisandro Dalcin SNESMSRegister("name",s,1,stability,NULL,NULL,alpha); 1829ce202e0SLisandro Dalcin .ve 1839ce202e0SLisandro Dalcin 18437e1895aSJed Brown Level: advanced 18537e1895aSJed Brown 186db781477SPatrick Sanan .seealso: `SNESMS` 18737e1895aSJed Brown @*/ 188*d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSRegister(SNESMSType name, PetscInt nstages, PetscInt nregisters, PetscReal stability, const PetscReal gamma[], const PetscReal delta[], const PetscReal betasub[]) 189*d71ae5a4SJacob Faibussowitsch { 19037e1895aSJed Brown SNESMSTableauLink link; 19137e1895aSJed Brown SNESMSTableau t; 19237e1895aSJed Brown 19337e1895aSJed Brown PetscFunctionBegin; 19437e1895aSJed Brown PetscValidCharPointer(name, 1); 19508401ef6SPierre Jolivet PetscCheck(nstages >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have at least one stage"); 1969ce202e0SLisandro Dalcin if (gamma || delta) { 19708401ef6SPierre Jolivet PetscCheck(nregisters == 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support for methods written in 3-register form"); 198064a246eSJacob Faibussowitsch PetscValidRealPointer(gamma, 5); 199064a246eSJacob Faibussowitsch PetscValidRealPointer(delta, 6); 2009ce202e0SLisandro Dalcin } else { 20108401ef6SPierre Jolivet PetscCheck(nregisters == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support for methods written in 1-register form"); 2029ce202e0SLisandro Dalcin } 203064a246eSJacob Faibussowitsch PetscValidRealPointer(betasub, 7); 20437e1895aSJed Brown 2059566063dSJacob Faibussowitsch PetscCall(SNESMSInitializePackage()); 2069566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 20737e1895aSJed Brown t = &link->tab; 2089566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &t->name)); 20937e1895aSJed Brown t->nstages = nstages; 21037e1895aSJed Brown t->nregisters = nregisters; 21137e1895aSJed Brown t->stability = stability; 2121aa26658SKarl Rupp 2139ce202e0SLisandro Dalcin if (gamma && delta) { 2149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nstages * nregisters, &t->gamma)); 2159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nstages, &t->delta)); 2169566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->gamma, gamma, nstages * nregisters)); 2179566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->delta, delta, nstages)); 2189ce202e0SLisandro Dalcin } 2199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nstages, &t->betasub)); 2209566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->betasub, betasub, nstages)); 22137e1895aSJed Brown 22237e1895aSJed Brown link->next = SNESMSTableauList; 22337e1895aSJed Brown SNESMSTableauList = link; 22437e1895aSJed Brown PetscFunctionReturn(0); 22537e1895aSJed Brown } 22637e1895aSJed Brown 22737e1895aSJed Brown /* 22837e1895aSJed Brown X - initial state, updated in-place. 22937e1895aSJed Brown F - residual, computed at the initial X on input 23037e1895aSJed Brown */ 231*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSStep_3Sstar(SNES snes, Vec X, Vec F) 232*d71ae5a4SJacob Faibussowitsch { 23337e1895aSJed Brown SNES_MS *ms = (SNES_MS *)snes->data; 23437e1895aSJed Brown SNESMSTableau t = ms->tableau; 23537e1895aSJed Brown const PetscReal *gamma = t->gamma, *delta = t->delta, *betasub = t->betasub; 23637e1895aSJed Brown Vec S1, S2, S3, Y; 237bf80fd75SLisandro Dalcin PetscInt i, nstages = t->nstages; 23837e1895aSJed Brown 23937e1895aSJed Brown PetscFunctionBegin; 24037e1895aSJed Brown Y = snes->work[0]; 24137e1895aSJed Brown S1 = X; 24237e1895aSJed Brown S2 = snes->work[1]; 24337e1895aSJed Brown S3 = snes->work[2]; 2449566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(S2)); 2459566063dSJacob Faibussowitsch PetscCall(VecCopy(X, S3)); 24637e1895aSJed Brown for (i = 0; i < nstages; i++) { 247da80777bSKarl Rupp Vec Ss[4]; 248da80777bSKarl Rupp PetscScalar scoeff[4]; 249da80777bSKarl Rupp 2509371c9d4SSatish Balay Ss[0] = S1; 2519371c9d4SSatish Balay Ss[1] = S2; 2529371c9d4SSatish Balay Ss[2] = S3; 2539371c9d4SSatish Balay Ss[3] = Y; 2541aa26658SKarl Rupp 255bf80fd75SLisandro Dalcin scoeff[0] = gamma[0 * nstages + i] - 1; 256da80777bSKarl Rupp scoeff[1] = gamma[1 * nstages + i]; 257da80777bSKarl Rupp scoeff[2] = gamma[2 * nstages + i]; 258da80777bSKarl Rupp scoeff[3] = -betasub[i] * ms->damping; 2591aa26658SKarl Rupp 2609566063dSJacob Faibussowitsch PetscCall(VecAXPY(S2, delta[i], S1)); 26148a46eb9SPierre Jolivet if (i > 0) PetscCall(SNESComputeFunction(snes, S1, F)); 2629566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, Y)); 2639566063dSJacob Faibussowitsch PetscCall(VecMAXPY(S1, 4, scoeff, Ss)); 26437e1895aSJed Brown } 26537e1895aSJed Brown PetscFunctionReturn(0); 26637e1895aSJed Brown } 26737e1895aSJed Brown 2689ce202e0SLisandro Dalcin /* 2699ce202e0SLisandro Dalcin X - initial state, updated in-place. 2709ce202e0SLisandro Dalcin F - residual, computed at the initial X on input 2719ce202e0SLisandro Dalcin */ 272*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSStep_Basic(SNES snes, Vec X, Vec F) 273*d71ae5a4SJacob Faibussowitsch { 2749ce202e0SLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 2759ce202e0SLisandro Dalcin SNESMSTableau tab = ms->tableau; 2769ce202e0SLisandro Dalcin const PetscReal *alpha = tab->betasub, h = ms->damping; 2779ce202e0SLisandro Dalcin PetscInt i, nstages = tab->nstages; 2789ce202e0SLisandro Dalcin Vec X0 = snes->work[0]; 2799ce202e0SLisandro Dalcin 2809ce202e0SLisandro Dalcin PetscFunctionBegin; 2819566063dSJacob Faibussowitsch PetscCall(VecCopy(X, X0)); 2829ce202e0SLisandro Dalcin for (i = 0; i < nstages; i++) { 28348a46eb9SPierre Jolivet if (i > 0) PetscCall(SNESComputeFunction(snes, X, F)); 2849566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, X)); 2859566063dSJacob Faibussowitsch PetscCall(VecAYPX(X, -alpha[i] * h, X0)); 2869ce202e0SLisandro Dalcin } 2879ce202e0SLisandro Dalcin PetscFunctionReturn(0); 2889ce202e0SLisandro Dalcin } 2899ce202e0SLisandro Dalcin 290*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSStep_Step(SNES snes, Vec X, Vec F) 291*d71ae5a4SJacob Faibussowitsch { 2929ce202e0SLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 2939ce202e0SLisandro Dalcin SNESMSTableau tab = ms->tableau; 2949ce202e0SLisandro Dalcin 2959ce202e0SLisandro Dalcin PetscFunctionBegin; 2969ce202e0SLisandro Dalcin if (tab->gamma && tab->delta) { 2979566063dSJacob Faibussowitsch PetscCall(SNESMSStep_3Sstar(snes, X, F)); 2989ce202e0SLisandro Dalcin } else { 2999566063dSJacob Faibussowitsch PetscCall(SNESMSStep_Basic(snes, X, F)); 3009ce202e0SLisandro Dalcin } 3019ce202e0SLisandro Dalcin PetscFunctionReturn(0); 3029ce202e0SLisandro Dalcin } 3039ce202e0SLisandro Dalcin 304*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSStep_Norms(SNES snes, PetscInt iter, Vec F) 305*d71ae5a4SJacob Faibussowitsch { 306bf80fd75SLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 307bf80fd75SLisandro Dalcin PetscReal fnorm; 308bf80fd75SLisandro Dalcin 309bf80fd75SLisandro Dalcin PetscFunctionBegin; 310bf80fd75SLisandro Dalcin if (ms->norms) { 3119566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 312bf80fd75SLisandro Dalcin SNESCheckFunctionNorm(snes, fnorm); 313bf80fd75SLisandro Dalcin /* Monitor convergence */ 3149566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 315bf80fd75SLisandro Dalcin snes->iter = iter; 316bf80fd75SLisandro Dalcin snes->norm = fnorm; 3179566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3189566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 3199566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 320bf80fd75SLisandro Dalcin /* Test for convergence */ 321dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 322bf80fd75SLisandro Dalcin } else if (iter > 0) { 3239566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 324bf80fd75SLisandro Dalcin snes->iter = iter; 3259566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 326bf80fd75SLisandro Dalcin } 327bf80fd75SLisandro Dalcin PetscFunctionReturn(0); 328bf80fd75SLisandro Dalcin } 329bf80fd75SLisandro Dalcin 330*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_MS(SNES snes) 331*d71ae5a4SJacob Faibussowitsch { 33237e1895aSJed Brown SNES_MS *ms = (SNES_MS *)snes->data; 33337e1895aSJed Brown Vec X = snes->vec_sol, F = snes->vec_func; 33437e1895aSJed Brown PetscInt i; 33537e1895aSJed Brown 33637e1895aSJed Brown PetscFunctionBegin; 3370b121fc5SBarry 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); 3389566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 339bf80fd75SLisandro Dalcin 34037e1895aSJed Brown snes->reason = SNES_CONVERGED_ITERATING; 3419566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 34237e1895aSJed Brown snes->iter = 0; 343bf80fd75SLisandro Dalcin snes->norm = 0; 3449566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 345bf80fd75SLisandro Dalcin 346e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 3479566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 3481aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 3491aa26658SKarl Rupp 3509566063dSJacob Faibussowitsch PetscCall(SNESMSStep_Norms(snes, 0, F)); 35137e1895aSJed Brown if (snes->reason) PetscFunctionReturn(0); 352bf80fd75SLisandro Dalcin 353bf80fd75SLisandro Dalcin for (i = 0; i < snes->max_its; i++) { 35437e1895aSJed Brown /* Call general purpose update function */ 355dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 356bf80fd75SLisandro Dalcin 357bf80fd75SLisandro Dalcin if (i == 0 && snes->jacobian) { 358bf80fd75SLisandro Dalcin /* This method does not require a Jacobian, but it is usually preconditioned by PBJacobi */ 3599566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, snes->vec_sol, snes->jacobian, snes->jacobian_pre)); 360bf80fd75SLisandro Dalcin SNESCheckJacobianDomainerror(snes); 3619566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 362bf80fd75SLisandro Dalcin } 363bf80fd75SLisandro Dalcin 3649566063dSJacob Faibussowitsch PetscCall(SNESMSStep_Step(snes, X, F)); 36537e1895aSJed Brown 36648a46eb9SPierre Jolivet if (i < snes->max_its - 1 || ms->norms) PetscCall(SNESComputeFunction(snes, X, F)); 36737e1895aSJed Brown 3689566063dSJacob Faibussowitsch PetscCall(SNESMSStep_Norms(snes, i + 1, F)); 36937e1895aSJed Brown if (snes->reason) PetscFunctionReturn(0); 37037e1895aSJed Brown } 37137e1895aSJed Brown if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; 37237e1895aSJed Brown PetscFunctionReturn(0); 37337e1895aSJed Brown } 37437e1895aSJed Brown 375*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_MS(SNES snes) 376*d71ae5a4SJacob Faibussowitsch { 3779ce202e0SLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 3789ce202e0SLisandro Dalcin SNESMSTableau tab = ms->tableau; 3799ce202e0SLisandro Dalcin PetscInt nwork = tab->nregisters; 38037e1895aSJed Brown 38137e1895aSJed Brown PetscFunctionBegin; 3829566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, nwork)); 3839566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 38437e1895aSJed Brown PetscFunctionReturn(0); 38537e1895aSJed Brown } 38637e1895aSJed Brown 387*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESReset_MS(SNES snes) 388*d71ae5a4SJacob Faibussowitsch { 38937e1895aSJed Brown PetscFunctionBegin; 39037e1895aSJed Brown PetscFunctionReturn(0); 39137e1895aSJed Brown } 39237e1895aSJed Brown 393*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_MS(SNES snes) 394*d71ae5a4SJacob Faibussowitsch { 39537e1895aSJed Brown PetscFunctionBegin; 3969566063dSJacob Faibussowitsch PetscCall(SNESReset_MS(snes)); 3979566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 3989566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetType_C", NULL)); 3999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetType_C", NULL)); 4009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetDamping_C", NULL)); 4019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetDamping_C", NULL)); 40237e1895aSJed Brown PetscFunctionReturn(0); 40337e1895aSJed Brown } 40437e1895aSJed Brown 405*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_MS(SNES snes, PetscViewer viewer) 406*d71ae5a4SJacob Faibussowitsch { 40737e1895aSJed Brown PetscBool iascii; 40837e1895aSJed Brown SNES_MS *ms = (SNES_MS *)snes->data; 40957715debSLisandro Dalcin SNESMSTableau tab = ms->tableau; 41037e1895aSJed Brown 41137e1895aSJed Brown PetscFunctionBegin; 4129566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 41348a46eb9SPierre Jolivet if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, " multi-stage method type: %s\n", tab->name)); 41437e1895aSJed Brown PetscFunctionReturn(0); 41537e1895aSJed Brown } 41637e1895aSJed Brown 417*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_MS(SNES snes, PetscOptionItems *PetscOptionsObject) 418*d71ae5a4SJacob Faibussowitsch { 419725b86d8SJed Brown SNES_MS *ms = (SNES_MS *)snes->data; 42037e1895aSJed Brown 42137e1895aSJed Brown PetscFunctionBegin; 422d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES MS options"); 42337e1895aSJed Brown { 42437e1895aSJed Brown SNESMSTableauLink link; 42537e1895aSJed Brown PetscInt count, choice; 42637e1895aSJed Brown PetscBool flg; 42737e1895aSJed Brown const char **namelist; 42857715debSLisandro Dalcin SNESMSType mstype; 42957715debSLisandro Dalcin PetscReal damping; 43037e1895aSJed Brown 4319566063dSJacob Faibussowitsch PetscCall(SNESMSGetType(snes, &mstype)); 4329371c9d4SSatish Balay for (link = SNESMSTableauList, count = 0; link; link = link->next, count++) 4339371c9d4SSatish Balay ; 4349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, (char ***)&namelist)); 43537e1895aSJed Brown for (link = SNESMSTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name; 4369566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-snes_ms_type", "Multistage smoother type", "SNESMSSetType", (const char *const *)namelist, count, mstype, &choice, &flg)); 4379566063dSJacob Faibussowitsch if (flg) PetscCall(SNESMSSetType(snes, namelist[choice])); 4389566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 4399566063dSJacob Faibussowitsch PetscCall(SNESMSGetDamping(snes, &damping)); 4409566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ms_damping", "Damping for multistage method", "SNESMSSetDamping", damping, &damping, &flg)); 4419566063dSJacob Faibussowitsch if (flg) PetscCall(SNESMSSetDamping(snes, damping)); 4429566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ms_norms", "Compute norms for monitoring", "none", ms->norms, &ms->norms, NULL)); 44337e1895aSJed Brown } 444d0609cedSBarry Smith PetscOptionsHeadEnd(); 44537e1895aSJed Brown PetscFunctionReturn(0); 44637e1895aSJed Brown } 44737e1895aSJed Brown 448*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSGetType_MS(SNES snes, SNESMSType *mstype) 449*d71ae5a4SJacob Faibussowitsch { 45057715debSLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 45157715debSLisandro Dalcin SNESMSTableau tab = ms->tableau; 45257715debSLisandro Dalcin 45357715debSLisandro Dalcin PetscFunctionBegin; 45457715debSLisandro Dalcin *mstype = tab->name; 45557715debSLisandro Dalcin PetscFunctionReturn(0); 45657715debSLisandro Dalcin } 45757715debSLisandro Dalcin 458*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSSetType_MS(SNES snes, SNESMSType mstype) 459*d71ae5a4SJacob Faibussowitsch { 46037e1895aSJed Brown SNES_MS *ms = (SNES_MS *)snes->data; 46137e1895aSJed Brown SNESMSTableauLink link; 46237e1895aSJed Brown PetscBool match; 46337e1895aSJed Brown 46437e1895aSJed Brown PetscFunctionBegin; 46537e1895aSJed Brown if (ms->tableau) { 4669566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ms->tableau->name, mstype, &match)); 46737e1895aSJed Brown if (match) PetscFunctionReturn(0); 46837e1895aSJed Brown } 46937e1895aSJed Brown for (link = SNESMSTableauList; link; link = link->next) { 4709566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name, mstype, &match)); 47137e1895aSJed Brown if (match) { 4729566063dSJacob Faibussowitsch if (snes->setupcalled) PetscCall(SNESReset_MS(snes)); 47337e1895aSJed Brown ms->tableau = &link->tab; 4749566063dSJacob Faibussowitsch if (snes->setupcalled) PetscCall(SNESSetUp_MS(snes)); 47537e1895aSJed Brown PetscFunctionReturn(0); 47637e1895aSJed Brown } 47737e1895aSJed Brown } 47898921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", mstype); 47937e1895aSJed Brown } 48037e1895aSJed Brown 48137e1895aSJed Brown /*@C 482f6dfbefdSBarry Smith SNESMSGetType - Get the type of multistage smoother `SNESMS` 48357715debSLisandro Dalcin 48457715debSLisandro Dalcin Not collective 48557715debSLisandro Dalcin 48657715debSLisandro Dalcin Input Parameter: 48757715debSLisandro Dalcin . snes - nonlinear solver context 48857715debSLisandro Dalcin 48957715debSLisandro Dalcin Output Parameter: 49057715debSLisandro Dalcin . mstype - type of multistage method 49157715debSLisandro Dalcin 492f6dfbefdSBarry Smith Level: advanced 49357715debSLisandro Dalcin 494f6dfbefdSBarry Smith .seealso: `SNESMS`, `SNESMSSetType()`, `SNESMSType`, `SNESMS` 49557715debSLisandro Dalcin @*/ 496*d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSGetType(SNES snes, SNESMSType *mstype) 497*d71ae5a4SJacob Faibussowitsch { 49857715debSLisandro Dalcin PetscFunctionBegin; 49957715debSLisandro Dalcin PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 50057715debSLisandro Dalcin PetscValidPointer(mstype, 2); 501cac4c232SBarry Smith PetscUseMethod(snes, "SNESMSGetType_C", (SNES, SNESMSType *), (snes, mstype)); 50257715debSLisandro Dalcin PetscFunctionReturn(0); 50357715debSLisandro Dalcin } 50457715debSLisandro Dalcin 50557715debSLisandro Dalcin /*@C 506f6dfbefdSBarry Smith SNESMSSetType - Set the type of multistage smoother `SNESMS` 50737e1895aSJed Brown 50837e1895aSJed Brown Logically collective 50937e1895aSJed Brown 510d8d19677SJose E. Roman Input Parameters: 51137e1895aSJed Brown + snes - nonlinear solver context 51237e1895aSJed Brown - mstype - type of multistage method 51337e1895aSJed Brown 514f6dfbefdSBarry Smith Level: advanced 51537e1895aSJed Brown 516f6dfbefdSBarry Smith .seealso: `SNESMS`, `SNESMSGetType()`, `SNESMSType`, `SNESMS` 51737e1895aSJed Brown @*/ 518*d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSSetType(SNES snes, SNESMSType mstype) 519*d71ae5a4SJacob Faibussowitsch { 52037e1895aSJed Brown PetscFunctionBegin; 52137e1895aSJed Brown PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 522dadcf809SJacob Faibussowitsch PetscValidCharPointer(mstype, 2); 523cac4c232SBarry Smith PetscTryMethod(snes, "SNESMSSetType_C", (SNES, SNESMSType), (snes, mstype)); 52457715debSLisandro Dalcin PetscFunctionReturn(0); 52557715debSLisandro Dalcin } 52657715debSLisandro Dalcin 527*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSGetDamping_MS(SNES snes, PetscReal *damping) 528*d71ae5a4SJacob Faibussowitsch { 52957715debSLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 53057715debSLisandro Dalcin 53157715debSLisandro Dalcin PetscFunctionBegin; 53257715debSLisandro Dalcin *damping = ms->damping; 53357715debSLisandro Dalcin PetscFunctionReturn(0); 53457715debSLisandro Dalcin } 53557715debSLisandro Dalcin 536*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSSetDamping_MS(SNES snes, PetscReal damping) 537*d71ae5a4SJacob Faibussowitsch { 53857715debSLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 53957715debSLisandro Dalcin 54057715debSLisandro Dalcin PetscFunctionBegin; 54157715debSLisandro Dalcin ms->damping = damping; 54257715debSLisandro Dalcin PetscFunctionReturn(0); 54357715debSLisandro Dalcin } 54457715debSLisandro Dalcin 54557715debSLisandro Dalcin /*@C 546f6dfbefdSBarry Smith SNESMSGetDamping - Get the damping parameter of `SNESMS` multistage scheme 54757715debSLisandro Dalcin 54857715debSLisandro Dalcin Not collective 54957715debSLisandro Dalcin 55057715debSLisandro Dalcin Input Parameter: 55157715debSLisandro Dalcin . snes - nonlinear solver context 55257715debSLisandro Dalcin 55357715debSLisandro Dalcin Output Parameter: 55457715debSLisandro Dalcin . damping - damping parameter 55557715debSLisandro Dalcin 55657715debSLisandro Dalcin Level: advanced 55757715debSLisandro Dalcin 558db781477SPatrick Sanan .seealso: `SNESMSSetDamping()`, `SNESMS` 55957715debSLisandro Dalcin @*/ 560*d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSGetDamping(SNES snes, PetscReal *damping) 561*d71ae5a4SJacob Faibussowitsch { 56257715debSLisandro Dalcin PetscFunctionBegin; 56357715debSLisandro Dalcin PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 564dadcf809SJacob Faibussowitsch PetscValidRealPointer(damping, 2); 565cac4c232SBarry Smith PetscUseMethod(snes, "SNESMSGetDamping_C", (SNES, PetscReal *), (snes, damping)); 56657715debSLisandro Dalcin PetscFunctionReturn(0); 56757715debSLisandro Dalcin } 56857715debSLisandro Dalcin 56957715debSLisandro Dalcin /*@C 570f6dfbefdSBarry Smith SNESMSSetDamping - Set the damping parameter for a `SNESMS` multistage scheme 57157715debSLisandro Dalcin 57257715debSLisandro Dalcin Logically collective 57357715debSLisandro Dalcin 574d8d19677SJose E. Roman Input Parameters: 57557715debSLisandro Dalcin + snes - nonlinear solver context 57657715debSLisandro Dalcin - damping - damping parameter 57757715debSLisandro Dalcin 57857715debSLisandro Dalcin Level: advanced 57957715debSLisandro Dalcin 580db781477SPatrick Sanan .seealso: `SNESMSGetDamping()`, `SNESMS` 58157715debSLisandro Dalcin @*/ 582*d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSSetDamping(SNES snes, PetscReal damping) 583*d71ae5a4SJacob Faibussowitsch { 58457715debSLisandro Dalcin PetscFunctionBegin; 58557715debSLisandro Dalcin PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 58657715debSLisandro Dalcin PetscValidLogicalCollectiveReal(snes, damping, 2); 587cac4c232SBarry Smith PetscTryMethod(snes, "SNESMSSetDamping_C", (SNES, PetscReal), (snes, damping)); 58837e1895aSJed Brown PetscFunctionReturn(0); 58937e1895aSJed Brown } 59037e1895aSJed Brown 59137e1895aSJed Brown /*MC 59237e1895aSJed Brown SNESMS - multi-stage smoothers 59337e1895aSJed Brown 594f6dfbefdSBarry Smith Options Database Keys: 59537e1895aSJed Brown + -snes_ms_type - type of multi-stage smoother 59637e1895aSJed Brown - -snes_ms_damping - damping for multi-stage method 59737e1895aSJed Brown 59837e1895aSJed Brown Notes: 59937e1895aSJed Brown These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for 6006c9de887SHong Zhang FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebyshev). 60137e1895aSJed Brown 60237e1895aSJed Brown Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds. 60337e1895aSJed Brown 604f6dfbefdSBarry Smith The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with `SNESMSRegister()`. 60537e1895aSJed Brown 60637e1895aSJed Brown References: 607606c0280SSatish Balay + * - Ketcheson (2010) Runge Kutta methods with minimum storage implementations (https://doi.org/10.1016/j.jcp.2009.11.006). 608606c0280SSatish 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). 609606c0280SSatish Balay . * - Pierce and Giles (1997) Preconditioned multigrid methods for compressible flow calculations on stretched meshes (https://doi.org/10.1006/jcph.1997.5772). 610606c0280SSatish 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). 61137e1895aSJed Brown 612f6dfbefdSBarry Smith Level: advanced 61337e1895aSJed Brown 614db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESMS`, `SNESFAS`, `KSPCHEBYSHEV` 61537e1895aSJed Brown 61637e1895aSJed Brown M*/ 617*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_MS(SNES snes) 618*d71ae5a4SJacob Faibussowitsch { 61937e1895aSJed Brown SNES_MS *ms; 62037e1895aSJed Brown 62137e1895aSJed Brown PetscFunctionBegin; 6229566063dSJacob Faibussowitsch PetscCall(SNESMSInitializePackage()); 62337e1895aSJed Brown 62437e1895aSJed Brown snes->ops->setup = SNESSetUp_MS; 62537e1895aSJed Brown snes->ops->solve = SNESSolve_MS; 62637e1895aSJed Brown snes->ops->destroy = SNESDestroy_MS; 62737e1895aSJed Brown snes->ops->setfromoptions = SNESSetFromOptions_MS; 62837e1895aSJed Brown snes->ops->view = SNESView_MS; 62937e1895aSJed Brown snes->ops->reset = SNESReset_MS; 63037e1895aSJed Brown 631efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 63237e1895aSJed Brown snes->usesksp = PETSC_TRUE; 63337e1895aSJed Brown 6344fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 6354fc747eaSLawrence Mitchell 6364dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ms)); 63737e1895aSJed Brown snes->data = (void *)ms; 63837e1895aSJed Brown ms->damping = 0.9; 63937e1895aSJed Brown ms->norms = PETSC_FALSE; 64037e1895aSJed Brown 6419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetType_C", SNESMSGetType_MS)); 6429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetType_C", SNESMSSetType_MS)); 6439566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetDamping_C", SNESMSGetDamping_MS)); 6449566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetDamping_C", SNESMSSetDamping_MS)); 64557715debSLisandro Dalcin 6469566063dSJacob Faibussowitsch PetscCall(SNESMSSetType(snes, SNESMSDefault)); 64737e1895aSJed Brown PetscFunctionReturn(0); 64837e1895aSJed Brown } 649