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