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