1*c4762a1bSJed Brown #include "finitevolume1d.h" 2*c4762a1bSJed Brown #include <petscdmda.h> 3*c4762a1bSJed Brown #include <petscdraw.h> 4*c4762a1bSJed Brown #include <petsc/private/tsimpl.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown #include <petsc/private/kernels/blockinvert.h> /* For the Kernel_*_gets_* stuff for BAIJ */ 7*c4762a1bSJed Brown const char *FVBCTypes[] = {"PERIODIC","OUTFLOW","FVBCType","FVBC_",0}; 8*c4762a1bSJed Brown 9*c4762a1bSJed Brown PETSC_STATIC_INLINE PetscReal Sgn(PetscReal a) { return (a<0) ? -1 : 1; } 10*c4762a1bSJed Brown PETSC_STATIC_INLINE PetscReal Abs(PetscReal a) { return (a<0) ? 0 : a; } 11*c4762a1bSJed Brown PETSC_STATIC_INLINE PetscReal Sqr(PetscReal a) { return a*a; } 12*c4762a1bSJed Brown //PETSC_STATIC_INLINE PetscReal MaxAbs(PetscReal a,PetscReal b) { return (PetscAbs(a) > PetscAbs(b)) ? a : b; } 13*c4762a1bSJed Brown PETSC_UNUSED PETSC_STATIC_INLINE PetscReal MinAbs(PetscReal a,PetscReal b) { return (PetscAbs(a) < PetscAbs(b)) ? a : b; } 14*c4762a1bSJed Brown PETSC_STATIC_INLINE PetscReal MinMod2(PetscReal a,PetscReal b) { return (a*b<0) ? 0 : Sgn(a)*PetscMin(PetscAbs(a),PetscAbs(b)); } 15*c4762a1bSJed Brown PETSC_STATIC_INLINE PetscReal MaxMod2(PetscReal a,PetscReal b) { return (a*b<0) ? 0 : Sgn(a)*PetscMax(PetscAbs(a),PetscAbs(b)); } 16*c4762a1bSJed Brown PETSC_STATIC_INLINE PetscReal MinMod3(PetscReal a,PetscReal b,PetscReal c) {return (a*b<0 || a*c<0) ? 0 : Sgn(a)*PetscMin(PetscAbs(a),PetscMin(PetscAbs(b),PetscAbs(c))); } 17*c4762a1bSJed Brown 18*c4762a1bSJed Brown //PETSC_STATIC_INLINE PetscReal RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax) { PetscReal range = xmax-xmin; return xmin +PetscFmodReal(range+PetscFmodReal(a,range),range); } 19*c4762a1bSJed Brown 20*c4762a1bSJed Brown 21*c4762a1bSJed Brown /* ----------------------- Lots of limiters, these could go in a separate library ------------------------- */ 22*c4762a1bSJed Brown void Limit_Upwind(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) 23*c4762a1bSJed Brown { 24*c4762a1bSJed Brown PetscInt i; 25*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0; 26*c4762a1bSJed Brown } 27*c4762a1bSJed Brown void Limit_LaxWendroff(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) 28*c4762a1bSJed Brown { 29*c4762a1bSJed Brown PetscInt i; 30*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jR[i]; 31*c4762a1bSJed Brown } 32*c4762a1bSJed Brown void Limit_BeamWarming(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) 33*c4762a1bSJed Brown { 34*c4762a1bSJed Brown PetscInt i; 35*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jL[i]; 36*c4762a1bSJed Brown } 37*c4762a1bSJed Brown void Limit_Fromm(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) 38*c4762a1bSJed Brown { 39*c4762a1bSJed Brown PetscInt i; 40*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]+jR[i]); 41*c4762a1bSJed Brown } 42*c4762a1bSJed Brown void Limit_Minmod(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) 43*c4762a1bSJed Brown { 44*c4762a1bSJed Brown PetscInt i; 45*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i],jR[i]); 46*c4762a1bSJed Brown } 47*c4762a1bSJed Brown void Limit_Superbee(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) 48*c4762a1bSJed Brown { 49*c4762a1bSJed Brown PetscInt i; 50*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i],2*jR[i]),MinMod2(2*jL[i],jR[i])); 51*c4762a1bSJed Brown } 52*c4762a1bSJed Brown void Limit_MC(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) 53*c4762a1bSJed Brown { 54*c4762a1bSJed Brown PetscInt i; 55*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],0.5*(jL[i]+jR[i]),2*jR[i]); 56*c4762a1bSJed Brown } 57*c4762a1bSJed Brown void Limit_VanLeer(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) 58*c4762a1bSJed Brown { /* phi = (t + abs(t)) / (1 + abs(t)) */ 59*c4762a1bSJed Brown PetscInt i; 60*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = (jL[i]*Abs(jR[i])+Abs(jL[i])*jR[i])/(Abs(jL[i])+Abs(jR[i])+1e-15); 61*c4762a1bSJed Brown } 62*c4762a1bSJed Brown void Limit_VanAlbada(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) /* differentiable */ 63*c4762a1bSJed Brown { /* phi = (t + t^2) / (1 + t^2) */ 64*c4762a1bSJed Brown PetscInt i; 65*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = (jL[i]*Sqr(jR[i])+Sqr(jL[i])*jR[i])/(Sqr(jL[i])+Sqr(jR[i])+1e-15); 66*c4762a1bSJed Brown } 67*c4762a1bSJed Brown void Limit_VanAlbadaTVD(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) 68*c4762a1bSJed Brown { /* phi = (t + t^2) / (1 + t^2) */ 69*c4762a1bSJed Brown PetscInt i; 70*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = (jL[i]*jR[i]<0) ? 0 : (jL[i]*Sqr(jR[i])+Sqr(jL[i])*jR[i])/(Sqr(jL[i])+Sqr(jR[i])+1e-15); 71*c4762a1bSJed Brown } 72*c4762a1bSJed Brown void Limit_Koren(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) /* differentiable */ 73*c4762a1bSJed Brown { /* phi = (t + 2*t^2) / (2 - t + 2*t^2) */ 74*c4762a1bSJed Brown PetscInt i; 75*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = ((jL[i]*Sqr(jR[i])+2*Sqr(jL[i])*jR[i])/(2*Sqr(jL[i])-jL[i]*jR[i]+2*Sqr(jR[i])+1e-15)); 76*c4762a1bSJed Brown } 77*c4762a1bSJed Brown void Limit_KorenSym(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) /* differentiable */ 78*c4762a1bSJed Brown { /* Symmetric version of above */ 79*c4762a1bSJed Brown PetscInt i; 80*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = (1.5*(jL[i]*Sqr(jR[i])+Sqr(jL[i])*jR[i])/(2*Sqr(jL[i])-jL[i]*jR[i]+2*Sqr(jR[i])+1e-15)); 81*c4762a1bSJed Brown } 82*c4762a1bSJed Brown void Limit_Koren3(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) 83*c4762a1bSJed Brown { /* Eq 11 of Cada-Torrilhon 2009 */ 84*c4762a1bSJed Brown PetscInt i; 85*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],(jL[i]+2*jR[i])/3,2*jR[i]); 86*c4762a1bSJed Brown } 87*c4762a1bSJed Brown static PetscReal CadaTorrilhonPhiHatR_Eq13(PetscReal L,PetscReal R) 88*c4762a1bSJed Brown { 89*c4762a1bSJed Brown return PetscMax(0,PetscMin((L+2*R)/3,PetscMax(-0.5*L,PetscMin(2*L,PetscMin((L+2*R)/3,1.6*R))))); 90*c4762a1bSJed Brown } 91*c4762a1bSJed Brown void Limit_CadaTorrilhon2(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) 92*c4762a1bSJed Brown { /* Cada-Torrilhon 2009, Eq 13 */ 93*c4762a1bSJed Brown PetscInt i; 94*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = CadaTorrilhonPhiHatR_Eq13(jL[i],jR[i]); 95*c4762a1bSJed Brown } 96*c4762a1bSJed Brown void Limit_CadaTorrilhon3R(PetscReal r,LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) 97*c4762a1bSJed Brown { /* Cada-Torrilhon 2009, Eq 22 */ 98*c4762a1bSJed Brown /* They recommend 0.001 < r < 1, but larger values are more accurate in smooth regions */ 99*c4762a1bSJed Brown const PetscReal eps = 1e-7,hx = info->hx; 100*c4762a1bSJed Brown PetscInt i; 101*c4762a1bSJed Brown for (i=0; i<info->m; i++) { 102*c4762a1bSJed Brown const PetscReal eta = (Sqr(jL[i])+Sqr(jR[i]))/Sqr(r*hx); 103*c4762a1bSJed Brown lmt[i] = ((eta < 1-eps) ? (jL[i]+2*jR[i])/3 : ((eta>1+eps) ? CadaTorrilhonPhiHatR_Eq13(jL[i],jR[i]) : 0.5*((1-(eta-1)/eps)*(jL[i]+2*jR[i])/3+(1+(eta+1)/eps)*CadaTorrilhonPhiHatR_Eq13(jL[i],jR[i])))); 104*c4762a1bSJed Brown } 105*c4762a1bSJed Brown } 106*c4762a1bSJed Brown void Limit_CadaTorrilhon3R0p1(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) 107*c4762a1bSJed Brown { 108*c4762a1bSJed Brown Limit_CadaTorrilhon3R(0.1,info,jL,jR,lmt); 109*c4762a1bSJed Brown } 110*c4762a1bSJed Brown void Limit_CadaTorrilhon3R1(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) 111*c4762a1bSJed Brown { 112*c4762a1bSJed Brown Limit_CadaTorrilhon3R(1,info,jL,jR,lmt); 113*c4762a1bSJed Brown } 114*c4762a1bSJed Brown void Limit_CadaTorrilhon3R10(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) 115*c4762a1bSJed Brown { 116*c4762a1bSJed Brown Limit_CadaTorrilhon3R(10,info,jL,jR,lmt); 117*c4762a1bSJed Brown } 118*c4762a1bSJed Brown void Limit_CadaTorrilhon3R100(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) 119*c4762a1bSJed Brown { 120*c4762a1bSJed Brown Limit_CadaTorrilhon3R(100,info,jL,jR,lmt); 121*c4762a1bSJed Brown } 122*c4762a1bSJed Brown 123*c4762a1bSJed Brown /* ----------------------- Limiters for split systems ------------------------- */ 124*c4762a1bSJed Brown 125*c4762a1bSJed Brown void Limit2_Upwind(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,const PetscInt sf,const PetscInt fs,PetscInt n,PetscScalar *lmt) 126*c4762a1bSJed Brown { 127*c4762a1bSJed Brown PetscInt i; 128*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0; 129*c4762a1bSJed Brown } 130*c4762a1bSJed Brown void Limit2_LaxWendroff(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,const PetscInt sf,const PetscInt fs,PetscInt n,PetscScalar *lmt) 131*c4762a1bSJed Brown { 132*c4762a1bSJed Brown PetscInt i; 133*c4762a1bSJed Brown if (n < sf-1) { /* slow components */ 134*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jR[i]/info->hxs; 135*c4762a1bSJed Brown } else if (n == sf-1) { /* slow component which is next to fast components */ 136*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jR[i]/(info->hxs/2.0+info->hxf/2.0); 137*c4762a1bSJed Brown } else if (n == sf) { /* fast component which is next to slow components */ 138*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jR[i]/info->hxf; 139*c4762a1bSJed Brown } else if (n > sf && n < fs-1) { /* fast components */ 140*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jR[i]/info->hxf; 141*c4762a1bSJed Brown } else if (n == fs-1) { /* fast component next to slow components */ 142*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jR[i]/(info->hxf/2.0+info->hxs/2.0); 143*c4762a1bSJed Brown } else if (n == fs) { /* slow component next to fast components */ 144*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jR[i]/info->hxs; 145*c4762a1bSJed Brown } else { /* slow components */ 146*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jR[i]/info->hxs; 147*c4762a1bSJed Brown } 148*c4762a1bSJed Brown } 149*c4762a1bSJed Brown void Limit2_BeamWarming(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,const PetscInt sf,const PetscInt fs,PetscInt n,PetscScalar *lmt) 150*c4762a1bSJed Brown { 151*c4762a1bSJed Brown PetscInt i; 152*c4762a1bSJed Brown if (n < sf-1) { 153*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jL[i]/info->hxs; 154*c4762a1bSJed Brown } else if (n == sf-1) { 155*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jL[i]/info->hxs; 156*c4762a1bSJed Brown } else if (n == sf) { 157*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jL[i]/(info->hxs/2.0+info->hxf/2.0); 158*c4762a1bSJed Brown } else if (n > sf && n < fs-1) { 159*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jL[i]/info->hxf; 160*c4762a1bSJed Brown } else if (n == fs-1) { 161*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jL[i]/info->hxf; 162*c4762a1bSJed Brown } else if (n == fs) { 163*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jL[i]/(info->hxf/2.0+info->hxs/2.0); 164*c4762a1bSJed Brown } else { 165*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jL[i]/info->hxs; 166*c4762a1bSJed Brown } 167*c4762a1bSJed Brown } 168*c4762a1bSJed Brown void Limit2_Fromm(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,const PetscInt sf,const PetscInt fs,PetscInt n,PetscScalar *lmt) 169*c4762a1bSJed Brown { 170*c4762a1bSJed Brown PetscInt i; 171*c4762a1bSJed Brown if (n < sf-1) { 172*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]+jR[i])/info->hxs; 173*c4762a1bSJed Brown } else if (n == sf-1) { 174*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]/info->hxs+jR[i]/(info->hxs/2.0+info->hxf/2.0)); 175*c4762a1bSJed Brown } else if (n == sf) { 176*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]/(info->hxs/2.0+info->hxf/2.0)+jR[i]/info->hxf); 177*c4762a1bSJed Brown } else if (n > sf && n < fs-1) { 178*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]+jR[i])/info->hxf; 179*c4762a1bSJed Brown } else if (n == fs-1) { 180*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]/info->hxf+jR[i]/(info->hxf/2.0+info->hxs/2.0)); 181*c4762a1bSJed Brown } else if (n == fs) { 182*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]/(info->hxf/2.0+info->hxs/2.0)+jR[i]/info->hxs); 183*c4762a1bSJed Brown } else { 184*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]+jR[i])/info->hxs; 185*c4762a1bSJed Brown } 186*c4762a1bSJed Brown } 187*c4762a1bSJed Brown void Limit2_Minmod(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,const PetscInt sf,const PetscInt fs,PetscInt n,PetscScalar *lmt) 188*c4762a1bSJed Brown { 189*c4762a1bSJed Brown PetscInt i; 190*c4762a1bSJed Brown if (n < sf-1) { 191*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i],jR[i])/info->hxs; 192*c4762a1bSJed Brown } else if (n == sf-1) { 193*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i]/info->hxs,jR[i]/(info->hxs/2.0+info->hxf/2.0)); 194*c4762a1bSJed Brown } else if (n == sf) { 195*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i]/(info->hxs/2.0+info->hxf/2.0),jR[i]/info->hxf); 196*c4762a1bSJed Brown } else if (n > sf && n < fs-1) { 197*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i],jR[i])/info->hxf; 198*c4762a1bSJed Brown } else if (n == fs-1) { 199*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i]/info->hxf,jR[i]/(info->hxf/2.0+info->hxs/2.0)); 200*c4762a1bSJed Brown } else if (n == fs) { 201*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i]/(info->hxf/2.0+info->hxs/2.0),jR[i]/info->hxs); 202*c4762a1bSJed Brown } else { 203*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i],jR[i])/info->hxs; 204*c4762a1bSJed Brown } 205*c4762a1bSJed Brown } 206*c4762a1bSJed Brown void Limit2_Superbee(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,const PetscInt sf,const PetscInt fs,PetscInt n,PetscScalar *lmt) 207*c4762a1bSJed Brown { 208*c4762a1bSJed Brown PetscInt i; 209*c4762a1bSJed Brown if (n < sf-1) { 210*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i],2*jR[i]),MinMod2(2*jL[i],jR[i]))/info->hxs; 211*c4762a1bSJed Brown } else if (n == sf-1) { 212*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i]/info->hxs,2*jR[i]/(info->hxs/2.0+info->hxf/2.0)),MinMod2(2*jL[i]/info->hxs,jR[i]/(info->hxs/2.0+info->hxf/2.0))); 213*c4762a1bSJed Brown } else if (n == sf) { 214*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i]/(info->hxs/2.0+info->hxf/2.0),2*jR[i]/info->hxf),MinMod2(2*jL[i]/(info->hxs/2.0+info->hxf/2.0),jR[i]/info->hxf)); 215*c4762a1bSJed Brown } else if (n > sf && n < fs-1) { 216*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(MinMod2(jL[i],2*jR[i]),MinMod2(2*jL[i],jR[i]))/info->hxf; 217*c4762a1bSJed Brown } else if (n == fs-1) { 218*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(MinMod2(jL[i]/info->hxf,2*jR[i]/(info->hxf/2.0+info->hxs/2.0)),MinMod2(2*jL[i]/info->hxf,jR[i]/(info->hxf/2.0+info->hxs/2.0))); 219*c4762a1bSJed Brown } else if (n == fs) { 220*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(MinMod2(jL[i]/(info->hxf/2.0+info->hxs/2.0),2*jR[i]/info->hxs),MinMod2(2*jL[i]/(info->hxf/2.0+info->hxs/2.0),jR[i]/info->hxs)); 221*c4762a1bSJed Brown } else { 222*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(MinMod2(jL[i],2*jR[i]),MinMod2(2*jL[i],jR[i]))/info->hxs; 223*c4762a1bSJed Brown } 224*c4762a1bSJed Brown } 225*c4762a1bSJed Brown void Limit2_MC(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,const PetscInt sf,const PetscInt fs,PetscInt n,PetscScalar *lmt) 226*c4762a1bSJed Brown { 227*c4762a1bSJed Brown PetscInt i; 228*c4762a1bSJed Brown if (n < sf-1) { 229*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],0.5*(jL[i]+jR[i]),2*jR[i])/info->hxs; 230*c4762a1bSJed Brown } else if (n == sf-1) { 231*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/info->hxs,0.5*(jL[i]/info->hxs+jR[i]/(info->hxf/2.0+info->hxs/2.0)),2*jR[i]/(info->hxf/2.0+info->hxs/2.0)); 232*c4762a1bSJed Brown } else if (n == sf) { 233*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/(info->hxs/2.0+info->hxf/2.0),0.5*(jL[i]/(info->hxs/2.0+info->hxf/2.0)+jR[i]/info->hxf),2*jR[i]/info->hxf); 234*c4762a1bSJed Brown } else if (n > sf && n < fs-1) { 235*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],0.5*(jL[i]+jR[i]),2*jR[i])/info->hxf; 236*c4762a1bSJed Brown } else if (n == fs-1) { 237*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/info->hxf,0.5*(jL[i]/info->hxf+jR[i]/(info->hxf/2.0+info->hxs/2.0)),2*jR[i]/(info->hxf/2.0+info->hxs/2.0)); 238*c4762a1bSJed Brown } else if (n == fs) { 239*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/(info->hxf/2.0+info->hxs/2.0),0.5*(jL[i]/(info->hxf/2.0+info->hxs/2.0)+jR[i]/info->hxs),2*jR[i]/info->hxs); 240*c4762a1bSJed Brown } else { 241*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],0.5*(jL[i]+jR[i]),2*jR[i])/info->hxs; 242*c4762a1bSJed Brown } 243*c4762a1bSJed Brown } 244*c4762a1bSJed Brown void Limit2_Koren3(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,const PetscInt sf,const PetscInt fs,PetscInt n,PetscScalar *lmt) 245*c4762a1bSJed Brown { /* Eq 11 of Cada-Torrilhon 2009 */ 246*c4762a1bSJed Brown PetscInt i; 247*c4762a1bSJed Brown if (n < sf-1) { 248*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],(jL[i]+2*jR[i])/3,2*jR[i])/info->hxs; 249*c4762a1bSJed Brown } else if (n == sf-1) { 250*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/info->hxs,(jL[i]/info->hxs+2*jR[i]/(info->hxf/2.0+info->hxs/2.0))/3,2*jR[i]/(info->hxf/2.0+info->hxs/2.0)); 251*c4762a1bSJed Brown } else if (n == sf) { 252*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/(info->hxs/2.0+info->hxf/2.0),(jL[i]/(info->hxs/2.0+info->hxf/2.0)+2*jR[i]/info->hxf)/3,2*jR[i]/info->hxf); 253*c4762a1bSJed Brown } else if (n > sf && n < fs-1) { 254*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],(jL[i]+2*jR[i])/3,2*jR[i])/info->hxf; 255*c4762a1bSJed Brown } else if (n == fs-1) { 256*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/info->hxf,(jL[i]/info->hxf+2*jR[i]/(info->hxf/2.0+info->hxs/2.0))/3,2*jR[i]/(info->hxf/2.0+info->hxs/2.0)); 257*c4762a1bSJed Brown } else if (n == fs) { 258*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/(info->hxf/2.0+info->hxs/2.0),(jL[i]/(info->hxf/2.0+info->hxs/2.0)+2*jR[i]/info->hxs)/3,2*jR[i]/info->hxs); 259*c4762a1bSJed Brown } else { 260*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],(jL[i]+2*jR[i])/3,2*jR[i])/info->hxs; 261*c4762a1bSJed Brown } 262*c4762a1bSJed Brown } 263*c4762a1bSJed Brown 264*c4762a1bSJed Brown /* ---- Three-way splitting ---- */ 265*c4762a1bSJed Brown void Limit3_Upwind(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,const PetscInt sm,const PetscInt mf,const PetscInt fm,const PetscInt ms,PetscInt n,PetscScalar *lmt) 266*c4762a1bSJed Brown { 267*c4762a1bSJed Brown PetscInt i; 268*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0; 269*c4762a1bSJed Brown } 270*c4762a1bSJed Brown void Limit3_LaxWendroff(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,const PetscInt sm,const PetscInt mf,const PetscInt fm,const PetscInt ms,PetscInt n,PetscScalar *lmt) 271*c4762a1bSJed Brown { 272*c4762a1bSJed Brown PetscInt i; 273*c4762a1bSJed Brown if (n < sm-1 || n > ms) { /* slow components */ 274*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jR[i]/info->hxs; 275*c4762a1bSJed Brown } else if (n == sm-1 || n == ms-1) { /* slow component which is next to medium components */ 276*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jR[i]/(info->hxs/2.0+info->hxm/2.0); 277*c4762a1bSJed Brown } else if (n < mf-1 || n > fm) { /* medium components */ 278*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jR[i]/info->hxm; 279*c4762a1bSJed Brown } else if (n == mf-1 || n == fm) { /* medium component next to fast components */ 280*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jR[i]/(info->hxm/2.0+info->hxf/2.0); 281*c4762a1bSJed Brown } else { /* fast components */ 282*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jR[i]/info->hxf; 283*c4762a1bSJed Brown } 284*c4762a1bSJed Brown } 285*c4762a1bSJed Brown void Limit3_BeamWarming(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,const PetscInt sm,const PetscInt mf,const PetscInt fm,const PetscInt ms,PetscInt n,PetscScalar *lmt) 286*c4762a1bSJed Brown { 287*c4762a1bSJed Brown PetscInt i; 288*c4762a1bSJed Brown if (n < sm || n > ms) { 289*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jL[i]/info->hxs; 290*c4762a1bSJed Brown } else if (n == sm || n == ms) { 291*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jL[i]/(info->hxs/2.0+info->hxf/2.0); 292*c4762a1bSJed Brown }else if (n < mf || n > fm) { 293*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jL[i]/info->hxm; 294*c4762a1bSJed Brown } else if (n == mf || n == fm) { 295*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jL[i]/(info->hxm/2.0+info->hxf/2.0); 296*c4762a1bSJed Brown } else { 297*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = jL[i]/info->hxf; 298*c4762a1bSJed Brown } 299*c4762a1bSJed Brown } 300*c4762a1bSJed Brown void Limit3_Fromm(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,const PetscInt sm,const PetscInt mf,const PetscInt fm,const PetscInt ms,PetscInt n,PetscScalar *lmt) 301*c4762a1bSJed Brown { 302*c4762a1bSJed Brown PetscInt i; 303*c4762a1bSJed Brown if (n < sm-1 || n > ms) { 304*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]+jR[i])/info->hxs; 305*c4762a1bSJed Brown } else if (n == sm-1) { 306*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]/info->hxs+jR[i]/(info->hxs/2.0+info->hxf/2.0)); 307*c4762a1bSJed Brown } else if (n == sm) { 308*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]/(info->hxs/2.0+info->hxm/2.0)+jR[i]/info->hxm); 309*c4762a1bSJed Brown } else if (n == ms-1) { 310*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]/info->hxm+jR[i]/(info->hxs/2.0+info->hxf/2.0)); 311*c4762a1bSJed Brown } else if (n == ms) { 312*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]/(info->hxm/2.0+info->hxs/2.0)+jR[i]/info->hxs); 313*c4762a1bSJed Brown } else if (n < mf-1 || n > fm) { 314*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]+jR[i])/info->hxm; 315*c4762a1bSJed Brown } else if (n == mf-1) { 316*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]/info->hxm+jR[i]/(info->hxm/2.0+info->hxf/2.0)); 317*c4762a1bSJed Brown } else if (n == mf) { 318*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]/(info->hxm/2.0+info->hxf/2.0)+jR[i]/info->hxf); 319*c4762a1bSJed Brown } else if (n == fm-1) { 320*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]/info->hxf+jR[i]/(info->hxf/2.0+info->hxm/2.0)); 321*c4762a1bSJed Brown } else if (n == fm) { 322*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]/(info->hxf/2.0+info->hxm/2.0)+jR[i]/info->hxm); 323*c4762a1bSJed Brown } else { 324*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]+jR[i])/info->hxf; 325*c4762a1bSJed Brown } 326*c4762a1bSJed Brown } 327*c4762a1bSJed Brown void Limit3_Minmod(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,const PetscInt sm,const PetscInt mf,const PetscInt fm,const PetscInt ms,PetscInt n,PetscScalar *lmt) 328*c4762a1bSJed Brown { 329*c4762a1bSJed Brown PetscInt i; 330*c4762a1bSJed Brown if (n < sm-1 || n > ms) { 331*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i],jR[i])/info->hxs; 332*c4762a1bSJed Brown } else if (n == sm-1) { 333*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i]/info->hxs,jR[i]/(info->hxs/2.0+info->hxf/2.0)); 334*c4762a1bSJed Brown } else if (n == sm) { 335*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i]/(info->hxs/2.0+info->hxf/2.0),jR[i]/info->hxf); 336*c4762a1bSJed Brown } else if (n == ms-1) { 337*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i]/info->hxm,jR[i]/(info->hxm/2.0+info->hxs/2.0)); 338*c4762a1bSJed Brown } else if (n == ms) { 339*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i]/(info->hxm/2.0+info->hxs/2.0),jR[i]/info->hxs); 340*c4762a1bSJed Brown } else if (n < mf-1 || n > fm) { 341*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i],jR[i])/info->hxm; 342*c4762a1bSJed Brown } else if (n == mf-1) { 343*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i]/info->hxm,jR[i]/(info->hxm/2.0+info->hxf/2.0)); 344*c4762a1bSJed Brown } else if (n == mf) { 345*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i]/(info->hxm/2.0+info->hxf/2.0),jR[i]/info->hxf); 346*c4762a1bSJed Brown } else if (n == fm-1) { 347*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i]/info->hxf,jR[i]/(info->hxf/2.0+info->hxm/2.0)); 348*c4762a1bSJed Brown } else if (n == fm) { 349*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i]/(info->hxf/2.0+info->hxm/2.0),jR[i]/info->hxm); 350*c4762a1bSJed Brown } else { 351*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i]+jR[i])/info->hxf; 352*c4762a1bSJed Brown } 353*c4762a1bSJed Brown } 354*c4762a1bSJed Brown void Limit3_Superbee(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,const PetscInt sm,const PetscInt mf,const PetscInt fm,const PetscInt ms,PetscInt n,PetscScalar *lmt) 355*c4762a1bSJed Brown { 356*c4762a1bSJed Brown PetscInt i; 357*c4762a1bSJed Brown if (n < sm-1 || n > ms) { 358*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i],2*jR[i]),MinMod2(2*jL[i],jR[i]))/info->hxs; 359*c4762a1bSJed Brown } else if (n == sm-1) { 360*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i]/info->hxs,2*jR[i]/(info->hxs/2.0+info->hxm/2.0)),MinMod2(2*jL[i]/info->hxs,jR[i]/(info->hxs/2.0+info->hxm/2.0))); 361*c4762a1bSJed Brown } else if (n == sm) { 362*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i]/(info->hxs/2.0+info->hxm/2.0),2*jR[i]/info->hxm),MinMod2(2*jL[i]/(info->hxs/2.0+info->hxm/2.0),jR[i]/info->hxm)); 363*c4762a1bSJed Brown } else if (n == ms-1) { 364*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(MinMod2(jL[i]/info->hxm,2*jR[i]/(info->hxm/2.0+info->hxs/2.0)),MinMod2(2*jL[i]/info->hxm,jR[i]/(info->hxm/2.0+info->hxs/2.0))); 365*c4762a1bSJed Brown } else if (n == ms) { 366*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(MinMod2(jL[i]/(info->hxm/2.0+info->hxs/2.0),2*jR[i]/info->hxs),MinMod2(2*jL[i]/(info->hxm/2.0+info->hxs/2.0),jR[i]/info->hxs)); 367*c4762a1bSJed Brown } else if (n < mf-1 || n > fm) { 368*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i],2*jR[i]),MinMod2(2*jL[i],jR[i]))/info->hxm; 369*c4762a1bSJed Brown } else if (n == mf-1) { 370*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i]/info->hxm,2*jR[i]/(info->hxm/2.0+info->hxf/2.0)),MinMod2(2*jL[i]/info->hxm,jR[i]/(info->hxm/2.0+info->hxf/2.0))); 371*c4762a1bSJed Brown } else if (n == mf) { 372*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i]/(info->hxm/2.0+info->hxf/2.0),2*jR[i]/info->hxf),MinMod2(2*jL[i]/(info->hxm/2.0+info->hxf/2.0),jR[i]/info->hxf)); 373*c4762a1bSJed Brown } else if (n == fm-1) { 374*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(MinMod2(jL[i]/info->hxf,2*jR[i]/(info->hxf/2.0+info->hxm/2.0)),MinMod2(2*jL[i]/info->hxf,jR[i]/(info->hxf/2.0+info->hxm/2.0))); 375*c4762a1bSJed Brown } else if (n == fm) { 376*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod2(MinMod2(jL[i]/(info->hxf/2.0+info->hxm/2.0),2*jR[i]/info->hxm),MinMod2(2*jL[i]/(info->hxf/2.0+info->hxm/2.0),jR[i]/info->hxm)); 377*c4762a1bSJed Brown } else { 378*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i],2*jR[i]),MinMod2(2*jL[i],jR[i]))/info->hxf; 379*c4762a1bSJed Brown } 380*c4762a1bSJed Brown } 381*c4762a1bSJed Brown void Limit3_MC(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,const PetscInt sm,const PetscInt mf,const PetscInt fm,const PetscInt ms,PetscInt n,PetscScalar *lmt) 382*c4762a1bSJed Brown { 383*c4762a1bSJed Brown PetscInt i; 384*c4762a1bSJed Brown if (n < sm-1 || n > ms) { 385*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],0.5*(jL[i]+jR[i]),2*jR[i])/info->hxs; 386*c4762a1bSJed Brown } else if (n == sm-1) { 387*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/info->hxs,0.5*(jL[i]/info->hxs+jR[i]/(info->hxs/2.0+info->hxm/2.0)),2*jR[i]/(info->hxs/2.0+info->hxm/2.0)); 388*c4762a1bSJed Brown } else if (n == sm) { 389*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/(info->hxs/2.0+info->hxm/2.0),0.5*(jL[i]/(info->hxs/2.0+info->hxm/2.0)+jR[i]/info->hxm),2*jR[i]/info->hxm); 390*c4762a1bSJed Brown } else if (n == ms-1) { 391*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/info->hxm,0.5*(jL[i]/info->hxm+jR[i]/(info->hxm/2.0+info->hxs/2.0)),2*jR[i]/(info->hxm/2.0+info->hxs/2.0)); 392*c4762a1bSJed Brown } else if (n == ms) { 393*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/(info->hxm/2.0+info->hxs/2.0),0.5*(jL[i]/(info->hxm/2.0+info->hxs/2.0)+jR[i]/info->hxs),2*jR[i]/info->hxs); 394*c4762a1bSJed Brown } else if (n < mf-1 || n > fm) { 395*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],0.5*(jL[i]+jR[i]),2*jR[i])/info->hxm; 396*c4762a1bSJed Brown } else if (n == mf-1) { 397*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/info->hxm,0.5*(jL[i]/info->hxm+jR[i]/(info->hxm/2.0+info->hxf/2.0)),2*jR[i]/(info->hxm/2.0+info->hxf/2.0)); 398*c4762a1bSJed Brown } else if (n == mf) { 399*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/(info->hxm/2.0+info->hxf/2.0),0.5*(jL[i]/(info->hxm/2.0+info->hxf/2.0)+jR[i]/info->hxf),2*jR[i]/info->hxf); 400*c4762a1bSJed Brown } else if (n == fm-1) { 401*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/info->hxf,0.5*(jL[i]/info->hxf+jR[i]/(info->hxf/2.0+info->hxm/2.0)),2*jR[i]/(info->hxf/2.0+info->hxm/2.0)); 402*c4762a1bSJed Brown } else if (n == fm) { 403*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/(info->hxf/2.0+info->hxm/2.0),0.5*(jL[i]/(info->hxf/2.0+info->hxm/2.0)+jR[i]/info->hxm),2*jR[i]/info->hxm); 404*c4762a1bSJed Brown } else { 405*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],0.5*(jL[i]+jR[i]),2*jR[i])/info->hxf; 406*c4762a1bSJed Brown } 407*c4762a1bSJed Brown } 408*c4762a1bSJed Brown void Limit3_Koren3(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,const PetscInt sm,const PetscInt mf,const PetscInt fm,const PetscInt ms,PetscInt n,PetscScalar *lmt) 409*c4762a1bSJed Brown { /* Eq 11 of Cada-Torrilhon 2009 */ 410*c4762a1bSJed Brown PetscInt i; 411*c4762a1bSJed Brown if (n < sm-1 || n > ms) { 412*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],(jL[i]+2*jR[i])/3,2*jR[i])/info->hxs; 413*c4762a1bSJed Brown } else if (n == sm-1) { 414*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/info->hxs,(jL[i]/info->hxs+2*jR[i]/(info->hxm/2.0+info->hxs/2.0))/3,2*jR[i]/(info->hxm/2.0+info->hxs/2.0)); 415*c4762a1bSJed Brown } else if (n == sm) { 416*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/(info->hxs/2.0+info->hxm/2.0),(jL[i]/(info->hxs/2.0+info->hxm/2.0)+2*jR[i]/info->hxm)/3,2*jR[i]/info->hxm); 417*c4762a1bSJed Brown } else if (n == ms-1) { 418*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/info->hxm,(jL[i]/info->hxm+2*jR[i]/(info->hxm/2.0+info->hxs/2.0))/3,2*jR[i]/(info->hxm/2.0+info->hxs/2.0)); 419*c4762a1bSJed Brown } else if (n == ms) { 420*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/(info->hxm/2.0+info->hxs/2.0),(jL[i]/(info->hxm/2.0+info->hxs/2.0)+2*jR[i]/info->hxs)/3,2*jR[i]/info->hxs); 421*c4762a1bSJed Brown } else if (n < mf-1 || n > fm) { 422*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],(jL[i]+2*jR[i])/3,2*jR[i])/info->hxm; 423*c4762a1bSJed Brown } else if (n == mf-1) { 424*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/info->hxm,(jL[i]/info->hxm+2*jR[i]/(info->hxm/2.0+info->hxf/2.0))/3,2*jR[i]/(info->hxm/2.0+info->hxf/2.0)); 425*c4762a1bSJed Brown } else if (n == mf) { 426*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/(info->hxm/2.0+info->hxf/2.0),(jL[i]/(info->hxm/2.0+info->hxf/2.0)+2*jR[i]/info->hxf)/3,2*jR[i]/info->hxf); 427*c4762a1bSJed Brown } else if (n == fm-1) { 428*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/info->hxf,(jL[i]/info->hxf+2*jR[i]/(info->hxf/2.0+info->hxm/2.0))/3,2*jR[i]/(info->hxf/2.0+info->hxm/2.0)); 429*c4762a1bSJed Brown } else if (n == fm) { 430*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i]/(info->hxf/2.0+info->hxm/2.0),(jL[i]/(info->hxf/2.0+info->hxm/2.0)+2*jR[i]/info->hxm)/3,2*jR[i]/info->hxm); 431*c4762a1bSJed Brown } else { 432*c4762a1bSJed Brown for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],(jL[i]+2*jR[i])/3,2*jR[i])/info->hxs; 433*c4762a1bSJed Brown } 434*c4762a1bSJed Brown } 435*c4762a1bSJed Brown PetscErrorCode RiemannListAdd(PetscFunctionList *flist,const char *name,RiemannFunction rsolve) 436*c4762a1bSJed Brown { 437*c4762a1bSJed Brown PetscErrorCode ierr; 438*c4762a1bSJed Brown 439*c4762a1bSJed Brown PetscFunctionBeginUser; 440*c4762a1bSJed Brown ierr = PetscFunctionListAdd(flist,name,rsolve);CHKERRQ(ierr); 441*c4762a1bSJed Brown PetscFunctionReturn(0); 442*c4762a1bSJed Brown } 443*c4762a1bSJed Brown 444*c4762a1bSJed Brown PetscErrorCode RiemannListFind(PetscFunctionList flist,const char *name,RiemannFunction *rsolve) 445*c4762a1bSJed Brown { 446*c4762a1bSJed Brown PetscErrorCode ierr; 447*c4762a1bSJed Brown 448*c4762a1bSJed Brown PetscFunctionBeginUser; 449*c4762a1bSJed Brown ierr = PetscFunctionListFind(flist,name,rsolve);CHKERRQ(ierr); 450*c4762a1bSJed Brown if (!*rsolve) SETERRQ1(PETSC_COMM_SELF,1,"Riemann solver \"%s\" could not be found",name); 451*c4762a1bSJed Brown PetscFunctionReturn(0); 452*c4762a1bSJed Brown } 453*c4762a1bSJed Brown 454*c4762a1bSJed Brown PetscErrorCode ReconstructListAdd(PetscFunctionList *flist,const char *name,ReconstructFunction r) 455*c4762a1bSJed Brown { 456*c4762a1bSJed Brown PetscErrorCode ierr; 457*c4762a1bSJed Brown 458*c4762a1bSJed Brown PetscFunctionBeginUser; 459*c4762a1bSJed Brown ierr = PetscFunctionListAdd(flist,name,r);CHKERRQ(ierr); 460*c4762a1bSJed Brown PetscFunctionReturn(0); 461*c4762a1bSJed Brown } 462*c4762a1bSJed Brown 463*c4762a1bSJed Brown PetscErrorCode ReconstructListFind(PetscFunctionList flist,const char *name,ReconstructFunction *r) 464*c4762a1bSJed Brown { 465*c4762a1bSJed Brown PetscErrorCode ierr; 466*c4762a1bSJed Brown 467*c4762a1bSJed Brown PetscFunctionBeginUser; 468*c4762a1bSJed Brown ierr = PetscFunctionListFind(flist,name,r);CHKERRQ(ierr); 469*c4762a1bSJed Brown if (!*r) SETERRQ1(PETSC_COMM_SELF,1,"Reconstruction \"%s\" could not be found",name); 470*c4762a1bSJed Brown PetscFunctionReturn(0); 471*c4762a1bSJed Brown } 472*c4762a1bSJed Brown 473*c4762a1bSJed Brown /* --------------------------------- Physics ----------------------------------- */ 474*c4762a1bSJed Brown 475*c4762a1bSJed Brown PetscErrorCode PhysicsDestroy_SimpleFree(void *vctx) 476*c4762a1bSJed Brown { 477*c4762a1bSJed Brown PetscErrorCode ierr; 478*c4762a1bSJed Brown 479*c4762a1bSJed Brown PetscFunctionBeginUser; 480*c4762a1bSJed Brown ierr = PetscFree(vctx);CHKERRQ(ierr); 481*c4762a1bSJed Brown PetscFunctionReturn(0); 482*c4762a1bSJed Brown } 483*c4762a1bSJed Brown 484*c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver ----------------------------------- */ 485*c4762a1bSJed Brown 486*c4762a1bSJed Brown PetscErrorCode FVRHSFunction(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 487*c4762a1bSJed Brown { 488*c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 489*c4762a1bSJed Brown PetscErrorCode ierr; 490*c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm; 491*c4762a1bSJed Brown PetscReal hx,cfl_idt = 0; 492*c4762a1bSJed Brown PetscScalar *x,*f,*slope; 493*c4762a1bSJed Brown Vec Xloc; 494*c4762a1bSJed Brown DM da; 495*c4762a1bSJed Brown 496*c4762a1bSJed Brown PetscFunctionBeginUser; 497*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 498*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); /* Xloc contains ghost points */ 499*c4762a1bSJed Brown ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); /* Mx is the number of center points */ 500*c4762a1bSJed Brown hx = (ctx->xmax-ctx->xmin)/Mx; 501*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); /* X is solution vector which does not contain ghost points */ 502*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 503*c4762a1bSJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); /* F is the right hand side function corresponds to center points */ 504*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 505*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr); 506*c4762a1bSJed Brown ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); /* contains ghost points */ 507*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 508*c4762a1bSJed Brown 509*c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 510*c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 511*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 512*c4762a1bSJed Brown } 513*c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 514*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 515*c4762a1bSJed Brown } 516*c4762a1bSJed Brown } 517*c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 518*c4762a1bSJed Brown struct _LimitInfo info; 519*c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 520*c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 521*c4762a1bSJed Brown ierr = (*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds,ctx->xmin+hx*i);CHKERRQ(ierr); 522*c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 523*c4762a1bSJed Brown ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 524*c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 525*c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 526*c4762a1bSJed Brown for (j=0; j<dof; j++) { 527*c4762a1bSJed Brown PetscScalar jmpL,jmpR; 528*c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 529*c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 530*c4762a1bSJed Brown for (k=0; k<dof; k++) { 531*c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 532*c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 533*c4762a1bSJed Brown } 534*c4762a1bSJed Brown } 535*c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 536*c4762a1bSJed Brown info.m = dof; 537*c4762a1bSJed Brown info.hx = hx; 538*c4762a1bSJed Brown (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope); 539*c4762a1bSJed Brown for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 540*c4762a1bSJed Brown for (j=0; j<dof; j++) { 541*c4762a1bSJed Brown PetscScalar tmp = 0; 542*c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 543*c4762a1bSJed Brown slope[i*dof+j] = tmp; 544*c4762a1bSJed Brown } 545*c4762a1bSJed Brown } 546*c4762a1bSJed Brown 547*c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 548*c4762a1bSJed Brown PetscReal maxspeed; 549*c4762a1bSJed Brown PetscScalar *uL,*uR; 550*c4762a1bSJed Brown uL = &ctx->uLR[0]; 551*c4762a1bSJed Brown uR = &ctx->uLR[dof]; 552*c4762a1bSJed Brown for (j=0; j<dof; j++) { 553*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 554*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 555*c4762a1bSJed Brown } 556*c4762a1bSJed Brown ierr = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr); 557*c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 558*c4762a1bSJed Brown if (i > xs) { 559*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hx; 560*c4762a1bSJed Brown } 561*c4762a1bSJed Brown if (i < xs+xm) { 562*c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hx; 563*c4762a1bSJed Brown } 564*c4762a1bSJed Brown } 565*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 566*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr); 567*c4762a1bSJed Brown ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 568*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 569*c4762a1bSJed Brown ierr = MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)da));CHKERRQ(ierr); 570*c4762a1bSJed Brown if (0) { 571*c4762a1bSJed Brown /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */ 572*c4762a1bSJed Brown PetscReal dt,tnow; 573*c4762a1bSJed Brown ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 574*c4762a1bSJed Brown ierr = TSGetTime(ts,&tnow);CHKERRQ(ierr); 575*c4762a1bSJed Brown if (dt > 0.5/ctx->cfl_idt) { 576*c4762a1bSJed Brown if (1) { 577*c4762a1bSJed Brown ierr = PetscPrintf(ctx->comm,"Stability constraint exceeded at t=%g, dt %g > %g\n",(double)tnow,(double)dt,(double)(0.5/ctx->cfl_idt));CHKERRQ(ierr); 578*c4762a1bSJed Brown } else SETERRQ2(PETSC_COMM_SELF,1,"Stability constraint exceeded, %g > %g",(double)dt,(double)(ctx->cfl/ctx->cfl_idt)); 579*c4762a1bSJed Brown } 580*c4762a1bSJed Brown } 581*c4762a1bSJed Brown PetscFunctionReturn(0); 582*c4762a1bSJed Brown } 583*c4762a1bSJed Brown 584*c4762a1bSJed Brown PetscErrorCode FVSample(FVCtx *ctx,DM da,PetscReal time,Vec U) 585*c4762a1bSJed Brown { 586*c4762a1bSJed Brown PetscErrorCode ierr; 587*c4762a1bSJed Brown PetscScalar *u,*uj; 588*c4762a1bSJed Brown PetscInt i,j,k,dof,xs,xm,Mx; 589*c4762a1bSJed Brown 590*c4762a1bSJed Brown PetscFunctionBeginUser; 591*c4762a1bSJed Brown if (!ctx->physics.sample) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Physics has not provided a sampling function"); 592*c4762a1bSJed Brown ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 593*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 594*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr); 595*c4762a1bSJed Brown ierr = PetscMalloc1(dof,&uj);CHKERRQ(ierr); 596*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 597*c4762a1bSJed Brown const PetscReal h = (ctx->xmax-ctx->xmin)/Mx,xi = ctx->xmin+h/2+i*h; 598*c4762a1bSJed Brown const PetscInt N = 200; 599*c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 600*c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 601*c4762a1bSJed Brown for (j=0; j<N+1; j++) { 602*c4762a1bSJed Brown PetscScalar xj = xi+h*(j-N/2)/(PetscReal)N; 603*c4762a1bSJed Brown ierr = (*ctx->physics.sample)(ctx->physics.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); 604*c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 605*c4762a1bSJed Brown } 606*c4762a1bSJed Brown } 607*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr); 608*c4762a1bSJed Brown ierr = PetscFree(uj);CHKERRQ(ierr); 609*c4762a1bSJed Brown PetscFunctionReturn(0); 610*c4762a1bSJed Brown } 611*c4762a1bSJed Brown 612*c4762a1bSJed Brown PetscErrorCode SolutionStatsView(DM da,Vec X,PetscViewer viewer) 613*c4762a1bSJed Brown { 614*c4762a1bSJed Brown PetscErrorCode ierr; 615*c4762a1bSJed Brown PetscReal xmin,xmax; 616*c4762a1bSJed Brown PetscScalar sum,tvsum,tvgsum; 617*c4762a1bSJed Brown const PetscScalar *x; 618*c4762a1bSJed Brown PetscInt imin,imax,Mx,i,j,xs,xm,dof; 619*c4762a1bSJed Brown Vec Xloc; 620*c4762a1bSJed Brown PetscBool iascii; 621*c4762a1bSJed Brown 622*c4762a1bSJed Brown PetscFunctionBeginUser; 623*c4762a1bSJed Brown ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 624*c4762a1bSJed Brown if (iascii) { 625*c4762a1bSJed Brown /* PETSc lacks a function to compute total variation norm (difficult in multiple dimensions), we do it here */ 626*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 627*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 628*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 629*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,Xloc,(void*)&x);CHKERRQ(ierr); 630*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 631*c4762a1bSJed Brown ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 632*c4762a1bSJed Brown tvsum = 0; 633*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 634*c4762a1bSJed Brown for (j=0; j<dof; j++) tvsum += PetscAbsScalar(x[i*dof+j]-x[(i-1)*dof+j]); 635*c4762a1bSJed Brown } 636*c4762a1bSJed Brown ierr = MPI_Allreduce(&tvsum,&tvgsum,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)da));CHKERRQ(ierr); 637*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,Xloc,(void*)&x);CHKERRQ(ierr); 638*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 639*c4762a1bSJed Brown ierr = VecMin(X,&imin,&xmin);CHKERRQ(ierr); 640*c4762a1bSJed Brown ierr = VecMax(X,&imax,&xmax);CHKERRQ(ierr); 641*c4762a1bSJed Brown ierr = VecSum(X,&sum);CHKERRQ(ierr); 642*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Solution range [%8.5f,%8.5f] with minimum at %D, mean %8.5f, ||x||_TV %8.5f\n",(double)xmin,(double)xmax,imin,(double)(sum/Mx),(double)(tvgsum/Mx));CHKERRQ(ierr); 643*c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type not supported"); 644*c4762a1bSJed Brown PetscFunctionReturn(0); 645*c4762a1bSJed Brown } 646