xref: /petsc/src/snes/impls/ms/ms.c (revision dcca6d9d80ebd869fe6029bd05a3aa9faafef49e)
1b45d2f2cSJed Brown #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 #undef __FUNCT__
3237e1895aSJed Brown #define __FUNCT__ "SNESMSRegisterAll"
3337e1895aSJed Brown /*@C
3437e1895aSJed Brown   SNESMSRegisterAll - Registers all of the multi-stage methods in SNESMS
3537e1895aSJed Brown 
3637e1895aSJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
3737e1895aSJed Brown 
3837e1895aSJed Brown   Level: advanced
3937e1895aSJed Brown 
4037e1895aSJed Brown .keywords: SNES, SNESMS, register, all
4137e1895aSJed Brown 
4237e1895aSJed Brown .seealso:  SNESMSRegisterDestroy()
4337e1895aSJed Brown @*/
4437e1895aSJed Brown PetscErrorCode SNESMSRegisterAll(void)
4537e1895aSJed Brown {
4637e1895aSJed Brown   PetscErrorCode ierr;
4737e1895aSJed Brown 
4837e1895aSJed Brown   PetscFunctionBegin;
4937e1895aSJed Brown   if (SNESMSRegisterAllCalled) PetscFunctionReturn(0);
5037e1895aSJed Brown   SNESMSRegisterAllCalled = PETSC_TRUE;
5137e1895aSJed Brown 
5237e1895aSJed Brown   {
5337e1895aSJed Brown     PetscReal gamma[3][1] = {{1.0},{0.0},{0.0}};
5437e1895aSJed Brown     PetscReal delta[1]    = {0.0};
5537e1895aSJed Brown     PetscReal betasub[1]  = {1.0};
5637e1895aSJed Brown     ierr = SNESMSRegister(SNESMSEULER,1,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
5737e1895aSJed Brown   }
5837e1895aSJed Brown 
5937e1895aSJed Brown   {
6037e1895aSJed Brown     PetscReal gamma[3][6] = {
6137e1895aSJed Brown       {0.0000000000000000E+00,  -7.0304722367110606E-01, -1.9836719667506464E-01, -1.6023843981863788E+00, 9.4483822882855284E-02,  -1.4204296130641869E-01},
6237e1895aSJed Brown       {1.0000000000000000E+00,  1.1111025767083920E+00,  5.6150921583923230E-01,  7.4151723494934041E-01,  3.1714538168600587E-01,  4.6479276238548706E-01},
6337e1895aSJed Brown       {0.0000000000000000E+00,  0.0000000000000000E+00,  0.0000000000000000E+00,  6.7968174970583317E-01,  -4.1755042846051737E-03, -1.9115668129923846E-01}
6437e1895aSJed Brown     };
6537e1895aSJed Brown     PetscReal delta[6]   = {1.0000000000000000E+00, 5.3275427433201750E-01, 6.0143544663985238E-01, 4.5874077053842177E-01, 2.7544386906104651E-01, 0.0000000000000000E+00};
6637e1895aSJed Brown     PetscReal betasub[6] = {8.4753115429481929E-01, 7.4018896368655618E-01, 6.5963574086583309E-03, 4.6747795645517759E-01, 1.3314545813643919E-01, 5.3260800028018784E-01};
6737e1895aSJed Brown     ierr = SNESMSRegister(SNESMSM62,6,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
6837e1895aSJed Brown   }
69a97cb6bcSJed Brown   {
70a97cb6bcSJed Brown     PetscReal gamma[3][4] = {{0,0,0,0},{0,0,0,0},{1,1,1,1}};
71a97cb6bcSJed Brown     PetscReal delta[4]    = {0,0,0,0};
72a97cb6bcSJed Brown     PetscReal betasub[4]  = {0.25, 0.5, 0.55, 1.0};
73a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSJAMESON83,4,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
74a97cb6bcSJed Brown   }
75a97cb6bcSJed Brown   {                             /* Van Leer, Tai, and Powell (1989) 2 stage, order 1 */
76a97cb6bcSJed Brown     PetscReal gamma[3][2] = {{0,0},{0,0},{1,1}};
77a97cb6bcSJed Brown     PetscReal delta[2]    = {0,0};
78a97cb6bcSJed Brown     PetscReal betasub[2]  = {0.3333,1.0};
79a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSVLTP21,2,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
80a97cb6bcSJed Brown   }
81a97cb6bcSJed Brown   {                             /* Van Leer, Tai, and Powell (1989) 3 stage, order 1 */
82a97cb6bcSJed Brown     PetscReal gamma[3][3] = {{0,0,0},{0,0,0},{1,1,1}};
83a97cb6bcSJed Brown     PetscReal delta[3]    = {0,0,0};
84a97cb6bcSJed Brown     PetscReal betasub[3]  = {0.1481,0.4000,1.0};
85a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSVLTP31,3,3,1.5,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
86a97cb6bcSJed Brown   }
87a97cb6bcSJed Brown   {                             /* Van Leer, Tai, and Powell (1989) 4 stage, order 1 */
88a97cb6bcSJed Brown     PetscReal gamma[3][4] = {{0,0,0,0},{0,0,0,0},{1,1,1,1}};
89a97cb6bcSJed Brown     PetscReal delta[4]    = {0,0,0,0};
90a97cb6bcSJed Brown     PetscReal betasub[4]  = {0.0833,0.2069,0.4265,1.0};
91a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSVLTP41,4,3,2.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
92a97cb6bcSJed Brown   }
93a97cb6bcSJed Brown   {                             /* Van Leer, Tai, and Powell (1989) 5 stage, order 1 */
94a97cb6bcSJed Brown     PetscReal gamma[3][5] = {{0,0,0,0,0},{0,0,0,0,0},{1,1,1,1,1}};
95a97cb6bcSJed Brown     PetscReal delta[5]    = {0,0,0,0,0};
96a97cb6bcSJed Brown     PetscReal betasub[5]  = {0.0533,0.1263,0.2375,0.4414,1.0};
97a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSVLTP51,5,3,2.5,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
98a97cb6bcSJed Brown   }
99a97cb6bcSJed Brown   {                             /* Van Leer, Tai, and Powell (1989) 6 stage, order 1 */
100a97cb6bcSJed Brown     PetscReal gamma[3][6] = {{0,0,0,0,0,0},{0,0,0,0,0,0},{1,1,1,1,1,1}};
101a97cb6bcSJed Brown     PetscReal delta[6]    = {0,0,0,0,0,0};
102a97cb6bcSJed Brown     PetscReal betasub[6]  = {0.0370,0.0851,0.1521,0.2562,0.4512,1.0};
103a97cb6bcSJed Brown     ierr = SNESMSRegister(SNESMSVLTP61,6,3,3.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr);
104a97cb6bcSJed Brown   }
10537e1895aSJed Brown   PetscFunctionReturn(0);
10637e1895aSJed Brown }
10737e1895aSJed Brown 
10837e1895aSJed Brown #undef __FUNCT__
10937e1895aSJed Brown #define __FUNCT__ "SNESMSRegisterDestroy"
11037e1895aSJed Brown /*@C
11137e1895aSJed Brown    SNESMSRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
11237e1895aSJed Brown 
11337e1895aSJed Brown    Not Collective
11437e1895aSJed Brown 
11537e1895aSJed Brown    Level: advanced
11637e1895aSJed Brown 
11737e1895aSJed Brown .keywords: TSRosW, register, destroy
1181c84c290SBarry Smith .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegister()
11937e1895aSJed Brown @*/
12037e1895aSJed Brown PetscErrorCode SNESMSRegisterDestroy(void)
12137e1895aSJed Brown {
12237e1895aSJed Brown   PetscErrorCode    ierr;
12337e1895aSJed Brown   SNESMSTableauLink link;
12437e1895aSJed Brown 
12537e1895aSJed Brown   PetscFunctionBegin;
12637e1895aSJed Brown   while ((link = SNESMSTableauList)) {
12737e1895aSJed Brown     SNESMSTableau t = &link->tab;
12837e1895aSJed Brown     SNESMSTableauList = link->next;
1291aa26658SKarl Rupp 
13037e1895aSJed Brown     ierr = PetscFree3(t->gamma,t->delta,t->betasub);CHKERRQ(ierr);
13137e1895aSJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
13237e1895aSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
13337e1895aSJed Brown   }
13437e1895aSJed Brown   SNESMSRegisterAllCalled = PETSC_FALSE;
13537e1895aSJed Brown   PetscFunctionReturn(0);
13637e1895aSJed Brown }
13737e1895aSJed Brown 
13837e1895aSJed Brown #undef __FUNCT__
13937e1895aSJed Brown #define __FUNCT__ "SNESMSInitializePackage"
14037e1895aSJed Brown /*@C
14137e1895aSJed Brown   SNESMSInitializePackage - This function initializes everything in the SNESMS package. It is called
14237e1895aSJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to SNESCreate_MS()
14337e1895aSJed Brown   when using static libraries.
14437e1895aSJed Brown 
14537e1895aSJed Brown   Level: developer
14637e1895aSJed Brown 
14737e1895aSJed Brown .keywords: SNES, SNESMS, initialize, package
14837e1895aSJed Brown .seealso: PetscInitialize()
14937e1895aSJed Brown @*/
150607a6623SBarry Smith PetscErrorCode SNESMSInitializePackage(void)
15137e1895aSJed Brown {
15237e1895aSJed Brown   PetscErrorCode ierr;
15337e1895aSJed Brown 
15437e1895aSJed Brown   PetscFunctionBegin;
15537e1895aSJed Brown   if (SNESMSPackageInitialized) PetscFunctionReturn(0);
15637e1895aSJed Brown   SNESMSPackageInitialized = PETSC_TRUE;
1571aa26658SKarl Rupp 
15837e1895aSJed Brown   ierr = SNESMSRegisterAll();CHKERRQ(ierr);
15937e1895aSJed Brown   ierr = PetscRegisterFinalize(SNESMSFinalizePackage);CHKERRQ(ierr);
16037e1895aSJed Brown   PetscFunctionReturn(0);
16137e1895aSJed Brown }
16237e1895aSJed Brown 
16337e1895aSJed Brown #undef __FUNCT__
16437e1895aSJed Brown #define __FUNCT__ "SNESMSFinalizePackage"
16537e1895aSJed Brown /*@C
16637e1895aSJed Brown   SNESMSFinalizePackage - This function destroys everything in the SNESMS package. It is
16737e1895aSJed Brown   called from PetscFinalize().
16837e1895aSJed Brown 
16937e1895aSJed Brown   Level: developer
17037e1895aSJed Brown 
17137e1895aSJed Brown .keywords: Petsc, destroy, package
17237e1895aSJed Brown .seealso: PetscFinalize()
17337e1895aSJed Brown @*/
17437e1895aSJed Brown PetscErrorCode SNESMSFinalizePackage(void)
17537e1895aSJed Brown {
17637e1895aSJed Brown   PetscErrorCode ierr;
17737e1895aSJed Brown 
17837e1895aSJed Brown   PetscFunctionBegin;
17937e1895aSJed Brown   SNESMSPackageInitialized = PETSC_FALSE;
1801aa26658SKarl Rupp 
18137e1895aSJed Brown   ierr = SNESMSRegisterDestroy();CHKERRQ(ierr);
18237e1895aSJed Brown   PetscFunctionReturn(0);
18337e1895aSJed Brown }
18437e1895aSJed Brown 
18537e1895aSJed Brown #undef __FUNCT__
18637e1895aSJed Brown #define __FUNCT__ "SNESMSRegister"
18737e1895aSJed Brown /*@C
18837e1895aSJed Brown    SNESMSRegister - register a multistage scheme
18937e1895aSJed Brown 
19037e1895aSJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
19137e1895aSJed Brown 
19237e1895aSJed Brown    Input Parameters:
19337e1895aSJed Brown +  name - identifier for method
19437e1895aSJed Brown .  nstages - number of stages
19537e1895aSJed Brown .  nregisters - number of registers used by low-storage implementation
19637e1895aSJed Brown .  gamma - coefficients, see Ketcheson's paper
19737e1895aSJed Brown .  delta - coefficients, see Ketcheson's paper
19837e1895aSJed Brown -  betasub - subdiagonal of Shu-Osher form
19937e1895aSJed Brown 
20037e1895aSJed Brown    Notes:
20137e1895aSJed Brown    The notation is described in Ketcheson (2010) Runge-Kutta methods with minimum storage implementations.
20237e1895aSJed Brown 
20337e1895aSJed Brown    Level: advanced
20437e1895aSJed Brown 
20537e1895aSJed Brown .keywords: SNES, register
20637e1895aSJed Brown 
20737e1895aSJed Brown .seealso: SNESMS
20837e1895aSJed Brown @*/
20919fd82e9SBarry Smith PetscErrorCode SNESMSRegister(SNESMSType name,PetscInt nstages,PetscInt nregisters,PetscReal stability,const PetscReal gamma[],const PetscReal delta[],const PetscReal betasub[])
21037e1895aSJed Brown {
21137e1895aSJed Brown   PetscErrorCode    ierr;
21237e1895aSJed Brown   SNESMSTableauLink link;
21337e1895aSJed Brown   SNESMSTableau     t;
21437e1895aSJed Brown 
21537e1895aSJed Brown   PetscFunctionBegin;
21637e1895aSJed Brown   PetscValidCharPointer(name,1);
21737e1895aSJed Brown   if (nstages < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have at least one stage");
21837e1895aSJed Brown   if (nregisters != 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only support for methods written in 3-register form");
21937e1895aSJed Brown   PetscValidPointer(gamma,4);
22037e1895aSJed Brown   PetscValidPointer(delta,5);
22137e1895aSJed Brown   PetscValidPointer(betasub,6);
22237e1895aSJed Brown 
22337e1895aSJed Brown   ierr          = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
22437e1895aSJed Brown   ierr          = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
22537e1895aSJed Brown   t             = &link->tab;
22637e1895aSJed Brown   ierr          = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
22737e1895aSJed Brown   t->nstages    = nstages;
22837e1895aSJed Brown   t->nregisters = nregisters;
22937e1895aSJed Brown   t->stability  = stability;
2301aa26658SKarl Rupp 
231*dcca6d9dSJed Brown   ierr = PetscMalloc3(nstages*nregisters,&t->gamma,nstages,&t->delta,nstages,&t->betasub);CHKERRQ(ierr);
23237e1895aSJed Brown   ierr = PetscMemcpy(t->gamma,gamma,nstages*nregisters*sizeof(PetscReal));CHKERRQ(ierr);
23337e1895aSJed Brown   ierr = PetscMemcpy(t->delta,delta,nstages*sizeof(PetscReal));CHKERRQ(ierr);
23437e1895aSJed Brown   ierr = PetscMemcpy(t->betasub,betasub,nstages*sizeof(PetscReal));CHKERRQ(ierr);
23537e1895aSJed Brown 
23637e1895aSJed Brown   link->next        = SNESMSTableauList;
23737e1895aSJed Brown   SNESMSTableauList = link;
23837e1895aSJed Brown   PetscFunctionReturn(0);
23937e1895aSJed Brown }
24037e1895aSJed Brown 
24137e1895aSJed Brown #undef __FUNCT__
24237e1895aSJed Brown #define __FUNCT__ "SNESMSStep_3Sstar"
24337e1895aSJed Brown /*
24437e1895aSJed Brown   X - initial state, updated in-place.
24537e1895aSJed Brown   F - residual, computed at the initial X on input
24637e1895aSJed Brown */
24737e1895aSJed Brown static PetscErrorCode SNESMSStep_3Sstar(SNES snes,Vec X,Vec F)
24837e1895aSJed Brown {
24937e1895aSJed Brown   PetscErrorCode  ierr;
25037e1895aSJed Brown   SNES_MS         *ms    = (SNES_MS*)snes->data;
25137e1895aSJed Brown   SNESMSTableau   t      = ms->tableau;
25237e1895aSJed Brown   const PetscReal *gamma = t->gamma,*delta = t->delta,*betasub = t->betasub;
25337e1895aSJed Brown   Vec             S1,S2,S3,Y;
25437e1895aSJed Brown   PetscInt        i,nstages = t->nstages;;
25537e1895aSJed Brown 
25637e1895aSJed Brown 
25737e1895aSJed Brown   PetscFunctionBegin;
25837e1895aSJed Brown   Y    = snes->work[0];
25937e1895aSJed Brown   S1   = X;
26037e1895aSJed Brown   S2   = snes->work[1];
26137e1895aSJed Brown   S3   = snes->work[2];
26237e1895aSJed Brown   ierr = VecZeroEntries(S2);CHKERRQ(ierr);
26337e1895aSJed Brown   ierr = VecCopy(X,S3);CHKERRQ(ierr);
26437e1895aSJed Brown   for (i=0; i<nstages; i++) {
265da80777bSKarl Rupp     Vec         Ss[4];
266da80777bSKarl Rupp     PetscScalar scoeff[4];
267da80777bSKarl Rupp 
268da80777bSKarl Rupp     Ss[0] = S1; Ss[1] = S2; Ss[2] = S3; Ss[3] = Y;
2691aa26658SKarl Rupp 
270da80777bSKarl Rupp     scoeff[0] = gamma[0*nstages+i]-(PetscReal)1.0;
271da80777bSKarl Rupp     scoeff[1] = gamma[1*nstages+i];
272da80777bSKarl Rupp     scoeff[2] = gamma[2*nstages+i];
273da80777bSKarl Rupp     scoeff[3] = -betasub[i]*ms->damping;
2741aa26658SKarl Rupp 
27537e1895aSJed Brown     ierr = VecAXPY(S2,delta[i],S1);CHKERRQ(ierr);
27637e1895aSJed Brown     if (i>0) {
27737e1895aSJed Brown       ierr = SNESComputeFunction(snes,S1,F);CHKERRQ(ierr);
27837e1895aSJed Brown       if (snes->domainerror) {
27937e1895aSJed Brown         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
28037e1895aSJed Brown         PetscFunctionReturn(0);
28137e1895aSJed Brown       }
28237e1895aSJed Brown     }
283d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
28437e1895aSJed Brown     ierr = VecMAXPY(S1,4,scoeff,Ss);CHKERRQ(ierr);
28537e1895aSJed Brown   }
28637e1895aSJed Brown   PetscFunctionReturn(0);
28737e1895aSJed Brown }
28837e1895aSJed Brown 
28937e1895aSJed Brown #undef __FUNCT__
29037e1895aSJed Brown #define __FUNCT__ "SNESSolve_MS"
29137e1895aSJed Brown static PetscErrorCode SNESSolve_MS(SNES snes)
29237e1895aSJed Brown {
29337e1895aSJed Brown   SNES_MS        *ms = (SNES_MS*)snes->data;
29437e1895aSJed Brown   Vec            X   = snes->vec_sol,F = snes->vec_func;
29537e1895aSJed Brown   PetscReal      fnorm;
29637e1895aSJed Brown   MatStructure   mstruct;
29737e1895aSJed Brown   PetscInt       i;
29837e1895aSJed Brown   PetscErrorCode ierr;
29937e1895aSJed Brown 
30037e1895aSJed Brown   PetscFunctionBegin;
30137e1895aSJed Brown   snes->reason = SNES_CONVERGED_ITERATING;
302e04113cfSBarry Smith   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
30337e1895aSJed Brown   snes->iter   = 0;
30437e1895aSJed Brown   snes->norm   = 0.;
305e04113cfSBarry Smith   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
306e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
30737e1895aSJed Brown     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
30837e1895aSJed Brown     if (snes->domainerror) {
30937e1895aSJed Brown       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
31037e1895aSJed Brown       PetscFunctionReturn(0);
31137e1895aSJed Brown     }
3121aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
3131aa26658SKarl Rupp 
31437e1895aSJed Brown   if (snes->jacobian) {         /* This method does not require a Jacobian, but it is usually preconditioned by PBJacobi */
31537e1895aSJed Brown     ierr = SNESComputeJacobian(snes,snes->vec_sol,&snes->jacobian,&snes->jacobian_pre,&mstruct);CHKERRQ(ierr);
31637e1895aSJed Brown     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,mstruct);CHKERRQ(ierr);
31737e1895aSJed Brown   }
31837e1895aSJed Brown   if (ms->norms) {
31937e1895aSJed Brown     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
320189a9710SBarry Smith     if (PetscIsInfOrNanReal(fnorm)) {
321189a9710SBarry Smith       snes->reason = SNES_DIVERGED_FNORM_NAN;
322189a9710SBarry Smith       PetscFunctionReturn(0);
323189a9710SBarry Smith     }
32437e1895aSJed Brown     /* Monitor convergence */
325e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
32637e1895aSJed Brown     snes->iter = 0;
32737e1895aSJed Brown     snes->norm = fnorm;
328e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
329a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
33037e1895aSJed Brown     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
33137e1895aSJed Brown 
33237e1895aSJed Brown     /* Test for convergence */
33337e1895aSJed Brown     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
33437e1895aSJed Brown     if (snes->reason) PetscFunctionReturn(0);
33537e1895aSJed Brown   }
33637e1895aSJed Brown 
33737e1895aSJed Brown   /* Call general purpose update function */
33837e1895aSJed Brown   if (snes->ops->update) {
33937e1895aSJed Brown     ierr = (*snes->ops->update)(snes,snes->iter);CHKERRQ(ierr);
34037e1895aSJed Brown   }
34137e1895aSJed Brown   for (i = 0; i < snes->max_its; i++) {
34237e1895aSJed Brown     ierr = SNESMSStep_3Sstar(snes,X,F);CHKERRQ(ierr);
34337e1895aSJed Brown 
34437e1895aSJed Brown     if (i+1 < snes->max_its || ms->norms) {
34537e1895aSJed Brown       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
34637e1895aSJed Brown       if (snes->domainerror) {
34737e1895aSJed Brown         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
34837e1895aSJed Brown         PetscFunctionReturn(0);
34937e1895aSJed Brown       }
35037e1895aSJed Brown     }
35137e1895aSJed Brown 
35237e1895aSJed Brown     if (ms->norms) {
35337e1895aSJed Brown       ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
354189a9710SBarry Smith       if (PetscIsInfOrNanReal(fnorm)) {
355189a9710SBarry Smith         snes->reason = SNES_DIVERGED_FNORM_NAN;
356189a9710SBarry Smith         PetscFunctionReturn(0);
357189a9710SBarry Smith       }
358189a9710SBarry Smith 
35937e1895aSJed Brown       /* Monitor convergence */
360e04113cfSBarry Smith       ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
36137e1895aSJed Brown       snes->iter = i+1;
36237e1895aSJed Brown       snes->norm = fnorm;
363e04113cfSBarry Smith       ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
364a71f0d7dSBarry Smith       ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
36537e1895aSJed Brown       ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
36637e1895aSJed Brown 
36737e1895aSJed Brown       /* Test for convergence */
36837e1895aSJed Brown       ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
36937e1895aSJed Brown       if (snes->reason) PetscFunctionReturn(0);
37037e1895aSJed Brown     }
37137e1895aSJed Brown 
37237e1895aSJed Brown     /* Call general purpose update function */
37337e1895aSJed Brown     if (snes->ops->update) {
37437e1895aSJed Brown       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
37537e1895aSJed Brown     }
37637e1895aSJed Brown   }
37737e1895aSJed Brown   if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
37837e1895aSJed Brown   PetscFunctionReturn(0);
37937e1895aSJed Brown }
38037e1895aSJed Brown 
38137e1895aSJed Brown #undef __FUNCT__
38237e1895aSJed Brown #define __FUNCT__ "SNESSetUp_MS"
38337e1895aSJed Brown static PetscErrorCode SNESSetUp_MS(SNES snes)
38437e1895aSJed Brown {
38537e1895aSJed Brown   SNES_MS        *ms = (SNES_MS*)snes->data;
38637e1895aSJed Brown   PetscErrorCode ierr;
38737e1895aSJed Brown 
38837e1895aSJed Brown   PetscFunctionBegin;
38937e1895aSJed Brown   if (!ms->tableau) {ierr = SNESMSSetType(snes,SNESMSDefault);CHKERRQ(ierr);}
390fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr);
3916cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
39237e1895aSJed Brown   PetscFunctionReturn(0);
39337e1895aSJed Brown }
39437e1895aSJed Brown 
39537e1895aSJed Brown #undef __FUNCT__
39637e1895aSJed Brown #define __FUNCT__ "SNESReset_MS"
39737e1895aSJed Brown static PetscErrorCode SNESReset_MS(SNES snes)
39837e1895aSJed Brown {
39937e1895aSJed Brown 
40037e1895aSJed Brown   PetscFunctionBegin;
40137e1895aSJed Brown   PetscFunctionReturn(0);
40237e1895aSJed Brown }
40337e1895aSJed Brown 
40437e1895aSJed Brown #undef __FUNCT__
40537e1895aSJed Brown #define __FUNCT__ "SNESDestroy_MS"
40637e1895aSJed Brown static PetscErrorCode SNESDestroy_MS(SNES snes)
40737e1895aSJed Brown {
40837e1895aSJed Brown   PetscErrorCode ierr;
40937e1895aSJed Brown 
41037e1895aSJed Brown   PetscFunctionBegin;
41137e1895aSJed Brown   ierr = PetscFree(snes->data);CHKERRQ(ierr);
412bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr);
41337e1895aSJed Brown   PetscFunctionReturn(0);
41437e1895aSJed Brown }
41537e1895aSJed Brown 
41637e1895aSJed Brown #undef __FUNCT__
41737e1895aSJed Brown #define __FUNCT__ "SNESView_MS"
41837e1895aSJed Brown static PetscErrorCode SNESView_MS(SNES snes,PetscViewer viewer)
41937e1895aSJed Brown {
42037e1895aSJed Brown   PetscBool      iascii;
42137e1895aSJed Brown   PetscErrorCode ierr;
42237e1895aSJed Brown   SNES_MS        *ms = (SNES_MS*)snes->data;
42337e1895aSJed Brown 
42437e1895aSJed Brown   PetscFunctionBegin;
425251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
42637e1895aSJed Brown   if (iascii) {
42737e1895aSJed Brown     SNESMSTableau tab = ms->tableau;
42837e1895aSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  multi-stage method type: %s\n",tab ? tab->name : "not yet set");CHKERRQ(ierr);
42937e1895aSJed Brown   }
43037e1895aSJed Brown   PetscFunctionReturn(0);
43137e1895aSJed Brown }
43237e1895aSJed Brown 
43337e1895aSJed Brown #undef __FUNCT__
43437e1895aSJed Brown #define __FUNCT__ "SNESSetFromOptions_MS"
43537e1895aSJed Brown static PetscErrorCode SNESSetFromOptions_MS(SNES snes)
43637e1895aSJed Brown {
437725b86d8SJed Brown   SNES_MS        *ms = (SNES_MS*)snes->data;
43837e1895aSJed Brown   PetscErrorCode ierr;
43937e1895aSJed Brown 
44037e1895aSJed Brown   PetscFunctionBegin;
44137e1895aSJed Brown   ierr = PetscOptionsHead("SNES MS options");CHKERRQ(ierr);
44237e1895aSJed Brown   {
44337e1895aSJed Brown     SNESMSTableauLink link;
44437e1895aSJed Brown     PetscInt          count,choice;
44537e1895aSJed Brown     PetscBool         flg;
44637e1895aSJed Brown     const char        **namelist;
44737e1895aSJed Brown     char              mstype[256];
44837e1895aSJed Brown 
4498caf3d72SBarry Smith     ierr = PetscStrncpy(mstype,SNESMSDefault,sizeof(mstype));CHKERRQ(ierr);
45037e1895aSJed Brown     for (link=SNESMSTableauList,count=0; link; link=link->next,count++) ;
45137e1895aSJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
45237e1895aSJed Brown     for (link=SNESMSTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
45337e1895aSJed Brown     ierr = PetscOptionsEList("-snes_ms_type","Multistage smoother type","SNESMSSetType",(const char*const*)namelist,count,mstype,&choice,&flg);CHKERRQ(ierr);
45437e1895aSJed Brown     ierr = SNESMSSetType(snes,flg ? namelist[choice] : mstype);CHKERRQ(ierr);
45537e1895aSJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
4560298fd71SBarry Smith     ierr = PetscOptionsReal("-snes_ms_damping","Damping for multistage method","SNESMSSetDamping",ms->damping,&ms->damping,NULL);CHKERRQ(ierr);
4570298fd71SBarry Smith     ierr = PetscOptionsBool("-snes_ms_norms","Compute norms for monitoring","none",ms->norms,&ms->norms,NULL);CHKERRQ(ierr);
45837e1895aSJed Brown   }
45937e1895aSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
46037e1895aSJed Brown   PetscFunctionReturn(0);
46137e1895aSJed Brown }
46237e1895aSJed Brown 
46337e1895aSJed Brown #undef __FUNCT__
46437e1895aSJed Brown #define __FUNCT__ "SNESMSSetType_MS"
46537e1895aSJed Brown PetscErrorCode  SNESMSSetType_MS(SNES snes,SNESMSType mstype)
46637e1895aSJed Brown {
46737e1895aSJed Brown   PetscErrorCode    ierr;
46837e1895aSJed Brown   SNES_MS           *ms = (SNES_MS*)snes->data;
46937e1895aSJed Brown   SNESMSTableauLink link;
47037e1895aSJed Brown   PetscBool         match;
47137e1895aSJed Brown 
47237e1895aSJed Brown   PetscFunctionBegin;
47337e1895aSJed Brown   if (ms->tableau) {
47437e1895aSJed Brown     ierr = PetscStrcmp(ms->tableau->name,mstype,&match);CHKERRQ(ierr);
47537e1895aSJed Brown     if (match) PetscFunctionReturn(0);
47637e1895aSJed Brown   }
47737e1895aSJed Brown   for (link = SNESMSTableauList; link; link=link->next) {
47837e1895aSJed Brown     ierr = PetscStrcmp(link->tab.name,mstype,&match);CHKERRQ(ierr);
47937e1895aSJed Brown     if (match) {
48037e1895aSJed Brown       ierr        = SNESReset_MS(snes);CHKERRQ(ierr);
48137e1895aSJed Brown       ms->tableau = &link->tab;
48237e1895aSJed Brown       PetscFunctionReturn(0);
48337e1895aSJed Brown     }
48437e1895aSJed Brown   }
485ce94432eSBarry Smith   SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mstype);
48637e1895aSJed Brown   PetscFunctionReturn(0);
48737e1895aSJed Brown }
48837e1895aSJed Brown 
48937e1895aSJed Brown #undef __FUNCT__
49037e1895aSJed Brown #define __FUNCT__ "SNESMSSetType"
49137e1895aSJed Brown /*@C
49237e1895aSJed Brown   SNESMSSetType - Set the type of multistage smoother
49337e1895aSJed Brown 
49437e1895aSJed Brown   Logically collective
49537e1895aSJed Brown 
49637e1895aSJed Brown   Input Parameter:
49737e1895aSJed Brown +  snes - nonlinear solver context
49837e1895aSJed Brown -  mstype - type of multistage method
49937e1895aSJed Brown 
50037e1895aSJed Brown   Level: beginner
50137e1895aSJed Brown 
50237e1895aSJed Brown .seealso: SNESMSGetType(), SNESMS
50337e1895aSJed Brown @*/
50419fd82e9SBarry Smith PetscErrorCode SNESMSSetType(SNES snes,SNESMSType rostype)
50537e1895aSJed Brown {
50637e1895aSJed Brown   PetscErrorCode ierr;
50737e1895aSJed Brown 
50837e1895aSJed Brown   PetscFunctionBegin;
50937e1895aSJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
51019fd82e9SBarry Smith   ierr = PetscTryMethod(snes,"SNESMSSetType_C",(SNES,SNESMSType),(snes,rostype));CHKERRQ(ierr);
51137e1895aSJed Brown   PetscFunctionReturn(0);
51237e1895aSJed Brown }
51337e1895aSJed Brown 
51437e1895aSJed Brown /* -------------------------------------------------------------------------- */
51537e1895aSJed Brown /*MC
51637e1895aSJed Brown       SNESMS - multi-stage smoothers
51737e1895aSJed Brown 
51837e1895aSJed Brown       Options Database:
51937e1895aSJed Brown 
52037e1895aSJed Brown +     -snes_ms_type - type of multi-stage smoother
52137e1895aSJed Brown -     -snes_ms_damping - damping for multi-stage method
52237e1895aSJed Brown 
52337e1895aSJed Brown       Notes:
52437e1895aSJed Brown       These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for
5256c9de887SHong Zhang       FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebyshev).
52637e1895aSJed Brown 
52737e1895aSJed Brown       Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds.
52837e1895aSJed Brown 
52937e1895aSJed Brown       The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with SNESMSRegister().
53037e1895aSJed Brown 
53137e1895aSJed Brown       References:
53237e1895aSJed Brown 
53337e1895aSJed Brown       Ketcheson (2010) Runge-Kutta methods with minimum storage implementations.
53437e1895aSJed Brown 
53537e1895aSJed Brown       Jameson (1983) Solution of the Euler equations for two dimensional transonic flow by a multigrid method.
53637e1895aSJed Brown 
53737e1895aSJed Brown       Pierce and Giles (1997) Preconditioned multigrid methods for compressible flow calculations on stretched meshes.
53837e1895aSJed Brown 
53937e1895aSJed Brown       Level: beginner
54037e1895aSJed Brown 
5416c9de887SHong Zhang .seealso:  SNESCreate(), SNES, SNESSetType(), SNESMS, SNESFAS, KSPCHEBYSHEV
54237e1895aSJed Brown 
54337e1895aSJed Brown M*/
54437e1895aSJed Brown #undef __FUNCT__
54537e1895aSJed Brown #define __FUNCT__ "SNESCreate_MS"
5468cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_MS(SNES snes)
54737e1895aSJed Brown {
54837e1895aSJed Brown   PetscErrorCode ierr;
54937e1895aSJed Brown   SNES_MS        *ms;
55037e1895aSJed Brown 
55137e1895aSJed Brown   PetscFunctionBegin;
552607a6623SBarry Smith   ierr = SNESMSInitializePackage();CHKERRQ(ierr);
55337e1895aSJed Brown 
55437e1895aSJed Brown   snes->ops->setup          = SNESSetUp_MS;
55537e1895aSJed Brown   snes->ops->solve          = SNESSolve_MS;
55637e1895aSJed Brown   snes->ops->destroy        = SNESDestroy_MS;
55737e1895aSJed Brown   snes->ops->setfromoptions = SNESSetFromOptions_MS;
55837e1895aSJed Brown   snes->ops->view           = SNESView_MS;
55937e1895aSJed Brown   snes->ops->reset          = SNESReset_MS;
56037e1895aSJed Brown 
56137e1895aSJed Brown   snes->usespc  = PETSC_FALSE;
56237e1895aSJed Brown   snes->usesksp = PETSC_TRUE;
56337e1895aSJed Brown 
56437e1895aSJed Brown   ierr        = PetscNewLog(snes,SNES_MS,&ms);CHKERRQ(ierr);
56537e1895aSJed Brown   snes->data  = (void*)ms;
56637e1895aSJed Brown   ms->damping = 0.9;
56737e1895aSJed Brown   ms->norms   = PETSC_FALSE;
56837e1895aSJed Brown 
569bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetType_C",SNESMSSetType_MS);CHKERRQ(ierr);
57037e1895aSJed Brown   PetscFunctionReturn(0);
57137e1895aSJed Brown }
57299e0435eSBarry Smith 
573