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 } SNES_MS; 2937e1895aSJed Brown 3037e1895aSJed Brown /*@C 31f6dfbefdSBarry Smith SNESMSRegisterAll - Registers all of the multi-stage methods in `SNESMS` 3237e1895aSJed Brown 33f6dfbefdSBarry Smith Logically Collective 3437e1895aSJed Brown 35f6dfbefdSBarry Smith Level: developer 3637e1895aSJed Brown 37f6dfbefdSBarry Smith .seealso: `SNES`, `SNESMS`, `SNESMSRegisterDestroy()` 3837e1895aSJed Brown @*/ 39d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSRegisterAll(void) 40d71ae5a4SJacob Faibussowitsch { 4137e1895aSJed Brown PetscFunctionBegin; 423ba16761SJacob Faibussowitsch if (SNESMSRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 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 } 903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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 @*/ 102d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSRegisterDestroy(void) 103d71ae5a4SJacob Faibussowitsch { 10437e1895aSJed Brown SNESMSTableauLink link; 10537e1895aSJed Brown 10637e1895aSJed Brown PetscFunctionBegin; 10737e1895aSJed Brown while ((link = SNESMSTableauList)) { 10837e1895aSJed Brown SNESMSTableau t = &link->tab; 10937e1895aSJed Brown SNESMSTableauList = link->next; 1101aa26658SKarl Rupp 1119566063dSJacob Faibussowitsch PetscCall(PetscFree(t->name)); 1129566063dSJacob Faibussowitsch PetscCall(PetscFree(t->gamma)); 1139566063dSJacob Faibussowitsch PetscCall(PetscFree(t->delta)); 1149566063dSJacob Faibussowitsch PetscCall(PetscFree(t->betasub)); 1159566063dSJacob Faibussowitsch PetscCall(PetscFree(link)); 11637e1895aSJed Brown } 11737e1895aSJed Brown SNESMSRegisterAllCalled = PETSC_FALSE; 1183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11937e1895aSJed Brown } 12037e1895aSJed Brown 12137e1895aSJed Brown /*@C 122f6dfbefdSBarry Smith SNESMSInitializePackage - This function initializes everything in the `SNESMS` package. It is called 123f6dfbefdSBarry Smith from `SNESInitializePackage()`. 12437e1895aSJed Brown 12537e1895aSJed Brown Level: developer 12637e1895aSJed Brown 127db781477SPatrick Sanan .seealso: `PetscInitialize()` 12837e1895aSJed Brown @*/ 129d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSInitializePackage(void) 130d71ae5a4SJacob Faibussowitsch { 13137e1895aSJed Brown PetscFunctionBegin; 1323ba16761SJacob Faibussowitsch if (SNESMSPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 13337e1895aSJed Brown SNESMSPackageInitialized = PETSC_TRUE; 1341aa26658SKarl Rupp 1359566063dSJacob Faibussowitsch PetscCall(SNESMSRegisterAll()); 1369566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(SNESMSFinalizePackage)); 1373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13837e1895aSJed Brown } 13937e1895aSJed Brown 14037e1895aSJed Brown /*@C 141f6dfbefdSBarry Smith SNESMSFinalizePackage - This function destroys everything in the `SNESMS` package. It is 142f6dfbefdSBarry Smith called from `PetscFinalize()`. 14337e1895aSJed Brown 14437e1895aSJed Brown Level: developer 14537e1895aSJed Brown 146db781477SPatrick Sanan .seealso: `PetscFinalize()` 14737e1895aSJed Brown @*/ 148d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSFinalizePackage(void) 149d71ae5a4SJacob Faibussowitsch { 15037e1895aSJed Brown PetscFunctionBegin; 15137e1895aSJed Brown SNESMSPackageInitialized = PETSC_FALSE; 1521aa26658SKarl Rupp 1539566063dSJacob Faibussowitsch PetscCall(SNESMSRegisterDestroy()); 1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15537e1895aSJed Brown } 15637e1895aSJed Brown 15737e1895aSJed Brown /*@C 158f6dfbefdSBarry Smith SNESMSRegister - register a multistage scheme for `SNESMS` 15937e1895aSJed Brown 160f6dfbefdSBarry Smith Logically Collective 16137e1895aSJed Brown 16237e1895aSJed Brown Input Parameters: 16337e1895aSJed Brown + name - identifier for method 16437e1895aSJed Brown . nstages - number of stages 16537e1895aSJed Brown . nregisters - number of registers used by low-storage implementation 1669ce202e0SLisandro Dalcin . stability - scaled stability region 16737e1895aSJed Brown . gamma - coefficients, see Ketcheson's paper 16837e1895aSJed Brown . delta - coefficients, see Ketcheson's paper 16937e1895aSJed Brown - betasub - subdiagonal of Shu-Osher form 17037e1895aSJed Brown 17120f4b53cSBarry Smith Level: advanced 17220f4b53cSBarry Smith 17337e1895aSJed Brown Notes: 17437e1895aSJed Brown The notation is described in Ketcheson (2010) Runge-Kutta methods with minimum storage implementations. 17537e1895aSJed Brown 1769ce202e0SLisandro Dalcin Many multistage schemes are of the form 1779ce202e0SLisandro Dalcin $ X_0 = X^{(n)} 1789ce202e0SLisandro Dalcin $ X_k = X_0 + \alpha_k * F(X_{k-1}), k = 1,\ldots,s 1799ce202e0SLisandro Dalcin $ X^{(n+1)} = X_s 1809ce202e0SLisandro Dalcin These methods can be registered with 1819ce202e0SLisandro Dalcin .vb 1829ce202e0SLisandro Dalcin SNESMSRegister("name",s,1,stability,NULL,NULL,alpha); 1839ce202e0SLisandro Dalcin .ve 1849ce202e0SLisandro Dalcin 185db781477SPatrick Sanan .seealso: `SNESMS` 18637e1895aSJed Brown @*/ 187d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSRegister(SNESMSType name, PetscInt nstages, PetscInt nregisters, PetscReal stability, const PetscReal gamma[], const PetscReal delta[], const PetscReal betasub[]) 188d71ae5a4SJacob Faibussowitsch { 18937e1895aSJed Brown SNESMSTableauLink link; 19037e1895aSJed Brown SNESMSTableau t; 19137e1895aSJed Brown 19237e1895aSJed Brown PetscFunctionBegin; 19337e1895aSJed Brown PetscValidCharPointer(name, 1); 19408401ef6SPierre Jolivet PetscCheck(nstages >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have at least one stage"); 1959ce202e0SLisandro Dalcin if (gamma || delta) { 19608401ef6SPierre Jolivet PetscCheck(nregisters == 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support for methods written in 3-register form"); 197064a246eSJacob Faibussowitsch PetscValidRealPointer(gamma, 5); 198064a246eSJacob Faibussowitsch PetscValidRealPointer(delta, 6); 1999ce202e0SLisandro Dalcin } else { 20008401ef6SPierre Jolivet PetscCheck(nregisters == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support for methods written in 1-register form"); 2019ce202e0SLisandro Dalcin } 202064a246eSJacob Faibussowitsch PetscValidRealPointer(betasub, 7); 20337e1895aSJed Brown 2049566063dSJacob Faibussowitsch PetscCall(SNESMSInitializePackage()); 2059566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 20637e1895aSJed Brown t = &link->tab; 2079566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &t->name)); 20837e1895aSJed Brown t->nstages = nstages; 20937e1895aSJed Brown t->nregisters = nregisters; 21037e1895aSJed Brown t->stability = stability; 2111aa26658SKarl Rupp 2129ce202e0SLisandro Dalcin if (gamma && delta) { 2139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nstages * nregisters, &t->gamma)); 2149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nstages, &t->delta)); 2159566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->gamma, gamma, nstages * nregisters)); 2169566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->delta, delta, nstages)); 2179ce202e0SLisandro Dalcin } 2189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nstages, &t->betasub)); 2199566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->betasub, betasub, nstages)); 22037e1895aSJed Brown 22137e1895aSJed Brown link->next = SNESMSTableauList; 22237e1895aSJed Brown SNESMSTableauList = link; 2233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22437e1895aSJed Brown } 22537e1895aSJed Brown 22637e1895aSJed Brown /* 22737e1895aSJed Brown X - initial state, updated in-place. 22837e1895aSJed Brown F - residual, computed at the initial X on input 22937e1895aSJed Brown */ 230d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSStep_3Sstar(SNES snes, Vec X, Vec F) 231d71ae5a4SJacob Faibussowitsch { 23237e1895aSJed Brown SNES_MS *ms = (SNES_MS *)snes->data; 23337e1895aSJed Brown SNESMSTableau t = ms->tableau; 2343e2ddb07SJacob Faibussowitsch const PetscInt nstages = t->nstages; 23537e1895aSJed Brown const PetscReal *gamma = t->gamma, *delta = t->delta, *betasub = t->betasub; 2363e2ddb07SJacob Faibussowitsch Vec S1 = X, S2 = snes->work[1], S3 = snes->work[2], Y = snes->work[0], S1copy = snes->work[3]; 23737e1895aSJed Brown 23837e1895aSJed Brown PetscFunctionBegin; 2399566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(S2)); 2409566063dSJacob Faibussowitsch PetscCall(VecCopy(X, S3)); 2413e2ddb07SJacob Faibussowitsch for (PetscInt i = 0; i < nstages; i++) { 2423e2ddb07SJacob Faibussowitsch Vec Ss[] = {S1copy, S2, S3, Y}; 2433e2ddb07SJacob Faibussowitsch const PetscScalar scoeff[] = {gamma[0 * nstages + i] - 1, gamma[1 * nstages + i], gamma[2 * nstages + i], -betasub[i] * ms->damping}; 2441aa26658SKarl Rupp 2459566063dSJacob Faibussowitsch PetscCall(VecAXPY(S2, delta[i], S1)); 24648a46eb9SPierre Jolivet if (i > 0) PetscCall(SNESComputeFunction(snes, S1, F)); 2479566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, Y)); 2483e2ddb07SJacob Faibussowitsch PetscCall(VecCopy(S1, S1copy)); 2499566063dSJacob Faibussowitsch PetscCall(VecMAXPY(S1, 4, scoeff, Ss)); 25037e1895aSJed Brown } 2513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 25237e1895aSJed Brown } 25337e1895aSJed Brown 2549ce202e0SLisandro Dalcin /* 2559ce202e0SLisandro Dalcin X - initial state, updated in-place. 2569ce202e0SLisandro Dalcin F - residual, computed at the initial X on input 2579ce202e0SLisandro Dalcin */ 258d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSStep_Basic(SNES snes, Vec X, Vec F) 259d71ae5a4SJacob Faibussowitsch { 2609ce202e0SLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 2619ce202e0SLisandro Dalcin SNESMSTableau tab = ms->tableau; 2629ce202e0SLisandro Dalcin const PetscReal *alpha = tab->betasub, h = ms->damping; 2639ce202e0SLisandro Dalcin PetscInt i, nstages = tab->nstages; 2649ce202e0SLisandro Dalcin Vec X0 = snes->work[0]; 2659ce202e0SLisandro Dalcin 2669ce202e0SLisandro Dalcin PetscFunctionBegin; 2679566063dSJacob Faibussowitsch PetscCall(VecCopy(X, X0)); 2689ce202e0SLisandro Dalcin for (i = 0; i < nstages; i++) { 26948a46eb9SPierre Jolivet if (i > 0) PetscCall(SNESComputeFunction(snes, X, F)); 2709566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, X)); 2719566063dSJacob Faibussowitsch PetscCall(VecAYPX(X, -alpha[i] * h, X0)); 2729ce202e0SLisandro Dalcin } 2733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2749ce202e0SLisandro Dalcin } 2759ce202e0SLisandro Dalcin 276d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSStep_Step(SNES snes, Vec X, Vec F) 277d71ae5a4SJacob Faibussowitsch { 2789ce202e0SLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 2799ce202e0SLisandro Dalcin SNESMSTableau tab = ms->tableau; 2809ce202e0SLisandro Dalcin 2819ce202e0SLisandro Dalcin PetscFunctionBegin; 2829ce202e0SLisandro Dalcin if (tab->gamma && tab->delta) { 2839566063dSJacob Faibussowitsch PetscCall(SNESMSStep_3Sstar(snes, X, F)); 2849ce202e0SLisandro Dalcin } else { 2859566063dSJacob Faibussowitsch PetscCall(SNESMSStep_Basic(snes, X, F)); 2869ce202e0SLisandro Dalcin } 2873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2889ce202e0SLisandro Dalcin } 2899ce202e0SLisandro Dalcin 290d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSStep_Norms(SNES snes, PetscInt iter, Vec F) 291d71ae5a4SJacob Faibussowitsch { 292bf80fd75SLisandro Dalcin PetscReal fnorm; 293bf80fd75SLisandro Dalcin 294bf80fd75SLisandro Dalcin PetscFunctionBegin; 295*7a14dab9SStefano Zampini if (SNESNeedNorm_Private(snes, iter)) { 2969566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 297bf80fd75SLisandro Dalcin SNESCheckFunctionNorm(snes, fnorm); 298bf80fd75SLisandro Dalcin /* Monitor convergence */ 2999566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 300bf80fd75SLisandro Dalcin snes->iter = iter; 301bf80fd75SLisandro Dalcin snes->norm = fnorm; 3029566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3039566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 304bf80fd75SLisandro Dalcin /* Test for convergence */ 3052d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, fnorm)); 3062d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 307bf80fd75SLisandro Dalcin } else if (iter > 0) { 3089566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 309bf80fd75SLisandro Dalcin snes->iter = iter; 3109566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 311bf80fd75SLisandro Dalcin } 3123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 313bf80fd75SLisandro Dalcin } 314bf80fd75SLisandro Dalcin 315d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_MS(SNES snes) 316d71ae5a4SJacob Faibussowitsch { 31737e1895aSJed Brown Vec X = snes->vec_sol, F = snes->vec_func; 31837e1895aSJed Brown PetscInt i; 31937e1895aSJed Brown 32037e1895aSJed Brown PetscFunctionBegin; 3210b121fc5SBarry 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); 3229566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 323bf80fd75SLisandro Dalcin 32437e1895aSJed Brown snes->reason = SNES_CONVERGED_ITERATING; 3259566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 32637e1895aSJed Brown snes->iter = 0; 327bf80fd75SLisandro Dalcin snes->norm = 0; 3289566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 329bf80fd75SLisandro Dalcin 330e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 3319566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 3321aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 3331aa26658SKarl Rupp 3349566063dSJacob Faibussowitsch PetscCall(SNESMSStep_Norms(snes, 0, F)); 3353ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 336bf80fd75SLisandro Dalcin 337bf80fd75SLisandro Dalcin for (i = 0; i < snes->max_its; i++) { 33837e1895aSJed Brown /* Call general purpose update function */ 339dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 340bf80fd75SLisandro Dalcin 341bf80fd75SLisandro Dalcin if (i == 0 && snes->jacobian) { 342bf80fd75SLisandro Dalcin /* This method does not require a Jacobian, but it is usually preconditioned by PBJacobi */ 3439566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, snes->vec_sol, snes->jacobian, snes->jacobian_pre)); 344bf80fd75SLisandro Dalcin SNESCheckJacobianDomainerror(snes); 3459566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 346bf80fd75SLisandro Dalcin } 347bf80fd75SLisandro Dalcin 3489566063dSJacob Faibussowitsch PetscCall(SNESMSStep_Step(snes, X, F)); 34937e1895aSJed Brown 350*7a14dab9SStefano Zampini if (i < snes->max_its - 1 || SNESNeedNorm_Private(snes, i + 1)) PetscCall(SNESComputeFunction(snes, X, F)); 35137e1895aSJed Brown 3529566063dSJacob Faibussowitsch PetscCall(SNESMSStep_Norms(snes, i + 1, F)); 3533ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 35437e1895aSJed Brown } 3553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35637e1895aSJed Brown } 35737e1895aSJed Brown 358d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_MS(SNES snes) 359d71ae5a4SJacob Faibussowitsch { 3609ce202e0SLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 3619ce202e0SLisandro Dalcin SNESMSTableau tab = ms->tableau; 3623e2ddb07SJacob Faibussowitsch PetscInt nwork = tab->nregisters + 1; // +1 because VecMAXPY() in SNESMSStep_3Sstar() 3633e2ddb07SJacob Faibussowitsch // needs an extra work vec 36437e1895aSJed Brown 36537e1895aSJed Brown PetscFunctionBegin; 3669566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, nwork)); 3679566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 3683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 36937e1895aSJed Brown } 37037e1895aSJed Brown 371d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESReset_MS(SNES snes) 372d71ae5a4SJacob Faibussowitsch { 37337e1895aSJed Brown PetscFunctionBegin; 3743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 37537e1895aSJed Brown } 37637e1895aSJed Brown 377d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_MS(SNES snes) 378d71ae5a4SJacob Faibussowitsch { 37937e1895aSJed Brown PetscFunctionBegin; 3809566063dSJacob Faibussowitsch PetscCall(SNESReset_MS(snes)); 3819566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 3829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetType_C", NULL)); 3839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetType_C", NULL)); 3849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetDamping_C", NULL)); 3859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetDamping_C", NULL)); 3863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38737e1895aSJed Brown } 38837e1895aSJed Brown 389d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_MS(SNES snes, PetscViewer viewer) 390d71ae5a4SJacob Faibussowitsch { 39137e1895aSJed Brown PetscBool iascii; 39237e1895aSJed Brown SNES_MS *ms = (SNES_MS *)snes->data; 39357715debSLisandro Dalcin SNESMSTableau tab = ms->tableau; 39437e1895aSJed Brown 39537e1895aSJed Brown PetscFunctionBegin; 3969566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 39748a46eb9SPierre Jolivet if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, " multi-stage method type: %s\n", tab->name)); 3983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 39937e1895aSJed Brown } 40037e1895aSJed Brown 401d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_MS(SNES snes, PetscOptionItems *PetscOptionsObject) 402d71ae5a4SJacob Faibussowitsch { 40337e1895aSJed Brown PetscFunctionBegin; 404d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES MS options"); 40537e1895aSJed Brown { 40637e1895aSJed Brown SNESMSTableauLink link; 40737e1895aSJed Brown PetscInt count, choice; 40837e1895aSJed Brown PetscBool flg; 40937e1895aSJed Brown const char **namelist; 41057715debSLisandro Dalcin SNESMSType mstype; 41157715debSLisandro Dalcin PetscReal damping; 412*7a14dab9SStefano Zampini PetscBool norms = PETSC_FALSE; 41337e1895aSJed Brown 4149566063dSJacob Faibussowitsch PetscCall(SNESMSGetType(snes, &mstype)); 4159371c9d4SSatish Balay for (link = SNESMSTableauList, count = 0; link; link = link->next, count++) 4169371c9d4SSatish Balay ; 4179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, (char ***)&namelist)); 41837e1895aSJed Brown for (link = SNESMSTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name; 4199566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-snes_ms_type", "Multistage smoother type", "SNESMSSetType", (const char *const *)namelist, count, mstype, &choice, &flg)); 4209566063dSJacob Faibussowitsch if (flg) PetscCall(SNESMSSetType(snes, namelist[choice])); 4219566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 4229566063dSJacob Faibussowitsch PetscCall(SNESMSGetDamping(snes, &damping)); 4239566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ms_damping", "Damping for multistage method", "SNESMSSetDamping", damping, &damping, &flg)); 4249566063dSJacob Faibussowitsch if (flg) PetscCall(SNESMSSetDamping(snes, damping)); 425*7a14dab9SStefano Zampini 426*7a14dab9SStefano Zampini /* deprecated option */ 427*7a14dab9SStefano Zampini PetscCall(PetscOptionsDeprecated("-snes_ms_norms", NULL, "3.20", "Use -snes_norm_schedule always")); 428*7a14dab9SStefano Zampini PetscCall(PetscOptionsBool("-snes_ms_norms", NULL, NULL, norms, &norms, NULL)); 429*7a14dab9SStefano Zampini if (norms) PetscCall(SNESSetNormSchedule(snes, SNES_NORM_ALWAYS)); 43037e1895aSJed Brown } 431d0609cedSBarry Smith PetscOptionsHeadEnd(); 4323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43337e1895aSJed Brown } 43437e1895aSJed Brown 435d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSGetType_MS(SNES snes, SNESMSType *mstype) 436d71ae5a4SJacob Faibussowitsch { 43757715debSLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 43857715debSLisandro Dalcin SNESMSTableau tab = ms->tableau; 43957715debSLisandro Dalcin 44057715debSLisandro Dalcin PetscFunctionBegin; 44157715debSLisandro Dalcin *mstype = tab->name; 4423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44357715debSLisandro Dalcin } 44457715debSLisandro Dalcin 445d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSSetType_MS(SNES snes, SNESMSType mstype) 446d71ae5a4SJacob Faibussowitsch { 44737e1895aSJed Brown SNES_MS *ms = (SNES_MS *)snes->data; 44837e1895aSJed Brown SNESMSTableauLink link; 44937e1895aSJed Brown PetscBool match; 45037e1895aSJed Brown 45137e1895aSJed Brown PetscFunctionBegin; 45237e1895aSJed Brown if (ms->tableau) { 4539566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ms->tableau->name, mstype, &match)); 4543ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 45537e1895aSJed Brown } 45637e1895aSJed Brown for (link = SNESMSTableauList; link; link = link->next) { 4579566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name, mstype, &match)); 45837e1895aSJed Brown if (match) { 4599566063dSJacob Faibussowitsch if (snes->setupcalled) PetscCall(SNESReset_MS(snes)); 46037e1895aSJed Brown ms->tableau = &link->tab; 4619566063dSJacob Faibussowitsch if (snes->setupcalled) PetscCall(SNESSetUp_MS(snes)); 4623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46337e1895aSJed Brown } 46437e1895aSJed Brown } 46598921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", mstype); 46637e1895aSJed Brown } 46737e1895aSJed Brown 46837e1895aSJed Brown /*@C 469f6dfbefdSBarry Smith SNESMSGetType - Get the type of multistage smoother `SNESMS` 47057715debSLisandro Dalcin 47120f4b53cSBarry Smith Not Collective 47257715debSLisandro Dalcin 47357715debSLisandro Dalcin Input Parameter: 47457715debSLisandro Dalcin . snes - nonlinear solver context 47557715debSLisandro Dalcin 47657715debSLisandro Dalcin Output Parameter: 47757715debSLisandro Dalcin . mstype - type of multistage method 47857715debSLisandro Dalcin 479f6dfbefdSBarry Smith Level: advanced 48057715debSLisandro Dalcin 481f6dfbefdSBarry Smith .seealso: `SNESMS`, `SNESMSSetType()`, `SNESMSType`, `SNESMS` 48257715debSLisandro Dalcin @*/ 483d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSGetType(SNES snes, SNESMSType *mstype) 484d71ae5a4SJacob Faibussowitsch { 48557715debSLisandro Dalcin PetscFunctionBegin; 48657715debSLisandro Dalcin PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 48757715debSLisandro Dalcin PetscValidPointer(mstype, 2); 488cac4c232SBarry Smith PetscUseMethod(snes, "SNESMSGetType_C", (SNES, SNESMSType *), (snes, mstype)); 4893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49057715debSLisandro Dalcin } 49157715debSLisandro Dalcin 49257715debSLisandro Dalcin /*@C 493f6dfbefdSBarry Smith SNESMSSetType - Set the type of multistage smoother `SNESMS` 49437e1895aSJed Brown 49520f4b53cSBarry Smith Logically Collective 49637e1895aSJed Brown 497d8d19677SJose E. Roman Input Parameters: 49837e1895aSJed Brown + snes - nonlinear solver context 49937e1895aSJed Brown - mstype - type of multistage method 50037e1895aSJed Brown 501f6dfbefdSBarry Smith Level: advanced 50237e1895aSJed Brown 503f6dfbefdSBarry Smith .seealso: `SNESMS`, `SNESMSGetType()`, `SNESMSType`, `SNESMS` 50437e1895aSJed Brown @*/ 505d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSSetType(SNES snes, SNESMSType mstype) 506d71ae5a4SJacob Faibussowitsch { 50737e1895aSJed Brown PetscFunctionBegin; 50837e1895aSJed Brown PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 509dadcf809SJacob Faibussowitsch PetscValidCharPointer(mstype, 2); 510cac4c232SBarry Smith PetscTryMethod(snes, "SNESMSSetType_C", (SNES, SNESMSType), (snes, mstype)); 5113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51257715debSLisandro Dalcin } 51357715debSLisandro Dalcin 514d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSGetDamping_MS(SNES snes, PetscReal *damping) 515d71ae5a4SJacob Faibussowitsch { 51657715debSLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 51757715debSLisandro Dalcin 51857715debSLisandro Dalcin PetscFunctionBegin; 51957715debSLisandro Dalcin *damping = ms->damping; 5203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 52157715debSLisandro Dalcin } 52257715debSLisandro Dalcin 523d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSSetDamping_MS(SNES snes, PetscReal damping) 524d71ae5a4SJacob Faibussowitsch { 52557715debSLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 52657715debSLisandro Dalcin 52757715debSLisandro Dalcin PetscFunctionBegin; 52857715debSLisandro Dalcin ms->damping = damping; 5293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53057715debSLisandro Dalcin } 53157715debSLisandro Dalcin 53257715debSLisandro Dalcin /*@C 533f6dfbefdSBarry Smith SNESMSGetDamping - Get the damping parameter of `SNESMS` multistage scheme 53457715debSLisandro Dalcin 53520f4b53cSBarry Smith Not Collective 53657715debSLisandro Dalcin 53757715debSLisandro Dalcin Input Parameter: 53857715debSLisandro Dalcin . snes - nonlinear solver context 53957715debSLisandro Dalcin 54057715debSLisandro Dalcin Output Parameter: 54157715debSLisandro Dalcin . damping - damping parameter 54257715debSLisandro Dalcin 54357715debSLisandro Dalcin Level: advanced 54457715debSLisandro Dalcin 545db781477SPatrick Sanan .seealso: `SNESMSSetDamping()`, `SNESMS` 54657715debSLisandro Dalcin @*/ 547d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSGetDamping(SNES snes, PetscReal *damping) 548d71ae5a4SJacob Faibussowitsch { 54957715debSLisandro Dalcin PetscFunctionBegin; 55057715debSLisandro Dalcin PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 551dadcf809SJacob Faibussowitsch PetscValidRealPointer(damping, 2); 552cac4c232SBarry Smith PetscUseMethod(snes, "SNESMSGetDamping_C", (SNES, PetscReal *), (snes, damping)); 5533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 55457715debSLisandro Dalcin } 55557715debSLisandro Dalcin 55657715debSLisandro Dalcin /*@C 557f6dfbefdSBarry Smith SNESMSSetDamping - Set the damping parameter for a `SNESMS` multistage scheme 55857715debSLisandro Dalcin 55920f4b53cSBarry Smith Logically Collective 56057715debSLisandro Dalcin 561d8d19677SJose E. Roman Input Parameters: 56257715debSLisandro Dalcin + snes - nonlinear solver context 56357715debSLisandro Dalcin - damping - damping parameter 56457715debSLisandro Dalcin 56557715debSLisandro Dalcin Level: advanced 56657715debSLisandro Dalcin 567db781477SPatrick Sanan .seealso: `SNESMSGetDamping()`, `SNESMS` 56857715debSLisandro Dalcin @*/ 569d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSSetDamping(SNES snes, PetscReal damping) 570d71ae5a4SJacob Faibussowitsch { 57157715debSLisandro Dalcin PetscFunctionBegin; 57257715debSLisandro Dalcin PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 57357715debSLisandro Dalcin PetscValidLogicalCollectiveReal(snes, damping, 2); 574cac4c232SBarry Smith PetscTryMethod(snes, "SNESMSSetDamping_C", (SNES, PetscReal), (snes, damping)); 5753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 57637e1895aSJed Brown } 57737e1895aSJed Brown 57837e1895aSJed Brown /*MC 57937e1895aSJed Brown SNESMS - multi-stage smoothers 58037e1895aSJed Brown 581f6dfbefdSBarry Smith Options Database Keys: 58237e1895aSJed Brown + -snes_ms_type - type of multi-stage smoother 58337e1895aSJed Brown - -snes_ms_damping - damping for multi-stage method 58437e1895aSJed Brown 58537e1895aSJed Brown Notes: 58637e1895aSJed Brown These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for 5876c9de887SHong Zhang FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebyshev). 58837e1895aSJed Brown 58937e1895aSJed Brown Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds. 59037e1895aSJed Brown 591f6dfbefdSBarry Smith The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with `SNESMSRegister()`. 59237e1895aSJed Brown 59337e1895aSJed Brown References: 594606c0280SSatish Balay + * - Ketcheson (2010) Runge Kutta methods with minimum storage implementations (https://doi.org/10.1016/j.jcp.2009.11.006). 595606c0280SSatish 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). 596606c0280SSatish Balay . * - Pierce and Giles (1997) Preconditioned multigrid methods for compressible flow calculations on stretched meshes (https://doi.org/10.1006/jcph.1997.5772). 597606c0280SSatish 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). 59837e1895aSJed Brown 599f6dfbefdSBarry Smith Level: advanced 60037e1895aSJed Brown 601db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESMS`, `SNESFAS`, `KSPCHEBYSHEV` 60237e1895aSJed Brown 60337e1895aSJed Brown M*/ 604d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_MS(SNES snes) 605d71ae5a4SJacob Faibussowitsch { 60637e1895aSJed Brown SNES_MS *ms; 60737e1895aSJed Brown 60837e1895aSJed Brown PetscFunctionBegin; 6099566063dSJacob Faibussowitsch PetscCall(SNESMSInitializePackage()); 61037e1895aSJed Brown 61137e1895aSJed Brown snes->ops->setup = SNESSetUp_MS; 61237e1895aSJed Brown snes->ops->solve = SNESSolve_MS; 61337e1895aSJed Brown snes->ops->destroy = SNESDestroy_MS; 61437e1895aSJed Brown snes->ops->setfromoptions = SNESSetFromOptions_MS; 61537e1895aSJed Brown snes->ops->view = SNESView_MS; 61637e1895aSJed Brown snes->ops->reset = SNESReset_MS; 61737e1895aSJed Brown 618efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 61937e1895aSJed Brown snes->usesksp = PETSC_TRUE; 62037e1895aSJed Brown 6214fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 6224fc747eaSLawrence Mitchell 6234dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ms)); 62437e1895aSJed Brown snes->data = (void *)ms; 62537e1895aSJed Brown ms->damping = 0.9; 62637e1895aSJed Brown 6279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetType_C", SNESMSGetType_MS)); 6289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetType_C", SNESMSSetType_MS)); 6299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetDamping_C", SNESMSGetDamping_MS)); 6309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetDamping_C", SNESMSSetDamping_MS)); 63157715debSLisandro Dalcin 6329566063dSJacob Faibussowitsch PetscCall(SNESMSSetType(snes, SNESMSDefault)); 6333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 63437e1895aSJed Brown } 635