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 } 71*3847c725SLisandro Dalcin 72*3847c725SLisandro Dalcin { /* Van Leer, Tai, and Powell (1989) 1 stage, order 1 */ 73*3847c725SLisandro Dalcin PetscReal gamma[3][1] = {{0},{0},{1}}; 74*3847c725SLisandro Dalcin PetscReal delta[1] = {0}; 75*3847c725SLisandro Dalcin PetscReal betasub[1] = {1.0}; 76*3847c725SLisandro Dalcin ierr = SNESMSRegister(SNESMSVLTP11,1,3,0.5,&gamma[0][0],delta,betasub);CHKERRQ(ierr); 77*3847c725SLisandro Dalcin } 78a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 2 stage, order 1 */ 79a97cb6bcSJed Brown PetscReal gamma[3][2] = {{0,0},{0,0},{1,1}}; 80a97cb6bcSJed Brown PetscReal delta[2] = {0,0}; 81a97cb6bcSJed Brown PetscReal betasub[2] = {0.3333,1.0}; 82a97cb6bcSJed Brown ierr = SNESMSRegister(SNESMSVLTP21,2,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr); 83a97cb6bcSJed Brown } 84a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 3 stage, order 1 */ 85a97cb6bcSJed Brown PetscReal gamma[3][3] = {{0,0,0},{0,0,0},{1,1,1}}; 86a97cb6bcSJed Brown PetscReal delta[3] = {0,0,0}; 87a97cb6bcSJed Brown PetscReal betasub[3] = {0.1481,0.4000,1.0}; 88a97cb6bcSJed Brown ierr = SNESMSRegister(SNESMSVLTP31,3,3,1.5,&gamma[0][0],delta,betasub);CHKERRQ(ierr); 89a97cb6bcSJed Brown } 90a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 4 stage, order 1 */ 91a97cb6bcSJed Brown PetscReal gamma[3][4] = {{0,0,0,0},{0,0,0,0},{1,1,1,1}}; 92a97cb6bcSJed Brown PetscReal delta[4] = {0,0,0,0}; 93a97cb6bcSJed Brown PetscReal betasub[4] = {0.0833,0.2069,0.4265,1.0}; 94a97cb6bcSJed Brown ierr = SNESMSRegister(SNESMSVLTP41,4,3,2.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr); 95a97cb6bcSJed Brown } 96a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 5 stage, order 1 */ 97a97cb6bcSJed Brown PetscReal gamma[3][5] = {{0,0,0,0,0},{0,0,0,0,0},{1,1,1,1,1}}; 98a97cb6bcSJed Brown PetscReal delta[5] = {0,0,0,0,0}; 99a97cb6bcSJed Brown PetscReal betasub[5] = {0.0533,0.1263,0.2375,0.4414,1.0}; 100a97cb6bcSJed Brown ierr = SNESMSRegister(SNESMSVLTP51,5,3,2.5,&gamma[0][0],delta,betasub);CHKERRQ(ierr); 101a97cb6bcSJed Brown } 102a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 6 stage, order 1 */ 103a97cb6bcSJed Brown PetscReal gamma[3][6] = {{0,0,0,0,0,0},{0,0,0,0,0,0},{1,1,1,1,1,1}}; 104a97cb6bcSJed Brown PetscReal delta[6] = {0,0,0,0,0,0}; 105a97cb6bcSJed Brown PetscReal betasub[6] = {0.0370,0.0851,0.1521,0.2562,0.4512,1.0}; 106a97cb6bcSJed Brown ierr = SNESMSRegister(SNESMSVLTP61,6,3,3.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr); 107a97cb6bcSJed Brown } 10837e1895aSJed Brown PetscFunctionReturn(0); 10937e1895aSJed Brown } 11037e1895aSJed Brown 11137e1895aSJed Brown /*@C 11257715debSLisandro Dalcin SNESMSRegisterDestroy - Frees the list of schemes that were registered by SNESMSRegister(). 11337e1895aSJed Brown 11437e1895aSJed Brown Not Collective 11537e1895aSJed Brown 11637e1895aSJed Brown Level: advanced 11737e1895aSJed Brown 11857715debSLisandro Dalcin .seealso: SNESMSRegister(), SNESMSRegisterAll() 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 /*@C 13937e1895aSJed Brown SNESMSInitializePackage - This function initializes everything in the SNESMS package. It is called 1408a690491SBarry Smith from SNESInitializePackage(). 14137e1895aSJed Brown 14237e1895aSJed Brown Level: developer 14337e1895aSJed Brown 14437e1895aSJed Brown .seealso: PetscInitialize() 14537e1895aSJed Brown @*/ 146607a6623SBarry Smith PetscErrorCode SNESMSInitializePackage(void) 14737e1895aSJed Brown { 14837e1895aSJed Brown PetscErrorCode ierr; 14937e1895aSJed Brown 15037e1895aSJed Brown PetscFunctionBegin; 15137e1895aSJed Brown if (SNESMSPackageInitialized) PetscFunctionReturn(0); 15237e1895aSJed Brown SNESMSPackageInitialized = PETSC_TRUE; 1531aa26658SKarl Rupp 15437e1895aSJed Brown ierr = SNESMSRegisterAll();CHKERRQ(ierr); 15537e1895aSJed Brown ierr = PetscRegisterFinalize(SNESMSFinalizePackage);CHKERRQ(ierr); 15637e1895aSJed Brown PetscFunctionReturn(0); 15737e1895aSJed Brown } 15837e1895aSJed Brown 15937e1895aSJed Brown /*@C 16037e1895aSJed Brown SNESMSFinalizePackage - This function destroys everything in the SNESMS package. It is 16137e1895aSJed Brown called from PetscFinalize(). 16237e1895aSJed Brown 16337e1895aSJed Brown Level: developer 16437e1895aSJed Brown 16537e1895aSJed Brown .seealso: PetscFinalize() 16637e1895aSJed Brown @*/ 16737e1895aSJed Brown PetscErrorCode SNESMSFinalizePackage(void) 16837e1895aSJed Brown { 16937e1895aSJed Brown PetscErrorCode ierr; 17037e1895aSJed Brown 17137e1895aSJed Brown PetscFunctionBegin; 17237e1895aSJed Brown SNESMSPackageInitialized = PETSC_FALSE; 1731aa26658SKarl Rupp 17437e1895aSJed Brown ierr = SNESMSRegisterDestroy();CHKERRQ(ierr); 17537e1895aSJed Brown PetscFunctionReturn(0); 17637e1895aSJed Brown } 17737e1895aSJed Brown 17837e1895aSJed Brown /*@C 17937e1895aSJed Brown SNESMSRegister - register a multistage scheme 18037e1895aSJed Brown 18137e1895aSJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 18237e1895aSJed Brown 18337e1895aSJed Brown Input Parameters: 18437e1895aSJed Brown + name - identifier for method 18537e1895aSJed Brown . nstages - number of stages 18637e1895aSJed Brown . nregisters - number of registers used by low-storage implementation 18737e1895aSJed Brown . gamma - coefficients, see Ketcheson's paper 18837e1895aSJed Brown . delta - coefficients, see Ketcheson's paper 18937e1895aSJed Brown - betasub - subdiagonal of Shu-Osher form 19037e1895aSJed Brown 19137e1895aSJed Brown Notes: 19237e1895aSJed Brown The notation is described in Ketcheson (2010) Runge-Kutta methods with minimum storage implementations. 19337e1895aSJed Brown 19437e1895aSJed Brown Level: advanced 19537e1895aSJed Brown 19637e1895aSJed Brown .seealso: SNESMS 19737e1895aSJed Brown @*/ 19819fd82e9SBarry Smith PetscErrorCode SNESMSRegister(SNESMSType name,PetscInt nstages,PetscInt nregisters,PetscReal stability,const PetscReal gamma[],const PetscReal delta[],const PetscReal betasub[]) 19937e1895aSJed Brown { 20037e1895aSJed Brown PetscErrorCode ierr; 20137e1895aSJed Brown SNESMSTableauLink link; 20237e1895aSJed Brown SNESMSTableau t; 20337e1895aSJed Brown 20437e1895aSJed Brown PetscFunctionBegin; 20537e1895aSJed Brown PetscValidCharPointer(name,1); 20637e1895aSJed Brown if (nstages < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have at least one stage"); 20737e1895aSJed Brown if (nregisters != 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only support for methods written in 3-register form"); 208534a8f05SLisandro Dalcin PetscValidRealPointer(gamma,4); 209534a8f05SLisandro Dalcin PetscValidRealPointer(delta,5); 210534a8f05SLisandro Dalcin PetscValidRealPointer(betasub,6); 21137e1895aSJed Brown 2129cc31a68SJed Brown ierr = SNESMSInitializePackage();CHKERRQ(ierr); 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); 221580bdb30SBarry Smith ierr = PetscArraycpy(t->gamma,gamma,nstages*nregisters);CHKERRQ(ierr); 222580bdb30SBarry Smith ierr = PetscArraycpy(t->delta,delta,nstages);CHKERRQ(ierr); 223580bdb30SBarry Smith ierr = PetscArraycpy(t->betasub,betasub,nstages);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; 241bf80fd75SLisandro Dalcin 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 257bf80fd75SLisandro Dalcin scoeff[0] = gamma[0*nstages+i] - 1; 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 272bf80fd75SLisandro Dalcin static PetscErrorCode SNESMSStep_Norms(SNES snes,PetscInt iter,Vec F) 273bf80fd75SLisandro Dalcin { 274bf80fd75SLisandro Dalcin SNES_MS *ms = (SNES_MS*)snes->data; 275bf80fd75SLisandro Dalcin PetscReal fnorm; 276bf80fd75SLisandro Dalcin PetscErrorCode ierr; 277bf80fd75SLisandro Dalcin 278bf80fd75SLisandro Dalcin PetscFunctionBegin; 279bf80fd75SLisandro Dalcin if (ms->norms) { 280bf80fd75SLisandro Dalcin ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 281bf80fd75SLisandro Dalcin SNESCheckFunctionNorm(snes,fnorm); 282bf80fd75SLisandro Dalcin /* Monitor convergence */ 283bf80fd75SLisandro Dalcin ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 284bf80fd75SLisandro Dalcin snes->iter = iter; 285bf80fd75SLisandro Dalcin snes->norm = fnorm; 286bf80fd75SLisandro Dalcin ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 287bf80fd75SLisandro Dalcin ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 288bf80fd75SLisandro Dalcin ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 289bf80fd75SLisandro Dalcin /* Test for convergence */ 290bf80fd75SLisandro Dalcin ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 291bf80fd75SLisandro Dalcin } else if (iter > 0) { 292bf80fd75SLisandro Dalcin ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 293bf80fd75SLisandro Dalcin snes->iter = iter; 294bf80fd75SLisandro Dalcin ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 295bf80fd75SLisandro Dalcin } 296bf80fd75SLisandro Dalcin PetscFunctionReturn(0); 297bf80fd75SLisandro Dalcin } 298bf80fd75SLisandro Dalcin 299bf80fd75SLisandro Dalcin 30037e1895aSJed Brown static PetscErrorCode SNESSolve_MS(SNES snes) 30137e1895aSJed Brown { 30237e1895aSJed Brown SNES_MS *ms = (SNES_MS*)snes->data; 30337e1895aSJed Brown Vec X = snes->vec_sol,F = snes->vec_func; 30437e1895aSJed Brown PetscInt i; 30537e1895aSJed Brown PetscErrorCode ierr; 30637e1895aSJed Brown 30737e1895aSJed Brown PetscFunctionBegin; 3086c4ed002SBarry 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); 309fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 310bf80fd75SLisandro Dalcin 31137e1895aSJed Brown snes->reason = SNES_CONVERGED_ITERATING; 312e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 31337e1895aSJed Brown snes->iter = 0; 314bf80fd75SLisandro Dalcin snes->norm = 0; 315e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 316bf80fd75SLisandro Dalcin 317e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 31837e1895aSJed Brown ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 3191aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 3201aa26658SKarl Rupp 321bf80fd75SLisandro Dalcin ierr = SNESMSStep_Norms(snes,0,F);CHKERRQ(ierr); 32237e1895aSJed Brown if (snes->reason) PetscFunctionReturn(0); 323bf80fd75SLisandro Dalcin 324bf80fd75SLisandro Dalcin for (i = 0; i < snes->max_its; i++) { 32537e1895aSJed Brown 32637e1895aSJed Brown /* Call general purpose update function */ 32737e1895aSJed Brown if (snes->ops->update) { 32837e1895aSJed Brown ierr = (*snes->ops->update)(snes,snes->iter);CHKERRQ(ierr); 32937e1895aSJed Brown } 330bf80fd75SLisandro Dalcin 331bf80fd75SLisandro Dalcin if (i == 0 && snes->jacobian) { 332bf80fd75SLisandro Dalcin /* This method does not require a Jacobian, but it is usually preconditioned by PBJacobi */ 333bf80fd75SLisandro Dalcin ierr = SNESComputeJacobian(snes,snes->vec_sol,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 334bf80fd75SLisandro Dalcin SNESCheckJacobianDomainerror(snes); 335bf80fd75SLisandro Dalcin ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 336bf80fd75SLisandro Dalcin } 337bf80fd75SLisandro Dalcin 33837e1895aSJed Brown ierr = SNESMSStep_3Sstar(snes,X,F);CHKERRQ(ierr); 33937e1895aSJed Brown 340bf80fd75SLisandro Dalcin if (i < snes->max_its-1 || ms->norms) { 34137e1895aSJed Brown ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 34237e1895aSJed Brown } 34337e1895aSJed Brown 344bf80fd75SLisandro Dalcin ierr = SNESMSStep_Norms(snes,i+1,F);CHKERRQ(ierr); 34537e1895aSJed Brown if (snes->reason) PetscFunctionReturn(0); 34637e1895aSJed Brown } 34737e1895aSJed Brown if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; 34837e1895aSJed Brown PetscFunctionReturn(0); 34937e1895aSJed Brown } 35037e1895aSJed Brown 35137e1895aSJed Brown static PetscErrorCode SNESSetUp_MS(SNES snes) 35237e1895aSJed Brown { 35337e1895aSJed Brown PetscErrorCode ierr; 35437e1895aSJed Brown 35537e1895aSJed Brown PetscFunctionBegin; 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 PetscFunctionBegin; 36437e1895aSJed Brown PetscFunctionReturn(0); 36537e1895aSJed Brown } 36637e1895aSJed Brown 36737e1895aSJed Brown static PetscErrorCode SNESDestroy_MS(SNES snes) 36837e1895aSJed Brown { 36937e1895aSJed Brown PetscErrorCode ierr; 37037e1895aSJed Brown 37137e1895aSJed Brown PetscFunctionBegin; 372bf80fd75SLisandro Dalcin ierr = SNESReset_MS(snes);CHKERRQ(ierr); 37337e1895aSJed Brown ierr = PetscFree(snes->data);CHKERRQ(ierr); 37457715debSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetType_C",NULL);CHKERRQ(ierr); 375bf80fd75SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetType_C",NULL);CHKERRQ(ierr); 37657715debSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetDamping_C",NULL);CHKERRQ(ierr); 37757715debSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetDamping_C",NULL);CHKERRQ(ierr); 37837e1895aSJed Brown PetscFunctionReturn(0); 37937e1895aSJed Brown } 38037e1895aSJed Brown 38137e1895aSJed Brown static PetscErrorCode SNESView_MS(SNES snes,PetscViewer viewer) 38237e1895aSJed Brown { 38337e1895aSJed Brown PetscBool iascii; 38437e1895aSJed Brown PetscErrorCode ierr; 38537e1895aSJed Brown SNES_MS *ms = (SNES_MS*)snes->data; 38657715debSLisandro Dalcin SNESMSTableau tab = ms->tableau; 38737e1895aSJed Brown 38837e1895aSJed Brown PetscFunctionBegin; 389251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 39037e1895aSJed Brown if (iascii) { 39157715debSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," multi-stage method type: %s\n",tab->name);CHKERRQ(ierr); 39237e1895aSJed Brown } 39337e1895aSJed Brown PetscFunctionReturn(0); 39437e1895aSJed Brown } 39537e1895aSJed Brown 3964416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_MS(PetscOptionItems *PetscOptionsObject,SNES snes) 39737e1895aSJed Brown { 398725b86d8SJed Brown SNES_MS *ms = (SNES_MS*)snes->data; 39937e1895aSJed Brown PetscErrorCode ierr; 40037e1895aSJed Brown 40137e1895aSJed Brown PetscFunctionBegin; 402e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES MS options");CHKERRQ(ierr); 40337e1895aSJed Brown { 40437e1895aSJed Brown SNESMSTableauLink link; 40537e1895aSJed Brown PetscInt count,choice; 40637e1895aSJed Brown PetscBool flg; 40737e1895aSJed Brown const char **namelist; 40857715debSLisandro Dalcin SNESMSType mstype; 40957715debSLisandro Dalcin PetscReal damping; 41037e1895aSJed Brown 41157715debSLisandro Dalcin ierr = SNESMSGetType(snes,&mstype);CHKERRQ(ierr); 41237e1895aSJed Brown for (link=SNESMSTableauList,count=0; link; link=link->next,count++) ; 413f489ac74SBarry Smith ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr); 41437e1895aSJed Brown for (link=SNESMSTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 41537e1895aSJed Brown ierr = PetscOptionsEList("-snes_ms_type","Multistage smoother type","SNESMSSetType",(const char*const*)namelist,count,mstype,&choice,&flg);CHKERRQ(ierr); 41657715debSLisandro Dalcin if (flg) {ierr = SNESMSSetType(snes,namelist[choice]);CHKERRQ(ierr);} 41737e1895aSJed Brown ierr = PetscFree(namelist);CHKERRQ(ierr); 41857715debSLisandro Dalcin ierr = SNESMSGetDamping(snes,&damping);CHKERRQ(ierr); 41957715debSLisandro Dalcin ierr = PetscOptionsReal("-snes_ms_damping","Damping for multistage method","SNESMSSetDamping",damping,&damping,&flg);CHKERRQ(ierr); 42057715debSLisandro Dalcin if (flg) {ierr = SNESMSSetDamping(snes,damping);CHKERRQ(ierr);} 4210298fd71SBarry Smith ierr = PetscOptionsBool("-snes_ms_norms","Compute norms for monitoring","none",ms->norms,&ms->norms,NULL);CHKERRQ(ierr); 42237e1895aSJed Brown } 42337e1895aSJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 42437e1895aSJed Brown PetscFunctionReturn(0); 42537e1895aSJed Brown } 42637e1895aSJed Brown 42757715debSLisandro Dalcin static PetscErrorCode SNESMSGetType_MS(SNES snes,SNESMSType *mstype) 42837e1895aSJed Brown { 42957715debSLisandro Dalcin SNES_MS *ms = (SNES_MS*)snes->data; 43057715debSLisandro Dalcin SNESMSTableau tab = ms->tableau; 43157715debSLisandro Dalcin 43257715debSLisandro Dalcin PetscFunctionBegin; 43357715debSLisandro Dalcin *mstype = tab->name; 43457715debSLisandro Dalcin PetscFunctionReturn(0); 43557715debSLisandro Dalcin } 43657715debSLisandro Dalcin 43757715debSLisandro Dalcin static PetscErrorCode SNESMSSetType_MS(SNES snes,SNESMSType mstype) 43857715debSLisandro Dalcin { 43937e1895aSJed Brown SNES_MS *ms = (SNES_MS*)snes->data; 44037e1895aSJed Brown SNESMSTableauLink link; 44137e1895aSJed Brown PetscBool match; 44257715debSLisandro Dalcin PetscErrorCode ierr; 44337e1895aSJed Brown 44437e1895aSJed Brown PetscFunctionBegin; 44537e1895aSJed Brown if (ms->tableau) { 44637e1895aSJed Brown ierr = PetscStrcmp(ms->tableau->name,mstype,&match);CHKERRQ(ierr); 44737e1895aSJed Brown if (match) PetscFunctionReturn(0); 44837e1895aSJed Brown } 44937e1895aSJed Brown for (link = SNESMSTableauList; link; link=link->next) { 45037e1895aSJed Brown ierr = PetscStrcmp(link->tab.name,mstype,&match);CHKERRQ(ierr); 45137e1895aSJed Brown if (match) { 45257715debSLisandro Dalcin if (snes->setupcalled) {ierr = SNESReset_MS(snes);CHKERRQ(ierr);} 45337e1895aSJed Brown ms->tableau = &link->tab; 45457715debSLisandro Dalcin if (snes->setupcalled) {ierr = SNESSetUp_MS(snes);CHKERRQ(ierr);} 45537e1895aSJed Brown PetscFunctionReturn(0); 45637e1895aSJed Brown } 45737e1895aSJed Brown } 458ce94432eSBarry Smith SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mstype); 45937e1895aSJed Brown PetscFunctionReturn(0); 46037e1895aSJed Brown } 46137e1895aSJed Brown 46237e1895aSJed Brown /*@C 46357715debSLisandro Dalcin SNESMSGetType - Get the type of multistage smoother 46457715debSLisandro Dalcin 46557715debSLisandro Dalcin Not collective 46657715debSLisandro Dalcin 46757715debSLisandro Dalcin Input Parameter: 46857715debSLisandro Dalcin . snes - nonlinear solver context 46957715debSLisandro Dalcin 47057715debSLisandro Dalcin Output Parameter: 47157715debSLisandro Dalcin . mstype - type of multistage method 47257715debSLisandro Dalcin 47357715debSLisandro Dalcin Level: beginner 47457715debSLisandro Dalcin 47557715debSLisandro Dalcin .seealso: SNESMSSetType(), SNESMSType, SNESMS 47657715debSLisandro Dalcin @*/ 47757715debSLisandro Dalcin PetscErrorCode SNESMSGetType(SNES snes,SNESMSType *mstype) 47857715debSLisandro Dalcin { 47957715debSLisandro Dalcin PetscErrorCode ierr; 48057715debSLisandro Dalcin 48157715debSLisandro Dalcin PetscFunctionBegin; 48257715debSLisandro Dalcin PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 48357715debSLisandro Dalcin PetscValidPointer(mstype,2); 48457715debSLisandro Dalcin ierr = PetscUseMethod(snes,"SNESMSGetType_C",(SNES,SNESMSType*),(snes,mstype));CHKERRQ(ierr); 48557715debSLisandro Dalcin PetscFunctionReturn(0); 48657715debSLisandro Dalcin } 48757715debSLisandro Dalcin 48857715debSLisandro Dalcin /*@C 48937e1895aSJed Brown SNESMSSetType - Set the type of multistage smoother 49037e1895aSJed Brown 49137e1895aSJed Brown Logically collective 49237e1895aSJed Brown 49337e1895aSJed Brown Input Parameter: 49437e1895aSJed Brown + snes - nonlinear solver context 49537e1895aSJed Brown - mstype - type of multistage method 49637e1895aSJed Brown 49737e1895aSJed Brown Level: beginner 49837e1895aSJed Brown 49957715debSLisandro Dalcin .seealso: SNESMSGetType(), SNESMSType, SNESMS 50037e1895aSJed Brown @*/ 50157715debSLisandro Dalcin PetscErrorCode SNESMSSetType(SNES snes,SNESMSType mstype) 50237e1895aSJed Brown { 50337e1895aSJed Brown PetscErrorCode ierr; 50437e1895aSJed Brown 50537e1895aSJed Brown PetscFunctionBegin; 50637e1895aSJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 50757715debSLisandro Dalcin PetscValidPointer(mstype,2); 50857715debSLisandro Dalcin ierr = PetscTryMethod(snes,"SNESMSSetType_C",(SNES,SNESMSType),(snes,mstype));CHKERRQ(ierr); 50957715debSLisandro Dalcin PetscFunctionReturn(0); 51057715debSLisandro Dalcin } 51157715debSLisandro Dalcin 51257715debSLisandro Dalcin static PetscErrorCode SNESMSGetDamping_MS(SNES snes,PetscReal *damping) 51357715debSLisandro Dalcin { 51457715debSLisandro Dalcin SNES_MS *ms = (SNES_MS*)snes->data; 51557715debSLisandro Dalcin 51657715debSLisandro Dalcin PetscFunctionBegin; 51757715debSLisandro Dalcin *damping = ms->damping; 51857715debSLisandro Dalcin PetscFunctionReturn(0); 51957715debSLisandro Dalcin } 52057715debSLisandro Dalcin 52157715debSLisandro Dalcin static PetscErrorCode SNESMSSetDamping_MS(SNES snes,PetscReal damping) 52257715debSLisandro Dalcin { 52357715debSLisandro Dalcin SNES_MS *ms = (SNES_MS*)snes->data; 52457715debSLisandro Dalcin 52557715debSLisandro Dalcin PetscFunctionBegin; 52657715debSLisandro Dalcin ms->damping = damping; 52757715debSLisandro Dalcin PetscFunctionReturn(0); 52857715debSLisandro Dalcin } 52957715debSLisandro Dalcin 53057715debSLisandro Dalcin /*@C 53157715debSLisandro Dalcin SNESMSGetDamping - Get the damping parameter 53257715debSLisandro Dalcin 53357715debSLisandro Dalcin Not collective 53457715debSLisandro Dalcin 53557715debSLisandro Dalcin Input Parameter: 53657715debSLisandro Dalcin . snes - nonlinear solver context 53757715debSLisandro Dalcin 53857715debSLisandro Dalcin Output Parameter: 53957715debSLisandro Dalcin . damping - damping parameter 54057715debSLisandro Dalcin 54157715debSLisandro Dalcin Level: advanced 54257715debSLisandro Dalcin 54357715debSLisandro Dalcin .seealso: SNESMSSetDamping(), SNESMS 54457715debSLisandro Dalcin @*/ 54557715debSLisandro Dalcin PetscErrorCode SNESMSGetDamping(SNES snes,PetscReal *damping) 54657715debSLisandro Dalcin { 54757715debSLisandro Dalcin PetscErrorCode ierr; 54857715debSLisandro Dalcin 54957715debSLisandro Dalcin PetscFunctionBegin; 55057715debSLisandro Dalcin PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 55157715debSLisandro Dalcin PetscValidPointer(damping,2); 55257715debSLisandro Dalcin ierr = PetscUseMethod(snes,"SNESMSGetDamping_C",(SNES,PetscReal*),(snes,damping));CHKERRQ(ierr); 55357715debSLisandro Dalcin PetscFunctionReturn(0); 55457715debSLisandro Dalcin } 55557715debSLisandro Dalcin 55657715debSLisandro Dalcin /*@C 55757715debSLisandro Dalcin SNESMSSetDamping - Set the damping parameter 55857715debSLisandro Dalcin 55957715debSLisandro Dalcin Logically collective 56057715debSLisandro Dalcin 56157715debSLisandro Dalcin Input Parameter: 56257715debSLisandro Dalcin + snes - nonlinear solver context 56357715debSLisandro Dalcin - damping - damping parameter 56457715debSLisandro Dalcin 56557715debSLisandro Dalcin Level: advanced 56657715debSLisandro Dalcin 56757715debSLisandro Dalcin .seealso: SNESMSGetDamping(), SNESMS 56857715debSLisandro Dalcin @*/ 56957715debSLisandro Dalcin PetscErrorCode SNESMSSetDamping(SNES snes,PetscReal damping) 57057715debSLisandro Dalcin { 57157715debSLisandro Dalcin PetscErrorCode ierr; 57257715debSLisandro Dalcin 57357715debSLisandro Dalcin PetscFunctionBegin; 57457715debSLisandro Dalcin PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 57557715debSLisandro Dalcin PetscValidLogicalCollectiveReal(snes,damping,2); 57657715debSLisandro Dalcin ierr = PetscTryMethod(snes,"SNESMSSetDamping_C",(SNES,PetscReal),(snes,damping));CHKERRQ(ierr); 57737e1895aSJed Brown PetscFunctionReturn(0); 57837e1895aSJed Brown } 57937e1895aSJed Brown 58037e1895aSJed Brown /* -------------------------------------------------------------------------- */ 58137e1895aSJed Brown /*MC 58237e1895aSJed Brown SNESMS - multi-stage smoothers 58337e1895aSJed Brown 58437e1895aSJed Brown Options Database: 58537e1895aSJed Brown 58637e1895aSJed Brown + -snes_ms_type - type of multi-stage smoother 58737e1895aSJed Brown - -snes_ms_damping - damping for multi-stage method 58837e1895aSJed Brown 58937e1895aSJed Brown Notes: 59037e1895aSJed Brown These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for 5916c9de887SHong Zhang FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebyshev). 59237e1895aSJed Brown 59337e1895aSJed Brown Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds. 59437e1895aSJed Brown 59537e1895aSJed Brown The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with SNESMSRegister(). 59637e1895aSJed Brown 59737e1895aSJed Brown References: 59857715debSLisandro Dalcin + 1. - Ketcheson (2010) Runge Kutta methods with minimum storage implementations (https://doi.org/10.1016/j.jcp.2009.11.006). 59957715debSLisandro Dalcin . 2. - Jameson (1983) Solution of the Euler equations for two dimensional transonic flow by a multigrid method (https://doi.org/10.1016/0096-3003(83)90019-X). 60057715debSLisandro Dalcin . 3. - Pierce and Giles (1997) Preconditioned multigrid methods for compressible flow calculations on stretched meshes (https://doi.org/10.1006/jcph.1997.5772). 60157715debSLisandro Dalcin - 4. - Van Leer, Tai, and Powell (1989) Design of optimally smoothing multi-stage schemes for the Euler equations (https://doi.org/10.2514/6.1989-1933). 60237e1895aSJed Brown 60337e1895aSJed Brown Level: beginner 60437e1895aSJed Brown 6056c9de887SHong Zhang .seealso: SNESCreate(), SNES, SNESSetType(), SNESMS, SNESFAS, KSPCHEBYSHEV 60637e1895aSJed Brown 60737e1895aSJed Brown M*/ 6088cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_MS(SNES snes) 60937e1895aSJed Brown { 61037e1895aSJed Brown PetscErrorCode ierr; 61137e1895aSJed Brown SNES_MS *ms; 61237e1895aSJed Brown 61337e1895aSJed Brown PetscFunctionBegin; 614607a6623SBarry Smith ierr = SNESMSInitializePackage();CHKERRQ(ierr); 61537e1895aSJed Brown 61637e1895aSJed Brown snes->ops->setup = SNESSetUp_MS; 61737e1895aSJed Brown snes->ops->solve = SNESSolve_MS; 61837e1895aSJed Brown snes->ops->destroy = SNESDestroy_MS; 61937e1895aSJed Brown snes->ops->setfromoptions = SNESSetFromOptions_MS; 62037e1895aSJed Brown snes->ops->view = SNESView_MS; 62137e1895aSJed Brown snes->ops->reset = SNESReset_MS; 62237e1895aSJed Brown 623efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 62437e1895aSJed Brown snes->usesksp = PETSC_TRUE; 62537e1895aSJed Brown 6264fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 6274fc747eaSLawrence Mitchell 628b00a9115SJed Brown ierr = PetscNewLog(snes,&ms);CHKERRQ(ierr); 62937e1895aSJed Brown snes->data = (void*)ms; 63037e1895aSJed Brown ms->damping = 0.9; 63137e1895aSJed Brown ms->norms = PETSC_FALSE; 63237e1895aSJed Brown 63357715debSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetType_C",SNESMSGetType_MS);CHKERRQ(ierr); 634bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetType_C",SNESMSSetType_MS);CHKERRQ(ierr); 63557715debSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetDamping_C",SNESMSGetDamping_MS);CHKERRQ(ierr); 63657715debSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetDamping_C",SNESMSSetDamping_MS);CHKERRQ(ierr); 63757715debSLisandro Dalcin 63857715debSLisandro Dalcin ierr = SNESMSSetType(snes,SNESMSDefault);CHKERRQ(ierr); 63937e1895aSJed Brown PetscFunctionReturn(0); 64037e1895aSJed Brown } 64199e0435eSBarry Smith 642