xref: /petsc/src/snes/impls/ms/ms.c (revision efd4aadf157bf1ba2d80c2be092fcf4247860003)
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 .keywords: SNES, SNESMS, register, all
3937e1895aSJed Brown 
4037e1895aSJed Brown .seealso:  SNESMSRegisterDestroy()
4137e1895aSJed Brown @*/
4237e1895aSJed Brown PetscErrorCode SNESMSRegisterAll(void)
4337e1895aSJed Brown {
4437e1895aSJed Brown   PetscErrorCode ierr;
4537e1895aSJed Brown 
4637e1895aSJed Brown   PetscFunctionBegin;
4737e1895aSJed Brown   if (SNESMSRegisterAllCalled) PetscFunctionReturn(0);
4837e1895aSJed Brown   SNESMSRegisterAllCalled = PETSC_TRUE;
4937e1895aSJed Brown 
5037e1895aSJed Brown   {
5137e1895aSJed Brown     PetscReal gamma[3][1] = {{1.0},{0.0},{0.0}};
5237e1895aSJed Brown     PetscReal delta[1]    = {0.0};
5337e1895aSJed Brown     PetscReal betasub[1]  = {1.0};
5437e1895aSJed Brown     ierr = SNESMSRegister(SNESMSEULER,1,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
5537e1895aSJed Brown   }
5637e1895aSJed Brown 
5737e1895aSJed Brown   {
5837e1895aSJed Brown     PetscReal gamma[3][6] = {
5937e1895aSJed Brown       {0.0000000000000000E+00,  -7.0304722367110606E-01, -1.9836719667506464E-01, -1.6023843981863788E+00, 9.4483822882855284E-02,  -1.4204296130641869E-01},
6037e1895aSJed Brown       {1.0000000000000000E+00,  1.1111025767083920E+00,  5.6150921583923230E-01,  7.4151723494934041E-01,  3.1714538168600587E-01,  4.6479276238548706E-01},
6137e1895aSJed Brown       {0.0000000000000000E+00,  0.0000000000000000E+00,  0.0000000000000000E+00,  6.7968174970583317E-01,  -4.1755042846051737E-03, -1.9115668129923846E-01}
6237e1895aSJed Brown     };
6337e1895aSJed Brown     PetscReal delta[6]   = {1.0000000000000000E+00, 5.3275427433201750E-01, 6.0143544663985238E-01, 4.5874077053842177E-01, 2.7544386906104651E-01, 0.0000000000000000E+00};
6437e1895aSJed Brown     PetscReal betasub[6] = {8.4753115429481929E-01, 7.4018896368655618E-01, 6.5963574086583309E-03, 4.6747795645517759E-01, 1.3314545813643919E-01, 5.3260800028018784E-01};
6537e1895aSJed Brown     ierr = SNESMSRegister(SNESMSM62,6,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
6637e1895aSJed Brown   }
67a97cb6bcSJed Brown   {
68a97cb6bcSJed Brown     PetscReal gamma[3][4] = {{0,0,0,0},{0,0,0,0},{1,1,1,1}};
69a97cb6bcSJed Brown     PetscReal delta[4]    = {0,0,0,0};
70a97cb6bcSJed Brown     PetscReal betasub[4]  = {0.25, 0.5, 0.55, 1.0};
71a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSJAMESON83,4,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
72a97cb6bcSJed Brown   }
73a97cb6bcSJed Brown   {                             /* Van Leer, Tai, and Powell (1989) 2 stage, order 1 */
74a97cb6bcSJed Brown     PetscReal gamma[3][2] = {{0,0},{0,0},{1,1}};
75a97cb6bcSJed Brown     PetscReal delta[2]    = {0,0};
76a97cb6bcSJed Brown     PetscReal betasub[2]  = {0.3333,1.0};
77a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSVLTP21,2,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
78a97cb6bcSJed Brown   }
79a97cb6bcSJed Brown   {                             /* Van Leer, Tai, and Powell (1989) 3 stage, order 1 */
80a97cb6bcSJed Brown     PetscReal gamma[3][3] = {{0,0,0},{0,0,0},{1,1,1}};
81a97cb6bcSJed Brown     PetscReal delta[3]    = {0,0,0};
82a97cb6bcSJed Brown     PetscReal betasub[3]  = {0.1481,0.4000,1.0};
83a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSVLTP31,3,3,1.5,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
84a97cb6bcSJed Brown   }
85a97cb6bcSJed Brown   {                             /* Van Leer, Tai, and Powell (1989) 4 stage, order 1 */
86a97cb6bcSJed Brown     PetscReal gamma[3][4] = {{0,0,0,0},{0,0,0,0},{1,1,1,1}};
87a97cb6bcSJed Brown     PetscReal delta[4]    = {0,0,0,0};
88a97cb6bcSJed Brown     PetscReal betasub[4]  = {0.0833,0.2069,0.4265,1.0};
89a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSVLTP41,4,3,2.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
90a97cb6bcSJed Brown   }
91a97cb6bcSJed Brown   {                             /* Van Leer, Tai, and Powell (1989) 5 stage, order 1 */
92a97cb6bcSJed Brown     PetscReal gamma[3][5] = {{0,0,0,0,0},{0,0,0,0,0},{1,1,1,1,1}};
93a97cb6bcSJed Brown     PetscReal delta[5]    = {0,0,0,0,0};
94a97cb6bcSJed Brown     PetscReal betasub[5]  = {0.0533,0.1263,0.2375,0.4414,1.0};
95a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSVLTP51,5,3,2.5,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
96a97cb6bcSJed Brown   }
97a97cb6bcSJed Brown   {                             /* Van Leer, Tai, and Powell (1989) 6 stage, order 1 */
98a97cb6bcSJed Brown     PetscReal gamma[3][6] = {{0,0,0,0,0,0},{0,0,0,0,0,0},{1,1,1,1,1,1}};
99a97cb6bcSJed Brown     PetscReal delta[6]    = {0,0,0,0,0,0};
100a97cb6bcSJed Brown     PetscReal betasub[6]  = {0.0370,0.0851,0.1521,0.2562,0.4512,1.0};
101a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSVLTP61,6,3,3.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
102a97cb6bcSJed Brown   }
10337e1895aSJed Brown   PetscFunctionReturn(0);
10437e1895aSJed Brown }
10537e1895aSJed Brown 
10637e1895aSJed Brown /*@C
10737e1895aSJed Brown    SNESMSRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
10837e1895aSJed Brown 
10937e1895aSJed Brown    Not Collective
11037e1895aSJed Brown 
11137e1895aSJed Brown    Level: advanced
11237e1895aSJed Brown 
11337e1895aSJed Brown .keywords: TSRosW, register, destroy
1141c84c290SBarry Smith .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegister()
11537e1895aSJed Brown @*/
11637e1895aSJed Brown PetscErrorCode SNESMSRegisterDestroy(void)
11737e1895aSJed Brown {
11837e1895aSJed Brown   PetscErrorCode    ierr;
11937e1895aSJed Brown   SNESMSTableauLink link;
12037e1895aSJed Brown 
12137e1895aSJed Brown   PetscFunctionBegin;
12237e1895aSJed Brown   while ((link = SNESMSTableauList)) {
12337e1895aSJed Brown     SNESMSTableau t = &link->tab;
12437e1895aSJed Brown     SNESMSTableauList = link->next;
1251aa26658SKarl Rupp 
12637e1895aSJed Brown     ierr = PetscFree3(t->gamma,t->delta,t->betasub);CHKERRQ(ierr);
12737e1895aSJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
12837e1895aSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
12937e1895aSJed Brown   }
13037e1895aSJed Brown   SNESMSRegisterAllCalled = PETSC_FALSE;
13137e1895aSJed Brown   PetscFunctionReturn(0);
13237e1895aSJed Brown }
13337e1895aSJed Brown 
13437e1895aSJed Brown /*@C
13537e1895aSJed Brown   SNESMSInitializePackage - This function initializes everything in the SNESMS package. It is called
13637e1895aSJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to SNESCreate_MS()
13737e1895aSJed Brown   when using static libraries.
13837e1895aSJed Brown 
13937e1895aSJed Brown   Level: developer
14037e1895aSJed Brown 
14137e1895aSJed Brown .keywords: SNES, SNESMS, initialize, package
14237e1895aSJed Brown .seealso: PetscInitialize()
14337e1895aSJed Brown @*/
144607a6623SBarry Smith PetscErrorCode SNESMSInitializePackage(void)
14537e1895aSJed Brown {
14637e1895aSJed Brown   PetscErrorCode ierr;
14737e1895aSJed Brown 
14837e1895aSJed Brown   PetscFunctionBegin;
14937e1895aSJed Brown   if (SNESMSPackageInitialized) PetscFunctionReturn(0);
15037e1895aSJed Brown   SNESMSPackageInitialized = PETSC_TRUE;
1511aa26658SKarl Rupp 
15237e1895aSJed Brown   ierr = SNESMSRegisterAll();CHKERRQ(ierr);
15337e1895aSJed Brown   ierr = PetscRegisterFinalize(SNESMSFinalizePackage);CHKERRQ(ierr);
15437e1895aSJed Brown   PetscFunctionReturn(0);
15537e1895aSJed Brown }
15637e1895aSJed Brown 
15737e1895aSJed Brown /*@C
15837e1895aSJed Brown   SNESMSFinalizePackage - This function destroys everything in the SNESMS package. It is
15937e1895aSJed Brown   called from PetscFinalize().
16037e1895aSJed Brown 
16137e1895aSJed Brown   Level: developer
16237e1895aSJed Brown 
16337e1895aSJed Brown .keywords: Petsc, destroy, package
16437e1895aSJed Brown .seealso: PetscFinalize()
16537e1895aSJed Brown @*/
16637e1895aSJed Brown PetscErrorCode SNESMSFinalizePackage(void)
16737e1895aSJed Brown {
16837e1895aSJed Brown   PetscErrorCode ierr;
16937e1895aSJed Brown 
17037e1895aSJed Brown   PetscFunctionBegin;
17137e1895aSJed Brown   SNESMSPackageInitialized = PETSC_FALSE;
1721aa26658SKarl Rupp 
17337e1895aSJed Brown   ierr = SNESMSRegisterDestroy();CHKERRQ(ierr);
17437e1895aSJed Brown   PetscFunctionReturn(0);
17537e1895aSJed Brown }
17637e1895aSJed Brown 
17737e1895aSJed Brown /*@C
17837e1895aSJed Brown    SNESMSRegister - register a multistage scheme
17937e1895aSJed Brown 
18037e1895aSJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
18137e1895aSJed Brown 
18237e1895aSJed Brown    Input Parameters:
18337e1895aSJed Brown +  name - identifier for method
18437e1895aSJed Brown .  nstages - number of stages
18537e1895aSJed Brown .  nregisters - number of registers used by low-storage implementation
18637e1895aSJed Brown .  gamma - coefficients, see Ketcheson's paper
18737e1895aSJed Brown .  delta - coefficients, see Ketcheson's paper
18837e1895aSJed Brown -  betasub - subdiagonal of Shu-Osher form
18937e1895aSJed Brown 
19037e1895aSJed Brown    Notes:
19137e1895aSJed Brown    The notation is described in Ketcheson (2010) Runge-Kutta methods with minimum storage implementations.
19237e1895aSJed Brown 
19337e1895aSJed Brown    Level: advanced
19437e1895aSJed Brown 
19537e1895aSJed Brown .keywords: SNES, register
19637e1895aSJed Brown 
19737e1895aSJed Brown .seealso: SNESMS
19837e1895aSJed Brown @*/
19919fd82e9SBarry Smith PetscErrorCode SNESMSRegister(SNESMSType name,PetscInt nstages,PetscInt nregisters,PetscReal stability,const PetscReal gamma[],const PetscReal delta[],const PetscReal betasub[])
20037e1895aSJed Brown {
20137e1895aSJed Brown   PetscErrorCode    ierr;
20237e1895aSJed Brown   SNESMSTableauLink link;
20337e1895aSJed Brown   SNESMSTableau     t;
20437e1895aSJed Brown 
20537e1895aSJed Brown   PetscFunctionBegin;
20637e1895aSJed Brown   PetscValidCharPointer(name,1);
20737e1895aSJed Brown   if (nstages < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have at least one stage");
20837e1895aSJed Brown   if (nregisters != 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only support for methods written in 3-register form");
20937e1895aSJed Brown   PetscValidPointer(gamma,4);
21037e1895aSJed Brown   PetscValidPointer(delta,5);
21137e1895aSJed Brown   PetscValidPointer(betasub,6);
21237e1895aSJed Brown 
213854ce69bSBarry Smith   ierr          = PetscNew(&link);CHKERRQ(ierr);
21437e1895aSJed Brown   t             = &link->tab;
21537e1895aSJed Brown   ierr          = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
21637e1895aSJed Brown   t->nstages    = nstages;
21737e1895aSJed Brown   t->nregisters = nregisters;
21837e1895aSJed Brown   t->stability  = stability;
2191aa26658SKarl Rupp 
220dcca6d9dSJed Brown   ierr = PetscMalloc3(nstages*nregisters,&t->gamma,nstages,&t->delta,nstages,&t->betasub);CHKERRQ(ierr);
22137e1895aSJed Brown   ierr = PetscMemcpy(t->gamma,gamma,nstages*nregisters*sizeof(PetscReal));CHKERRQ(ierr);
22237e1895aSJed Brown   ierr = PetscMemcpy(t->delta,delta,nstages*sizeof(PetscReal));CHKERRQ(ierr);
22337e1895aSJed Brown   ierr = PetscMemcpy(t->betasub,betasub,nstages*sizeof(PetscReal));CHKERRQ(ierr);
22437e1895aSJed Brown 
22537e1895aSJed Brown   link->next        = SNESMSTableauList;
22637e1895aSJed Brown   SNESMSTableauList = link;
22737e1895aSJed Brown   PetscFunctionReturn(0);
22837e1895aSJed Brown }
22937e1895aSJed Brown 
23037e1895aSJed Brown /*
23137e1895aSJed Brown   X - initial state, updated in-place.
23237e1895aSJed Brown   F - residual, computed at the initial X on input
23337e1895aSJed Brown */
23437e1895aSJed Brown static PetscErrorCode SNESMSStep_3Sstar(SNES snes,Vec X,Vec F)
23537e1895aSJed Brown {
23637e1895aSJed Brown   PetscErrorCode  ierr;
23737e1895aSJed Brown   SNES_MS         *ms    = (SNES_MS*)snes->data;
23837e1895aSJed Brown   SNESMSTableau   t      = ms->tableau;
23937e1895aSJed Brown   const PetscReal *gamma = t->gamma,*delta = t->delta,*betasub = t->betasub;
24037e1895aSJed Brown   Vec             S1,S2,S3,Y;
24137e1895aSJed Brown   PetscInt        i,nstages = t->nstages;;
24237e1895aSJed Brown 
24337e1895aSJed Brown 
24437e1895aSJed Brown   PetscFunctionBegin;
24537e1895aSJed Brown   Y    = snes->work[0];
24637e1895aSJed Brown   S1   = X;
24737e1895aSJed Brown   S2   = snes->work[1];
24837e1895aSJed Brown   S3   = snes->work[2];
24937e1895aSJed Brown   ierr = VecZeroEntries(S2);CHKERRQ(ierr);
25037e1895aSJed Brown   ierr = VecCopy(X,S3);CHKERRQ(ierr);
25137e1895aSJed Brown   for (i=0; i<nstages; i++) {
252da80777bSKarl Rupp     Vec         Ss[4];
253da80777bSKarl Rupp     PetscScalar scoeff[4];
254da80777bSKarl Rupp 
255da80777bSKarl Rupp     Ss[0] = S1; Ss[1] = S2; Ss[2] = S3; Ss[3] = Y;
2561aa26658SKarl Rupp 
257da80777bSKarl Rupp     scoeff[0] = gamma[0*nstages+i]-(PetscReal)1.0;
258da80777bSKarl Rupp     scoeff[1] = gamma[1*nstages+i];
259da80777bSKarl Rupp     scoeff[2] = gamma[2*nstages+i];
260da80777bSKarl Rupp     scoeff[3] = -betasub[i]*ms->damping;
2611aa26658SKarl Rupp 
26237e1895aSJed Brown     ierr = VecAXPY(S2,delta[i],S1);CHKERRQ(ierr);
26337e1895aSJed Brown     if (i>0) {
26437e1895aSJed Brown       ierr = SNESComputeFunction(snes,S1,F);CHKERRQ(ierr);
26537e1895aSJed Brown     }
266d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
26737e1895aSJed Brown     ierr = VecMAXPY(S1,4,scoeff,Ss);CHKERRQ(ierr);
26837e1895aSJed Brown   }
26937e1895aSJed Brown   PetscFunctionReturn(0);
27037e1895aSJed Brown }
27137e1895aSJed Brown 
27237e1895aSJed Brown static PetscErrorCode SNESSolve_MS(SNES snes)
27337e1895aSJed Brown {
27437e1895aSJed Brown   SNES_MS        *ms = (SNES_MS*)snes->data;
27537e1895aSJed Brown   Vec            X   = snes->vec_sol,F = snes->vec_func;
27637e1895aSJed Brown   PetscReal      fnorm;
27737e1895aSJed Brown   PetscInt       i;
27837e1895aSJed Brown   PetscErrorCode ierr;
27937e1895aSJed Brown 
28037e1895aSJed Brown   PetscFunctionBegin;
2816c4ed002SBarry 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);
282c579b300SPatrick Farrell 
283fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
28437e1895aSJed Brown   snes->reason = SNES_CONVERGED_ITERATING;
285e04113cfSBarry Smith   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
28637e1895aSJed Brown   snes->iter   = 0;
28737e1895aSJed Brown   snes->norm   = 0.;
288e04113cfSBarry Smith   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
289e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
29037e1895aSJed Brown     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
2911aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
2921aa26658SKarl Rupp 
29337e1895aSJed Brown   if (snes->jacobian) {         /* This method does not require a Jacobian, but it is usually preconditioned by PBJacobi */
294d1e9a80fSBarry Smith     ierr = SNESComputeJacobian(snes,snes->vec_sol,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
29537e1895aSJed Brown   }
29637e1895aSJed Brown   if (ms->norms) {
29737e1895aSJed Brown     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
298422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
29937e1895aSJed Brown     /* Monitor convergence */
300e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
30137e1895aSJed Brown     snes->iter = 0;
30237e1895aSJed Brown     snes->norm = fnorm;
303e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
304a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
30537e1895aSJed Brown     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
30637e1895aSJed Brown 
30737e1895aSJed Brown     /* Test for convergence */
30837e1895aSJed Brown     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
30937e1895aSJed Brown     if (snes->reason) PetscFunctionReturn(0);
31037e1895aSJed Brown   }
31137e1895aSJed Brown 
31237e1895aSJed Brown   /* Call general purpose update function */
31337e1895aSJed Brown   if (snes->ops->update) {
31437e1895aSJed Brown     ierr = (*snes->ops->update)(snes,snes->iter);CHKERRQ(ierr);
31537e1895aSJed Brown   }
31637e1895aSJed Brown   for (i = 0; i < snes->max_its; i++) {
31737e1895aSJed Brown     ierr = SNESMSStep_3Sstar(snes,X,F);CHKERRQ(ierr);
31837e1895aSJed Brown 
31937e1895aSJed Brown     if (i+1 < snes->max_its || ms->norms) {
32037e1895aSJed Brown       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
32137e1895aSJed Brown     }
32237e1895aSJed Brown 
32337e1895aSJed Brown     if (ms->norms) {
32437e1895aSJed Brown       ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
325422a814eSBarry Smith       SNESCheckFunctionNorm(snes,fnorm);
326189a9710SBarry Smith 
32737e1895aSJed Brown       /* Monitor convergence */
328e04113cfSBarry Smith       ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
32937e1895aSJed Brown       snes->iter = i+1;
33037e1895aSJed Brown       snes->norm = fnorm;
331e04113cfSBarry Smith       ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
332a71f0d7dSBarry Smith       ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
33337e1895aSJed Brown       ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
33437e1895aSJed Brown 
33537e1895aSJed Brown       /* Test for convergence */
33637e1895aSJed Brown       ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
33737e1895aSJed Brown       if (snes->reason) PetscFunctionReturn(0);
33837e1895aSJed Brown     }
33937e1895aSJed Brown 
34037e1895aSJed Brown     /* Call general purpose update function */
34137e1895aSJed Brown     if (snes->ops->update) {
34237e1895aSJed Brown       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
34337e1895aSJed Brown     }
34437e1895aSJed Brown   }
34537e1895aSJed Brown   if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
34637e1895aSJed Brown   PetscFunctionReturn(0);
34737e1895aSJed Brown }
34837e1895aSJed Brown 
34937e1895aSJed Brown static PetscErrorCode SNESSetUp_MS(SNES snes)
35037e1895aSJed Brown {
35137e1895aSJed Brown   SNES_MS        *ms = (SNES_MS*)snes->data;
35237e1895aSJed Brown   PetscErrorCode ierr;
35337e1895aSJed Brown 
35437e1895aSJed Brown   PetscFunctionBegin;
35537e1895aSJed Brown   if (!ms->tableau) {ierr = SNESMSSetType(snes,SNESMSDefault);CHKERRQ(ierr);}
356fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr);
3576cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
35837e1895aSJed Brown   PetscFunctionReturn(0);
35937e1895aSJed Brown }
36037e1895aSJed Brown 
36137e1895aSJed Brown static PetscErrorCode SNESReset_MS(SNES snes)
36237e1895aSJed Brown {
36337e1895aSJed Brown 
36437e1895aSJed Brown   PetscFunctionBegin;
36537e1895aSJed Brown   PetscFunctionReturn(0);
36637e1895aSJed Brown }
36737e1895aSJed Brown 
36837e1895aSJed Brown static PetscErrorCode SNESDestroy_MS(SNES snes)
36937e1895aSJed Brown {
37037e1895aSJed Brown   PetscErrorCode ierr;
37137e1895aSJed Brown 
37237e1895aSJed Brown   PetscFunctionBegin;
37337e1895aSJed Brown   ierr = PetscFree(snes->data);CHKERRQ(ierr);
374bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr);
37537e1895aSJed Brown   PetscFunctionReturn(0);
37637e1895aSJed Brown }
37737e1895aSJed Brown 
37837e1895aSJed Brown static PetscErrorCode SNESView_MS(SNES snes,PetscViewer viewer)
37937e1895aSJed Brown {
38037e1895aSJed Brown   PetscBool      iascii;
38137e1895aSJed Brown   PetscErrorCode ierr;
38237e1895aSJed Brown   SNES_MS        *ms = (SNES_MS*)snes->data;
38337e1895aSJed Brown 
38437e1895aSJed Brown   PetscFunctionBegin;
385251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
38637e1895aSJed Brown   if (iascii) {
38737e1895aSJed Brown     SNESMSTableau tab = ms->tableau;
38837e1895aSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  multi-stage method type: %s\n",tab ? tab->name : "not yet set");CHKERRQ(ierr);
38937e1895aSJed Brown   }
39037e1895aSJed Brown   PetscFunctionReturn(0);
39137e1895aSJed Brown }
39237e1895aSJed Brown 
3934416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_MS(PetscOptionItems *PetscOptionsObject,SNES snes)
39437e1895aSJed Brown {
395725b86d8SJed Brown   SNES_MS        *ms = (SNES_MS*)snes->data;
39637e1895aSJed Brown   PetscErrorCode ierr;
39737e1895aSJed Brown 
39837e1895aSJed Brown   PetscFunctionBegin;
399e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES MS options");CHKERRQ(ierr);
40037e1895aSJed Brown   {
40137e1895aSJed Brown     SNESMSTableauLink link;
40237e1895aSJed Brown     PetscInt          count,choice;
40337e1895aSJed Brown     PetscBool         flg;
40437e1895aSJed Brown     const char        **namelist;
40537e1895aSJed Brown     char              mstype[256];
40637e1895aSJed Brown 
4078caf3d72SBarry Smith     ierr = PetscStrncpy(mstype,SNESMSDefault,sizeof(mstype));CHKERRQ(ierr);
40837e1895aSJed Brown     for (link=SNESMSTableauList,count=0; link; link=link->next,count++) ;
409785e854fSJed Brown     ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr);
41037e1895aSJed Brown     for (link=SNESMSTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
41137e1895aSJed Brown     ierr = PetscOptionsEList("-snes_ms_type","Multistage smoother type","SNESMSSetType",(const char*const*)namelist,count,mstype,&choice,&flg);CHKERRQ(ierr);
41237e1895aSJed Brown     ierr = SNESMSSetType(snes,flg ? namelist[choice] : mstype);CHKERRQ(ierr);
41337e1895aSJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
4140298fd71SBarry Smith     ierr = PetscOptionsReal("-snes_ms_damping","Damping for multistage method","SNESMSSetDamping",ms->damping,&ms->damping,NULL);CHKERRQ(ierr);
4150298fd71SBarry Smith     ierr = PetscOptionsBool("-snes_ms_norms","Compute norms for monitoring","none",ms->norms,&ms->norms,NULL);CHKERRQ(ierr);
41637e1895aSJed Brown   }
41737e1895aSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
41837e1895aSJed Brown   PetscFunctionReturn(0);
41937e1895aSJed Brown }
42037e1895aSJed Brown 
42137e1895aSJed Brown PetscErrorCode  SNESMSSetType_MS(SNES snes,SNESMSType mstype)
42237e1895aSJed Brown {
42337e1895aSJed Brown   PetscErrorCode    ierr;
42437e1895aSJed Brown   SNES_MS           *ms = (SNES_MS*)snes->data;
42537e1895aSJed Brown   SNESMSTableauLink link;
42637e1895aSJed Brown   PetscBool         match;
42737e1895aSJed Brown 
42837e1895aSJed Brown   PetscFunctionBegin;
42937e1895aSJed Brown   if (ms->tableau) {
43037e1895aSJed Brown     ierr = PetscStrcmp(ms->tableau->name,mstype,&match);CHKERRQ(ierr);
43137e1895aSJed Brown     if (match) PetscFunctionReturn(0);
43237e1895aSJed Brown   }
43337e1895aSJed Brown   for (link = SNESMSTableauList; link; link=link->next) {
43437e1895aSJed Brown     ierr = PetscStrcmp(link->tab.name,mstype,&match);CHKERRQ(ierr);
43537e1895aSJed Brown     if (match) {
43637e1895aSJed Brown       ierr        = SNESReset_MS(snes);CHKERRQ(ierr);
43737e1895aSJed Brown       ms->tableau = &link->tab;
43837e1895aSJed Brown       PetscFunctionReturn(0);
43937e1895aSJed Brown     }
44037e1895aSJed Brown   }
441ce94432eSBarry Smith   SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mstype);
44237e1895aSJed Brown   PetscFunctionReturn(0);
44337e1895aSJed Brown }
44437e1895aSJed Brown 
44537e1895aSJed Brown /*@C
44637e1895aSJed Brown   SNESMSSetType - Set the type of multistage smoother
44737e1895aSJed Brown 
44837e1895aSJed Brown   Logically collective
44937e1895aSJed Brown 
45037e1895aSJed Brown   Input Parameter:
45137e1895aSJed Brown +  snes - nonlinear solver context
45237e1895aSJed Brown -  mstype - type of multistage method
45337e1895aSJed Brown 
45437e1895aSJed Brown   Level: beginner
45537e1895aSJed Brown 
45637e1895aSJed Brown .seealso: SNESMSGetType(), SNESMS
45737e1895aSJed Brown @*/
45819fd82e9SBarry Smith PetscErrorCode SNESMSSetType(SNES snes,SNESMSType rostype)
45937e1895aSJed Brown {
46037e1895aSJed Brown   PetscErrorCode ierr;
46137e1895aSJed Brown 
46237e1895aSJed Brown   PetscFunctionBegin;
46337e1895aSJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
46419fd82e9SBarry Smith   ierr = PetscTryMethod(snes,"SNESMSSetType_C",(SNES,SNESMSType),(snes,rostype));CHKERRQ(ierr);
46537e1895aSJed Brown   PetscFunctionReturn(0);
46637e1895aSJed Brown }
46737e1895aSJed Brown 
46837e1895aSJed Brown /* -------------------------------------------------------------------------- */
46937e1895aSJed Brown /*MC
47037e1895aSJed Brown       SNESMS - multi-stage smoothers
47137e1895aSJed Brown 
47237e1895aSJed Brown       Options Database:
47337e1895aSJed Brown 
47437e1895aSJed Brown +     -snes_ms_type - type of multi-stage smoother
47537e1895aSJed Brown -     -snes_ms_damping - damping for multi-stage method
47637e1895aSJed Brown 
47737e1895aSJed Brown       Notes:
47837e1895aSJed Brown       These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for
4796c9de887SHong Zhang       FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebyshev).
48037e1895aSJed Brown 
48137e1895aSJed Brown       Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds.
48237e1895aSJed Brown 
48337e1895aSJed Brown       The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with SNESMSRegister().
48437e1895aSJed Brown 
48537e1895aSJed Brown       References:
48696a0c994SBarry Smith +   1. -   Ketcheson (2010) Runge Kutta methods with minimum storage implementations.
48796a0c994SBarry Smith .   2. -   Jameson (1983) Solution of the Euler equations for two dimensional transonic flow by a multigrid method.
48896a0c994SBarry Smith -   3. -   Pierce and Giles (1997) Preconditioned multigrid methods for compressible flow calculations on stretched meshes.
48937e1895aSJed Brown 
49037e1895aSJed Brown       Level: beginner
49137e1895aSJed Brown 
4926c9de887SHong Zhang .seealso:  SNESCreate(), SNES, SNESSetType(), SNESMS, SNESFAS, KSPCHEBYSHEV
49337e1895aSJed Brown 
49437e1895aSJed Brown M*/
4958cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_MS(SNES snes)
49637e1895aSJed Brown {
49737e1895aSJed Brown   PetscErrorCode ierr;
49837e1895aSJed Brown   SNES_MS        *ms;
49937e1895aSJed Brown 
50037e1895aSJed Brown   PetscFunctionBegin;
501607a6623SBarry Smith   ierr = SNESMSInitializePackage();CHKERRQ(ierr);
50237e1895aSJed Brown 
50337e1895aSJed Brown   snes->ops->setup          = SNESSetUp_MS;
50437e1895aSJed Brown   snes->ops->solve          = SNESSolve_MS;
50537e1895aSJed Brown   snes->ops->destroy        = SNESDestroy_MS;
50637e1895aSJed Brown   snes->ops->setfromoptions = SNESSetFromOptions_MS;
50737e1895aSJed Brown   snes->ops->view           = SNESView_MS;
50837e1895aSJed Brown   snes->ops->reset          = SNESReset_MS;
50937e1895aSJed Brown 
510*efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
51137e1895aSJed Brown   snes->usesksp = PETSC_TRUE;
51237e1895aSJed Brown 
5134fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
5144fc747eaSLawrence Mitchell 
515b00a9115SJed Brown   ierr        = PetscNewLog(snes,&ms);CHKERRQ(ierr);
51637e1895aSJed Brown   snes->data  = (void*)ms;
51737e1895aSJed Brown   ms->damping = 0.9;
51837e1895aSJed Brown   ms->norms   = PETSC_FALSE;
51937e1895aSJed Brown 
520bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetType_C",SNESMSSetType_MS);CHKERRQ(ierr);
52137e1895aSJed Brown   PetscFunctionReturn(0);
52237e1895aSJed Brown }
52399e0435eSBarry Smith 
524