xref: /petsc/src/snes/impls/ms/ms.c (revision d8d19677bbccf95218448bee62e6b87f4513e133)
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
3237e1895aSJed Brown   SNESMSRegisterAll - Registers all of the multi-stage methods in SNESMS
3337e1895aSJed Brown 
3437e1895aSJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
3537e1895aSJed Brown 
3637e1895aSJed Brown   Level: advanced
3737e1895aSJed Brown 
3837e1895aSJed Brown .seealso:  SNESMSRegisterDestroy()
3937e1895aSJed Brown @*/
4037e1895aSJed Brown PetscErrorCode SNESMSRegisterAll(void)
4137e1895aSJed Brown {
4237e1895aSJed Brown   PetscErrorCode ierr;
4337e1895aSJed Brown 
4437e1895aSJed Brown   PetscFunctionBegin;
4537e1895aSJed Brown   if (SNESMSRegisterAllCalled) PetscFunctionReturn(0);
4637e1895aSJed Brown   SNESMSRegisterAllCalled = PETSC_TRUE;
4737e1895aSJed Brown 
4837e1895aSJed Brown   {
499ce202e0SLisandro Dalcin     PetscReal alpha[1] = {1.0};
509ce202e0SLisandro Dalcin     ierr = SNESMSRegister(SNESMSEULER,1,1,1.0,NULL,NULL,alpha);CHKERRQ(ierr);
5137e1895aSJed Brown   }
5237e1895aSJed Brown 
5337e1895aSJed Brown   {
5437e1895aSJed Brown     PetscReal gamma[3][6] = {
5537e1895aSJed Brown       {0.0000000000000000E+00, -7.0304722367110606E-01, -1.9836719667506464E-01, -1.6023843981863788E+00,  9.4483822882855284E-02, -1.4204296130641869E-01},
5637e1895aSJed Brown       {1.0000000000000000E+00,  1.1111025767083920E+00,  5.6150921583923230E-01,  7.4151723494934041E-01,  3.1714538168600587E-01,  4.6479276238548706E-01},
5737e1895aSJed Brown       {0.0000000000000000E+00,  0.0000000000000000E+00,  0.0000000000000000E+00,  6.7968174970583317E-01, -4.1755042846051737E-03, -1.9115668129923846E-01}
5837e1895aSJed Brown     };
5937e1895aSJed Brown     PetscReal delta[6]   = {1.0000000000000000E+00, 5.3275427433201750E-01, 6.0143544663985238E-01, 4.5874077053842177E-01, 2.7544386906104651E-01, 0.0000000000000000E+00};
6037e1895aSJed Brown     PetscReal betasub[6] = {8.4753115429481929E-01, 7.4018896368655618E-01, 6.5963574086583309E-03, 4.6747795645517759E-01, 1.3314545813643919E-01, 5.3260800028018784E-01};
6137e1895aSJed Brown     ierr = SNESMSRegister(SNESMSM62,6,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
6237e1895aSJed Brown   }
639ce202e0SLisandro Dalcin 
649ce202e0SLisandro Dalcin   { /* Jameson (1983) */
659ce202e0SLisandro Dalcin     PetscReal alpha[4] = {0.25, 0.5, 0.55, 1.0};
669ce202e0SLisandro Dalcin     ierr = SNESMSRegister(SNESMSJAMESON83,4,1,1.0,NULL,NULL,alpha);CHKERRQ(ierr);
67a97cb6bcSJed Brown   }
683847c725SLisandro Dalcin 
693847c725SLisandro Dalcin   { /* Van Leer, Tai, and Powell (1989) 1 stage, order 1 */
709ce202e0SLisandro Dalcin     PetscReal alpha[1]  = {1.0};
719ce202e0SLisandro Dalcin     ierr = SNESMSRegister(SNESMSVLTP11,1,1,0.5,NULL,NULL,alpha);CHKERRQ(ierr);
723847c725SLisandro Dalcin   }
73a97cb6bcSJed Brown   { /* Van Leer, Tai, and Powell (1989) 2 stage, order 1 */
749ce202e0SLisandro Dalcin     PetscReal alpha[2] = {0.3333, 1.0};
759ce202e0SLisandro Dalcin     ierr = SNESMSRegister(SNESMSVLTP21,2,1,1.0,NULL,NULL,alpha);CHKERRQ(ierr);
76a97cb6bcSJed Brown   }
77a97cb6bcSJed Brown   { /* Van Leer, Tai, and Powell (1989) 3 stage, order 1 */
789ce202e0SLisandro Dalcin     PetscReal alpha[3] = {0.1481, 0.4000, 1.0};
799ce202e0SLisandro Dalcin     ierr = SNESMSRegister(SNESMSVLTP31,3,1,1.5,NULL,NULL,alpha);CHKERRQ(ierr);
80a97cb6bcSJed Brown   }
81a97cb6bcSJed Brown   { /* Van Leer, Tai, and Powell (1989) 4 stage, order 1 */
829ce202e0SLisandro Dalcin     PetscReal alpha[4] = {0.0833, 0.2069, 0.4265, 1.0};
839ce202e0SLisandro Dalcin     ierr = SNESMSRegister(SNESMSVLTP41,4,1,2.0,NULL,NULL,alpha);CHKERRQ(ierr);
84a97cb6bcSJed Brown   }
85a97cb6bcSJed Brown   { /* Van Leer, Tai, and Powell (1989) 5 stage, order 1 */
869ce202e0SLisandro Dalcin     PetscReal alpha[5] = {0.0533, 0.1263, 0.2375, 0.4414,1.0};
879ce202e0SLisandro Dalcin     ierr = SNESMSRegister(SNESMSVLTP51,5,1,2.5,NULL,NULL,alpha);CHKERRQ(ierr);
88a97cb6bcSJed Brown   }
89a97cb6bcSJed Brown   { /* Van Leer, Tai, and Powell (1989) 6 stage, order 1 */
909ce202e0SLisandro Dalcin     PetscReal alpha[6] = {0.0370, 0.0851, 0.1521, 0.2562, 0.4512, 1.0};
919ce202e0SLisandro Dalcin     ierr = SNESMSRegister(SNESMSVLTP61,6,1,3.0,NULL,NULL,alpha);CHKERRQ(ierr);
92a97cb6bcSJed Brown   }
9337e1895aSJed Brown   PetscFunctionReturn(0);
9437e1895aSJed Brown }
9537e1895aSJed Brown 
9637e1895aSJed Brown /*@C
9757715debSLisandro Dalcin    SNESMSRegisterDestroy - Frees the list of schemes that were registered by SNESMSRegister().
9837e1895aSJed Brown 
9937e1895aSJed Brown    Not Collective
10037e1895aSJed Brown 
10137e1895aSJed Brown    Level: advanced
10237e1895aSJed Brown 
10357715debSLisandro Dalcin .seealso: SNESMSRegister(), SNESMSRegisterAll()
10437e1895aSJed Brown @*/
10537e1895aSJed Brown PetscErrorCode SNESMSRegisterDestroy(void)
10637e1895aSJed Brown {
10737e1895aSJed Brown   PetscErrorCode    ierr;
10837e1895aSJed Brown   SNESMSTableauLink link;
10937e1895aSJed Brown 
11037e1895aSJed Brown   PetscFunctionBegin;
11137e1895aSJed Brown   while ((link = SNESMSTableauList)) {
11237e1895aSJed Brown     SNESMSTableau t = &link->tab;
11337e1895aSJed Brown     SNESMSTableauList = link->next;
1141aa26658SKarl Rupp 
11537e1895aSJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
1169ce202e0SLisandro Dalcin     ierr = PetscFree(t->gamma);CHKERRQ(ierr);
1179ce202e0SLisandro Dalcin     ierr = PetscFree(t->delta);CHKERRQ(ierr);
1189ce202e0SLisandro Dalcin     ierr = PetscFree(t->betasub);CHKERRQ(ierr);
11937e1895aSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
12037e1895aSJed Brown   }
12137e1895aSJed Brown   SNESMSRegisterAllCalled = PETSC_FALSE;
12237e1895aSJed Brown   PetscFunctionReturn(0);
12337e1895aSJed Brown }
12437e1895aSJed Brown 
12537e1895aSJed Brown /*@C
12637e1895aSJed Brown   SNESMSInitializePackage - This function initializes everything in the SNESMS package. It is called
1278a690491SBarry Smith   from SNESInitializePackage().
12837e1895aSJed Brown 
12937e1895aSJed Brown   Level: developer
13037e1895aSJed Brown 
13137e1895aSJed Brown .seealso: PetscInitialize()
13237e1895aSJed Brown @*/
133607a6623SBarry Smith PetscErrorCode SNESMSInitializePackage(void)
13437e1895aSJed Brown {
13537e1895aSJed Brown   PetscErrorCode ierr;
13637e1895aSJed Brown 
13737e1895aSJed Brown   PetscFunctionBegin;
13837e1895aSJed Brown   if (SNESMSPackageInitialized) PetscFunctionReturn(0);
13937e1895aSJed Brown   SNESMSPackageInitialized = PETSC_TRUE;
1401aa26658SKarl Rupp 
14137e1895aSJed Brown   ierr = SNESMSRegisterAll();CHKERRQ(ierr);
14237e1895aSJed Brown   ierr = PetscRegisterFinalize(SNESMSFinalizePackage);CHKERRQ(ierr);
14337e1895aSJed Brown   PetscFunctionReturn(0);
14437e1895aSJed Brown }
14537e1895aSJed Brown 
14637e1895aSJed Brown /*@C
14737e1895aSJed Brown   SNESMSFinalizePackage - This function destroys everything in the SNESMS package. It is
14837e1895aSJed Brown   called from PetscFinalize().
14937e1895aSJed Brown 
15037e1895aSJed Brown   Level: developer
15137e1895aSJed Brown 
15237e1895aSJed Brown .seealso: PetscFinalize()
15337e1895aSJed Brown @*/
15437e1895aSJed Brown PetscErrorCode SNESMSFinalizePackage(void)
15537e1895aSJed Brown {
15637e1895aSJed Brown   PetscErrorCode ierr;
15737e1895aSJed Brown 
15837e1895aSJed Brown   PetscFunctionBegin;
15937e1895aSJed Brown   SNESMSPackageInitialized = PETSC_FALSE;
1601aa26658SKarl Rupp 
16137e1895aSJed Brown   ierr = SNESMSRegisterDestroy();CHKERRQ(ierr);
16237e1895aSJed Brown   PetscFunctionReturn(0);
16337e1895aSJed Brown }
16437e1895aSJed Brown 
16537e1895aSJed Brown /*@C
16637e1895aSJed Brown    SNESMSRegister - register a multistage scheme
16737e1895aSJed Brown 
16837e1895aSJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
16937e1895aSJed Brown 
17037e1895aSJed Brown    Input Parameters:
17137e1895aSJed Brown +  name - identifier for method
17237e1895aSJed Brown .  nstages - number of stages
17337e1895aSJed Brown .  nregisters - number of registers used by low-storage implementation
1749ce202e0SLisandro Dalcin .  stability - scaled stability region
17537e1895aSJed Brown .  gamma - coefficients, see Ketcheson's paper
17637e1895aSJed Brown .  delta - coefficients, see Ketcheson's paper
17737e1895aSJed Brown -  betasub - subdiagonal of Shu-Osher form
17837e1895aSJed Brown 
17937e1895aSJed Brown    Notes:
18037e1895aSJed Brown    The notation is described in Ketcheson (2010) Runge-Kutta methods with minimum storage implementations.
18137e1895aSJed Brown 
1829ce202e0SLisandro Dalcin    Many multistage schemes are of the form
1839ce202e0SLisandro Dalcin    $ X_0 = X^{(n)}
1849ce202e0SLisandro Dalcin    $ X_k = X_0 + \alpha_k * F(X_{k-1}), k = 1,\ldots,s
1859ce202e0SLisandro Dalcin    $ X^{(n+1)} = X_s
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 
19137e1895aSJed Brown    Level: advanced
19237e1895aSJed Brown 
19337e1895aSJed Brown .seealso: SNESMS
19437e1895aSJed Brown @*/
19519fd82e9SBarry Smith PetscErrorCode SNESMSRegister(SNESMSType name,PetscInt nstages,PetscInt nregisters,PetscReal stability,const PetscReal gamma[],const PetscReal delta[],const PetscReal betasub[])
19637e1895aSJed Brown {
19737e1895aSJed Brown   PetscErrorCode    ierr;
19837e1895aSJed Brown   SNESMSTableauLink link;
19937e1895aSJed Brown   SNESMSTableau     t;
20037e1895aSJed Brown 
20137e1895aSJed Brown   PetscFunctionBegin;
20237e1895aSJed Brown   PetscValidCharPointer(name,1);
20337e1895aSJed Brown   if (nstages < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have at least one stage");
2049ce202e0SLisandro Dalcin   if (gamma || delta) {
20537e1895aSJed Brown     if (nregisters != 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only support for methods written in 3-register form");
206064a246eSJacob Faibussowitsch     PetscValidRealPointer(gamma,5);
207064a246eSJacob Faibussowitsch     PetscValidRealPointer(delta,6);
2089ce202e0SLisandro Dalcin   } else {
2099ce202e0SLisandro Dalcin     if (nregisters != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only support for methods written in 1-register form");
2109ce202e0SLisandro Dalcin   }
211064a246eSJacob Faibussowitsch   PetscValidRealPointer(betasub,7);
21237e1895aSJed Brown 
2139cc31a68SJed Brown   ierr          = SNESMSInitializePackage();CHKERRQ(ierr);
214854ce69bSBarry Smith   ierr          = PetscNew(&link);CHKERRQ(ierr);
21537e1895aSJed Brown   t             = &link->tab;
21637e1895aSJed Brown   ierr          = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
21737e1895aSJed Brown   t->nstages    = nstages;
21837e1895aSJed Brown   t->nregisters = nregisters;
21937e1895aSJed Brown   t->stability  = stability;
2201aa26658SKarl Rupp 
2219ce202e0SLisandro Dalcin   if (gamma && delta) {
2229ce202e0SLisandro Dalcin     ierr = PetscMalloc1(nstages*nregisters,&t->gamma);CHKERRQ(ierr);
2239ce202e0SLisandro Dalcin     ierr = PetscMalloc1(nstages,&t->delta);CHKERRQ(ierr);
224580bdb30SBarry Smith     ierr = PetscArraycpy(t->gamma,gamma,nstages*nregisters);CHKERRQ(ierr);
225580bdb30SBarry Smith     ierr = PetscArraycpy(t->delta,delta,nstages);CHKERRQ(ierr);
2269ce202e0SLisandro Dalcin   }
2279ce202e0SLisandro Dalcin   ierr = PetscMalloc1(nstages,&t->betasub);CHKERRQ(ierr);
228580bdb30SBarry Smith   ierr = PetscArraycpy(t->betasub,betasub,nstages);CHKERRQ(ierr);
22937e1895aSJed Brown 
23037e1895aSJed Brown   link->next        = SNESMSTableauList;
23137e1895aSJed Brown   SNESMSTableauList = link;
23237e1895aSJed Brown   PetscFunctionReturn(0);
23337e1895aSJed Brown }
23437e1895aSJed Brown 
23537e1895aSJed Brown /*
23637e1895aSJed Brown   X - initial state, updated in-place.
23737e1895aSJed Brown   F - residual, computed at the initial X on input
23837e1895aSJed Brown */
23937e1895aSJed Brown static PetscErrorCode SNESMSStep_3Sstar(SNES snes,Vec X,Vec F)
24037e1895aSJed Brown {
24137e1895aSJed Brown   PetscErrorCode  ierr;
24237e1895aSJed Brown   SNES_MS         *ms    = (SNES_MS*)snes->data;
24337e1895aSJed Brown   SNESMSTableau   t      = ms->tableau;
24437e1895aSJed Brown   const PetscReal *gamma = t->gamma,*delta = t->delta,*betasub = t->betasub;
24537e1895aSJed Brown   Vec             S1,S2,S3,Y;
246bf80fd75SLisandro Dalcin   PetscInt        i,nstages = t->nstages;
24737e1895aSJed Brown 
24837e1895aSJed Brown   PetscFunctionBegin;
24937e1895aSJed Brown   Y    = snes->work[0];
25037e1895aSJed Brown   S1   = X;
25137e1895aSJed Brown   S2   = snes->work[1];
25237e1895aSJed Brown   S3   = snes->work[2];
25337e1895aSJed Brown   ierr = VecZeroEntries(S2);CHKERRQ(ierr);
25437e1895aSJed Brown   ierr = VecCopy(X,S3);CHKERRQ(ierr);
25537e1895aSJed Brown   for (i = 0; i < nstages; i++) {
256da80777bSKarl Rupp     Vec         Ss[4];
257da80777bSKarl Rupp     PetscScalar scoeff[4];
258da80777bSKarl Rupp 
259da80777bSKarl Rupp     Ss[0] = S1; Ss[1] = S2; Ss[2] = S3; Ss[3] = Y;
2601aa26658SKarl Rupp 
261bf80fd75SLisandro Dalcin     scoeff[0] = gamma[0*nstages+i] - 1;
262da80777bSKarl Rupp     scoeff[1] = gamma[1*nstages+i];
263da80777bSKarl Rupp     scoeff[2] = gamma[2*nstages+i];
264da80777bSKarl Rupp     scoeff[3] = -betasub[i]*ms->damping;
2651aa26658SKarl Rupp 
26637e1895aSJed Brown     ierr = VecAXPY(S2,delta[i],S1);CHKERRQ(ierr);
26737e1895aSJed Brown     if (i > 0) {
26837e1895aSJed Brown       ierr = SNESComputeFunction(snes,S1,F);CHKERRQ(ierr);
26937e1895aSJed Brown     }
270d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
27137e1895aSJed Brown     ierr = VecMAXPY(S1,4,scoeff,Ss);CHKERRQ(ierr);
27237e1895aSJed Brown   }
27337e1895aSJed Brown   PetscFunctionReturn(0);
27437e1895aSJed Brown }
27537e1895aSJed Brown 
2769ce202e0SLisandro Dalcin /*
2779ce202e0SLisandro Dalcin   X - initial state, updated in-place.
2789ce202e0SLisandro Dalcin   F - residual, computed at the initial X on input
2799ce202e0SLisandro Dalcin */
2809ce202e0SLisandro Dalcin static PetscErrorCode SNESMSStep_Basic(SNES snes,Vec X,Vec F)
2819ce202e0SLisandro Dalcin {
2829ce202e0SLisandro Dalcin   PetscErrorCode  ierr;
2839ce202e0SLisandro Dalcin   SNES_MS         *ms    = (SNES_MS*)snes->data;
2849ce202e0SLisandro Dalcin   SNESMSTableau   tab    = ms->tableau;
2859ce202e0SLisandro Dalcin   const PetscReal *alpha = tab->betasub, h = ms->damping;
2869ce202e0SLisandro Dalcin   PetscInt        i,nstages = tab->nstages;
2879ce202e0SLisandro Dalcin   Vec             X0 = snes->work[0];
2889ce202e0SLisandro Dalcin 
2899ce202e0SLisandro Dalcin   PetscFunctionBegin;
2909ce202e0SLisandro Dalcin   ierr = VecCopy(X,X0);CHKERRQ(ierr);
2919ce202e0SLisandro Dalcin   for (i = 0; i < nstages; i++) {
2929ce202e0SLisandro Dalcin     if (i > 0) {
2939ce202e0SLisandro Dalcin       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
2949ce202e0SLisandro Dalcin     }
2959ce202e0SLisandro Dalcin     ierr = KSPSolve(snes->ksp,F,X);CHKERRQ(ierr);
2969ce202e0SLisandro Dalcin     ierr = VecAYPX(X,-alpha[i]*h,X0);CHKERRQ(ierr);
2979ce202e0SLisandro Dalcin   }
2989ce202e0SLisandro Dalcin   PetscFunctionReturn(0);
2999ce202e0SLisandro Dalcin }
3009ce202e0SLisandro Dalcin 
3019ce202e0SLisandro Dalcin static PetscErrorCode SNESMSStep_Step(SNES snes,Vec X,Vec F)
3029ce202e0SLisandro Dalcin {
3039ce202e0SLisandro Dalcin   SNES_MS        *ms = (SNES_MS*)snes->data;
3049ce202e0SLisandro Dalcin   SNESMSTableau  tab = ms->tableau;
3059ce202e0SLisandro Dalcin   PetscErrorCode ierr;
3069ce202e0SLisandro Dalcin 
3079ce202e0SLisandro Dalcin   PetscFunctionBegin;
3089ce202e0SLisandro Dalcin   if (tab->gamma && tab->delta) {
3099ce202e0SLisandro Dalcin     ierr = SNESMSStep_3Sstar(snes,X,F);CHKERRQ(ierr);
3109ce202e0SLisandro Dalcin   } else {
3119ce202e0SLisandro Dalcin     ierr = SNESMSStep_Basic(snes,X,F);CHKERRQ(ierr);
3129ce202e0SLisandro Dalcin   }
3139ce202e0SLisandro Dalcin   PetscFunctionReturn(0);
3149ce202e0SLisandro Dalcin }
3159ce202e0SLisandro Dalcin 
316bf80fd75SLisandro Dalcin static PetscErrorCode SNESMSStep_Norms(SNES snes,PetscInt iter,Vec F)
317bf80fd75SLisandro Dalcin {
318bf80fd75SLisandro Dalcin   SNES_MS        *ms = (SNES_MS*)snes->data;
319bf80fd75SLisandro Dalcin   PetscReal      fnorm;
320bf80fd75SLisandro Dalcin   PetscErrorCode ierr;
321bf80fd75SLisandro Dalcin 
322bf80fd75SLisandro Dalcin   PetscFunctionBegin;
323bf80fd75SLisandro Dalcin   if (ms->norms) {
324bf80fd75SLisandro Dalcin     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
325bf80fd75SLisandro Dalcin     SNESCheckFunctionNorm(snes,fnorm);
326bf80fd75SLisandro Dalcin     /* Monitor convergence */
327bf80fd75SLisandro Dalcin     ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
328bf80fd75SLisandro Dalcin     snes->iter = iter;
329bf80fd75SLisandro Dalcin     snes->norm = fnorm;
330bf80fd75SLisandro Dalcin     ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
331bf80fd75SLisandro Dalcin     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
332bf80fd75SLisandro Dalcin     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
333bf80fd75SLisandro Dalcin     /* Test for convergence */
334bf80fd75SLisandro Dalcin     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
335bf80fd75SLisandro Dalcin   } else if (iter > 0) {
336bf80fd75SLisandro Dalcin     ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
337bf80fd75SLisandro Dalcin     snes->iter = iter;
338bf80fd75SLisandro Dalcin     ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
339bf80fd75SLisandro Dalcin   }
340bf80fd75SLisandro Dalcin   PetscFunctionReturn(0);
341bf80fd75SLisandro Dalcin }
342bf80fd75SLisandro Dalcin 
34337e1895aSJed Brown static PetscErrorCode SNESSolve_MS(SNES snes)
34437e1895aSJed Brown {
34537e1895aSJed Brown   SNES_MS        *ms = (SNES_MS*)snes->data;
34637e1895aSJed Brown   Vec            X   = snes->vec_sol,F = snes->vec_func;
34737e1895aSJed Brown   PetscInt       i;
34837e1895aSJed Brown   PetscErrorCode ierr;
34937e1895aSJed Brown 
35037e1895aSJed Brown   PetscFunctionBegin;
3516c4ed002SBarry Smith   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
352fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
353bf80fd75SLisandro Dalcin 
35437e1895aSJed Brown   snes->reason = SNES_CONVERGED_ITERATING;
355e04113cfSBarry Smith   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
35637e1895aSJed Brown   snes->iter   = 0;
357bf80fd75SLisandro Dalcin   snes->norm   = 0;
358e04113cfSBarry Smith   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
359bf80fd75SLisandro Dalcin 
360e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
36137e1895aSJed Brown     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
3621aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
3631aa26658SKarl Rupp 
364bf80fd75SLisandro Dalcin   ierr = SNESMSStep_Norms(snes,0,F);CHKERRQ(ierr);
36537e1895aSJed Brown   if (snes->reason) PetscFunctionReturn(0);
366bf80fd75SLisandro Dalcin 
367bf80fd75SLisandro Dalcin   for (i = 0; i < snes->max_its; i++) {
36837e1895aSJed Brown 
36937e1895aSJed Brown     /* Call general purpose update function */
37037e1895aSJed Brown     if (snes->ops->update) {
37137e1895aSJed Brown       ierr = (*snes->ops->update)(snes,snes->iter);CHKERRQ(ierr);
37237e1895aSJed Brown     }
373bf80fd75SLisandro Dalcin 
374bf80fd75SLisandro Dalcin     if (i == 0 && snes->jacobian) {
375bf80fd75SLisandro Dalcin       /* This method does not require a Jacobian, but it is usually preconditioned by PBJacobi */
376bf80fd75SLisandro Dalcin       ierr = SNESComputeJacobian(snes,snes->vec_sol,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
377bf80fd75SLisandro Dalcin       SNESCheckJacobianDomainerror(snes);
378bf80fd75SLisandro Dalcin       ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
379bf80fd75SLisandro Dalcin     }
380bf80fd75SLisandro Dalcin 
3819ce202e0SLisandro Dalcin     ierr = SNESMSStep_Step(snes,X,F);CHKERRQ(ierr);
38237e1895aSJed Brown 
383bf80fd75SLisandro Dalcin     if (i < snes->max_its-1 || ms->norms) {
38437e1895aSJed Brown       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
38537e1895aSJed Brown     }
38637e1895aSJed Brown 
387bf80fd75SLisandro Dalcin     ierr = SNESMSStep_Norms(snes,i+1,F);CHKERRQ(ierr);
38837e1895aSJed Brown     if (snes->reason) PetscFunctionReturn(0);
38937e1895aSJed Brown   }
39037e1895aSJed Brown   if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
39137e1895aSJed Brown   PetscFunctionReturn(0);
39237e1895aSJed Brown }
39337e1895aSJed Brown 
39437e1895aSJed Brown static PetscErrorCode SNESSetUp_MS(SNES snes)
39537e1895aSJed Brown {
3969ce202e0SLisandro Dalcin   SNES_MS        *ms   = (SNES_MS*)snes->data;
3979ce202e0SLisandro Dalcin   SNESMSTableau  tab   = ms->tableau;
3989ce202e0SLisandro Dalcin   PetscInt       nwork = tab->nregisters;
39937e1895aSJed Brown   PetscErrorCode ierr;
40037e1895aSJed Brown 
40137e1895aSJed Brown   PetscFunctionBegin;
4029ce202e0SLisandro Dalcin   ierr = SNESSetWorkVecs(snes,nwork);CHKERRQ(ierr);
4036cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
40437e1895aSJed Brown   PetscFunctionReturn(0);
40537e1895aSJed Brown }
40637e1895aSJed Brown 
40737e1895aSJed Brown static PetscErrorCode SNESReset_MS(SNES snes)
40837e1895aSJed Brown {
40937e1895aSJed Brown   PetscFunctionBegin;
41037e1895aSJed Brown   PetscFunctionReturn(0);
41137e1895aSJed Brown }
41237e1895aSJed Brown 
41337e1895aSJed Brown static PetscErrorCode SNESDestroy_MS(SNES snes)
41437e1895aSJed Brown {
41537e1895aSJed Brown   PetscErrorCode ierr;
41637e1895aSJed Brown 
41737e1895aSJed Brown   PetscFunctionBegin;
418bf80fd75SLisandro Dalcin   ierr = SNESReset_MS(snes);CHKERRQ(ierr);
41937e1895aSJed Brown   ierr = PetscFree(snes->data);CHKERRQ(ierr);
42057715debSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetType_C",NULL);CHKERRQ(ierr);
421bf80fd75SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetType_C",NULL);CHKERRQ(ierr);
42257715debSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetDamping_C",NULL);CHKERRQ(ierr);
42357715debSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetDamping_C",NULL);CHKERRQ(ierr);
42437e1895aSJed Brown   PetscFunctionReturn(0);
42537e1895aSJed Brown }
42637e1895aSJed Brown 
42737e1895aSJed Brown static PetscErrorCode SNESView_MS(SNES snes,PetscViewer viewer)
42837e1895aSJed Brown {
42937e1895aSJed Brown   PetscBool      iascii;
43037e1895aSJed Brown   PetscErrorCode ierr;
43137e1895aSJed Brown   SNES_MS        *ms = (SNES_MS*)snes->data;
43257715debSLisandro Dalcin   SNESMSTableau  tab = ms->tableau;
43337e1895aSJed Brown 
43437e1895aSJed Brown   PetscFunctionBegin;
435251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
43637e1895aSJed Brown   if (iascii) {
43757715debSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  multi-stage method type: %s\n",tab->name);CHKERRQ(ierr);
43837e1895aSJed Brown   }
43937e1895aSJed Brown   PetscFunctionReturn(0);
44037e1895aSJed Brown }
44137e1895aSJed Brown 
4424416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_MS(PetscOptionItems *PetscOptionsObject,SNES snes)
44337e1895aSJed Brown {
444725b86d8SJed Brown   SNES_MS        *ms = (SNES_MS*)snes->data;
44537e1895aSJed Brown   PetscErrorCode ierr;
44637e1895aSJed Brown 
44737e1895aSJed Brown   PetscFunctionBegin;
448e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES MS options");CHKERRQ(ierr);
44937e1895aSJed Brown   {
45037e1895aSJed Brown     SNESMSTableauLink link;
45137e1895aSJed Brown     PetscInt          count,choice;
45237e1895aSJed Brown     PetscBool         flg;
45337e1895aSJed Brown     const char        **namelist;
45457715debSLisandro Dalcin     SNESMSType        mstype;
45557715debSLisandro Dalcin     PetscReal         damping;
45637e1895aSJed Brown 
45757715debSLisandro Dalcin     ierr = SNESMSGetType(snes,&mstype);CHKERRQ(ierr);
45837e1895aSJed Brown     for (link=SNESMSTableauList,count=0; link; link=link->next,count++) ;
459f489ac74SBarry Smith     ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr);
46037e1895aSJed Brown     for (link=SNESMSTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
46137e1895aSJed Brown     ierr = PetscOptionsEList("-snes_ms_type","Multistage smoother type","SNESMSSetType",(const char*const*)namelist,count,mstype,&choice,&flg);CHKERRQ(ierr);
46257715debSLisandro Dalcin     if (flg) {ierr = SNESMSSetType(snes,namelist[choice]);CHKERRQ(ierr);}
46337e1895aSJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
46457715debSLisandro Dalcin     ierr = SNESMSGetDamping(snes,&damping);CHKERRQ(ierr);
46557715debSLisandro Dalcin     ierr = PetscOptionsReal("-snes_ms_damping","Damping for multistage method","SNESMSSetDamping",damping,&damping,&flg);CHKERRQ(ierr);
46657715debSLisandro Dalcin     if (flg) {ierr = SNESMSSetDamping(snes,damping);CHKERRQ(ierr);}
4670298fd71SBarry Smith     ierr = PetscOptionsBool("-snes_ms_norms","Compute norms for monitoring","none",ms->norms,&ms->norms,NULL);CHKERRQ(ierr);
46837e1895aSJed Brown   }
46937e1895aSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
47037e1895aSJed Brown   PetscFunctionReturn(0);
47137e1895aSJed Brown }
47237e1895aSJed Brown 
47357715debSLisandro Dalcin static PetscErrorCode SNESMSGetType_MS(SNES snes,SNESMSType *mstype)
47437e1895aSJed Brown {
47557715debSLisandro Dalcin   SNES_MS        *ms = (SNES_MS*)snes->data;
47657715debSLisandro Dalcin   SNESMSTableau  tab = ms->tableau;
47757715debSLisandro Dalcin 
47857715debSLisandro Dalcin   PetscFunctionBegin;
47957715debSLisandro Dalcin   *mstype = tab->name;
48057715debSLisandro Dalcin   PetscFunctionReturn(0);
48157715debSLisandro Dalcin }
48257715debSLisandro Dalcin 
48357715debSLisandro Dalcin static PetscErrorCode SNESMSSetType_MS(SNES snes,SNESMSType mstype)
48457715debSLisandro Dalcin {
48537e1895aSJed Brown   SNES_MS           *ms = (SNES_MS*)snes->data;
48637e1895aSJed Brown   SNESMSTableauLink link;
48737e1895aSJed Brown   PetscBool         match;
48857715debSLisandro Dalcin   PetscErrorCode    ierr;
48937e1895aSJed Brown 
49037e1895aSJed Brown   PetscFunctionBegin;
49137e1895aSJed Brown   if (ms->tableau) {
49237e1895aSJed Brown     ierr = PetscStrcmp(ms->tableau->name,mstype,&match);CHKERRQ(ierr);
49337e1895aSJed Brown     if (match) PetscFunctionReturn(0);
49437e1895aSJed Brown   }
49537e1895aSJed Brown   for (link = SNESMSTableauList; link; link=link->next) {
49637e1895aSJed Brown     ierr = PetscStrcmp(link->tab.name,mstype,&match);CHKERRQ(ierr);
49737e1895aSJed Brown     if (match) {
49857715debSLisandro Dalcin       if (snes->setupcalled)  {ierr = SNESReset_MS(snes);CHKERRQ(ierr);}
49937e1895aSJed Brown       ms->tableau = &link->tab;
50057715debSLisandro Dalcin       if (snes->setupcalled)  {ierr = SNESSetUp_MS(snes);CHKERRQ(ierr);}
50137e1895aSJed Brown       PetscFunctionReturn(0);
50237e1895aSJed Brown     }
50337e1895aSJed Brown   }
504ce94432eSBarry Smith   SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mstype);
50537e1895aSJed Brown }
50637e1895aSJed Brown 
50737e1895aSJed Brown /*@C
50857715debSLisandro Dalcin   SNESMSGetType - Get the type of multistage smoother
50957715debSLisandro Dalcin 
51057715debSLisandro Dalcin   Not collective
51157715debSLisandro Dalcin 
51257715debSLisandro Dalcin   Input Parameter:
51357715debSLisandro Dalcin .  snes - nonlinear solver context
51457715debSLisandro Dalcin 
51557715debSLisandro Dalcin   Output Parameter:
51657715debSLisandro Dalcin .  mstype - type of multistage method
51757715debSLisandro Dalcin 
51857715debSLisandro Dalcin   Level: beginner
51957715debSLisandro Dalcin 
52057715debSLisandro Dalcin .seealso: SNESMSSetType(), SNESMSType, SNESMS
52157715debSLisandro Dalcin @*/
52257715debSLisandro Dalcin PetscErrorCode SNESMSGetType(SNES snes,SNESMSType *mstype)
52357715debSLisandro Dalcin {
52457715debSLisandro Dalcin   PetscErrorCode ierr;
52557715debSLisandro Dalcin 
52657715debSLisandro Dalcin   PetscFunctionBegin;
52757715debSLisandro Dalcin   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
52857715debSLisandro Dalcin   PetscValidPointer(mstype,2);
52957715debSLisandro Dalcin   ierr = PetscUseMethod(snes,"SNESMSGetType_C",(SNES,SNESMSType*),(snes,mstype));CHKERRQ(ierr);
53057715debSLisandro Dalcin   PetscFunctionReturn(0);
53157715debSLisandro Dalcin }
53257715debSLisandro Dalcin 
53357715debSLisandro Dalcin /*@C
53437e1895aSJed Brown   SNESMSSetType - Set the type of multistage smoother
53537e1895aSJed Brown 
53637e1895aSJed Brown   Logically collective
53737e1895aSJed Brown 
538*d8d19677SJose E. Roman   Input Parameters:
53937e1895aSJed Brown +  snes - nonlinear solver context
54037e1895aSJed Brown -  mstype - type of multistage method
54137e1895aSJed Brown 
54237e1895aSJed Brown   Level: beginner
54337e1895aSJed Brown 
54457715debSLisandro Dalcin .seealso: SNESMSGetType(), SNESMSType, SNESMS
54537e1895aSJed Brown @*/
54657715debSLisandro Dalcin PetscErrorCode SNESMSSetType(SNES snes,SNESMSType mstype)
54737e1895aSJed Brown {
54837e1895aSJed Brown   PetscErrorCode ierr;
54937e1895aSJed Brown 
55037e1895aSJed Brown   PetscFunctionBegin;
55137e1895aSJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
55257715debSLisandro Dalcin   PetscValidPointer(mstype,2);
55357715debSLisandro Dalcin   ierr = PetscTryMethod(snes,"SNESMSSetType_C",(SNES,SNESMSType),(snes,mstype));CHKERRQ(ierr);
55457715debSLisandro Dalcin   PetscFunctionReturn(0);
55557715debSLisandro Dalcin }
55657715debSLisandro Dalcin 
55757715debSLisandro Dalcin static PetscErrorCode SNESMSGetDamping_MS(SNES snes,PetscReal *damping)
55857715debSLisandro Dalcin {
55957715debSLisandro Dalcin   SNES_MS        *ms = (SNES_MS*)snes->data;
56057715debSLisandro Dalcin 
56157715debSLisandro Dalcin   PetscFunctionBegin;
56257715debSLisandro Dalcin   *damping = ms->damping;
56357715debSLisandro Dalcin   PetscFunctionReturn(0);
56457715debSLisandro Dalcin }
56557715debSLisandro Dalcin 
56657715debSLisandro Dalcin static PetscErrorCode SNESMSSetDamping_MS(SNES snes,PetscReal damping)
56757715debSLisandro Dalcin {
56857715debSLisandro Dalcin   SNES_MS           *ms = (SNES_MS*)snes->data;
56957715debSLisandro Dalcin 
57057715debSLisandro Dalcin   PetscFunctionBegin;
57157715debSLisandro Dalcin   ms->damping = damping;
57257715debSLisandro Dalcin   PetscFunctionReturn(0);
57357715debSLisandro Dalcin }
57457715debSLisandro Dalcin 
57557715debSLisandro Dalcin /*@C
57657715debSLisandro Dalcin   SNESMSGetDamping - Get the damping parameter
57757715debSLisandro Dalcin 
57857715debSLisandro Dalcin   Not collective
57957715debSLisandro Dalcin 
58057715debSLisandro Dalcin   Input Parameter:
58157715debSLisandro Dalcin .  snes - nonlinear solver context
58257715debSLisandro Dalcin 
58357715debSLisandro Dalcin   Output Parameter:
58457715debSLisandro Dalcin .  damping - damping parameter
58557715debSLisandro Dalcin 
58657715debSLisandro Dalcin   Level: advanced
58757715debSLisandro Dalcin 
58857715debSLisandro Dalcin .seealso: SNESMSSetDamping(), SNESMS
58957715debSLisandro Dalcin @*/
59057715debSLisandro Dalcin PetscErrorCode SNESMSGetDamping(SNES snes,PetscReal *damping)
59157715debSLisandro Dalcin {
59257715debSLisandro Dalcin   PetscErrorCode ierr;
59357715debSLisandro Dalcin 
59457715debSLisandro Dalcin   PetscFunctionBegin;
59557715debSLisandro Dalcin   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
59657715debSLisandro Dalcin   PetscValidPointer(damping,2);
59757715debSLisandro Dalcin   ierr = PetscUseMethod(snes,"SNESMSGetDamping_C",(SNES,PetscReal*),(snes,damping));CHKERRQ(ierr);
59857715debSLisandro Dalcin   PetscFunctionReturn(0);
59957715debSLisandro Dalcin }
60057715debSLisandro Dalcin 
60157715debSLisandro Dalcin /*@C
60257715debSLisandro Dalcin   SNESMSSetDamping - Set the damping parameter
60357715debSLisandro Dalcin 
60457715debSLisandro Dalcin   Logically collective
60557715debSLisandro Dalcin 
606*d8d19677SJose E. Roman   Input Parameters:
60757715debSLisandro Dalcin +  snes - nonlinear solver context
60857715debSLisandro Dalcin -  damping - damping parameter
60957715debSLisandro Dalcin 
61057715debSLisandro Dalcin   Level: advanced
61157715debSLisandro Dalcin 
61257715debSLisandro Dalcin .seealso: SNESMSGetDamping(), SNESMS
61357715debSLisandro Dalcin @*/
61457715debSLisandro Dalcin PetscErrorCode SNESMSSetDamping(SNES snes,PetscReal damping)
61557715debSLisandro Dalcin {
61657715debSLisandro Dalcin   PetscErrorCode ierr;
61757715debSLisandro Dalcin 
61857715debSLisandro Dalcin   PetscFunctionBegin;
61957715debSLisandro Dalcin   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
62057715debSLisandro Dalcin   PetscValidLogicalCollectiveReal(snes,damping,2);
62157715debSLisandro Dalcin   ierr = PetscTryMethod(snes,"SNESMSSetDamping_C",(SNES,PetscReal),(snes,damping));CHKERRQ(ierr);
62237e1895aSJed Brown   PetscFunctionReturn(0);
62337e1895aSJed Brown }
62437e1895aSJed Brown 
62537e1895aSJed Brown /* -------------------------------------------------------------------------- */
62637e1895aSJed Brown /*MC
62737e1895aSJed Brown       SNESMS - multi-stage smoothers
62837e1895aSJed Brown 
62937e1895aSJed Brown       Options Database:
63037e1895aSJed Brown 
63137e1895aSJed Brown +     -snes_ms_type - type of multi-stage smoother
63237e1895aSJed Brown -     -snes_ms_damping - damping for multi-stage method
63337e1895aSJed Brown 
63437e1895aSJed Brown       Notes:
63537e1895aSJed Brown       These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for
6366c9de887SHong Zhang       FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebyshev).
63737e1895aSJed Brown 
63837e1895aSJed Brown       Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds.
63937e1895aSJed Brown 
64037e1895aSJed Brown       The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with SNESMSRegister().
64137e1895aSJed Brown 
64237e1895aSJed Brown       References:
64357715debSLisandro Dalcin +   1. -   Ketcheson (2010) Runge Kutta methods with minimum storage implementations (https://doi.org/10.1016/j.jcp.2009.11.006).
64457715debSLisandro Dalcin .   2. -   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).
64557715debSLisandro Dalcin .   3. -   Pierce and Giles (1997) Preconditioned multigrid methods for compressible flow calculations on stretched meshes (https://doi.org/10.1006/jcph.1997.5772).
64657715debSLisandro Dalcin -   4. -   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).
64737e1895aSJed Brown 
64837e1895aSJed Brown       Level: beginner
64937e1895aSJed Brown 
6506c9de887SHong Zhang .seealso:  SNESCreate(), SNES, SNESSetType(), SNESMS, SNESFAS, KSPCHEBYSHEV
65137e1895aSJed Brown 
65237e1895aSJed Brown M*/
6538cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_MS(SNES snes)
65437e1895aSJed Brown {
65537e1895aSJed Brown   PetscErrorCode ierr;
65637e1895aSJed Brown   SNES_MS        *ms;
65737e1895aSJed Brown 
65837e1895aSJed Brown   PetscFunctionBegin;
659607a6623SBarry Smith   ierr = SNESMSInitializePackage();CHKERRQ(ierr);
66037e1895aSJed Brown 
66137e1895aSJed Brown   snes->ops->setup          = SNESSetUp_MS;
66237e1895aSJed Brown   snes->ops->solve          = SNESSolve_MS;
66337e1895aSJed Brown   snes->ops->destroy        = SNESDestroy_MS;
66437e1895aSJed Brown   snes->ops->setfromoptions = SNESSetFromOptions_MS;
66537e1895aSJed Brown   snes->ops->view           = SNESView_MS;
66637e1895aSJed Brown   snes->ops->reset          = SNESReset_MS;
66737e1895aSJed Brown 
668efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
66937e1895aSJed Brown   snes->usesksp = PETSC_TRUE;
67037e1895aSJed Brown 
6714fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
6724fc747eaSLawrence Mitchell 
673b00a9115SJed Brown   ierr        = PetscNewLog(snes,&ms);CHKERRQ(ierr);
67437e1895aSJed Brown   snes->data  = (void*)ms;
67537e1895aSJed Brown   ms->damping = 0.9;
67637e1895aSJed Brown   ms->norms   = PETSC_FALSE;
67737e1895aSJed Brown 
67857715debSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetType_C",SNESMSGetType_MS);CHKERRQ(ierr);
679bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetType_C",SNESMSSetType_MS);CHKERRQ(ierr);
68057715debSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetDamping_C",SNESMSGetDamping_MS);CHKERRQ(ierr);
68157715debSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetDamping_C",SNESMSSetDamping_MS);CHKERRQ(ierr);
68257715debSLisandro Dalcin 
68357715debSLisandro Dalcin   ierr = SNESMSSetType(snes,SNESMSDefault);CHKERRQ(ierr);
68437e1895aSJed Brown   PetscFunctionReturn(0);
68537e1895aSJed Brown }
686