1*37e1895aSJed Brown #include <private/snesimpl.h> /*I "petscsnes.h" I*/ 2*37e1895aSJed Brown 3*37e1895aSJed Brown static const SNESMSType SNESMSDefault = SNESMSM62; 4*37e1895aSJed Brown static PetscBool SNESMSRegisterAllCalled; 5*37e1895aSJed Brown static PetscBool SNESMSPackageInitialized; 6*37e1895aSJed Brown 7*37e1895aSJed Brown typedef struct _SNESMSTableau *SNESMSTableau; 8*37e1895aSJed Brown struct _SNESMSTableau { 9*37e1895aSJed Brown char *name; 10*37e1895aSJed Brown PetscInt nstages; /* Number of stages */ 11*37e1895aSJed Brown PetscInt nregisters; /* Number of registers */ 12*37e1895aSJed Brown PetscReal stability; /* Scaled stability region */ 13*37e1895aSJed Brown PetscReal *gamma; /* Coefficients of 3S* method */ 14*37e1895aSJed Brown PetscReal *delta; /* Coefficients of 3S* method */ 15*37e1895aSJed Brown PetscReal *betasub; /* Subdiagonal of beta in Shu-Osher form */ 16*37e1895aSJed Brown }; 17*37e1895aSJed Brown 18*37e1895aSJed Brown typedef struct _SNESMSTableauLink *SNESMSTableauLink; 19*37e1895aSJed Brown struct _SNESMSTableauLink { 20*37e1895aSJed Brown struct _SNESMSTableau tab; 21*37e1895aSJed Brown SNESMSTableauLink next; 22*37e1895aSJed Brown }; 23*37e1895aSJed Brown static SNESMSTableauLink SNESMSTableauList; 24*37e1895aSJed Brown 25*37e1895aSJed Brown typedef struct { 26*37e1895aSJed Brown SNESMSTableau tableau; /* Tableau in low-storage form */ 27*37e1895aSJed Brown PetscReal damping; /* Damping parameter, like length of (pseudo) time step */ 28*37e1895aSJed Brown PetscBool norms; /* Compute norms, usually only for monitoring purposes */ 29*37e1895aSJed Brown } SNES_MS; 30*37e1895aSJed Brown 31*37e1895aSJed Brown #undef __FUNCT__ 32*37e1895aSJed Brown #define __FUNCT__ "SNESMSRegisterAll" 33*37e1895aSJed Brown /*@C 34*37e1895aSJed Brown SNESMSRegisterAll - Registers all of the multi-stage methods in SNESMS 35*37e1895aSJed Brown 36*37e1895aSJed Brown Not Collective, but should be called by all processes which will need the schemes to be registered 37*37e1895aSJed Brown 38*37e1895aSJed Brown Level: advanced 39*37e1895aSJed Brown 40*37e1895aSJed Brown .keywords: SNES, SNESMS, register, all 41*37e1895aSJed Brown 42*37e1895aSJed Brown .seealso: SNESMSRegisterDestroy() 43*37e1895aSJed Brown @*/ 44*37e1895aSJed Brown PetscErrorCode SNESMSRegisterAll(void) 45*37e1895aSJed Brown { 46*37e1895aSJed Brown PetscErrorCode ierr; 47*37e1895aSJed Brown 48*37e1895aSJed Brown PetscFunctionBegin; 49*37e1895aSJed Brown if (SNESMSRegisterAllCalled) PetscFunctionReturn(0); 50*37e1895aSJed Brown SNESMSRegisterAllCalled = PETSC_TRUE; 51*37e1895aSJed Brown 52*37e1895aSJed Brown { 53*37e1895aSJed Brown PetscReal gamma[3][1] = {{1.0},{0.0},{0.0}}; 54*37e1895aSJed Brown PetscReal delta[1] = {0.0}; 55*37e1895aSJed Brown PetscReal betasub[1] = {1.0}; 56*37e1895aSJed Brown ierr = SNESMSRegister(SNESMSEULER,1,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr); 57*37e1895aSJed Brown } 58*37e1895aSJed Brown 59*37e1895aSJed Brown { 60*37e1895aSJed Brown PetscReal gamma[3][6] = { 61*37e1895aSJed Brown {0.0000000000000000E+00, -7.0304722367110606E-01, -1.9836719667506464E-01, -1.6023843981863788E+00, 9.4483822882855284E-02, -1.4204296130641869E-01}, 62*37e1895aSJed Brown {1.0000000000000000E+00, 1.1111025767083920E+00, 5.6150921583923230E-01, 7.4151723494934041E-01, 3.1714538168600587E-01, 4.6479276238548706E-01}, 63*37e1895aSJed Brown {0.0000000000000000E+00, 0.0000000000000000E+00, 0.0000000000000000E+00, 6.7968174970583317E-01, -4.1755042846051737E-03, -1.9115668129923846E-01} 64*37e1895aSJed Brown }; 65*37e1895aSJed Brown PetscReal delta[6] = {1.0000000000000000E+00, 5.3275427433201750E-01, 6.0143544663985238E-01, 4.5874077053842177E-01, 2.7544386906104651E-01, 0.0000000000000000E+00}; 66*37e1895aSJed Brown PetscReal betasub[6] = {8.4753115429481929E-01, 7.4018896368655618E-01, 6.5963574086583309E-03, 4.6747795645517759E-01, 1.3314545813643919E-01, 5.3260800028018784E-01}; 67*37e1895aSJed Brown ierr = SNESMSRegister(SNESMSM62,6,3,1.0,&gamma[0][0],delta,betasub);CHKERRQ(ierr); 68*37e1895aSJed Brown } 69*37e1895aSJed Brown PetscFunctionReturn(0); 70*37e1895aSJed Brown } 71*37e1895aSJed Brown 72*37e1895aSJed Brown #undef __FUNCT__ 73*37e1895aSJed Brown #define __FUNCT__ "SNESMSRegisterDestroy" 74*37e1895aSJed Brown /*@C 75*37e1895aSJed Brown SNESMSRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister(). 76*37e1895aSJed Brown 77*37e1895aSJed Brown Not Collective 78*37e1895aSJed Brown 79*37e1895aSJed Brown Level: advanced 80*37e1895aSJed Brown 81*37e1895aSJed Brown .keywords: TSRosW, register, destroy 82*37e1895aSJed Brown .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic() 83*37e1895aSJed Brown @*/ 84*37e1895aSJed Brown PetscErrorCode SNESMSRegisterDestroy(void) 85*37e1895aSJed Brown { 86*37e1895aSJed Brown PetscErrorCode ierr; 87*37e1895aSJed Brown SNESMSTableauLink link; 88*37e1895aSJed Brown 89*37e1895aSJed Brown PetscFunctionBegin; 90*37e1895aSJed Brown while ((link = SNESMSTableauList)) { 91*37e1895aSJed Brown SNESMSTableau t = &link->tab; 92*37e1895aSJed Brown SNESMSTableauList = link->next; 93*37e1895aSJed Brown ierr = PetscFree3(t->gamma,t->delta,t->betasub);CHKERRQ(ierr); 94*37e1895aSJed Brown ierr = PetscFree(t->name);CHKERRQ(ierr); 95*37e1895aSJed Brown ierr = PetscFree(link);CHKERRQ(ierr); 96*37e1895aSJed Brown } 97*37e1895aSJed Brown SNESMSRegisterAllCalled = PETSC_FALSE; 98*37e1895aSJed Brown PetscFunctionReturn(0); 99*37e1895aSJed Brown } 100*37e1895aSJed Brown 101*37e1895aSJed Brown #undef __FUNCT__ 102*37e1895aSJed Brown #define __FUNCT__ "SNESMSInitializePackage" 103*37e1895aSJed Brown /*@C 104*37e1895aSJed Brown SNESMSInitializePackage - This function initializes everything in the SNESMS package. It is called 105*37e1895aSJed Brown from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to SNESCreate_MS() 106*37e1895aSJed Brown when using static libraries. 107*37e1895aSJed Brown 108*37e1895aSJed Brown Input Parameter: 109*37e1895aSJed Brown path - The dynamic library path, or PETSC_NULL 110*37e1895aSJed Brown 111*37e1895aSJed Brown Level: developer 112*37e1895aSJed Brown 113*37e1895aSJed Brown .keywords: SNES, SNESMS, initialize, package 114*37e1895aSJed Brown .seealso: PetscInitialize() 115*37e1895aSJed Brown @*/ 116*37e1895aSJed Brown PetscErrorCode SNESMSInitializePackage(const char path[]) 117*37e1895aSJed Brown { 118*37e1895aSJed Brown PetscErrorCode ierr; 119*37e1895aSJed Brown 120*37e1895aSJed Brown PetscFunctionBegin; 121*37e1895aSJed Brown if (SNESMSPackageInitialized) PetscFunctionReturn(0); 122*37e1895aSJed Brown SNESMSPackageInitialized = PETSC_TRUE; 123*37e1895aSJed Brown ierr = SNESMSRegisterAll();CHKERRQ(ierr); 124*37e1895aSJed Brown ierr = PetscRegisterFinalize(SNESMSFinalizePackage);CHKERRQ(ierr); 125*37e1895aSJed Brown PetscFunctionReturn(0); 126*37e1895aSJed Brown } 127*37e1895aSJed Brown 128*37e1895aSJed Brown #undef __FUNCT__ 129*37e1895aSJed Brown #define __FUNCT__ "SNESMSFinalizePackage" 130*37e1895aSJed Brown /*@C 131*37e1895aSJed Brown SNESMSFinalizePackage - This function destroys everything in the SNESMS package. It is 132*37e1895aSJed Brown called from PetscFinalize(). 133*37e1895aSJed Brown 134*37e1895aSJed Brown Level: developer 135*37e1895aSJed Brown 136*37e1895aSJed Brown .keywords: Petsc, destroy, package 137*37e1895aSJed Brown .seealso: PetscFinalize() 138*37e1895aSJed Brown @*/ 139*37e1895aSJed Brown PetscErrorCode SNESMSFinalizePackage(void) 140*37e1895aSJed Brown { 141*37e1895aSJed Brown PetscErrorCode ierr; 142*37e1895aSJed Brown 143*37e1895aSJed Brown PetscFunctionBegin; 144*37e1895aSJed Brown SNESMSPackageInitialized = PETSC_FALSE; 145*37e1895aSJed Brown ierr = SNESMSRegisterDestroy();CHKERRQ(ierr); 146*37e1895aSJed Brown PetscFunctionReturn(0); 147*37e1895aSJed Brown } 148*37e1895aSJed Brown 149*37e1895aSJed Brown #undef __FUNCT__ 150*37e1895aSJed Brown #define __FUNCT__ "SNESMSRegister" 151*37e1895aSJed Brown /*@C 152*37e1895aSJed Brown SNESMSRegister - register a multistage scheme 153*37e1895aSJed Brown 154*37e1895aSJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 155*37e1895aSJed Brown 156*37e1895aSJed Brown Input Parameters: 157*37e1895aSJed Brown + name - identifier for method 158*37e1895aSJed Brown . nstages - number of stages 159*37e1895aSJed Brown . nregisters - number of registers used by low-storage implementation 160*37e1895aSJed Brown . gamma - coefficients, see Ketcheson's paper 161*37e1895aSJed Brown . delta - coefficients, see Ketcheson's paper 162*37e1895aSJed Brown - betasub - subdiagonal of Shu-Osher form 163*37e1895aSJed Brown 164*37e1895aSJed Brown Notes: 165*37e1895aSJed Brown The notation is described in Ketcheson (2010) Runge-Kutta methods with minimum storage implementations. 166*37e1895aSJed Brown 167*37e1895aSJed Brown Level: advanced 168*37e1895aSJed Brown 169*37e1895aSJed Brown .keywords: SNES, register 170*37e1895aSJed Brown 171*37e1895aSJed Brown .seealso: SNESMS 172*37e1895aSJed Brown @*/ 173*37e1895aSJed Brown PetscErrorCode SNESMSRegister(const SNESMSType name,PetscInt nstages,PetscInt nregisters,PetscReal stability,const PetscReal gamma[],const PetscReal delta[],const PetscReal betasub[]) 174*37e1895aSJed Brown { 175*37e1895aSJed Brown PetscErrorCode ierr; 176*37e1895aSJed Brown SNESMSTableauLink link; 177*37e1895aSJed Brown SNESMSTableau t; 178*37e1895aSJed Brown 179*37e1895aSJed Brown PetscFunctionBegin; 180*37e1895aSJed Brown PetscValidCharPointer(name,1); 181*37e1895aSJed Brown if (nstages < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have at least one stage"); 182*37e1895aSJed Brown if (nregisters != 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only support for methods written in 3-register form"); 183*37e1895aSJed Brown PetscValidPointer(gamma,4); 184*37e1895aSJed Brown PetscValidPointer(delta,5); 185*37e1895aSJed Brown PetscValidPointer(betasub,6); 186*37e1895aSJed Brown 187*37e1895aSJed Brown ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 188*37e1895aSJed Brown ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr); 189*37e1895aSJed Brown t = &link->tab; 190*37e1895aSJed Brown ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 191*37e1895aSJed Brown t->nstages = nstages; 192*37e1895aSJed Brown t->nregisters = nregisters; 193*37e1895aSJed Brown t->stability = stability; 194*37e1895aSJed Brown ierr = PetscMalloc3(nstages*nregisters,PetscReal,&t->gamma,nstages,PetscReal,&t->delta,nstages,PetscReal,&t->betasub);CHKERRQ(ierr); 195*37e1895aSJed Brown ierr = PetscMemcpy(t->gamma,gamma,nstages*nregisters*sizeof(PetscReal));CHKERRQ(ierr); 196*37e1895aSJed Brown ierr = PetscMemcpy(t->delta,delta,nstages*sizeof(PetscReal));CHKERRQ(ierr); 197*37e1895aSJed Brown ierr = PetscMemcpy(t->betasub,betasub,nstages*sizeof(PetscReal));CHKERRQ(ierr); 198*37e1895aSJed Brown 199*37e1895aSJed Brown link->next = SNESMSTableauList; 200*37e1895aSJed Brown SNESMSTableauList = link; 201*37e1895aSJed Brown PetscFunctionReturn(0); 202*37e1895aSJed Brown } 203*37e1895aSJed Brown 204*37e1895aSJed Brown #undef __FUNCT__ 205*37e1895aSJed Brown #define __FUNCT__ "SNESMSStep_3Sstar" 206*37e1895aSJed Brown /* 207*37e1895aSJed Brown X - initial state, updated in-place. 208*37e1895aSJed Brown F - residual, computed at the initial X on input 209*37e1895aSJed Brown */ 210*37e1895aSJed Brown static PetscErrorCode SNESMSStep_3Sstar(SNES snes,Vec X,Vec F) 211*37e1895aSJed Brown { 212*37e1895aSJed Brown PetscErrorCode ierr; 213*37e1895aSJed Brown SNES_MS *ms = (SNES_MS*)snes->data; 214*37e1895aSJed Brown SNESMSTableau t = ms->tableau; 215*37e1895aSJed Brown const PetscReal *gamma = t->gamma,*delta = t->delta,*betasub = t->betasub; 216*37e1895aSJed Brown Vec S1,S2,S3,Y; 217*37e1895aSJed Brown PetscInt i,nstages = t->nstages;; 218*37e1895aSJed Brown 219*37e1895aSJed Brown 220*37e1895aSJed Brown PetscFunctionBegin; 221*37e1895aSJed Brown Y = snes->work[0]; 222*37e1895aSJed Brown S1 = X; 223*37e1895aSJed Brown S2 = snes->work[1]; 224*37e1895aSJed Brown S3 = snes->work[2]; 225*37e1895aSJed Brown ierr = VecZeroEntries(S2);CHKERRQ(ierr); 226*37e1895aSJed Brown ierr = VecCopy(X,S3);CHKERRQ(ierr); 227*37e1895aSJed Brown for (i=0; i<nstages; i++) { 228*37e1895aSJed Brown Vec Ss[4] = {S1,S2,S3,Y}; 229*37e1895aSJed Brown PetscReal scoeff[4] = {gamma[0*nstages+i]-1.0,gamma[1*nstages+i],gamma[2*nstages+i],-betasub[i]*ms->damping}; 230*37e1895aSJed Brown ierr = VecAXPY(S2,delta[i],S1);CHKERRQ(ierr); 231*37e1895aSJed Brown if (i>0) { 232*37e1895aSJed Brown ierr = SNESComputeFunction(snes,S1,F);CHKERRQ(ierr); 233*37e1895aSJed Brown if (snes->domainerror) { 234*37e1895aSJed Brown snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 235*37e1895aSJed Brown PetscFunctionReturn(0); 236*37e1895aSJed Brown } 237*37e1895aSJed Brown } 238*37e1895aSJed Brown ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr); 239*37e1895aSJed Brown ierr = VecMAXPY(S1,4,scoeff,Ss);CHKERRQ(ierr); 240*37e1895aSJed Brown } 241*37e1895aSJed Brown PetscFunctionReturn(0); 242*37e1895aSJed Brown } 243*37e1895aSJed Brown 244*37e1895aSJed Brown #undef __FUNCT__ 245*37e1895aSJed Brown #define __FUNCT__ "SNESSolve_MS" 246*37e1895aSJed Brown static PetscErrorCode SNESSolve_MS(SNES snes) 247*37e1895aSJed Brown { 248*37e1895aSJed Brown SNES_MS *ms = (SNES_MS*)snes->data; 249*37e1895aSJed Brown Vec X = snes->vec_sol,F = snes->vec_func; 250*37e1895aSJed Brown PetscReal fnorm; 251*37e1895aSJed Brown MatStructure mstruct; 252*37e1895aSJed Brown PetscInt i; 253*37e1895aSJed Brown PetscErrorCode ierr; 254*37e1895aSJed Brown 255*37e1895aSJed Brown PetscFunctionBegin; 256*37e1895aSJed Brown snes->reason = SNES_CONVERGED_ITERATING; 257*37e1895aSJed Brown ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 258*37e1895aSJed Brown snes->iter = 0; 259*37e1895aSJed Brown snes->norm = 0.; 260*37e1895aSJed Brown ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 261*37e1895aSJed Brown 262*37e1895aSJed Brown ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 263*37e1895aSJed Brown if (snes->domainerror) { 264*37e1895aSJed Brown snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 265*37e1895aSJed Brown PetscFunctionReturn(0); 266*37e1895aSJed Brown } 267*37e1895aSJed Brown if (snes->jacobian) { /* This method does not require a Jacobian, but it is usually preconditioned by PBJacobi */ 268*37e1895aSJed Brown ierr = SNESComputeJacobian(snes,snes->vec_sol,&snes->jacobian,&snes->jacobian_pre,&mstruct);CHKERRQ(ierr); 269*37e1895aSJed Brown ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,mstruct);CHKERRQ(ierr); 270*37e1895aSJed Brown } 271*37e1895aSJed Brown if (ms->norms) { 272*37e1895aSJed Brown ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 273*37e1895aSJed Brown if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 274*37e1895aSJed Brown /* Monitor convergence */ 275*37e1895aSJed Brown ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 276*37e1895aSJed Brown snes->iter = 0; 277*37e1895aSJed Brown snes->norm = fnorm; 278*37e1895aSJed Brown ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 279*37e1895aSJed Brown SNESLogConvHistory(snes,snes->norm,0); 280*37e1895aSJed Brown ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 281*37e1895aSJed Brown 282*37e1895aSJed Brown /* set parameter for default relative tolerance convergence test */ 283*37e1895aSJed Brown snes->ttol = fnorm*snes->rtol; 284*37e1895aSJed Brown /* Test for convergence */ 285*37e1895aSJed Brown ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 286*37e1895aSJed Brown if (snes->reason) PetscFunctionReturn(0); 287*37e1895aSJed Brown } 288*37e1895aSJed Brown 289*37e1895aSJed Brown /* Call general purpose update function */ 290*37e1895aSJed Brown if (snes->ops->update) { 291*37e1895aSJed Brown ierr = (*snes->ops->update)(snes,snes->iter);CHKERRQ(ierr); 292*37e1895aSJed Brown } 293*37e1895aSJed Brown for (i = 0; i < snes->max_its; i++) { 294*37e1895aSJed Brown ierr = SNESMSStep_3Sstar(snes,X,F);CHKERRQ(ierr); 295*37e1895aSJed Brown 296*37e1895aSJed Brown if (i+1 < snes->max_its || ms->norms) { 297*37e1895aSJed Brown ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 298*37e1895aSJed Brown if (snes->domainerror) { 299*37e1895aSJed Brown snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 300*37e1895aSJed Brown PetscFunctionReturn(0); 301*37e1895aSJed Brown } 302*37e1895aSJed Brown } 303*37e1895aSJed Brown 304*37e1895aSJed Brown if (ms->norms) { 305*37e1895aSJed Brown ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 306*37e1895aSJed Brown if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 307*37e1895aSJed Brown /* Monitor convergence */ 308*37e1895aSJed Brown ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 309*37e1895aSJed Brown snes->iter = i+1; 310*37e1895aSJed Brown snes->norm = fnorm; 311*37e1895aSJed Brown ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 312*37e1895aSJed Brown SNESLogConvHistory(snes,snes->norm,0); 313*37e1895aSJed Brown ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 314*37e1895aSJed Brown 315*37e1895aSJed Brown /* Test for convergence */ 316*37e1895aSJed Brown ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 317*37e1895aSJed Brown if (snes->reason) PetscFunctionReturn(0); 318*37e1895aSJed Brown } 319*37e1895aSJed Brown 320*37e1895aSJed Brown /* Call general purpose update function */ 321*37e1895aSJed Brown if (snes->ops->update) { 322*37e1895aSJed Brown ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 323*37e1895aSJed Brown } 324*37e1895aSJed Brown } 325*37e1895aSJed Brown if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; 326*37e1895aSJed Brown PetscFunctionReturn(0); 327*37e1895aSJed Brown } 328*37e1895aSJed Brown 329*37e1895aSJed Brown #undef __FUNCT__ 330*37e1895aSJed Brown #define __FUNCT__ "SNESSetUp_MS" 331*37e1895aSJed Brown static PetscErrorCode SNESSetUp_MS(SNES snes) 332*37e1895aSJed Brown { 333*37e1895aSJed Brown SNES_MS * ms = (SNES_MS *)snes->data; 334*37e1895aSJed Brown PetscErrorCode ierr; 335*37e1895aSJed Brown 336*37e1895aSJed Brown PetscFunctionBegin; 337*37e1895aSJed Brown if (!ms->tableau) {ierr = SNESMSSetType(snes,SNESMSDefault);CHKERRQ(ierr);} 338*37e1895aSJed Brown ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr); 339*37e1895aSJed Brown PetscFunctionReturn(0); 340*37e1895aSJed Brown } 341*37e1895aSJed Brown 342*37e1895aSJed Brown #undef __FUNCT__ 343*37e1895aSJed Brown #define __FUNCT__ "SNESReset_MS" 344*37e1895aSJed Brown static PetscErrorCode SNESReset_MS(SNES snes) 345*37e1895aSJed Brown { 346*37e1895aSJed Brown 347*37e1895aSJed Brown PetscFunctionBegin; 348*37e1895aSJed Brown PetscFunctionReturn(0); 349*37e1895aSJed Brown } 350*37e1895aSJed Brown 351*37e1895aSJed Brown #undef __FUNCT__ 352*37e1895aSJed Brown #define __FUNCT__ "SNESDestroy_MS" 353*37e1895aSJed Brown static PetscErrorCode SNESDestroy_MS(SNES snes) 354*37e1895aSJed Brown { 355*37e1895aSJed Brown PetscErrorCode ierr; 356*37e1895aSJed Brown 357*37e1895aSJed Brown PetscFunctionBegin; 358*37e1895aSJed Brown ierr = PetscFree(snes->data);CHKERRQ(ierr); 359*37e1895aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr); 360*37e1895aSJed Brown PetscFunctionReturn(0); 361*37e1895aSJed Brown } 362*37e1895aSJed Brown 363*37e1895aSJed Brown #undef __FUNCT__ 364*37e1895aSJed Brown #define __FUNCT__ "SNESView_MS" 365*37e1895aSJed Brown static PetscErrorCode SNESView_MS(SNES snes,PetscViewer viewer) 366*37e1895aSJed Brown { 367*37e1895aSJed Brown PetscBool iascii; 368*37e1895aSJed Brown PetscErrorCode ierr; 369*37e1895aSJed Brown SNES_MS *ms = (SNES_MS*)snes->data; 370*37e1895aSJed Brown 371*37e1895aSJed Brown PetscFunctionBegin; 372*37e1895aSJed Brown ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 373*37e1895aSJed Brown if (iascii) { 374*37e1895aSJed Brown SNESMSTableau tab = ms->tableau; 375*37e1895aSJed Brown ierr = PetscViewerASCIIPrintf(viewer," multi-stage method type: %s\n",tab?tab->name:"not yet set");CHKERRQ(ierr); 376*37e1895aSJed Brown } 377*37e1895aSJed Brown PetscFunctionReturn(0); 378*37e1895aSJed Brown } 379*37e1895aSJed Brown 380*37e1895aSJed Brown #undef __FUNCT__ 381*37e1895aSJed Brown #define __FUNCT__ "SNESSetFromOptions_MS" 382*37e1895aSJed Brown static PetscErrorCode SNESSetFromOptions_MS(SNES snes) 383*37e1895aSJed Brown { 384*37e1895aSJed Brown SNES_MS *ms = (SNES_MS*)snes->data;; 385*37e1895aSJed Brown PetscErrorCode ierr; 386*37e1895aSJed Brown 387*37e1895aSJed Brown PetscFunctionBegin; 388*37e1895aSJed Brown ierr = PetscOptionsHead("SNES MS options");CHKERRQ(ierr); 389*37e1895aSJed Brown { 390*37e1895aSJed Brown SNESMSTableauLink link; 391*37e1895aSJed Brown PetscInt count,choice; 392*37e1895aSJed Brown PetscBool flg; 393*37e1895aSJed Brown const char **namelist; 394*37e1895aSJed Brown char mstype[256]; 395*37e1895aSJed Brown 396*37e1895aSJed Brown ierr = PetscStrncpy(mstype,SNESMSDefault,sizeof mstype);CHKERRQ(ierr); 397*37e1895aSJed Brown for (link=SNESMSTableauList,count=0; link; link=link->next,count++) ; 398*37e1895aSJed Brown ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr); 399*37e1895aSJed Brown for (link=SNESMSTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 400*37e1895aSJed Brown ierr = PetscOptionsEList("-snes_ms_type","Multistage smoother type","SNESMSSetType",(const char*const*)namelist,count,mstype,&choice,&flg);CHKERRQ(ierr); 401*37e1895aSJed Brown ierr = SNESMSSetType(snes,flg ? namelist[choice] : mstype);CHKERRQ(ierr); 402*37e1895aSJed Brown ierr = PetscFree(namelist);CHKERRQ(ierr); 403*37e1895aSJed Brown ierr = PetscOptionsReal("-snes_ms_damping","Damping for multistage method","SNESMSSetDamping",ms->damping,&ms->damping,PETSC_NULL);CHKERRQ(ierr); 404*37e1895aSJed Brown ierr = PetscOptionsBool("-snes_ms_norms","Compute norms for monitoring","none",ms->norms,&ms->norms,PETSC_NULL);CHKERRQ(ierr); 405*37e1895aSJed Brown } 406*37e1895aSJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 407*37e1895aSJed Brown PetscFunctionReturn(0); 408*37e1895aSJed Brown } 409*37e1895aSJed Brown 410*37e1895aSJed Brown EXTERN_C_BEGIN 411*37e1895aSJed Brown #undef __FUNCT__ 412*37e1895aSJed Brown #define __FUNCT__ "SNESMSSetType_MS" 413*37e1895aSJed Brown PetscErrorCode SNESMSSetType_MS(SNES snes,SNESMSType mstype) 414*37e1895aSJed Brown { 415*37e1895aSJed Brown PetscErrorCode ierr; 416*37e1895aSJed Brown SNES_MS *ms = (SNES_MS*)snes->data; 417*37e1895aSJed Brown SNESMSTableauLink link; 418*37e1895aSJed Brown PetscBool match; 419*37e1895aSJed Brown 420*37e1895aSJed Brown PetscFunctionBegin; 421*37e1895aSJed Brown if (ms->tableau) { 422*37e1895aSJed Brown ierr = PetscStrcmp(ms->tableau->name,mstype,&match);CHKERRQ(ierr); 423*37e1895aSJed Brown if (match) PetscFunctionReturn(0); 424*37e1895aSJed Brown } 425*37e1895aSJed Brown for (link = SNESMSTableauList; link; link=link->next) { 426*37e1895aSJed Brown ierr = PetscStrcmp(link->tab.name,mstype,&match);CHKERRQ(ierr); 427*37e1895aSJed Brown if (match) { 428*37e1895aSJed Brown ierr = SNESReset_MS(snes);CHKERRQ(ierr); 429*37e1895aSJed Brown ms->tableau = &link->tab; 430*37e1895aSJed Brown PetscFunctionReturn(0); 431*37e1895aSJed Brown } 432*37e1895aSJed Brown } 433*37e1895aSJed Brown SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mstype); 434*37e1895aSJed Brown PetscFunctionReturn(0); 435*37e1895aSJed Brown } 436*37e1895aSJed Brown EXTERN_C_END 437*37e1895aSJed Brown 438*37e1895aSJed Brown 439*37e1895aSJed Brown #undef __FUNCT__ 440*37e1895aSJed Brown #define __FUNCT__ "SNESMSSetType" 441*37e1895aSJed Brown /*@C 442*37e1895aSJed Brown SNESMSSetType - Set the type of multistage smoother 443*37e1895aSJed Brown 444*37e1895aSJed Brown Logically collective 445*37e1895aSJed Brown 446*37e1895aSJed Brown Input Parameter: 447*37e1895aSJed Brown + snes - nonlinear solver context 448*37e1895aSJed Brown - mstype - type of multistage method 449*37e1895aSJed Brown 450*37e1895aSJed Brown Level: beginner 451*37e1895aSJed Brown 452*37e1895aSJed Brown .seealso: SNESMSGetType(), SNESMS 453*37e1895aSJed Brown @*/ 454*37e1895aSJed Brown PetscErrorCode SNESMSSetType(SNES snes,const SNESMSType rostype) 455*37e1895aSJed Brown { 456*37e1895aSJed Brown PetscErrorCode ierr; 457*37e1895aSJed Brown 458*37e1895aSJed Brown PetscFunctionBegin; 459*37e1895aSJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 460*37e1895aSJed Brown ierr = PetscTryMethod(snes,"SNESMSSetType_C",(SNES,const SNESMSType),(snes,rostype));CHKERRQ(ierr); 461*37e1895aSJed Brown PetscFunctionReturn(0); 462*37e1895aSJed Brown } 463*37e1895aSJed Brown 464*37e1895aSJed Brown /* -------------------------------------------------------------------------- */ 465*37e1895aSJed Brown /*MC 466*37e1895aSJed Brown SNESMS - multi-stage smoothers 467*37e1895aSJed Brown 468*37e1895aSJed Brown Options Database: 469*37e1895aSJed Brown 470*37e1895aSJed Brown + -snes_ms_type - type of multi-stage smoother 471*37e1895aSJed Brown - -snes_ms_damping - damping for multi-stage method 472*37e1895aSJed Brown 473*37e1895aSJed Brown Notes: 474*37e1895aSJed Brown These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for 475*37e1895aSJed Brown FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebychev). 476*37e1895aSJed Brown 477*37e1895aSJed Brown Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds. 478*37e1895aSJed Brown 479*37e1895aSJed Brown The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with SNESMSRegister(). 480*37e1895aSJed Brown 481*37e1895aSJed Brown References: 482*37e1895aSJed Brown 483*37e1895aSJed Brown Ketcheson (2010) Runge-Kutta methods with minimum storage implementations. 484*37e1895aSJed Brown 485*37e1895aSJed Brown Jameson (1983) Solution of the Euler equations for two dimensional transonic flow by a multigrid method. 486*37e1895aSJed Brown 487*37e1895aSJed Brown Pierce and Giles (1997) Preconditioned multigrid methods for compressible flow calculations on stretched meshes. 488*37e1895aSJed Brown 489*37e1895aSJed Brown Level: beginner 490*37e1895aSJed Brown 491*37e1895aSJed Brown .seealso: SNESCreate(), SNES, SNESSetType(), SNESMS, SNESFAS, KSPCHEBYCHEV 492*37e1895aSJed Brown 493*37e1895aSJed Brown M*/ 494*37e1895aSJed Brown EXTERN_C_BEGIN 495*37e1895aSJed Brown #undef __FUNCT__ 496*37e1895aSJed Brown #define __FUNCT__ "SNESCreate_MS" 497*37e1895aSJed Brown PetscErrorCode SNESCreate_MS(SNES snes) 498*37e1895aSJed Brown { 499*37e1895aSJed Brown PetscErrorCode ierr; 500*37e1895aSJed Brown SNES_MS *ms; 501*37e1895aSJed Brown 502*37e1895aSJed Brown PetscFunctionBegin; 503*37e1895aSJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 504*37e1895aSJed Brown ierr = SNESMSInitializePackage(PETSC_NULL);CHKERRQ(ierr); 505*37e1895aSJed Brown #endif 506*37e1895aSJed Brown 507*37e1895aSJed Brown snes->ops->setup = SNESSetUp_MS; 508*37e1895aSJed Brown snes->ops->solve = SNESSolve_MS; 509*37e1895aSJed Brown snes->ops->destroy = SNESDestroy_MS; 510*37e1895aSJed Brown snes->ops->setfromoptions = SNESSetFromOptions_MS; 511*37e1895aSJed Brown snes->ops->view = SNESView_MS; 512*37e1895aSJed Brown snes->ops->reset = SNESReset_MS; 513*37e1895aSJed Brown 514*37e1895aSJed Brown snes->usespc = PETSC_FALSE; 515*37e1895aSJed Brown snes->usesksp = PETSC_TRUE; 516*37e1895aSJed Brown 517*37e1895aSJed Brown ierr = PetscNewLog(snes,SNES_MS,&ms);CHKERRQ(ierr); 518*37e1895aSJed Brown snes->data = (void*)ms; 519*37e1895aSJed Brown ms->damping = 0.9; 520*37e1895aSJed Brown ms->norms = PETSC_FALSE; 521*37e1895aSJed Brown 522*37e1895aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESMSSetType_C","SNESMSSetType_MS",SNESMSSetType_MS);CHKERRQ(ierr); 523*37e1895aSJed Brown PetscFunctionReturn(0); 524*37e1895aSJed Brown } 525*37e1895aSJed Brown EXTERN_C_END 526