xref: /petsc/src/snes/impls/ms/ms.c (revision 580bdb303e1ee3b1222b2042810b4c26340259c6)
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   {
4937e1895aSJed Brown     PetscReal gamma[3][1] = {{1.0},{0.0},{0.0}};
5037e1895aSJed Brown     PetscReal delta[1]    = {0.0};
5137e1895aSJed Brown     PetscReal betasub[1]  = {1.0};
5237e1895aSJed Brown     ierr = SNESMSRegister(SNESMSEULER,1,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
5337e1895aSJed Brown   }
5437e1895aSJed Brown 
5537e1895aSJed Brown   {
5637e1895aSJed Brown     PetscReal gamma[3][6] = {
5737e1895aSJed Brown       {0.0000000000000000E+00,  -7.0304722367110606E-01, -1.9836719667506464E-01, -1.6023843981863788E+00, 9.4483822882855284E-02,  -1.4204296130641869E-01},
5837e1895aSJed Brown       {1.0000000000000000E+00,  1.1111025767083920E+00,  5.6150921583923230E-01,  7.4151723494934041E-01,  3.1714538168600587E-01,  4.6479276238548706E-01},
5937e1895aSJed Brown       {0.0000000000000000E+00,  0.0000000000000000E+00,  0.0000000000000000E+00,  6.7968174970583317E-01,  -4.1755042846051737E-03, -1.9115668129923846E-01}
6037e1895aSJed Brown     };
6137e1895aSJed Brown     PetscReal delta[6]   = {1.0000000000000000E+00, 5.3275427433201750E-01, 6.0143544663985238E-01, 4.5874077053842177E-01, 2.7544386906104651E-01, 0.0000000000000000E+00};
6237e1895aSJed Brown     PetscReal betasub[6] = {8.4753115429481929E-01, 7.4018896368655618E-01, 6.5963574086583309E-03, 4.6747795645517759E-01, 1.3314545813643919E-01, 5.3260800028018784E-01};
6337e1895aSJed Brown     ierr = SNESMSRegister(SNESMSM62,6,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
6437e1895aSJed Brown   }
65a97cb6bcSJed Brown   {
66a97cb6bcSJed Brown     PetscReal gamma[3][4] = {{0,0,0,0},{0,0,0,0},{1,1,1,1}};
67a97cb6bcSJed Brown     PetscReal delta[4]    = {0,0,0,0};
68a97cb6bcSJed Brown     PetscReal betasub[4]  = {0.25, 0.5, 0.55, 1.0};
69a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSJAMESON83,4,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
70a97cb6bcSJed Brown   }
71a97cb6bcSJed Brown   {                             /* Van Leer, Tai, and Powell (1989) 2 stage, order 1 */
72a97cb6bcSJed Brown     PetscReal gamma[3][2] = {{0,0},{0,0},{1,1}};
73a97cb6bcSJed Brown     PetscReal delta[2]    = {0,0};
74a97cb6bcSJed Brown     PetscReal betasub[2]  = {0.3333,1.0};
75a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSVLTP21,2,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
76a97cb6bcSJed Brown   }
77a97cb6bcSJed Brown   {                             /* Van Leer, Tai, and Powell (1989) 3 stage, order 1 */
78a97cb6bcSJed Brown     PetscReal gamma[3][3] = {{0,0,0},{0,0,0},{1,1,1}};
79a97cb6bcSJed Brown     PetscReal delta[3]    = {0,0,0};
80a97cb6bcSJed Brown     PetscReal betasub[3]  = {0.1481,0.4000,1.0};
81a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSVLTP31,3,3,1.5,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
82a97cb6bcSJed Brown   }
83a97cb6bcSJed Brown   {                             /* Van Leer, Tai, and Powell (1989) 4 stage, order 1 */
84a97cb6bcSJed Brown     PetscReal gamma[3][4] = {{0,0,0,0},{0,0,0,0},{1,1,1,1}};
85a97cb6bcSJed Brown     PetscReal delta[4]    = {0,0,0,0};
86a97cb6bcSJed Brown     PetscReal betasub[4]  = {0.0833,0.2069,0.4265,1.0};
87a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSVLTP41,4,3,2.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
88a97cb6bcSJed Brown   }
89a97cb6bcSJed Brown   {                             /* Van Leer, Tai, and Powell (1989) 5 stage, order 1 */
90a97cb6bcSJed Brown     PetscReal gamma[3][5] = {{0,0,0,0,0},{0,0,0,0,0},{1,1,1,1,1}};
91a97cb6bcSJed Brown     PetscReal delta[5]    = {0,0,0,0,0};
92a97cb6bcSJed Brown     PetscReal betasub[5]  = {0.0533,0.1263,0.2375,0.4414,1.0};
93a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSVLTP51,5,3,2.5,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
94a97cb6bcSJed Brown   }
95a97cb6bcSJed Brown   {                             /* Van Leer, Tai, and Powell (1989) 6 stage, order 1 */
96a97cb6bcSJed Brown     PetscReal gamma[3][6] = {{0,0,0,0,0,0},{0,0,0,0,0,0},{1,1,1,1,1,1}};
97a97cb6bcSJed Brown     PetscReal delta[6]    = {0,0,0,0,0,0};
98a97cb6bcSJed Brown     PetscReal betasub[6]  = {0.0370,0.0851,0.1521,0.2562,0.4512,1.0};
99a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSVLTP61,6,3,3.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
100a97cb6bcSJed Brown   }
10137e1895aSJed Brown   PetscFunctionReturn(0);
10237e1895aSJed Brown }
10337e1895aSJed Brown 
10437e1895aSJed Brown /*@C
10537e1895aSJed Brown    SNESMSRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
10637e1895aSJed Brown 
10737e1895aSJed Brown    Not Collective
10837e1895aSJed Brown 
10937e1895aSJed Brown    Level: advanced
11037e1895aSJed Brown 
1111c84c290SBarry Smith .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegister()
11237e1895aSJed Brown @*/
11337e1895aSJed Brown PetscErrorCode SNESMSRegisterDestroy(void)
11437e1895aSJed Brown {
11537e1895aSJed Brown   PetscErrorCode    ierr;
11637e1895aSJed Brown   SNESMSTableauLink link;
11737e1895aSJed Brown 
11837e1895aSJed Brown   PetscFunctionBegin;
11937e1895aSJed Brown   while ((link = SNESMSTableauList)) {
12037e1895aSJed Brown     SNESMSTableau t = &link->tab;
12137e1895aSJed Brown     SNESMSTableauList = link->next;
1221aa26658SKarl Rupp 
12337e1895aSJed Brown     ierr = PetscFree3(t->gamma,t->delta,t->betasub);CHKERRQ(ierr);
12437e1895aSJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
12537e1895aSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
12637e1895aSJed Brown   }
12737e1895aSJed Brown   SNESMSRegisterAllCalled = PETSC_FALSE;
12837e1895aSJed Brown   PetscFunctionReturn(0);
12937e1895aSJed Brown }
13037e1895aSJed Brown 
13137e1895aSJed Brown /*@C
13237e1895aSJed Brown   SNESMSInitializePackage - This function initializes everything in the SNESMS package. It is called
1338a690491SBarry Smith   from SNESInitializePackage().
13437e1895aSJed Brown 
13537e1895aSJed Brown   Level: developer
13637e1895aSJed Brown 
13737e1895aSJed Brown .seealso: PetscInitialize()
13837e1895aSJed Brown @*/
139607a6623SBarry Smith PetscErrorCode SNESMSInitializePackage(void)
14037e1895aSJed Brown {
14137e1895aSJed Brown   PetscErrorCode ierr;
14237e1895aSJed Brown 
14337e1895aSJed Brown   PetscFunctionBegin;
14437e1895aSJed Brown   if (SNESMSPackageInitialized) PetscFunctionReturn(0);
14537e1895aSJed Brown   SNESMSPackageInitialized = PETSC_TRUE;
1461aa26658SKarl Rupp 
14737e1895aSJed Brown   ierr = SNESMSRegisterAll();CHKERRQ(ierr);
14837e1895aSJed Brown   ierr = PetscRegisterFinalize(SNESMSFinalizePackage);CHKERRQ(ierr);
14937e1895aSJed Brown   PetscFunctionReturn(0);
15037e1895aSJed Brown }
15137e1895aSJed Brown 
15237e1895aSJed Brown /*@C
15337e1895aSJed Brown   SNESMSFinalizePackage - This function destroys everything in the SNESMS package. It is
15437e1895aSJed Brown   called from PetscFinalize().
15537e1895aSJed Brown 
15637e1895aSJed Brown   Level: developer
15737e1895aSJed Brown 
15837e1895aSJed Brown .seealso: PetscFinalize()
15937e1895aSJed Brown @*/
16037e1895aSJed Brown PetscErrorCode SNESMSFinalizePackage(void)
16137e1895aSJed Brown {
16237e1895aSJed Brown   PetscErrorCode ierr;
16337e1895aSJed Brown 
16437e1895aSJed Brown   PetscFunctionBegin;
16537e1895aSJed Brown   SNESMSPackageInitialized = PETSC_FALSE;
1661aa26658SKarl Rupp 
16737e1895aSJed Brown   ierr = SNESMSRegisterDestroy();CHKERRQ(ierr);
16837e1895aSJed Brown   PetscFunctionReturn(0);
16937e1895aSJed Brown }
17037e1895aSJed Brown 
17137e1895aSJed Brown /*@C
17237e1895aSJed Brown    SNESMSRegister - register a multistage scheme
17337e1895aSJed Brown 
17437e1895aSJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
17537e1895aSJed Brown 
17637e1895aSJed Brown    Input Parameters:
17737e1895aSJed Brown +  name - identifier for method
17837e1895aSJed Brown .  nstages - number of stages
17937e1895aSJed Brown .  nregisters - number of registers used by low-storage implementation
18037e1895aSJed Brown .  gamma - coefficients, see Ketcheson's paper
18137e1895aSJed Brown .  delta - coefficients, see Ketcheson's paper
18237e1895aSJed Brown -  betasub - subdiagonal of Shu-Osher form
18337e1895aSJed Brown 
18437e1895aSJed Brown    Notes:
18537e1895aSJed Brown    The notation is described in Ketcheson (2010) Runge-Kutta methods with minimum storage implementations.
18637e1895aSJed Brown 
18737e1895aSJed Brown    Level: advanced
18837e1895aSJed Brown 
18937e1895aSJed Brown .seealso: SNESMS
19037e1895aSJed Brown @*/
19119fd82e9SBarry Smith PetscErrorCode SNESMSRegister(SNESMSType name,PetscInt nstages,PetscInt nregisters,PetscReal stability,const PetscReal gamma[],const PetscReal delta[],const PetscReal betasub[])
19237e1895aSJed Brown {
19337e1895aSJed Brown   PetscErrorCode    ierr;
19437e1895aSJed Brown   SNESMSTableauLink link;
19537e1895aSJed Brown   SNESMSTableau     t;
19637e1895aSJed Brown 
19737e1895aSJed Brown   PetscFunctionBegin;
19837e1895aSJed Brown   PetscValidCharPointer(name,1);
19937e1895aSJed Brown   if (nstages < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have at least one stage");
20037e1895aSJed Brown   if (nregisters != 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only support for methods written in 3-register form");
20137e1895aSJed Brown   PetscValidPointer(gamma,4);
20237e1895aSJed Brown   PetscValidPointer(delta,5);
20337e1895aSJed Brown   PetscValidPointer(betasub,6);
20437e1895aSJed Brown 
2059cc31a68SJed Brown   ierr          = SNESMSInitializePackage();CHKERRQ(ierr);
206854ce69bSBarry Smith   ierr          = PetscNew(&link);CHKERRQ(ierr);
20737e1895aSJed Brown   t             = &link->tab;
20837e1895aSJed Brown   ierr          = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
20937e1895aSJed Brown   t->nstages    = nstages;
21037e1895aSJed Brown   t->nregisters = nregisters;
21137e1895aSJed Brown   t->stability  = stability;
2121aa26658SKarl Rupp 
213dcca6d9dSJed Brown   ierr = PetscMalloc3(nstages*nregisters,&t->gamma,nstages,&t->delta,nstages,&t->betasub);CHKERRQ(ierr);
214*580bdb30SBarry Smith   ierr = PetscArraycpy(t->gamma,gamma,nstages*nregisters);CHKERRQ(ierr);
215*580bdb30SBarry Smith   ierr = PetscArraycpy(t->delta,delta,nstages);CHKERRQ(ierr);
216*580bdb30SBarry Smith   ierr = PetscArraycpy(t->betasub,betasub,nstages);CHKERRQ(ierr);
21737e1895aSJed Brown 
21837e1895aSJed Brown   link->next        = SNESMSTableauList;
21937e1895aSJed Brown   SNESMSTableauList = link;
22037e1895aSJed Brown   PetscFunctionReturn(0);
22137e1895aSJed Brown }
22237e1895aSJed Brown 
22337e1895aSJed Brown /*
22437e1895aSJed Brown   X - initial state, updated in-place.
22537e1895aSJed Brown   F - residual, computed at the initial X on input
22637e1895aSJed Brown */
22737e1895aSJed Brown static PetscErrorCode SNESMSStep_3Sstar(SNES snes,Vec X,Vec F)
22837e1895aSJed Brown {
22937e1895aSJed Brown   PetscErrorCode  ierr;
23037e1895aSJed Brown   SNES_MS         *ms    = (SNES_MS*)snes->data;
23137e1895aSJed Brown   SNESMSTableau   t      = ms->tableau;
23237e1895aSJed Brown   const PetscReal *gamma = t->gamma,*delta = t->delta,*betasub = t->betasub;
23337e1895aSJed Brown   Vec             S1,S2,S3,Y;
23437e1895aSJed Brown   PetscInt        i,nstages = t->nstages;;
23537e1895aSJed Brown 
23637e1895aSJed Brown 
23737e1895aSJed Brown   PetscFunctionBegin;
23837e1895aSJed Brown   Y    = snes->work[0];
23937e1895aSJed Brown   S1   = X;
24037e1895aSJed Brown   S2   = snes->work[1];
24137e1895aSJed Brown   S3   = snes->work[2];
24237e1895aSJed Brown   ierr = VecZeroEntries(S2);CHKERRQ(ierr);
24337e1895aSJed Brown   ierr = VecCopy(X,S3);CHKERRQ(ierr);
24437e1895aSJed Brown   for (i=0; i<nstages; i++) {
245da80777bSKarl Rupp     Vec         Ss[4];
246da80777bSKarl Rupp     PetscScalar scoeff[4];
247da80777bSKarl Rupp 
248da80777bSKarl Rupp     Ss[0] = S1; Ss[1] = S2; Ss[2] = S3; Ss[3] = Y;
2491aa26658SKarl Rupp 
250da80777bSKarl Rupp     scoeff[0] = gamma[0*nstages+i]-(PetscReal)1.0;
251da80777bSKarl Rupp     scoeff[1] = gamma[1*nstages+i];
252da80777bSKarl Rupp     scoeff[2] = gamma[2*nstages+i];
253da80777bSKarl Rupp     scoeff[3] = -betasub[i]*ms->damping;
2541aa26658SKarl Rupp 
25537e1895aSJed Brown     ierr = VecAXPY(S2,delta[i],S1);CHKERRQ(ierr);
25637e1895aSJed Brown     if (i>0) {
25737e1895aSJed Brown       ierr = SNESComputeFunction(snes,S1,F);CHKERRQ(ierr);
25837e1895aSJed Brown     }
259d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
26037e1895aSJed Brown     ierr = VecMAXPY(S1,4,scoeff,Ss);CHKERRQ(ierr);
26137e1895aSJed Brown   }
26237e1895aSJed Brown   PetscFunctionReturn(0);
26337e1895aSJed Brown }
26437e1895aSJed Brown 
26537e1895aSJed Brown static PetscErrorCode SNESSolve_MS(SNES snes)
26637e1895aSJed Brown {
26737e1895aSJed Brown   SNES_MS        *ms = (SNES_MS*)snes->data;
26837e1895aSJed Brown   Vec            X   = snes->vec_sol,F = snes->vec_func;
26937e1895aSJed Brown   PetscReal      fnorm;
27037e1895aSJed Brown   PetscInt       i;
27137e1895aSJed Brown   PetscErrorCode ierr;
27237e1895aSJed Brown 
27337e1895aSJed Brown   PetscFunctionBegin;
2746c4ed002SBarry 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);
275c579b300SPatrick Farrell 
276fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
27737e1895aSJed Brown   snes->reason = SNES_CONVERGED_ITERATING;
278e04113cfSBarry Smith   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
27937e1895aSJed Brown   snes->iter   = 0;
28037e1895aSJed Brown   snes->norm   = 0.;
281e04113cfSBarry Smith   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
282e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
28337e1895aSJed Brown     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
2841aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
2851aa26658SKarl Rupp 
28637e1895aSJed Brown   if (snes->jacobian) {         /* This method does not require a Jacobian, but it is usually preconditioned by PBJacobi */
287d1e9a80fSBarry Smith     ierr = SNESComputeJacobian(snes,snes->vec_sol,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
28807b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
28937e1895aSJed Brown   }
29037e1895aSJed Brown   if (ms->norms) {
29137e1895aSJed Brown     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
292422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
29337e1895aSJed Brown     /* Monitor convergence */
294e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
29537e1895aSJed Brown     snes->iter = 0;
29637e1895aSJed Brown     snes->norm = fnorm;
297e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
298a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
29937e1895aSJed Brown     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
30037e1895aSJed Brown 
30137e1895aSJed Brown     /* Test for convergence */
30237e1895aSJed Brown     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
30337e1895aSJed Brown     if (snes->reason) PetscFunctionReturn(0);
30437e1895aSJed Brown   }
30537e1895aSJed Brown 
30637e1895aSJed Brown   /* Call general purpose update function */
30737e1895aSJed Brown   if (snes->ops->update) {
30837e1895aSJed Brown     ierr = (*snes->ops->update)(snes,snes->iter);CHKERRQ(ierr);
30937e1895aSJed Brown   }
31037e1895aSJed Brown   for (i = 0; i < snes->max_its; i++) {
31137e1895aSJed Brown     ierr = SNESMSStep_3Sstar(snes,X,F);CHKERRQ(ierr);
31237e1895aSJed Brown 
31337e1895aSJed Brown     if (i+1 < snes->max_its || ms->norms) {
31437e1895aSJed Brown       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
31537e1895aSJed Brown     }
31637e1895aSJed Brown 
31737e1895aSJed Brown     if (ms->norms) {
31837e1895aSJed Brown       ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
319422a814eSBarry Smith       SNESCheckFunctionNorm(snes,fnorm);
320189a9710SBarry Smith 
32137e1895aSJed Brown       /* Monitor convergence */
322e04113cfSBarry Smith       ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
32337e1895aSJed Brown       snes->iter = i+1;
32437e1895aSJed Brown       snes->norm = fnorm;
325e04113cfSBarry Smith       ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
326a71f0d7dSBarry Smith       ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
32737e1895aSJed Brown       ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
32837e1895aSJed Brown 
32937e1895aSJed Brown       /* Test for convergence */
33037e1895aSJed Brown       ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
33137e1895aSJed Brown       if (snes->reason) PetscFunctionReturn(0);
33237e1895aSJed Brown     }
33337e1895aSJed Brown 
33437e1895aSJed Brown     /* Call general purpose update function */
33537e1895aSJed Brown     if (snes->ops->update) {
33637e1895aSJed Brown       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
33737e1895aSJed Brown     }
33837e1895aSJed Brown   }
33937e1895aSJed Brown   if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
34037e1895aSJed Brown   PetscFunctionReturn(0);
34137e1895aSJed Brown }
34237e1895aSJed Brown 
34337e1895aSJed Brown static PetscErrorCode SNESSetUp_MS(SNES snes)
34437e1895aSJed Brown {
34537e1895aSJed Brown   SNES_MS        *ms = (SNES_MS*)snes->data;
34637e1895aSJed Brown   PetscErrorCode ierr;
34737e1895aSJed Brown 
34837e1895aSJed Brown   PetscFunctionBegin;
34937e1895aSJed Brown   if (!ms->tableau) {ierr = SNESMSSetType(snes,SNESMSDefault);CHKERRQ(ierr);}
350fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr);
3516cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
35237e1895aSJed Brown   PetscFunctionReturn(0);
35337e1895aSJed Brown }
35437e1895aSJed Brown 
35537e1895aSJed Brown static PetscErrorCode SNESReset_MS(SNES snes)
35637e1895aSJed Brown {
35737e1895aSJed Brown 
35837e1895aSJed Brown   PetscFunctionBegin;
35937e1895aSJed Brown   PetscFunctionReturn(0);
36037e1895aSJed Brown }
36137e1895aSJed Brown 
36237e1895aSJed Brown static PetscErrorCode SNESDestroy_MS(SNES snes)
36337e1895aSJed Brown {
36437e1895aSJed Brown   PetscErrorCode ierr;
36537e1895aSJed Brown 
36637e1895aSJed Brown   PetscFunctionBegin;
36737e1895aSJed Brown   ierr = PetscFree(snes->data);CHKERRQ(ierr);
368bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr);
36937e1895aSJed Brown   PetscFunctionReturn(0);
37037e1895aSJed Brown }
37137e1895aSJed Brown 
37237e1895aSJed Brown static PetscErrorCode SNESView_MS(SNES snes,PetscViewer viewer)
37337e1895aSJed Brown {
37437e1895aSJed Brown   PetscBool      iascii;
37537e1895aSJed Brown   PetscErrorCode ierr;
37637e1895aSJed Brown   SNES_MS        *ms = (SNES_MS*)snes->data;
37737e1895aSJed Brown 
37837e1895aSJed Brown   PetscFunctionBegin;
379251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
38037e1895aSJed Brown   if (iascii) {
38137e1895aSJed Brown     SNESMSTableau tab = ms->tableau;
38237e1895aSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  multi-stage method type: %s\n",tab ? tab->name : "not yet set");CHKERRQ(ierr);
38337e1895aSJed Brown   }
38437e1895aSJed Brown   PetscFunctionReturn(0);
38537e1895aSJed Brown }
38637e1895aSJed Brown 
3874416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_MS(PetscOptionItems *PetscOptionsObject,SNES snes)
38837e1895aSJed Brown {
389725b86d8SJed Brown   SNES_MS        *ms = (SNES_MS*)snes->data;
39037e1895aSJed Brown   PetscErrorCode ierr;
39137e1895aSJed Brown 
39237e1895aSJed Brown   PetscFunctionBegin;
393e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES MS options");CHKERRQ(ierr);
39437e1895aSJed Brown   {
39537e1895aSJed Brown     SNESMSTableauLink link;
39637e1895aSJed Brown     PetscInt          count,choice;
39737e1895aSJed Brown     PetscBool         flg;
39837e1895aSJed Brown     const char        **namelist;
39937e1895aSJed Brown     char              mstype[256];
40037e1895aSJed Brown 
4018caf3d72SBarry Smith     ierr = PetscStrncpy(mstype,SNESMSDefault,sizeof(mstype));CHKERRQ(ierr);
40237e1895aSJed Brown     for (link=SNESMSTableauList,count=0; link; link=link->next,count++) ;
403f489ac74SBarry Smith     ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr);
40437e1895aSJed Brown     for (link=SNESMSTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
40537e1895aSJed Brown     ierr = PetscOptionsEList("-snes_ms_type","Multistage smoother type","SNESMSSetType",(const char*const*)namelist,count,mstype,&choice,&flg);CHKERRQ(ierr);
40637e1895aSJed Brown     ierr = SNESMSSetType(snes,flg ? namelist[choice] : mstype);CHKERRQ(ierr);
40737e1895aSJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
4080298fd71SBarry Smith     ierr = PetscOptionsReal("-snes_ms_damping","Damping for multistage method","SNESMSSetDamping",ms->damping,&ms->damping,NULL);CHKERRQ(ierr);
4090298fd71SBarry Smith     ierr = PetscOptionsBool("-snes_ms_norms","Compute norms for monitoring","none",ms->norms,&ms->norms,NULL);CHKERRQ(ierr);
41037e1895aSJed Brown   }
41137e1895aSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
41237e1895aSJed Brown   PetscFunctionReturn(0);
41337e1895aSJed Brown }
41437e1895aSJed Brown 
41537e1895aSJed Brown PetscErrorCode  SNESMSSetType_MS(SNES snes,SNESMSType mstype)
41637e1895aSJed Brown {
41737e1895aSJed Brown   PetscErrorCode    ierr;
41837e1895aSJed Brown   SNES_MS           *ms = (SNES_MS*)snes->data;
41937e1895aSJed Brown   SNESMSTableauLink link;
42037e1895aSJed Brown   PetscBool         match;
42137e1895aSJed Brown 
42237e1895aSJed Brown   PetscFunctionBegin;
42337e1895aSJed Brown   if (ms->tableau) {
42437e1895aSJed Brown     ierr = PetscStrcmp(ms->tableau->name,mstype,&match);CHKERRQ(ierr);
42537e1895aSJed Brown     if (match) PetscFunctionReturn(0);
42637e1895aSJed Brown   }
42737e1895aSJed Brown   for (link = SNESMSTableauList; link; link=link->next) {
42837e1895aSJed Brown     ierr = PetscStrcmp(link->tab.name,mstype,&match);CHKERRQ(ierr);
42937e1895aSJed Brown     if (match) {
43037e1895aSJed Brown       ierr        = SNESReset_MS(snes);CHKERRQ(ierr);
43137e1895aSJed Brown       ms->tableau = &link->tab;
43237e1895aSJed Brown       PetscFunctionReturn(0);
43337e1895aSJed Brown     }
43437e1895aSJed Brown   }
435ce94432eSBarry Smith   SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mstype);
43637e1895aSJed Brown   PetscFunctionReturn(0);
43737e1895aSJed Brown }
43837e1895aSJed Brown 
43937e1895aSJed Brown /*@C
44037e1895aSJed Brown   SNESMSSetType - Set the type of multistage smoother
44137e1895aSJed Brown 
44237e1895aSJed Brown   Logically collective
44337e1895aSJed Brown 
44437e1895aSJed Brown   Input Parameter:
44537e1895aSJed Brown +  snes - nonlinear solver context
44637e1895aSJed Brown -  mstype - type of multistage method
44737e1895aSJed Brown 
44837e1895aSJed Brown   Level: beginner
44937e1895aSJed Brown 
45037e1895aSJed Brown .seealso: SNESMSGetType(), SNESMS
45137e1895aSJed Brown @*/
45219fd82e9SBarry Smith PetscErrorCode SNESMSSetType(SNES snes,SNESMSType rostype)
45337e1895aSJed Brown {
45437e1895aSJed Brown   PetscErrorCode ierr;
45537e1895aSJed Brown 
45637e1895aSJed Brown   PetscFunctionBegin;
45737e1895aSJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
45819fd82e9SBarry Smith   ierr = PetscTryMethod(snes,"SNESMSSetType_C",(SNES,SNESMSType),(snes,rostype));CHKERRQ(ierr);
45937e1895aSJed Brown   PetscFunctionReturn(0);
46037e1895aSJed Brown }
46137e1895aSJed Brown 
46237e1895aSJed Brown /* -------------------------------------------------------------------------- */
46337e1895aSJed Brown /*MC
46437e1895aSJed Brown       SNESMS - multi-stage smoothers
46537e1895aSJed Brown 
46637e1895aSJed Brown       Options Database:
46737e1895aSJed Brown 
46837e1895aSJed Brown +     -snes_ms_type - type of multi-stage smoother
46937e1895aSJed Brown -     -snes_ms_damping - damping for multi-stage method
47037e1895aSJed Brown 
47137e1895aSJed Brown       Notes:
47237e1895aSJed Brown       These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for
4736c9de887SHong Zhang       FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebyshev).
47437e1895aSJed Brown 
47537e1895aSJed Brown       Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds.
47637e1895aSJed Brown 
47737e1895aSJed Brown       The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with SNESMSRegister().
47837e1895aSJed Brown 
47937e1895aSJed Brown       References:
48096a0c994SBarry Smith +   1. -   Ketcheson (2010) Runge Kutta methods with minimum storage implementations.
48196a0c994SBarry Smith .   2. -   Jameson (1983) Solution of the Euler equations for two dimensional transonic flow by a multigrid method.
48296a0c994SBarry Smith -   3. -   Pierce and Giles (1997) Preconditioned multigrid methods for compressible flow calculations on stretched meshes.
48337e1895aSJed Brown 
48437e1895aSJed Brown       Level: beginner
48537e1895aSJed Brown 
4866c9de887SHong Zhang .seealso:  SNESCreate(), SNES, SNESSetType(), SNESMS, SNESFAS, KSPCHEBYSHEV
48737e1895aSJed Brown 
48837e1895aSJed Brown M*/
4898cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_MS(SNES snes)
49037e1895aSJed Brown {
49137e1895aSJed Brown   PetscErrorCode ierr;
49237e1895aSJed Brown   SNES_MS        *ms;
49337e1895aSJed Brown 
49437e1895aSJed Brown   PetscFunctionBegin;
495607a6623SBarry Smith   ierr = SNESMSInitializePackage();CHKERRQ(ierr);
49637e1895aSJed Brown 
49737e1895aSJed Brown   snes->ops->setup          = SNESSetUp_MS;
49837e1895aSJed Brown   snes->ops->solve          = SNESSolve_MS;
49937e1895aSJed Brown   snes->ops->destroy        = SNESDestroy_MS;
50037e1895aSJed Brown   snes->ops->setfromoptions = SNESSetFromOptions_MS;
50137e1895aSJed Brown   snes->ops->view           = SNESView_MS;
50237e1895aSJed Brown   snes->ops->reset          = SNESReset_MS;
50337e1895aSJed Brown 
504efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
50537e1895aSJed Brown   snes->usesksp = PETSC_TRUE;
50637e1895aSJed Brown 
5074fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
5084fc747eaSLawrence Mitchell 
509b00a9115SJed Brown   ierr        = PetscNewLog(snes,&ms);CHKERRQ(ierr);
51037e1895aSJed Brown   snes->data  = (void*)ms;
51137e1895aSJed Brown   ms->damping = 0.9;
51237e1895aSJed Brown   ms->norms   = PETSC_FALSE;
51337e1895aSJed Brown 
514bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetType_C",SNESMSSetType_MS);CHKERRQ(ierr);
51537e1895aSJed Brown   PetscFunctionReturn(0);
51637e1895aSJed Brown }
51799e0435eSBarry Smith 
518