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 37420bcc1bSBarry Smith .seealso: [](ch_snes), `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 96420bcc1bSBarry Smith Logically Collective 9737e1895aSJed Brown 98f6dfbefdSBarry Smith Level: developer 9937e1895aSJed Brown 100420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESMS`, `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 127420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESMS`, `SNESMSRegister()`, `SNESMSRegisterAll()`, `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 146420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESMS`, `SNESMSRegister()`, `SNESMSRegisterAll()`, `SNESMSInitializePackage()`, `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 160cc4c1da9SBarry Smith Logically Collective, No Fortran Support 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 167420bcc1bSBarry Smith . gamma - coefficients, see Ketcheson's paper {cite}`ketcheson2010runge` 168420bcc1bSBarry Smith . delta - coefficients, see Ketcheson's paper {cite}`ketcheson2010runge` 16937e1895aSJed Brown - betasub - subdiagonal of Shu-Osher form 17037e1895aSJed Brown 17120f4b53cSBarry Smith Level: advanced 17220f4b53cSBarry Smith 17337e1895aSJed Brown Notes: 174420bcc1bSBarry Smith The notation is described in {cite}`ketcheson2010runge` Ketcheson (2010) Runge-Kutta methods with minimum storage implementations. 17537e1895aSJed Brown 1769ce202e0SLisandro Dalcin Many multistage schemes are of the form 177420bcc1bSBarry Smith 178420bcc1bSBarry Smith $$ 179420bcc1bSBarry Smith \begin{align*} 180420bcc1bSBarry Smith X_0 = X^{(n)} \\ 181420bcc1bSBarry Smith X_k = X_0 + \alpha_k * F(X_{k-1}), k = 1,\ldots,s \\ 182420bcc1bSBarry Smith X^{(n+1)} = X_s 183420bcc1bSBarry Smith \end{align*} 184420bcc1bSBarry Smith $$ 185420bcc1bSBarry Smith 1869ce202e0SLisandro Dalcin These methods can be registered with 1879ce202e0SLisandro Dalcin .vb 1889ce202e0SLisandro Dalcin SNESMSRegister("name",s,1,stability,NULL,NULL,alpha); 1899ce202e0SLisandro Dalcin .ve 1909ce202e0SLisandro Dalcin 191420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESMS` 19237e1895aSJed Brown @*/ 193d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSRegister(SNESMSType name, PetscInt nstages, PetscInt nregisters, PetscReal stability, const PetscReal gamma[], const PetscReal delta[], const PetscReal betasub[]) 194d71ae5a4SJacob Faibussowitsch { 19537e1895aSJed Brown SNESMSTableauLink link; 19637e1895aSJed Brown SNESMSTableau t; 19737e1895aSJed Brown 19837e1895aSJed Brown PetscFunctionBegin; 1994f572ea9SToby Isaac PetscAssertPointer(name, 1); 20008401ef6SPierre Jolivet PetscCheck(nstages >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have at least one stage"); 2019ce202e0SLisandro Dalcin if (gamma || delta) { 20208401ef6SPierre Jolivet PetscCheck(nregisters == 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support for methods written in 3-register form"); 2034f572ea9SToby Isaac PetscAssertPointer(gamma, 5); 2044f572ea9SToby Isaac PetscAssertPointer(delta, 6); 2059ce202e0SLisandro Dalcin } else { 20608401ef6SPierre Jolivet PetscCheck(nregisters == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support for methods written in 1-register form"); 2079ce202e0SLisandro Dalcin } 2084f572ea9SToby Isaac PetscAssertPointer(betasub, 7); 20937e1895aSJed Brown 2109566063dSJacob Faibussowitsch PetscCall(SNESMSInitializePackage()); 2119566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 21237e1895aSJed Brown t = &link->tab; 2139566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &t->name)); 21437e1895aSJed Brown t->nstages = nstages; 21537e1895aSJed Brown t->nregisters = nregisters; 21637e1895aSJed Brown t->stability = stability; 2171aa26658SKarl Rupp 2189ce202e0SLisandro Dalcin if (gamma && delta) { 2199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nstages * nregisters, &t->gamma)); 2209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nstages, &t->delta)); 2219566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->gamma, gamma, nstages * nregisters)); 2229566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->delta, delta, nstages)); 2239ce202e0SLisandro Dalcin } 2249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nstages, &t->betasub)); 2259566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->betasub, betasub, nstages)); 22637e1895aSJed Brown 22737e1895aSJed Brown link->next = SNESMSTableauList; 22837e1895aSJed Brown SNESMSTableauList = link; 2293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 23037e1895aSJed Brown } 23137e1895aSJed Brown 23237e1895aSJed Brown /* 23337e1895aSJed Brown X - initial state, updated in-place. 23437e1895aSJed Brown F - residual, computed at the initial X on input 23537e1895aSJed Brown */ 236d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSStep_3Sstar(SNES snes, Vec X, Vec F) 237d71ae5a4SJacob Faibussowitsch { 23837e1895aSJed Brown SNES_MS *ms = (SNES_MS *)snes->data; 23937e1895aSJed Brown SNESMSTableau t = ms->tableau; 2403e2ddb07SJacob Faibussowitsch const PetscInt nstages = t->nstages; 24137e1895aSJed Brown const PetscReal *gamma = t->gamma, *delta = t->delta, *betasub = t->betasub; 2423e2ddb07SJacob Faibussowitsch Vec S1 = X, S2 = snes->work[1], S3 = snes->work[2], Y = snes->work[0], S1copy = snes->work[3]; 24337e1895aSJed Brown 24437e1895aSJed Brown PetscFunctionBegin; 2459566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(S2)); 2469566063dSJacob Faibussowitsch PetscCall(VecCopy(X, S3)); 2473e2ddb07SJacob Faibussowitsch for (PetscInt i = 0; i < nstages; i++) { 2483e2ddb07SJacob Faibussowitsch Vec Ss[] = {S1copy, S2, S3, Y}; 2493e2ddb07SJacob Faibussowitsch const PetscScalar scoeff[] = {gamma[0 * nstages + i] - 1, gamma[1 * nstages + i], gamma[2 * nstages + i], -betasub[i] * ms->damping}; 2501aa26658SKarl Rupp 2519566063dSJacob Faibussowitsch PetscCall(VecAXPY(S2, delta[i], S1)); 25248a46eb9SPierre Jolivet if (i > 0) PetscCall(SNESComputeFunction(snes, S1, F)); 2539566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, Y)); 2543e2ddb07SJacob Faibussowitsch PetscCall(VecCopy(S1, S1copy)); 2559566063dSJacob Faibussowitsch PetscCall(VecMAXPY(S1, 4, scoeff, Ss)); 25637e1895aSJed Brown } 2573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 25837e1895aSJed Brown } 25937e1895aSJed Brown 2609ce202e0SLisandro Dalcin /* 2619ce202e0SLisandro Dalcin X - initial state, updated in-place. 2629ce202e0SLisandro Dalcin F - residual, computed at the initial X on input 2639ce202e0SLisandro Dalcin */ 264d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSStep_Basic(SNES snes, Vec X, Vec F) 265d71ae5a4SJacob Faibussowitsch { 2669ce202e0SLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 2679ce202e0SLisandro Dalcin SNESMSTableau tab = ms->tableau; 2689ce202e0SLisandro Dalcin const PetscReal *alpha = tab->betasub, h = ms->damping; 2699ce202e0SLisandro Dalcin PetscInt i, nstages = tab->nstages; 2709ce202e0SLisandro Dalcin Vec X0 = snes->work[0]; 2719ce202e0SLisandro Dalcin 2729ce202e0SLisandro Dalcin PetscFunctionBegin; 2739566063dSJacob Faibussowitsch PetscCall(VecCopy(X, X0)); 2749ce202e0SLisandro Dalcin for (i = 0; i < nstages; i++) { 27548a46eb9SPierre Jolivet if (i > 0) PetscCall(SNESComputeFunction(snes, X, F)); 2769566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, X)); 2779566063dSJacob Faibussowitsch PetscCall(VecAYPX(X, -alpha[i] * h, X0)); 2789ce202e0SLisandro Dalcin } 2793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2809ce202e0SLisandro Dalcin } 2819ce202e0SLisandro Dalcin 282d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSStep_Step(SNES snes, Vec X, Vec F) 283d71ae5a4SJacob Faibussowitsch { 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 } 2933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2949ce202e0SLisandro Dalcin } 2959ce202e0SLisandro Dalcin 296d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSStep_Norms(SNES snes, PetscInt iter, Vec F) 297d71ae5a4SJacob Faibussowitsch { 298bf80fd75SLisandro Dalcin PetscReal fnorm; 299bf80fd75SLisandro Dalcin 300bf80fd75SLisandro Dalcin PetscFunctionBegin; 3017a14dab9SStefano Zampini if (SNESNeedNorm_Private(snes, iter)) { 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)); 310bf80fd75SLisandro Dalcin /* Test for convergence */ 3112d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, fnorm)); 3122d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 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 } 3183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 319bf80fd75SLisandro Dalcin } 320bf80fd75SLisandro Dalcin 321d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_MS(SNES snes) 322d71ae5a4SJacob Faibussowitsch { 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)); 3413ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 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 3567a14dab9SStefano Zampini if (i < snes->max_its - 1 || SNESNeedNorm_Private(snes, i + 1)) PetscCall(SNESComputeFunction(snes, X, F)); 35737e1895aSJed Brown 3589566063dSJacob Faibussowitsch PetscCall(SNESMSStep_Norms(snes, i + 1, F)); 3593ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 36037e1895aSJed Brown } 3613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 36237e1895aSJed Brown } 36337e1895aSJed Brown 364d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_MS(SNES snes) 365d71ae5a4SJacob Faibussowitsch { 3669ce202e0SLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 3679ce202e0SLisandro Dalcin SNESMSTableau tab = ms->tableau; 3683e2ddb07SJacob Faibussowitsch PetscInt nwork = tab->nregisters + 1; // +1 because VecMAXPY() in SNESMSStep_3Sstar() 3693e2ddb07SJacob Faibussowitsch // needs an extra work vec 37037e1895aSJed Brown 37137e1895aSJed Brown PetscFunctionBegin; 3729566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, nwork)); 3739566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 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(PetscFree(snes->data)); 3819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetType_C", NULL)); 3829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetType_C", NULL)); 3839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetDamping_C", NULL)); 3849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetDamping_C", NULL)); 3853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38637e1895aSJed Brown } 38737e1895aSJed Brown 388d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_MS(SNES snes, PetscViewer viewer) 389d71ae5a4SJacob Faibussowitsch { 390*9f196a02SMartin Diehl PetscBool isascii; 39137e1895aSJed Brown SNES_MS *ms = (SNES_MS *)snes->data; 39257715debSLisandro Dalcin SNESMSTableau tab = ms->tableau; 39337e1895aSJed Brown 39437e1895aSJed Brown PetscFunctionBegin; 395*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 396*9f196a02SMartin Diehl if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, " multi-stage method type: %s\n", tab->name)); 3973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 39837e1895aSJed Brown } 39937e1895aSJed Brown 400ce78bad3SBarry Smith static PetscErrorCode SNESSetFromOptions_MS(SNES snes, PetscOptionItems PetscOptionsObject) 401d71ae5a4SJacob Faibussowitsch { 40237e1895aSJed Brown PetscFunctionBegin; 403d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES MS options"); 40437e1895aSJed Brown { 40537e1895aSJed Brown SNESMSTableauLink link; 40637e1895aSJed Brown PetscInt count, choice; 40737e1895aSJed Brown PetscBool flg; 40837e1895aSJed Brown const char **namelist; 40957715debSLisandro Dalcin SNESMSType mstype; 41057715debSLisandro Dalcin PetscReal damping; 4117a14dab9SStefano Zampini PetscBool norms = PETSC_FALSE; 41237e1895aSJed Brown 4139566063dSJacob Faibussowitsch PetscCall(SNESMSGetType(snes, &mstype)); 414fbccb6d4SPierre Jolivet for (link = SNESMSTableauList, count = 0; link; link = link->next, count++); 4159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, (char ***)&namelist)); 41637e1895aSJed Brown for (link = SNESMSTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name; 4179566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-snes_ms_type", "Multistage smoother type", "SNESMSSetType", (const char *const *)namelist, count, mstype, &choice, &flg)); 4189566063dSJacob Faibussowitsch if (flg) PetscCall(SNESMSSetType(snes, namelist[choice])); 4199566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 4209566063dSJacob Faibussowitsch PetscCall(SNESMSGetDamping(snes, &damping)); 4219566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ms_damping", "Damping for multistage method", "SNESMSSetDamping", damping, &damping, &flg)); 4229566063dSJacob Faibussowitsch if (flg) PetscCall(SNESMSSetDamping(snes, damping)); 4237a14dab9SStefano Zampini 4247a14dab9SStefano Zampini /* deprecated option */ 4257a14dab9SStefano Zampini PetscCall(PetscOptionsDeprecated("-snes_ms_norms", NULL, "3.20", "Use -snes_norm_schedule always")); 4267a14dab9SStefano Zampini PetscCall(PetscOptionsBool("-snes_ms_norms", NULL, NULL, norms, &norms, NULL)); 4277a14dab9SStefano Zampini if (norms) PetscCall(SNESSetNormSchedule(snes, SNES_NORM_ALWAYS)); 42837e1895aSJed Brown } 429d0609cedSBarry Smith PetscOptionsHeadEnd(); 4303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43137e1895aSJed Brown } 43237e1895aSJed Brown 433d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSGetType_MS(SNES snes, SNESMSType *mstype) 434d71ae5a4SJacob Faibussowitsch { 43557715debSLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 43657715debSLisandro Dalcin SNESMSTableau tab = ms->tableau; 43757715debSLisandro Dalcin 43857715debSLisandro Dalcin PetscFunctionBegin; 43957715debSLisandro Dalcin *mstype = tab->name; 4403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44157715debSLisandro Dalcin } 44257715debSLisandro Dalcin 443d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSSetType_MS(SNES snes, SNESMSType mstype) 444d71ae5a4SJacob Faibussowitsch { 44537e1895aSJed Brown SNES_MS *ms = (SNES_MS *)snes->data; 44637e1895aSJed Brown SNESMSTableauLink link; 44737e1895aSJed Brown PetscBool match; 44837e1895aSJed Brown 44937e1895aSJed Brown PetscFunctionBegin; 45037e1895aSJed Brown if (ms->tableau) { 4519566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ms->tableau->name, mstype, &match)); 4523ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 45337e1895aSJed Brown } 45437e1895aSJed Brown for (link = SNESMSTableauList; link; link = link->next) { 4559566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name, mstype, &match)); 45637e1895aSJed Brown if (match) { 45737e1895aSJed Brown ms->tableau = &link->tab; 4589566063dSJacob Faibussowitsch if (snes->setupcalled) PetscCall(SNESSetUp_MS(snes)); 4593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46037e1895aSJed Brown } 46137e1895aSJed Brown } 46298921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", mstype); 46337e1895aSJed Brown } 46437e1895aSJed Brown 465cc4c1da9SBarry Smith /*@ 466f6dfbefdSBarry Smith SNESMSGetType - Get the type of multistage smoother `SNESMS` 46757715debSLisandro Dalcin 46820f4b53cSBarry Smith Not Collective 46957715debSLisandro Dalcin 47057715debSLisandro Dalcin Input Parameter: 47157715debSLisandro Dalcin . snes - nonlinear solver context 47257715debSLisandro Dalcin 47357715debSLisandro Dalcin Output Parameter: 47457715debSLisandro Dalcin . mstype - type of multistage method 47557715debSLisandro Dalcin 476f6dfbefdSBarry Smith Level: advanced 47757715debSLisandro Dalcin 478420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESMS`, `SNESMSSetType()`, `SNESMSType` 47957715debSLisandro Dalcin @*/ 480d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSGetType(SNES snes, SNESMSType *mstype) 481d71ae5a4SJacob Faibussowitsch { 48257715debSLisandro Dalcin PetscFunctionBegin; 48357715debSLisandro Dalcin PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 4844f572ea9SToby Isaac PetscAssertPointer(mstype, 2); 485cac4c232SBarry Smith PetscUseMethod(snes, "SNESMSGetType_C", (SNES, SNESMSType *), (snes, mstype)); 4863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 48757715debSLisandro Dalcin } 48857715debSLisandro Dalcin 489cc4c1da9SBarry Smith /*@ 490f6dfbefdSBarry Smith SNESMSSetType - Set the type of multistage smoother `SNESMS` 49137e1895aSJed Brown 49220f4b53cSBarry Smith Logically Collective 49337e1895aSJed Brown 494d8d19677SJose E. Roman Input Parameters: 49537e1895aSJed Brown + snes - nonlinear solver context 49637e1895aSJed Brown - mstype - type of multistage method 49737e1895aSJed Brown 498f6dfbefdSBarry Smith Level: advanced 49937e1895aSJed Brown 500420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESMS`, `SNESMSGetType()`, `SNESMSType` 50137e1895aSJed Brown @*/ 502d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSSetType(SNES snes, SNESMSType mstype) 503d71ae5a4SJacob Faibussowitsch { 50437e1895aSJed Brown PetscFunctionBegin; 50537e1895aSJed Brown PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 5064f572ea9SToby Isaac PetscAssertPointer(mstype, 2); 507cac4c232SBarry Smith PetscTryMethod(snes, "SNESMSSetType_C", (SNES, SNESMSType), (snes, mstype)); 5083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 50957715debSLisandro Dalcin } 51057715debSLisandro Dalcin 511d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSGetDamping_MS(SNES snes, PetscReal *damping) 512d71ae5a4SJacob Faibussowitsch { 51357715debSLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 51457715debSLisandro Dalcin 51557715debSLisandro Dalcin PetscFunctionBegin; 51657715debSLisandro Dalcin *damping = ms->damping; 5173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51857715debSLisandro Dalcin } 51957715debSLisandro Dalcin 520d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESMSSetDamping_MS(SNES snes, PetscReal damping) 521d71ae5a4SJacob Faibussowitsch { 52257715debSLisandro Dalcin SNES_MS *ms = (SNES_MS *)snes->data; 52357715debSLisandro Dalcin 52457715debSLisandro Dalcin PetscFunctionBegin; 52557715debSLisandro Dalcin ms->damping = damping; 5263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 52757715debSLisandro Dalcin } 52857715debSLisandro Dalcin 5295d83a8b1SBarry Smith /*@ 530f6dfbefdSBarry Smith SNESMSGetDamping - Get the damping parameter of `SNESMS` multistage scheme 53157715debSLisandro Dalcin 53220f4b53cSBarry Smith Not Collective 53357715debSLisandro Dalcin 53457715debSLisandro Dalcin Input Parameter: 53557715debSLisandro Dalcin . snes - nonlinear solver context 53657715debSLisandro Dalcin 53757715debSLisandro Dalcin Output Parameter: 53857715debSLisandro Dalcin . damping - damping parameter 53957715debSLisandro Dalcin 54057715debSLisandro Dalcin Level: advanced 54157715debSLisandro Dalcin 542420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESMSSetDamping()`, `SNESMS` 54357715debSLisandro Dalcin @*/ 544d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSGetDamping(SNES snes, PetscReal *damping) 545d71ae5a4SJacob Faibussowitsch { 54657715debSLisandro Dalcin PetscFunctionBegin; 54757715debSLisandro Dalcin PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 5484f572ea9SToby Isaac PetscAssertPointer(damping, 2); 549cac4c232SBarry Smith PetscUseMethod(snes, "SNESMSGetDamping_C", (SNES, PetscReal *), (snes, damping)); 5503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 55157715debSLisandro Dalcin } 55257715debSLisandro Dalcin 5535d83a8b1SBarry Smith /*@ 554f6dfbefdSBarry Smith SNESMSSetDamping - Set the damping parameter for a `SNESMS` multistage scheme 55557715debSLisandro Dalcin 55620f4b53cSBarry Smith Logically Collective 55757715debSLisandro Dalcin 558d8d19677SJose E. Roman Input Parameters: 55957715debSLisandro Dalcin + snes - nonlinear solver context 56057715debSLisandro Dalcin - damping - damping parameter 56157715debSLisandro Dalcin 56257715debSLisandro Dalcin Level: advanced 56357715debSLisandro Dalcin 564420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESMSGetDamping()`, `SNESMS` 56557715debSLisandro Dalcin @*/ 566d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMSSetDamping(SNES snes, PetscReal damping) 567d71ae5a4SJacob Faibussowitsch { 56857715debSLisandro Dalcin PetscFunctionBegin; 56957715debSLisandro Dalcin PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 57057715debSLisandro Dalcin PetscValidLogicalCollectiveReal(snes, damping, 2); 571cac4c232SBarry Smith PetscTryMethod(snes, "SNESMSSetDamping_C", (SNES, PetscReal), (snes, damping)); 5723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 57337e1895aSJed Brown } 57437e1895aSJed Brown 57537e1895aSJed Brown /*MC 57637e1895aSJed Brown SNESMS - multi-stage smoothers 57737e1895aSJed Brown 578f6dfbefdSBarry Smith Options Database Keys: 57937e1895aSJed Brown + -snes_ms_type - type of multi-stage smoother 58037e1895aSJed Brown - -snes_ms_damping - damping for multi-stage method 58137e1895aSJed Brown 582420bcc1bSBarry Smith Level: advanced 583420bcc1bSBarry Smith 58437e1895aSJed Brown Notes: 58537e1895aSJed Brown These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for 5866c9de887SHong Zhang FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebyshev). 58737e1895aSJed Brown 58837e1895aSJed Brown Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds. 58937e1895aSJed Brown 590f6dfbefdSBarry Smith The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with `SNESMSRegister()`. 59137e1895aSJed Brown 5921d27aa22SBarry Smith See {cite}`ketcheson2010runge`, {cite}`jameson1983`, {cite}`pierce1997preconditioned`, and {cite}`va1981design` 59337e1895aSJed Brown 594420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESMS`, `SNESFAS`, `KSPCHEBYSHEV`, `SNESMSSetDamping()`, `SNESMSGetDamping()`, 595420bcc1bSBarry Smith `SNESMSSetType()`, `SNESMSGetType()` 59637e1895aSJed Brown M*/ 597420bcc1bSBarry Smith 598d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_MS(SNES snes) 599d71ae5a4SJacob Faibussowitsch { 60037e1895aSJed Brown SNES_MS *ms; 60137e1895aSJed Brown 60237e1895aSJed Brown PetscFunctionBegin; 6039566063dSJacob Faibussowitsch PetscCall(SNESMSInitializePackage()); 60437e1895aSJed Brown 60537e1895aSJed Brown snes->ops->setup = SNESSetUp_MS; 60637e1895aSJed Brown snes->ops->solve = SNESSolve_MS; 60737e1895aSJed Brown snes->ops->destroy = SNESDestroy_MS; 60837e1895aSJed Brown snes->ops->setfromoptions = SNESSetFromOptions_MS; 60937e1895aSJed Brown snes->ops->view = SNESView_MS; 61037e1895aSJed Brown 611efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 61237e1895aSJed Brown snes->usesksp = PETSC_TRUE; 61337e1895aSJed Brown 6144fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 6154fc747eaSLawrence Mitchell 61677e5a1f9SBarry Smith PetscCall(SNESParametersInitialize(snes)); 61777e5a1f9SBarry Smith 6184dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ms)); 61937e1895aSJed Brown snes->data = (void *)ms; 62037e1895aSJed Brown ms->damping = 0.9; 62137e1895aSJed Brown 6229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetType_C", SNESMSGetType_MS)); 6239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetType_C", SNESMSSetType_MS)); 6249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetDamping_C", SNESMSGetDamping_MS)); 6259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetDamping_C", SNESMSSetDamping_MS)); 62657715debSLisandro Dalcin 6279566063dSJacob Faibussowitsch PetscCall(SNESMSSetType(snes, SNESMSDefault)); 6283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 62937e1895aSJed Brown } 630