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 PetscFunctionBegin; 4337e1895aSJed Brown if (SNESMSRegisterAllCalled) PetscFunctionReturn(0); 4437e1895aSJed Brown SNESMSRegisterAllCalled = PETSC_TRUE; 4537e1895aSJed Brown 4637e1895aSJed Brown { 479ce202e0SLisandro Dalcin PetscReal alpha[1] = {1.0}; 489566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSEULER,1,1,1.0,NULL,NULL,alpha)); 4937e1895aSJed Brown } 5037e1895aSJed Brown 5137e1895aSJed Brown { 5237e1895aSJed Brown PetscReal gamma[3][6] = { 5337e1895aSJed Brown {0.0000000000000000E+00, -7.0304722367110606E-01, -1.9836719667506464E-01, -1.6023843981863788E+00, 9.4483822882855284E-02, -1.4204296130641869E-01}, 5437e1895aSJed Brown {1.0000000000000000E+00, 1.1111025767083920E+00, 5.6150921583923230E-01, 7.4151723494934041E-01, 3.1714538168600587E-01, 4.6479276238548706E-01}, 5537e1895aSJed Brown {0.0000000000000000E+00, 0.0000000000000000E+00, 0.0000000000000000E+00, 6.7968174970583317E-01, -4.1755042846051737E-03, -1.9115668129923846E-01} 5637e1895aSJed Brown }; 5737e1895aSJed Brown PetscReal delta[6] = {1.0000000000000000E+00, 5.3275427433201750E-01, 6.0143544663985238E-01, 4.5874077053842177E-01, 2.7544386906104651E-01, 0.0000000000000000E+00}; 5837e1895aSJed Brown PetscReal betasub[6] = {8.4753115429481929E-01, 7.4018896368655618E-01, 6.5963574086583309E-03, 4.6747795645517759E-01, 1.3314545813643919E-01, 5.3260800028018784E-01}; 599566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSM62,6,3,1.0,&gamma[0][0],delta,betasub)); 6037e1895aSJed Brown } 619ce202e0SLisandro Dalcin 629ce202e0SLisandro Dalcin { /* Jameson (1983) */ 639ce202e0SLisandro Dalcin PetscReal alpha[4] = {0.25, 0.5, 0.55, 1.0}; 649566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSJAMESON83,4,1,1.0,NULL,NULL,alpha)); 65a97cb6bcSJed Brown } 663847c725SLisandro Dalcin 673847c725SLisandro Dalcin { /* Van Leer, Tai, and Powell (1989) 1 stage, order 1 */ 689ce202e0SLisandro Dalcin PetscReal alpha[1] = {1.0}; 699566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP11,1,1,0.5,NULL,NULL,alpha)); 703847c725SLisandro Dalcin } 71a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 2 stage, order 1 */ 729ce202e0SLisandro Dalcin PetscReal alpha[2] = {0.3333, 1.0}; 739566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP21,2,1,1.0,NULL,NULL,alpha)); 74a97cb6bcSJed Brown } 75a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 3 stage, order 1 */ 769ce202e0SLisandro Dalcin PetscReal alpha[3] = {0.1481, 0.4000, 1.0}; 779566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP31,3,1,1.5,NULL,NULL,alpha)); 78a97cb6bcSJed Brown } 79a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 4 stage, order 1 */ 809ce202e0SLisandro Dalcin PetscReal alpha[4] = {0.0833, 0.2069, 0.4265, 1.0}; 819566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP41,4,1,2.0,NULL,NULL,alpha)); 82a97cb6bcSJed Brown } 83a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 5 stage, order 1 */ 849ce202e0SLisandro Dalcin PetscReal alpha[5] = {0.0533, 0.1263, 0.2375, 0.4414,1.0}; 859566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP51,5,1,2.5,NULL,NULL,alpha)); 86a97cb6bcSJed Brown } 87a97cb6bcSJed Brown { /* Van Leer, Tai, and Powell (1989) 6 stage, order 1 */ 889ce202e0SLisandro Dalcin PetscReal alpha[6] = {0.0370, 0.0851, 0.1521, 0.2562, 0.4512, 1.0}; 899566063dSJacob Faibussowitsch PetscCall(SNESMSRegister(SNESMSVLTP61,6,1,3.0,NULL,NULL,alpha)); 90a97cb6bcSJed Brown } 9137e1895aSJed Brown PetscFunctionReturn(0); 9237e1895aSJed Brown } 9337e1895aSJed Brown 9437e1895aSJed Brown /*@C 9557715debSLisandro Dalcin SNESMSRegisterDestroy - Frees the list of schemes that were registered by SNESMSRegister(). 9637e1895aSJed Brown 9737e1895aSJed Brown Not Collective 9837e1895aSJed Brown 9937e1895aSJed Brown Level: advanced 10037e1895aSJed Brown 10157715debSLisandro Dalcin .seealso: SNESMSRegister(), SNESMSRegisterAll() 10237e1895aSJed Brown @*/ 10337e1895aSJed Brown PetscErrorCode SNESMSRegisterDestroy(void) 10437e1895aSJed Brown { 10537e1895aSJed Brown SNESMSTableauLink link; 10637e1895aSJed Brown 10737e1895aSJed Brown PetscFunctionBegin; 10837e1895aSJed Brown while ((link = SNESMSTableauList)) { 10937e1895aSJed Brown SNESMSTableau t = &link->tab; 11037e1895aSJed Brown SNESMSTableauList = link->next; 1111aa26658SKarl Rupp 1129566063dSJacob Faibussowitsch PetscCall(PetscFree(t->name)); 1139566063dSJacob Faibussowitsch PetscCall(PetscFree(t->gamma)); 1149566063dSJacob Faibussowitsch PetscCall(PetscFree(t->delta)); 1159566063dSJacob Faibussowitsch PetscCall(PetscFree(t->betasub)); 1169566063dSJacob Faibussowitsch PetscCall(PetscFree(link)); 11737e1895aSJed Brown } 11837e1895aSJed Brown SNESMSRegisterAllCalled = PETSC_FALSE; 11937e1895aSJed Brown PetscFunctionReturn(0); 12037e1895aSJed Brown } 12137e1895aSJed Brown 12237e1895aSJed Brown /*@C 12337e1895aSJed Brown SNESMSInitializePackage - This function initializes everything in the SNESMS package. It is called 1248a690491SBarry Smith from SNESInitializePackage(). 12537e1895aSJed Brown 12637e1895aSJed Brown Level: developer 12737e1895aSJed Brown 12837e1895aSJed Brown .seealso: PetscInitialize() 12937e1895aSJed Brown @*/ 130607a6623SBarry Smith PetscErrorCode SNESMSInitializePackage(void) 13137e1895aSJed Brown { 13237e1895aSJed Brown PetscFunctionBegin; 13337e1895aSJed Brown if (SNESMSPackageInitialized) PetscFunctionReturn(0); 13437e1895aSJed Brown SNESMSPackageInitialized = PETSC_TRUE; 1351aa26658SKarl Rupp 1369566063dSJacob Faibussowitsch PetscCall(SNESMSRegisterAll()); 1379566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(SNESMSFinalizePackage)); 13837e1895aSJed Brown PetscFunctionReturn(0); 13937e1895aSJed Brown } 14037e1895aSJed Brown 14137e1895aSJed Brown /*@C 14237e1895aSJed Brown SNESMSFinalizePackage - This function destroys everything in the SNESMS package. It is 14337e1895aSJed Brown called from PetscFinalize(). 14437e1895aSJed Brown 14537e1895aSJed Brown Level: developer 14637e1895aSJed Brown 14737e1895aSJed Brown .seealso: PetscFinalize() 14837e1895aSJed Brown @*/ 14937e1895aSJed Brown PetscErrorCode SNESMSFinalizePackage(void) 15037e1895aSJed Brown { 15137e1895aSJed Brown PetscFunctionBegin; 15237e1895aSJed Brown SNESMSPackageInitialized = PETSC_FALSE; 1531aa26658SKarl Rupp 1549566063dSJacob Faibussowitsch PetscCall(SNESMSRegisterDestroy()); 15537e1895aSJed Brown PetscFunctionReturn(0); 15637e1895aSJed Brown } 15737e1895aSJed Brown 15837e1895aSJed Brown /*@C 15937e1895aSJed Brown SNESMSRegister - register a multistage scheme 16037e1895aSJed Brown 16137e1895aSJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 16237e1895aSJed Brown 16337e1895aSJed Brown Input Parameters: 16437e1895aSJed Brown + name - identifier for method 16537e1895aSJed Brown . nstages - number of stages 16637e1895aSJed Brown . nregisters - number of registers used by low-storage implementation 1679ce202e0SLisandro Dalcin . stability - scaled stability region 16837e1895aSJed Brown . gamma - coefficients, see Ketcheson's paper 16937e1895aSJed Brown . delta - coefficients, see Ketcheson's paper 17037e1895aSJed Brown - betasub - subdiagonal of Shu-Osher form 17137e1895aSJed Brown 17237e1895aSJed Brown Notes: 17337e1895aSJed Brown The notation is described in Ketcheson (2010) Runge-Kutta methods with minimum storage implementations. 17437e1895aSJed Brown 1759ce202e0SLisandro Dalcin Many multistage schemes are of the form 1769ce202e0SLisandro Dalcin $ X_0 = X^{(n)} 1779ce202e0SLisandro Dalcin $ X_k = X_0 + \alpha_k * F(X_{k-1}), k = 1,\ldots,s 1789ce202e0SLisandro Dalcin $ X^{(n+1)} = X_s 1799ce202e0SLisandro Dalcin These methods can be registered with 1809ce202e0SLisandro Dalcin .vb 1819ce202e0SLisandro Dalcin SNESMSRegister("name",s,1,stability,NULL,NULL,alpha); 1829ce202e0SLisandro Dalcin .ve 1839ce202e0SLisandro Dalcin 18437e1895aSJed Brown Level: advanced 18537e1895aSJed Brown 18637e1895aSJed Brown .seealso: SNESMS 18737e1895aSJed Brown @*/ 18819fd82e9SBarry Smith PetscErrorCode SNESMSRegister(SNESMSType name,PetscInt nstages,PetscInt nregisters,PetscReal stability,const PetscReal gamma[],const PetscReal delta[],const PetscReal betasub[]) 18937e1895aSJed Brown { 19037e1895aSJed Brown SNESMSTableauLink link; 19137e1895aSJed Brown SNESMSTableau t; 19237e1895aSJed Brown 19337e1895aSJed Brown PetscFunctionBegin; 19437e1895aSJed Brown PetscValidCharPointer(name,1); 1952c71b3e2SJacob Faibussowitsch PetscCheckFalse(nstages < 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have at least one stage"); 1969ce202e0SLisandro Dalcin if (gamma || delta) { 1972c71b3e2SJacob Faibussowitsch PetscCheckFalse(nregisters != 3,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only support for methods written in 3-register form"); 198064a246eSJacob Faibussowitsch PetscValidRealPointer(gamma,5); 199064a246eSJacob Faibussowitsch PetscValidRealPointer(delta,6); 2009ce202e0SLisandro Dalcin } else { 2012c71b3e2SJacob Faibussowitsch PetscCheckFalse(nregisters != 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only support for methods written in 1-register form"); 2029ce202e0SLisandro Dalcin } 203064a246eSJacob Faibussowitsch PetscValidRealPointer(betasub,7); 20437e1895aSJed Brown 2059566063dSJacob Faibussowitsch PetscCall(SNESMSInitializePackage()); 2069566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 20737e1895aSJed Brown t = &link->tab; 2089566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name,&t->name)); 20937e1895aSJed Brown t->nstages = nstages; 21037e1895aSJed Brown t->nregisters = nregisters; 21137e1895aSJed Brown t->stability = stability; 2121aa26658SKarl Rupp 2139ce202e0SLisandro Dalcin if (gamma && delta) { 2149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nstages*nregisters,&t->gamma)); 2159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nstages,&t->delta)); 2169566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->gamma,gamma,nstages*nregisters)); 2179566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->delta,delta,nstages)); 2189ce202e0SLisandro Dalcin } 2199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nstages,&t->betasub)); 2209566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->betasub,betasub,nstages)); 22137e1895aSJed Brown 22237e1895aSJed Brown link->next = SNESMSTableauList; 22337e1895aSJed Brown SNESMSTableauList = link; 22437e1895aSJed Brown PetscFunctionReturn(0); 22537e1895aSJed Brown } 22637e1895aSJed Brown 22737e1895aSJed Brown /* 22837e1895aSJed Brown X - initial state, updated in-place. 22937e1895aSJed Brown F - residual, computed at the initial X on input 23037e1895aSJed Brown */ 23137e1895aSJed Brown static PetscErrorCode SNESMSStep_3Sstar(SNES snes,Vec X,Vec F) 23237e1895aSJed Brown { 23337e1895aSJed Brown SNES_MS *ms = (SNES_MS*)snes->data; 23437e1895aSJed Brown SNESMSTableau t = ms->tableau; 23537e1895aSJed Brown const PetscReal *gamma = t->gamma,*delta = t->delta,*betasub = t->betasub; 23637e1895aSJed Brown Vec S1,S2,S3,Y; 237bf80fd75SLisandro Dalcin PetscInt i,nstages = t->nstages; 23837e1895aSJed Brown 23937e1895aSJed Brown PetscFunctionBegin; 24037e1895aSJed Brown Y = snes->work[0]; 24137e1895aSJed Brown S1 = X; 24237e1895aSJed Brown S2 = snes->work[1]; 24337e1895aSJed Brown S3 = snes->work[2]; 2449566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(S2)); 2459566063dSJacob Faibussowitsch PetscCall(VecCopy(X,S3)); 24637e1895aSJed Brown for (i = 0; i < nstages; i++) { 247da80777bSKarl Rupp Vec Ss[4]; 248da80777bSKarl Rupp PetscScalar scoeff[4]; 249da80777bSKarl Rupp 250da80777bSKarl Rupp Ss[0] = S1; Ss[1] = S2; Ss[2] = S3; Ss[3] = Y; 2511aa26658SKarl Rupp 252bf80fd75SLisandro Dalcin scoeff[0] = gamma[0*nstages+i] - 1; 253da80777bSKarl Rupp scoeff[1] = gamma[1*nstages+i]; 254da80777bSKarl Rupp scoeff[2] = gamma[2*nstages+i]; 255da80777bSKarl Rupp scoeff[3] = -betasub[i]*ms->damping; 2561aa26658SKarl Rupp 2579566063dSJacob Faibussowitsch PetscCall(VecAXPY(S2,delta[i],S1)); 25837e1895aSJed Brown if (i > 0) { 2599566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,S1,F)); 26037e1895aSJed Brown } 2619566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp,F,Y)); 2629566063dSJacob Faibussowitsch PetscCall(VecMAXPY(S1,4,scoeff,Ss)); 26337e1895aSJed Brown } 26437e1895aSJed Brown PetscFunctionReturn(0); 26537e1895aSJed Brown } 26637e1895aSJed Brown 2679ce202e0SLisandro Dalcin /* 2689ce202e0SLisandro Dalcin X - initial state, updated in-place. 2699ce202e0SLisandro Dalcin F - residual, computed at the initial X on input 2709ce202e0SLisandro Dalcin */ 2719ce202e0SLisandro Dalcin static PetscErrorCode SNESMSStep_Basic(SNES snes,Vec X,Vec F) 2729ce202e0SLisandro Dalcin { 2739ce202e0SLisandro Dalcin SNES_MS *ms = (SNES_MS*)snes->data; 2749ce202e0SLisandro Dalcin SNESMSTableau tab = ms->tableau; 2759ce202e0SLisandro Dalcin const PetscReal *alpha = tab->betasub, h = ms->damping; 2769ce202e0SLisandro Dalcin PetscInt i,nstages = tab->nstages; 2779ce202e0SLisandro Dalcin Vec X0 = snes->work[0]; 2789ce202e0SLisandro Dalcin 2799ce202e0SLisandro Dalcin PetscFunctionBegin; 2809566063dSJacob Faibussowitsch PetscCall(VecCopy(X,X0)); 2819ce202e0SLisandro Dalcin for (i = 0; i < nstages; i++) { 2829ce202e0SLisandro Dalcin if (i > 0) { 2839566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,X,F)); 2849ce202e0SLisandro Dalcin } 2859566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp,F,X)); 2869566063dSJacob Faibussowitsch PetscCall(VecAYPX(X,-alpha[i]*h,X0)); 2879ce202e0SLisandro Dalcin } 2889ce202e0SLisandro Dalcin PetscFunctionReturn(0); 2899ce202e0SLisandro Dalcin } 2909ce202e0SLisandro Dalcin 2919ce202e0SLisandro Dalcin static PetscErrorCode SNESMSStep_Step(SNES snes,Vec X,Vec F) 2929ce202e0SLisandro Dalcin { 2939ce202e0SLisandro Dalcin SNES_MS *ms = (SNES_MS*)snes->data; 2949ce202e0SLisandro Dalcin SNESMSTableau tab = ms->tableau; 2959ce202e0SLisandro Dalcin 2969ce202e0SLisandro Dalcin PetscFunctionBegin; 2979ce202e0SLisandro Dalcin if (tab->gamma && tab->delta) { 2989566063dSJacob Faibussowitsch PetscCall(SNESMSStep_3Sstar(snes,X,F)); 2999ce202e0SLisandro Dalcin } else { 3009566063dSJacob Faibussowitsch PetscCall(SNESMSStep_Basic(snes,X,F)); 3019ce202e0SLisandro Dalcin } 3029ce202e0SLisandro Dalcin PetscFunctionReturn(0); 3039ce202e0SLisandro Dalcin } 3049ce202e0SLisandro Dalcin 305bf80fd75SLisandro Dalcin static PetscErrorCode SNESMSStep_Norms(SNES snes,PetscInt iter,Vec F) 306bf80fd75SLisandro Dalcin { 307bf80fd75SLisandro Dalcin SNES_MS *ms = (SNES_MS*)snes->data; 308bf80fd75SLisandro Dalcin PetscReal fnorm; 309bf80fd75SLisandro Dalcin 310bf80fd75SLisandro Dalcin PetscFunctionBegin; 311bf80fd75SLisandro Dalcin if (ms->norms) { 3129566063dSJacob Faibussowitsch PetscCall(VecNorm(F,NORM_2,&fnorm)); /* fnorm <- ||F|| */ 313bf80fd75SLisandro Dalcin SNESCheckFunctionNorm(snes,fnorm); 314bf80fd75SLisandro Dalcin /* Monitor convergence */ 3159566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 316bf80fd75SLisandro Dalcin snes->iter = iter; 317bf80fd75SLisandro Dalcin snes->norm = fnorm; 3189566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3199566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0)); 3209566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes,snes->iter,snes->norm)); 321bf80fd75SLisandro Dalcin /* Test for convergence */ 3229566063dSJacob Faibussowitsch PetscCall((*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP)); 323bf80fd75SLisandro Dalcin } else if (iter > 0) { 3249566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 325bf80fd75SLisandro Dalcin snes->iter = iter; 3269566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 327bf80fd75SLisandro Dalcin } 328bf80fd75SLisandro Dalcin PetscFunctionReturn(0); 329bf80fd75SLisandro Dalcin } 330bf80fd75SLisandro Dalcin 33137e1895aSJed Brown static PetscErrorCode SNESSolve_MS(SNES snes) 33237e1895aSJed Brown { 33337e1895aSJed Brown SNES_MS *ms = (SNES_MS*)snes->data; 33437e1895aSJed Brown Vec X = snes->vec_sol,F = snes->vec_func; 33537e1895aSJed Brown PetscInt i; 33637e1895aSJed Brown 33737e1895aSJed Brown PetscFunctionBegin; 3382c71b3e2SJacob Faibussowitsch PetscCheckFalse(snes->xl || snes->xu || snes->ops->computevariablebounds,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 3399566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation,&SNEScite)); 340bf80fd75SLisandro Dalcin 34137e1895aSJed Brown snes->reason = SNES_CONVERGED_ITERATING; 3429566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 34337e1895aSJed Brown snes->iter = 0; 344bf80fd75SLisandro Dalcin snes->norm = 0; 3459566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 346bf80fd75SLisandro Dalcin 347e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 3489566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,X,F)); 3491aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 3501aa26658SKarl Rupp 3519566063dSJacob Faibussowitsch PetscCall(SNESMSStep_Norms(snes,0,F)); 35237e1895aSJed Brown if (snes->reason) PetscFunctionReturn(0); 353bf80fd75SLisandro Dalcin 354bf80fd75SLisandro Dalcin for (i = 0; i < snes->max_its; i++) { 35537e1895aSJed Brown 35637e1895aSJed Brown /* Call general purpose update function */ 35737e1895aSJed Brown if (snes->ops->update) { 3589566063dSJacob Faibussowitsch PetscCall((*snes->ops->update)(snes,snes->iter)); 35937e1895aSJed Brown } 360bf80fd75SLisandro Dalcin 361bf80fd75SLisandro Dalcin if (i == 0 && snes->jacobian) { 362bf80fd75SLisandro Dalcin /* This method does not require a Jacobian, but it is usually preconditioned by PBJacobi */ 3639566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes,snes->vec_sol,snes->jacobian,snes->jacobian_pre)); 364bf80fd75SLisandro Dalcin SNESCheckJacobianDomainerror(snes); 3659566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre)); 366bf80fd75SLisandro Dalcin } 367bf80fd75SLisandro Dalcin 3689566063dSJacob Faibussowitsch PetscCall(SNESMSStep_Step(snes,X,F)); 36937e1895aSJed Brown 370bf80fd75SLisandro Dalcin if (i < snes->max_its-1 || ms->norms) { 3719566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,X,F)); 37237e1895aSJed Brown } 37337e1895aSJed Brown 3749566063dSJacob Faibussowitsch PetscCall(SNESMSStep_Norms(snes,i+1,F)); 37537e1895aSJed Brown if (snes->reason) PetscFunctionReturn(0); 37637e1895aSJed Brown } 37737e1895aSJed Brown if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; 37837e1895aSJed Brown PetscFunctionReturn(0); 37937e1895aSJed Brown } 38037e1895aSJed Brown 38137e1895aSJed Brown static PetscErrorCode SNESSetUp_MS(SNES snes) 38237e1895aSJed Brown { 3839ce202e0SLisandro Dalcin SNES_MS *ms = (SNES_MS*)snes->data; 3849ce202e0SLisandro Dalcin SNESMSTableau tab = ms->tableau; 3859ce202e0SLisandro Dalcin PetscInt nwork = tab->nregisters; 38637e1895aSJed Brown 38737e1895aSJed Brown PetscFunctionBegin; 3889566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes,nwork)); 3899566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 39037e1895aSJed Brown PetscFunctionReturn(0); 39137e1895aSJed Brown } 39237e1895aSJed Brown 39337e1895aSJed Brown static PetscErrorCode SNESReset_MS(SNES snes) 39437e1895aSJed Brown { 39537e1895aSJed Brown PetscFunctionBegin; 39637e1895aSJed Brown PetscFunctionReturn(0); 39737e1895aSJed Brown } 39837e1895aSJed Brown 39937e1895aSJed Brown static PetscErrorCode SNESDestroy_MS(SNES snes) 40037e1895aSJed Brown { 40137e1895aSJed Brown PetscFunctionBegin; 4029566063dSJacob Faibussowitsch PetscCall(SNESReset_MS(snes)); 4039566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 4049566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetType_C",NULL)); 4059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetType_C",NULL)); 4069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetDamping_C",NULL)); 4079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetDamping_C",NULL)); 40837e1895aSJed Brown PetscFunctionReturn(0); 40937e1895aSJed Brown } 41037e1895aSJed Brown 41137e1895aSJed Brown static PetscErrorCode SNESView_MS(SNES snes,PetscViewer viewer) 41237e1895aSJed Brown { 41337e1895aSJed Brown PetscBool iascii; 41437e1895aSJed Brown SNES_MS *ms = (SNES_MS*)snes->data; 41557715debSLisandro Dalcin SNESMSTableau tab = ms->tableau; 41637e1895aSJed Brown 41737e1895aSJed Brown PetscFunctionBegin; 4189566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 41937e1895aSJed Brown if (iascii) { 4209566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," multi-stage method type: %s\n",tab->name)); 42137e1895aSJed Brown } 42237e1895aSJed Brown PetscFunctionReturn(0); 42337e1895aSJed Brown } 42437e1895aSJed Brown 4254416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_MS(PetscOptionItems *PetscOptionsObject,SNES snes) 42637e1895aSJed Brown { 427725b86d8SJed Brown SNES_MS *ms = (SNES_MS*)snes->data; 42837e1895aSJed Brown 42937e1895aSJed Brown PetscFunctionBegin; 4309566063dSJacob Faibussowitsch PetscCall(PetscOptionsHead(PetscOptionsObject,"SNES MS options")); 43137e1895aSJed Brown { 43237e1895aSJed Brown SNESMSTableauLink link; 43337e1895aSJed Brown PetscInt count,choice; 43437e1895aSJed Brown PetscBool flg; 43537e1895aSJed Brown const char **namelist; 43657715debSLisandro Dalcin SNESMSType mstype; 43757715debSLisandro Dalcin PetscReal damping; 43837e1895aSJed Brown 4399566063dSJacob Faibussowitsch PetscCall(SNESMSGetType(snes,&mstype)); 44037e1895aSJed Brown for (link=SNESMSTableauList,count=0; link; link=link->next,count++) ; 4419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count,(char***)&namelist)); 44237e1895aSJed Brown for (link=SNESMSTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 4439566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-snes_ms_type","Multistage smoother type","SNESMSSetType",(const char*const*)namelist,count,mstype,&choice,&flg)); 4449566063dSJacob Faibussowitsch if (flg) PetscCall(SNESMSSetType(snes,namelist[choice])); 4459566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 4469566063dSJacob Faibussowitsch PetscCall(SNESMSGetDamping(snes,&damping)); 4479566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ms_damping","Damping for multistage method","SNESMSSetDamping",damping,&damping,&flg)); 4489566063dSJacob Faibussowitsch if (flg) PetscCall(SNESMSSetDamping(snes,damping)); 4499566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ms_norms","Compute norms for monitoring","none",ms->norms,&ms->norms,NULL)); 45037e1895aSJed Brown } 4519566063dSJacob Faibussowitsch PetscCall(PetscOptionsTail()); 45237e1895aSJed Brown PetscFunctionReturn(0); 45337e1895aSJed Brown } 45437e1895aSJed Brown 45557715debSLisandro Dalcin static PetscErrorCode SNESMSGetType_MS(SNES snes,SNESMSType *mstype) 45637e1895aSJed Brown { 45757715debSLisandro Dalcin SNES_MS *ms = (SNES_MS*)snes->data; 45857715debSLisandro Dalcin SNESMSTableau tab = ms->tableau; 45957715debSLisandro Dalcin 46057715debSLisandro Dalcin PetscFunctionBegin; 46157715debSLisandro Dalcin *mstype = tab->name; 46257715debSLisandro Dalcin PetscFunctionReturn(0); 46357715debSLisandro Dalcin } 46457715debSLisandro Dalcin 46557715debSLisandro Dalcin static PetscErrorCode SNESMSSetType_MS(SNES snes,SNESMSType mstype) 46657715debSLisandro Dalcin { 46737e1895aSJed Brown SNES_MS *ms = (SNES_MS*)snes->data; 46837e1895aSJed Brown SNESMSTableauLink link; 46937e1895aSJed Brown PetscBool match; 47037e1895aSJed Brown 47137e1895aSJed Brown PetscFunctionBegin; 47237e1895aSJed Brown if (ms->tableau) { 4739566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ms->tableau->name,mstype,&match)); 47437e1895aSJed Brown if (match) PetscFunctionReturn(0); 47537e1895aSJed Brown } 47637e1895aSJed Brown for (link = SNESMSTableauList; link; link=link->next) { 4779566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name,mstype,&match)); 47837e1895aSJed Brown if (match) { 4799566063dSJacob Faibussowitsch if (snes->setupcalled) PetscCall(SNESReset_MS(snes)); 48037e1895aSJed Brown ms->tableau = &link->tab; 4819566063dSJacob Faibussowitsch if (snes->setupcalled) PetscCall(SNESSetUp_MS(snes)); 48237e1895aSJed Brown PetscFunctionReturn(0); 48337e1895aSJed Brown } 48437e1895aSJed Brown } 48598921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mstype); 48637e1895aSJed Brown } 48737e1895aSJed Brown 48837e1895aSJed Brown /*@C 48957715debSLisandro Dalcin SNESMSGetType - Get the type of multistage smoother 49057715debSLisandro Dalcin 49157715debSLisandro Dalcin Not collective 49257715debSLisandro Dalcin 49357715debSLisandro Dalcin Input Parameter: 49457715debSLisandro Dalcin . snes - nonlinear solver context 49557715debSLisandro Dalcin 49657715debSLisandro Dalcin Output Parameter: 49757715debSLisandro Dalcin . mstype - type of multistage method 49857715debSLisandro Dalcin 49957715debSLisandro Dalcin Level: beginner 50057715debSLisandro Dalcin 50157715debSLisandro Dalcin .seealso: SNESMSSetType(), SNESMSType, SNESMS 50257715debSLisandro Dalcin @*/ 50357715debSLisandro Dalcin PetscErrorCode SNESMSGetType(SNES snes,SNESMSType *mstype) 50457715debSLisandro Dalcin { 50557715debSLisandro Dalcin PetscFunctionBegin; 50657715debSLisandro Dalcin PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 50757715debSLisandro Dalcin PetscValidPointer(mstype,2); 508*cac4c232SBarry Smith PetscUseMethod(snes,"SNESMSGetType_C",(SNES,SNESMSType*),(snes,mstype)); 50957715debSLisandro Dalcin PetscFunctionReturn(0); 51057715debSLisandro Dalcin } 51157715debSLisandro Dalcin 51257715debSLisandro Dalcin /*@C 51337e1895aSJed Brown SNESMSSetType - Set the type of multistage smoother 51437e1895aSJed Brown 51537e1895aSJed Brown Logically collective 51637e1895aSJed Brown 517d8d19677SJose E. Roman Input Parameters: 51837e1895aSJed Brown + snes - nonlinear solver context 51937e1895aSJed Brown - mstype - type of multistage method 52037e1895aSJed Brown 52137e1895aSJed Brown Level: beginner 52237e1895aSJed Brown 52357715debSLisandro Dalcin .seealso: SNESMSGetType(), SNESMSType, SNESMS 52437e1895aSJed Brown @*/ 52557715debSLisandro Dalcin PetscErrorCode SNESMSSetType(SNES snes,SNESMSType mstype) 52637e1895aSJed Brown { 52737e1895aSJed Brown PetscFunctionBegin; 52837e1895aSJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 529dadcf809SJacob Faibussowitsch PetscValidCharPointer(mstype,2); 530*cac4c232SBarry Smith PetscTryMethod(snes,"SNESMSSetType_C",(SNES,SNESMSType),(snes,mstype)); 53157715debSLisandro Dalcin PetscFunctionReturn(0); 53257715debSLisandro Dalcin } 53357715debSLisandro Dalcin 53457715debSLisandro Dalcin static PetscErrorCode SNESMSGetDamping_MS(SNES snes,PetscReal *damping) 53557715debSLisandro Dalcin { 53657715debSLisandro Dalcin SNES_MS *ms = (SNES_MS*)snes->data; 53757715debSLisandro Dalcin 53857715debSLisandro Dalcin PetscFunctionBegin; 53957715debSLisandro Dalcin *damping = ms->damping; 54057715debSLisandro Dalcin PetscFunctionReturn(0); 54157715debSLisandro Dalcin } 54257715debSLisandro Dalcin 54357715debSLisandro Dalcin static PetscErrorCode SNESMSSetDamping_MS(SNES snes,PetscReal damping) 54457715debSLisandro Dalcin { 54557715debSLisandro Dalcin SNES_MS *ms = (SNES_MS*)snes->data; 54657715debSLisandro Dalcin 54757715debSLisandro Dalcin PetscFunctionBegin; 54857715debSLisandro Dalcin ms->damping = damping; 54957715debSLisandro Dalcin PetscFunctionReturn(0); 55057715debSLisandro Dalcin } 55157715debSLisandro Dalcin 55257715debSLisandro Dalcin /*@C 55357715debSLisandro Dalcin SNESMSGetDamping - Get the damping parameter 55457715debSLisandro Dalcin 55557715debSLisandro Dalcin Not collective 55657715debSLisandro Dalcin 55757715debSLisandro Dalcin Input Parameter: 55857715debSLisandro Dalcin . snes - nonlinear solver context 55957715debSLisandro Dalcin 56057715debSLisandro Dalcin Output Parameter: 56157715debSLisandro Dalcin . damping - damping parameter 56257715debSLisandro Dalcin 56357715debSLisandro Dalcin Level: advanced 56457715debSLisandro Dalcin 56557715debSLisandro Dalcin .seealso: SNESMSSetDamping(), SNESMS 56657715debSLisandro Dalcin @*/ 56757715debSLisandro Dalcin PetscErrorCode SNESMSGetDamping(SNES snes,PetscReal *damping) 56857715debSLisandro Dalcin { 56957715debSLisandro Dalcin PetscFunctionBegin; 57057715debSLisandro Dalcin PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 571dadcf809SJacob Faibussowitsch PetscValidRealPointer(damping,2); 572*cac4c232SBarry Smith PetscUseMethod(snes,"SNESMSGetDamping_C",(SNES,PetscReal*),(snes,damping)); 57357715debSLisandro Dalcin PetscFunctionReturn(0); 57457715debSLisandro Dalcin } 57557715debSLisandro Dalcin 57657715debSLisandro Dalcin /*@C 57757715debSLisandro Dalcin SNESMSSetDamping - Set the damping parameter 57857715debSLisandro Dalcin 57957715debSLisandro Dalcin Logically collective 58057715debSLisandro Dalcin 581d8d19677SJose E. Roman Input Parameters: 58257715debSLisandro Dalcin + snes - nonlinear solver context 58357715debSLisandro Dalcin - damping - damping parameter 58457715debSLisandro Dalcin 58557715debSLisandro Dalcin Level: advanced 58657715debSLisandro Dalcin 58757715debSLisandro Dalcin .seealso: SNESMSGetDamping(), SNESMS 58857715debSLisandro Dalcin @*/ 58957715debSLisandro Dalcin PetscErrorCode SNESMSSetDamping(SNES snes,PetscReal damping) 59057715debSLisandro Dalcin { 59157715debSLisandro Dalcin PetscFunctionBegin; 59257715debSLisandro Dalcin PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 59357715debSLisandro Dalcin PetscValidLogicalCollectiveReal(snes,damping,2); 594*cac4c232SBarry Smith PetscTryMethod(snes,"SNESMSSetDamping_C",(SNES,PetscReal),(snes,damping)); 59537e1895aSJed Brown PetscFunctionReturn(0); 59637e1895aSJed Brown } 59737e1895aSJed Brown 59837e1895aSJed Brown /* -------------------------------------------------------------------------- */ 59937e1895aSJed Brown /*MC 60037e1895aSJed Brown SNESMS - multi-stage smoothers 60137e1895aSJed Brown 60237e1895aSJed Brown Options Database: 60337e1895aSJed Brown 60437e1895aSJed Brown + -snes_ms_type - type of multi-stage smoother 60537e1895aSJed Brown - -snes_ms_damping - damping for multi-stage method 60637e1895aSJed Brown 60737e1895aSJed Brown Notes: 60837e1895aSJed Brown These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for 6096c9de887SHong Zhang FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebyshev). 61037e1895aSJed Brown 61137e1895aSJed Brown Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds. 61237e1895aSJed Brown 61337e1895aSJed Brown The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with SNESMSRegister(). 61437e1895aSJed Brown 61537e1895aSJed Brown References: 616606c0280SSatish Balay + * - Ketcheson (2010) Runge Kutta methods with minimum storage implementations (https://doi.org/10.1016/j.jcp.2009.11.006). 617606c0280SSatish Balay . * - 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). 618606c0280SSatish Balay . * - Pierce and Giles (1997) Preconditioned multigrid methods for compressible flow calculations on stretched meshes (https://doi.org/10.1006/jcph.1997.5772). 619606c0280SSatish Balay - * - 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). 62037e1895aSJed Brown 62137e1895aSJed Brown Level: beginner 62237e1895aSJed Brown 6236c9de887SHong Zhang .seealso: SNESCreate(), SNES, SNESSetType(), SNESMS, SNESFAS, KSPCHEBYSHEV 62437e1895aSJed Brown 62537e1895aSJed Brown M*/ 6268cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_MS(SNES snes) 62737e1895aSJed Brown { 62837e1895aSJed Brown SNES_MS *ms; 62937e1895aSJed Brown 63037e1895aSJed Brown PetscFunctionBegin; 6319566063dSJacob Faibussowitsch PetscCall(SNESMSInitializePackage()); 63237e1895aSJed Brown 63337e1895aSJed Brown snes->ops->setup = SNESSetUp_MS; 63437e1895aSJed Brown snes->ops->solve = SNESSolve_MS; 63537e1895aSJed Brown snes->ops->destroy = SNESDestroy_MS; 63637e1895aSJed Brown snes->ops->setfromoptions = SNESSetFromOptions_MS; 63737e1895aSJed Brown snes->ops->view = SNESView_MS; 63837e1895aSJed Brown snes->ops->reset = SNESReset_MS; 63937e1895aSJed Brown 640efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 64137e1895aSJed Brown snes->usesksp = PETSC_TRUE; 64237e1895aSJed Brown 6434fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 6444fc747eaSLawrence Mitchell 6459566063dSJacob Faibussowitsch PetscCall(PetscNewLog(snes,&ms)); 64637e1895aSJed Brown snes->data = (void*)ms; 64737e1895aSJed Brown ms->damping = 0.9; 64837e1895aSJed Brown ms->norms = PETSC_FALSE; 64937e1895aSJed Brown 6509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetType_C",SNESMSGetType_MS)); 6519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetType_C",SNESMSSetType_MS)); 6529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESMSGetDamping_C",SNESMSGetDamping_MS)); 6539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESMSSetDamping_C",SNESMSSetDamping_MS)); 65457715debSLisandro Dalcin 6559566063dSJacob Faibussowitsch PetscCall(SNESMSSetType(snes,SNESMSDefault)); 65637e1895aSJed Brown PetscFunctionReturn(0); 65737e1895aSJed Brown } 658