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 105*57715debSLisandro Dalcin SNESMSRegisterDestroy - Frees the list of schemes that were registered by SNESMSRegister(). 10637e1895aSJed Brown 10737e1895aSJed Brown Not Collective 10837e1895aSJed Brown 10937e1895aSJed Brown Level: advanced 11037e1895aSJed Brown 111*57715debSLisandro Dalcin .seealso: SNESMSRegister(), SNESMSRegisterAll() 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"); 201534a8f05SLisandro Dalcin PetscValidRealPointer(gamma,4); 202534a8f05SLisandro Dalcin PetscValidRealPointer(delta,5); 203534a8f05SLisandro Dalcin PetscValidRealPointer(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); 214580bdb30SBarry Smith ierr = PetscArraycpy(t->gamma,gamma,nstages*nregisters);CHKERRQ(ierr); 215580bdb30SBarry Smith ierr = PetscArraycpy(t->delta,delta,nstages);CHKERRQ(ierr); 216580bdb30SBarry 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; 234bf80fd75SLisandro Dalcin 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 250bf80fd75SLisandro Dalcin scoeff[0] = gamma[0*nstages+i] - 1; 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 265bf80fd75SLisandro Dalcin static PetscErrorCode SNESMSStep_Norms(SNES snes,PetscInt iter,Vec F) 266bf80fd75SLisandro Dalcin { 267bf80fd75SLisandro Dalcin SNES_MS *ms = (SNES_MS*)snes->data; 268bf80fd75SLisandro Dalcin PetscReal fnorm; 269bf80fd75SLisandro Dalcin PetscErrorCode ierr; 270bf80fd75SLisandro Dalcin 271bf80fd75SLisandro Dalcin PetscFunctionBegin; 272bf80fd75SLisandro Dalcin if (ms->norms) { 273bf80fd75SLisandro Dalcin ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 274bf80fd75SLisandro Dalcin SNESCheckFunctionNorm(snes,fnorm); 275bf80fd75SLisandro Dalcin /* Monitor convergence */ 276bf80fd75SLisandro Dalcin ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 277bf80fd75SLisandro Dalcin snes->iter = iter; 278bf80fd75SLisandro Dalcin snes->norm = fnorm; 279bf80fd75SLisandro Dalcin ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 280bf80fd75SLisandro Dalcin ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 281bf80fd75SLisandro Dalcin ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 282bf80fd75SLisandro Dalcin /* Test for convergence */ 283bf80fd75SLisandro Dalcin ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 284bf80fd75SLisandro Dalcin } else if (iter > 0) { 285bf80fd75SLisandro Dalcin ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 286bf80fd75SLisandro Dalcin snes->iter = iter; 287bf80fd75SLisandro Dalcin ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 288bf80fd75SLisandro Dalcin } 289bf80fd75SLisandro Dalcin PetscFunctionReturn(0); 290bf80fd75SLisandro Dalcin } 291bf80fd75SLisandro Dalcin 292bf80fd75SLisandro Dalcin 29337e1895aSJed Brown static PetscErrorCode SNESSolve_MS(SNES snes) 29437e1895aSJed Brown { 29537e1895aSJed Brown SNES_MS *ms = (SNES_MS*)snes->data; 29637e1895aSJed Brown Vec X = snes->vec_sol,F = snes->vec_func; 29737e1895aSJed Brown PetscInt i; 29837e1895aSJed Brown PetscErrorCode ierr; 29937e1895aSJed Brown 30037e1895aSJed Brown PetscFunctionBegin; 3016c4ed002SBarry 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); 302fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 303bf80fd75SLisandro Dalcin 30437e1895aSJed Brown snes->reason = SNES_CONVERGED_ITERATING; 305e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 30637e1895aSJed Brown snes->iter = 0; 307bf80fd75SLisandro Dalcin snes->norm = 0; 308e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 309bf80fd75SLisandro Dalcin 310e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 31137e1895aSJed Brown ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 3121aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 3131aa26658SKarl Rupp 314bf80fd75SLisandro Dalcin ierr = SNESMSStep_Norms(snes,0,F);CHKERRQ(ierr); 31537e1895aSJed Brown if (snes->reason) PetscFunctionReturn(0); 316bf80fd75SLisandro Dalcin 317bf80fd75SLisandro Dalcin for (i = 0; i < snes->max_its; i++) { 31837e1895aSJed Brown 31937e1895aSJed Brown /* Call general purpose update function */ 32037e1895aSJed Brown if (snes->ops->update) { 32137e1895aSJed Brown ierr = (*snes->ops->update)(snes,snes->iter);CHKERRQ(ierr); 32237e1895aSJed Brown } 323bf80fd75SLisandro Dalcin 324bf80fd75SLisandro Dalcin if (i == 0 && snes->jacobian) { 325bf80fd75SLisandro Dalcin /* This method does not require a Jacobian, but it is usually preconditioned by PBJacobi */ 326bf80fd75SLisandro Dalcin ierr = SNESComputeJacobian(snes,snes->vec_sol,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 327bf80fd75SLisandro Dalcin SNESCheckJacobianDomainerror(snes); 328bf80fd75SLisandro Dalcin ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 329bf80fd75SLisandro Dalcin } 330bf80fd75SLisandro Dalcin 33137e1895aSJed Brown ierr = SNESMSStep_3Sstar(snes,X,F);CHKERRQ(ierr); 33237e1895aSJed Brown 333bf80fd75SLisandro Dalcin if (i < snes->max_its-1 || ms->norms) { 33437e1895aSJed Brown ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 33537e1895aSJed Brown } 33637e1895aSJed Brown 337bf80fd75SLisandro Dalcin ierr = SNESMSStep_Norms(snes,i+1,F);CHKERRQ(ierr); 33837e1895aSJed Brown if (snes->reason) PetscFunctionReturn(0); 33937e1895aSJed Brown } 34037e1895aSJed Brown if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; 34137e1895aSJed Brown PetscFunctionReturn(0); 34237e1895aSJed Brown } 34337e1895aSJed Brown 34437e1895aSJed Brown static PetscErrorCode SNESSetUp_MS(SNES snes) 34537e1895aSJed Brown { 34637e1895aSJed Brown PetscErrorCode ierr; 34737e1895aSJed Brown 34837e1895aSJed Brown PetscFunctionBegin; 349fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr); 3506cab3a1bSJed Brown ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 35137e1895aSJed Brown PetscFunctionReturn(0); 35237e1895aSJed Brown } 35337e1895aSJed Brown 35437e1895aSJed Brown static PetscErrorCode SNESReset_MS(SNES snes) 35537e1895aSJed Brown { 35637e1895aSJed Brown PetscFunctionBegin; 35737e1895aSJed Brown PetscFunctionReturn(0); 35837e1895aSJed Brown } 35937e1895aSJed Brown 36037e1895aSJed Brown static PetscErrorCode SNESDestroy_MS(SNES snes) 36137e1895aSJed Brown { 36237e1895aSJed Brown PetscErrorCode ierr; 36337e1895aSJed Brown 36437e1895aSJed Brown PetscFunctionBegin; 365bf80fd75SLisandro Dalcin ierr = SNESReset_MS(snes);CHKERRQ(ierr); 36637e1895aSJed Brown ierr = PetscFree(snes->data);CHKERRQ(ierr); 367*57715debSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetType_C",NULL);CHKERRQ(ierr); 368bf80fd75SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetType_C",NULL);CHKERRQ(ierr); 369*57715debSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetDamping_C",NULL);CHKERRQ(ierr); 370*57715debSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetDamping_C",NULL);CHKERRQ(ierr); 37137e1895aSJed Brown PetscFunctionReturn(0); 37237e1895aSJed Brown } 37337e1895aSJed Brown 37437e1895aSJed Brown static PetscErrorCode SNESView_MS(SNES snes,PetscViewer viewer) 37537e1895aSJed Brown { 37637e1895aSJed Brown PetscBool iascii; 37737e1895aSJed Brown PetscErrorCode ierr; 37837e1895aSJed Brown SNES_MS *ms = (SNES_MS*)snes->data; 379*57715debSLisandro Dalcin SNESMSTableau tab = ms->tableau; 38037e1895aSJed Brown 38137e1895aSJed Brown PetscFunctionBegin; 382251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 38337e1895aSJed Brown if (iascii) { 384*57715debSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," multi-stage method type: %s\n",tab->name);CHKERRQ(ierr); 38537e1895aSJed Brown } 38637e1895aSJed Brown PetscFunctionReturn(0); 38737e1895aSJed Brown } 38837e1895aSJed Brown 3894416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_MS(PetscOptionItems *PetscOptionsObject,SNES snes) 39037e1895aSJed Brown { 391725b86d8SJed Brown SNES_MS *ms = (SNES_MS*)snes->data; 39237e1895aSJed Brown PetscErrorCode ierr; 39337e1895aSJed Brown 39437e1895aSJed Brown PetscFunctionBegin; 395e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES MS options");CHKERRQ(ierr); 39637e1895aSJed Brown { 39737e1895aSJed Brown SNESMSTableauLink link; 39837e1895aSJed Brown PetscInt count,choice; 39937e1895aSJed Brown PetscBool flg; 40037e1895aSJed Brown const char **namelist; 401*57715debSLisandro Dalcin SNESMSType mstype; 402*57715debSLisandro Dalcin PetscReal damping; 40337e1895aSJed Brown 404*57715debSLisandro Dalcin ierr = SNESMSGetType(snes,&mstype);CHKERRQ(ierr); 40537e1895aSJed Brown for (link=SNESMSTableauList,count=0; link; link=link->next,count++) ; 406f489ac74SBarry Smith ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr); 40737e1895aSJed Brown for (link=SNESMSTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 40837e1895aSJed Brown ierr = PetscOptionsEList("-snes_ms_type","Multistage smoother type","SNESMSSetType",(const char*const*)namelist,count,mstype,&choice,&flg);CHKERRQ(ierr); 409*57715debSLisandro Dalcin if (flg) {ierr = SNESMSSetType(snes,namelist[choice]);CHKERRQ(ierr);} 41037e1895aSJed Brown ierr = PetscFree(namelist);CHKERRQ(ierr); 411*57715debSLisandro Dalcin ierr = SNESMSGetDamping(snes,&damping);CHKERRQ(ierr); 412*57715debSLisandro Dalcin ierr = PetscOptionsReal("-snes_ms_damping","Damping for multistage method","SNESMSSetDamping",damping,&damping,&flg);CHKERRQ(ierr); 413*57715debSLisandro Dalcin if (flg) {ierr = SNESMSSetDamping(snes,damping);CHKERRQ(ierr);} 4140298fd71SBarry Smith ierr = PetscOptionsBool("-snes_ms_norms","Compute norms for monitoring","none",ms->norms,&ms->norms,NULL);CHKERRQ(ierr); 41537e1895aSJed Brown } 41637e1895aSJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 41737e1895aSJed Brown PetscFunctionReturn(0); 41837e1895aSJed Brown } 41937e1895aSJed Brown 420*57715debSLisandro Dalcin static PetscErrorCode SNESMSGetType_MS(SNES snes,SNESMSType *mstype) 42137e1895aSJed Brown { 422*57715debSLisandro Dalcin SNES_MS *ms = (SNES_MS*)snes->data; 423*57715debSLisandro Dalcin SNESMSTableau tab = ms->tableau; 424*57715debSLisandro Dalcin 425*57715debSLisandro Dalcin PetscFunctionBegin; 426*57715debSLisandro Dalcin *mstype = tab->name; 427*57715debSLisandro Dalcin PetscFunctionReturn(0); 428*57715debSLisandro Dalcin } 429*57715debSLisandro Dalcin 430*57715debSLisandro Dalcin static PetscErrorCode SNESMSSetType_MS(SNES snes,SNESMSType mstype) 431*57715debSLisandro Dalcin { 43237e1895aSJed Brown SNES_MS *ms = (SNES_MS*)snes->data; 43337e1895aSJed Brown SNESMSTableauLink link; 43437e1895aSJed Brown PetscBool match; 435*57715debSLisandro Dalcin PetscErrorCode ierr; 43637e1895aSJed Brown 43737e1895aSJed Brown PetscFunctionBegin; 43837e1895aSJed Brown if (ms->tableau) { 43937e1895aSJed Brown ierr = PetscStrcmp(ms->tableau->name,mstype,&match);CHKERRQ(ierr); 44037e1895aSJed Brown if (match) PetscFunctionReturn(0); 44137e1895aSJed Brown } 44237e1895aSJed Brown for (link = SNESMSTableauList; link; link=link->next) { 44337e1895aSJed Brown ierr = PetscStrcmp(link->tab.name,mstype,&match);CHKERRQ(ierr); 44437e1895aSJed Brown if (match) { 445*57715debSLisandro Dalcin if (snes->setupcalled) {ierr = SNESReset_MS(snes);CHKERRQ(ierr);} 44637e1895aSJed Brown ms->tableau = &link->tab; 447*57715debSLisandro Dalcin if (snes->setupcalled) {ierr = SNESSetUp_MS(snes);CHKERRQ(ierr);} 44837e1895aSJed Brown PetscFunctionReturn(0); 44937e1895aSJed Brown } 45037e1895aSJed Brown } 451ce94432eSBarry Smith SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mstype); 45237e1895aSJed Brown PetscFunctionReturn(0); 45337e1895aSJed Brown } 45437e1895aSJed Brown 45537e1895aSJed Brown /*@C 456*57715debSLisandro Dalcin SNESMSGetType - Get the type of multistage smoother 457*57715debSLisandro Dalcin 458*57715debSLisandro Dalcin Not collective 459*57715debSLisandro Dalcin 460*57715debSLisandro Dalcin Input Parameter: 461*57715debSLisandro Dalcin . snes - nonlinear solver context 462*57715debSLisandro Dalcin 463*57715debSLisandro Dalcin Output Parameter: 464*57715debSLisandro Dalcin . mstype - type of multistage method 465*57715debSLisandro Dalcin 466*57715debSLisandro Dalcin Level: beginner 467*57715debSLisandro Dalcin 468*57715debSLisandro Dalcin .seealso: SNESMSSetType(), SNESMSType, SNESMS 469*57715debSLisandro Dalcin @*/ 470*57715debSLisandro Dalcin PetscErrorCode SNESMSGetType(SNES snes,SNESMSType *mstype) 471*57715debSLisandro Dalcin { 472*57715debSLisandro Dalcin PetscErrorCode ierr; 473*57715debSLisandro Dalcin 474*57715debSLisandro Dalcin PetscFunctionBegin; 475*57715debSLisandro Dalcin PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 476*57715debSLisandro Dalcin PetscValidPointer(mstype,2); 477*57715debSLisandro Dalcin ierr = PetscUseMethod(snes,"SNESMSGetType_C",(SNES,SNESMSType*),(snes,mstype));CHKERRQ(ierr); 478*57715debSLisandro Dalcin PetscFunctionReturn(0); 479*57715debSLisandro Dalcin } 480*57715debSLisandro Dalcin 481*57715debSLisandro Dalcin /*@C 48237e1895aSJed Brown SNESMSSetType - Set the type of multistage smoother 48337e1895aSJed Brown 48437e1895aSJed Brown Logically collective 48537e1895aSJed Brown 48637e1895aSJed Brown Input Parameter: 48737e1895aSJed Brown + snes - nonlinear solver context 48837e1895aSJed Brown - mstype - type of multistage method 48937e1895aSJed Brown 49037e1895aSJed Brown Level: beginner 49137e1895aSJed Brown 492*57715debSLisandro Dalcin .seealso: SNESMSGetType(), SNESMSType, SNESMS 49337e1895aSJed Brown @*/ 494*57715debSLisandro Dalcin PetscErrorCode SNESMSSetType(SNES snes,SNESMSType mstype) 49537e1895aSJed Brown { 49637e1895aSJed Brown PetscErrorCode ierr; 49737e1895aSJed Brown 49837e1895aSJed Brown PetscFunctionBegin; 49937e1895aSJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 500*57715debSLisandro Dalcin PetscValidPointer(mstype,2); 501*57715debSLisandro Dalcin ierr = PetscTryMethod(snes,"SNESMSSetType_C",(SNES,SNESMSType),(snes,mstype));CHKERRQ(ierr); 502*57715debSLisandro Dalcin PetscFunctionReturn(0); 503*57715debSLisandro Dalcin } 504*57715debSLisandro Dalcin 505*57715debSLisandro Dalcin static PetscErrorCode SNESMSGetDamping_MS(SNES snes,PetscReal *damping) 506*57715debSLisandro Dalcin { 507*57715debSLisandro Dalcin SNES_MS *ms = (SNES_MS*)snes->data; 508*57715debSLisandro Dalcin 509*57715debSLisandro Dalcin PetscFunctionBegin; 510*57715debSLisandro Dalcin *damping = ms->damping; 511*57715debSLisandro Dalcin PetscFunctionReturn(0); 512*57715debSLisandro Dalcin } 513*57715debSLisandro Dalcin 514*57715debSLisandro Dalcin static PetscErrorCode SNESMSSetDamping_MS(SNES snes,PetscReal damping) 515*57715debSLisandro Dalcin { 516*57715debSLisandro Dalcin SNES_MS *ms = (SNES_MS*)snes->data; 517*57715debSLisandro Dalcin 518*57715debSLisandro Dalcin PetscFunctionBegin; 519*57715debSLisandro Dalcin ms->damping = damping; 520*57715debSLisandro Dalcin PetscFunctionReturn(0); 521*57715debSLisandro Dalcin } 522*57715debSLisandro Dalcin 523*57715debSLisandro Dalcin /*@C 524*57715debSLisandro Dalcin SNESMSGetDamping - Get the damping parameter 525*57715debSLisandro Dalcin 526*57715debSLisandro Dalcin Not collective 527*57715debSLisandro Dalcin 528*57715debSLisandro Dalcin Input Parameter: 529*57715debSLisandro Dalcin . snes - nonlinear solver context 530*57715debSLisandro Dalcin 531*57715debSLisandro Dalcin Output Parameter: 532*57715debSLisandro Dalcin . damping - damping parameter 533*57715debSLisandro Dalcin 534*57715debSLisandro Dalcin Level: advanced 535*57715debSLisandro Dalcin 536*57715debSLisandro Dalcin .seealso: SNESMSSetDamping(), SNESMS 537*57715debSLisandro Dalcin @*/ 538*57715debSLisandro Dalcin PetscErrorCode SNESMSGetDamping(SNES snes,PetscReal *damping) 539*57715debSLisandro Dalcin { 540*57715debSLisandro Dalcin PetscErrorCode ierr; 541*57715debSLisandro Dalcin 542*57715debSLisandro Dalcin PetscFunctionBegin; 543*57715debSLisandro Dalcin PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 544*57715debSLisandro Dalcin PetscValidPointer(damping,2); 545*57715debSLisandro Dalcin ierr = PetscUseMethod(snes,"SNESMSGetDamping_C",(SNES,PetscReal*),(snes,damping));CHKERRQ(ierr); 546*57715debSLisandro Dalcin PetscFunctionReturn(0); 547*57715debSLisandro Dalcin } 548*57715debSLisandro Dalcin 549*57715debSLisandro Dalcin /*@C 550*57715debSLisandro Dalcin SNESMSSetDamping - Set the damping parameter 551*57715debSLisandro Dalcin 552*57715debSLisandro Dalcin Logically collective 553*57715debSLisandro Dalcin 554*57715debSLisandro Dalcin Input Parameter: 555*57715debSLisandro Dalcin + snes - nonlinear solver context 556*57715debSLisandro Dalcin - damping - damping parameter 557*57715debSLisandro Dalcin 558*57715debSLisandro Dalcin Level: advanced 559*57715debSLisandro Dalcin 560*57715debSLisandro Dalcin .seealso: SNESMSGetDamping(), SNESMS 561*57715debSLisandro Dalcin @*/ 562*57715debSLisandro Dalcin PetscErrorCode SNESMSSetDamping(SNES snes,PetscReal damping) 563*57715debSLisandro Dalcin { 564*57715debSLisandro Dalcin PetscErrorCode ierr; 565*57715debSLisandro Dalcin 566*57715debSLisandro Dalcin PetscFunctionBegin; 567*57715debSLisandro Dalcin PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 568*57715debSLisandro Dalcin PetscValidLogicalCollectiveReal(snes,damping,2); 569*57715debSLisandro Dalcin ierr = PetscTryMethod(snes,"SNESMSSetDamping_C",(SNES,PetscReal),(snes,damping));CHKERRQ(ierr); 57037e1895aSJed Brown PetscFunctionReturn(0); 57137e1895aSJed Brown } 57237e1895aSJed Brown 57337e1895aSJed Brown /* -------------------------------------------------------------------------- */ 57437e1895aSJed Brown /*MC 57537e1895aSJed Brown SNESMS - multi-stage smoothers 57637e1895aSJed Brown 57737e1895aSJed Brown Options Database: 57837e1895aSJed Brown 57937e1895aSJed Brown + -snes_ms_type - type of multi-stage smoother 58037e1895aSJed Brown - -snes_ms_damping - damping for multi-stage method 58137e1895aSJed Brown 58237e1895aSJed Brown Notes: 58337e1895aSJed Brown These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for 5846c9de887SHong Zhang FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebyshev). 58537e1895aSJed Brown 58637e1895aSJed Brown Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds. 58737e1895aSJed Brown 58837e1895aSJed Brown The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with SNESMSRegister(). 58937e1895aSJed Brown 59037e1895aSJed Brown References: 591*57715debSLisandro Dalcin + 1. - Ketcheson (2010) Runge Kutta methods with minimum storage implementations (https://doi.org/10.1016/j.jcp.2009.11.006). 592*57715debSLisandro 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). 593*57715debSLisandro Dalcin . 3. - Pierce and Giles (1997) Preconditioned multigrid methods for compressible flow calculations on stretched meshes (https://doi.org/10.1006/jcph.1997.5772). 594*57715debSLisandro 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). 59537e1895aSJed Brown 59637e1895aSJed Brown Level: beginner 59737e1895aSJed Brown 5986c9de887SHong Zhang .seealso: SNESCreate(), SNES, SNESSetType(), SNESMS, SNESFAS, KSPCHEBYSHEV 59937e1895aSJed Brown 60037e1895aSJed Brown M*/ 6018cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_MS(SNES snes) 60237e1895aSJed Brown { 60337e1895aSJed Brown PetscErrorCode ierr; 60437e1895aSJed Brown SNES_MS *ms; 60537e1895aSJed Brown 60637e1895aSJed Brown PetscFunctionBegin; 607607a6623SBarry Smith ierr = SNESMSInitializePackage();CHKERRQ(ierr); 60837e1895aSJed Brown 60937e1895aSJed Brown snes->ops->setup = SNESSetUp_MS; 61037e1895aSJed Brown snes->ops->solve = SNESSolve_MS; 61137e1895aSJed Brown snes->ops->destroy = SNESDestroy_MS; 61237e1895aSJed Brown snes->ops->setfromoptions = SNESSetFromOptions_MS; 61337e1895aSJed Brown snes->ops->view = SNESView_MS; 61437e1895aSJed Brown snes->ops->reset = SNESReset_MS; 61537e1895aSJed Brown 616efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 61737e1895aSJed Brown snes->usesksp = PETSC_TRUE; 61837e1895aSJed Brown 6194fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 6204fc747eaSLawrence Mitchell 621b00a9115SJed Brown ierr = PetscNewLog(snes,&ms);CHKERRQ(ierr); 62237e1895aSJed Brown snes->data = (void*)ms; 62337e1895aSJed Brown ms->damping = 0.9; 62437e1895aSJed Brown ms->norms = PETSC_FALSE; 62537e1895aSJed Brown 626*57715debSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetType_C",SNESMSGetType_MS);CHKERRQ(ierr); 627bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetType_C",SNESMSSetType_MS);CHKERRQ(ierr); 628*57715debSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetDamping_C",SNESMSGetDamping_MS);CHKERRQ(ierr); 629*57715debSLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetDamping_C",SNESMSSetDamping_MS);CHKERRQ(ierr); 630*57715debSLisandro Dalcin 631*57715debSLisandro Dalcin ierr = SNESMSSetType(snes,SNESMSDefault);CHKERRQ(ierr); 63237e1895aSJed Brown PetscFunctionReturn(0); 63337e1895aSJed Brown } 63499e0435eSBarry Smith 635