xref: /petsc/src/ts/tutorials/multirate/ex4.c (revision 84ff8c8b54fd7c9cb88641c01dfe6357ec5f72d0)
1*84ff8c8bSHong Zhang /*
2*84ff8c8bSHong Zhang   Note:
3*84ff8c8bSHong Zhang     -hratio is the ratio between mesh size of corse grids and fine grids
4*84ff8c8bSHong Zhang */
5*84ff8c8bSHong Zhang 
6*84ff8c8bSHong Zhang static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n"
7*84ff8c8bSHong Zhang   "  advect      - Constant coefficient scalar advection\n"
8*84ff8c8bSHong Zhang   "                u_t       + (a*u)_x               = 0\n"
9*84ff8c8bSHong Zhang   "  shallow     - 1D Shallow water equations (Saint Venant System)\n"
10*84ff8c8bSHong Zhang   "                h_t + (q)_x = 0 \n"
11*84ff8c8bSHong Zhang   "                q_t + (\frac{q^2}{h} + g/2*h^2)_x = 0  \n"
12*84ff8c8bSHong Zhang   "                where, h(x,t) denotes the height of the water and q(x,t) the momentum.\n"
13*84ff8c8bSHong Zhang   "  for this toy problem, we choose different meshsizes for different sub-domains (slow-fast-slow), say\n"
14*84ff8c8bSHong Zhang   "                hxs  = hratio*hxf \n"
15*84ff8c8bSHong Zhang   "  where hxs and hxf are the grid spacings for coarse and fine grids respectively.\n"
16*84ff8c8bSHong Zhang   "  exact       - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n"
17*84ff8c8bSHong Zhang   "                the states across shocks and rarefactions\n"
18*84ff8c8bSHong Zhang   "  simulation  - use reference solution which is generated by smaller time step size to be true solution,\n"
19*84ff8c8bSHong Zhang   "                also the reference solution should be generated by user and stored in a binary file.\n"
20*84ff8c8bSHong Zhang   "  characteristic - Limit the characteristic variables, this is usually preferred (default)\n"
21*84ff8c8bSHong Zhang   "  bc_type     - Boundary condition for the problem, options are: periodic, outflow, inflow "
22*84ff8c8bSHong Zhang   "Several problem descriptions (initial data, physics specific features, boundary data) can be chosen with -initial N\n\n"
23*84ff8c8bSHong Zhang   "The problem size should be set with -da_grid_x M\n\n";
24*84ff8c8bSHong Zhang 
25*84ff8c8bSHong Zhang /*
26*84ff8c8bSHong Zhang   Example:
27*84ff8c8bSHong Zhang      ./ex4 -da_grid_x 40 -initial 1 -hratio 1 -limit mc -ts_dt 0.01 -ts_max_time 7.0 -ts_type mprk -ts_mprk_type 2a22 -ts_monitor_draw_solution -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0
28*84ff8c8bSHong Zhang      ./ex4 -da_grid_x 40 -initial 2 -hratio 1 -limit mc -ts_dt 0.01 -ts_max_time 2.5 -ts_type mprk -ts_mprk_type 2a22 -ts_monitor_draw_solution -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0
29*84ff8c8bSHong Zhang      ./ex4 -da_grid_x 40 -initial 3 -hratio 1 -limit mc -ts_dt 0.01 -ts_max_time 4.0 -ts_type mprk -ts_mprk_type 2a22 -ts_monitor_draw_solution -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0
30*84ff8c8bSHong Zhang      ./ex4 -da_grid_x 40 -initial 4 -hratio 1 -limit koren3 -ts_dt 0.01 -ts_max_time 4.0 -ts_type mprk -ts_mprk_type 2a22 -ts_monitor_draw_solution -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0
31*84ff8c8bSHong Zhang      ./ex4 -da_grid_x 40 -initial 5 -hratio 1 -limit mc -ts_dt 0.01 -ts_max_time 5.0 -ts_type mprk -ts_mprk_type 2a22 -ts_monitor_draw_solution -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0
32*84ff8c8bSHong Zhang 
33*84ff8c8bSHong Zhang   Contributed by: Aidan Hamilton <aidan@udel.edu>
34*84ff8c8bSHong Zhang */
35*84ff8c8bSHong Zhang 
36*84ff8c8bSHong Zhang #include <petscts.h>
37*84ff8c8bSHong Zhang #include <petscdm.h>
38*84ff8c8bSHong Zhang #include <petscdmda.h>
39*84ff8c8bSHong Zhang #include <petscdraw.h>
40*84ff8c8bSHong Zhang #include "finitevolume1d.h"
41*84ff8c8bSHong Zhang #include <petsc/private/kernels/blockinvert.h>
42*84ff8c8bSHong Zhang 
43*84ff8c8bSHong Zhang PETSC_STATIC_INLINE PetscReal RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax) { PetscReal range = xmax-xmin; return xmin +PetscFmodReal(range+PetscFmodReal(a,range),range); }
44*84ff8c8bSHong Zhang PETSC_STATIC_INLINE PetscReal MaxAbs(PetscReal a,PetscReal b) { return (PetscAbs(a) > PetscAbs(b)) ? a : b; }
45*84ff8c8bSHong Zhang /* --------------------------------- Advection ----------------------------------- */
46*84ff8c8bSHong Zhang typedef struct {
47*84ff8c8bSHong Zhang   PetscReal a;                  /* advective velocity */
48*84ff8c8bSHong Zhang } AdvectCtx;
49*84ff8c8bSHong Zhang 
50*84ff8c8bSHong Zhang static PetscErrorCode PhysicsRiemann_Advect(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
51*84ff8c8bSHong Zhang {
52*84ff8c8bSHong Zhang   AdvectCtx *ctx = (AdvectCtx*)vctx;
53*84ff8c8bSHong Zhang   PetscReal speed;
54*84ff8c8bSHong Zhang 
55*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
56*84ff8c8bSHong Zhang   speed     = ctx->a;
57*84ff8c8bSHong Zhang   flux[0]   = PetscMax(0,speed)*uL[0] + PetscMin(0,speed)*uR[0];
58*84ff8c8bSHong Zhang   *maxspeed = speed;
59*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
60*84ff8c8bSHong Zhang }
61*84ff8c8bSHong Zhang 
62*84ff8c8bSHong Zhang static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
63*84ff8c8bSHong Zhang {
64*84ff8c8bSHong Zhang   AdvectCtx *ctx = (AdvectCtx*)vctx;
65*84ff8c8bSHong Zhang 
66*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
67*84ff8c8bSHong Zhang   X[0]      = 1.;
68*84ff8c8bSHong Zhang   Xi[0]     = 1.;
69*84ff8c8bSHong Zhang   speeds[0] = ctx->a;
70*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
71*84ff8c8bSHong Zhang }
72*84ff8c8bSHong Zhang 
73*84ff8c8bSHong Zhang static PetscErrorCode PhysicsSample_Advect(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
74*84ff8c8bSHong Zhang {
75*84ff8c8bSHong Zhang   AdvectCtx *ctx = (AdvectCtx*)vctx;
76*84ff8c8bSHong Zhang   PetscReal a    = ctx->a,x0;
77*84ff8c8bSHong Zhang 
78*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
79*84ff8c8bSHong Zhang   switch (bctype) {
80*84ff8c8bSHong Zhang     case FVBC_OUTFLOW:   x0 = x-a*t; break;
81*84ff8c8bSHong Zhang     case FVBC_PERIODIC: x0 = RangeMod(x-a*t,xmin,xmax); break;
82*84ff8c8bSHong Zhang     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown BCType");
83*84ff8c8bSHong Zhang   }
84*84ff8c8bSHong Zhang   switch (initial) {
85*84ff8c8bSHong Zhang     case 0: u[0] = (x0 < 0) ? 1 : -1; break;
86*84ff8c8bSHong Zhang     case 1: u[0] = (x0 < 0) ? -1 : 1; break;
87*84ff8c8bSHong Zhang     case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break;
88*84ff8c8bSHong Zhang     case 3: u[0] = PetscSinReal(2*PETSC_PI*x0); break;
89*84ff8c8bSHong Zhang     case 4: u[0] = PetscAbs(x0); break;
90*84ff8c8bSHong Zhang     case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2*PETSC_PI*x0)); break;
91*84ff8c8bSHong Zhang     case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2-x0 : 0)); break;
92*84ff8c8bSHong Zhang     case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI*x0),10.0);break;
93*84ff8c8bSHong Zhang     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition");
94*84ff8c8bSHong Zhang   }
95*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
96*84ff8c8bSHong Zhang }
97*84ff8c8bSHong Zhang 
98*84ff8c8bSHong Zhang static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx)
99*84ff8c8bSHong Zhang {
100*84ff8c8bSHong Zhang   PetscErrorCode ierr;
101*84ff8c8bSHong Zhang   AdvectCtx      *user;
102*84ff8c8bSHong Zhang 
103*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
104*84ff8c8bSHong Zhang   ierr = PetscNew(&user);CHKERRQ(ierr);
105*84ff8c8bSHong Zhang   ctx->physics2.sample2         = PhysicsSample_Advect;
106*84ff8c8bSHong Zhang   ctx->physics2.riemann2        = PhysicsRiemann_Advect;
107*84ff8c8bSHong Zhang   ctx->physics2.characteristic2 = PhysicsCharacteristic_Advect;
108*84ff8c8bSHong Zhang   ctx->physics2.destroy         = PhysicsDestroy_SimpleFree;
109*84ff8c8bSHong Zhang   ctx->physics2.user            = user;
110*84ff8c8bSHong Zhang   ctx->physics2.dof             = 1;
111*84ff8c8bSHong Zhang 
112*84ff8c8bSHong Zhang   ierr = PetscStrallocpy("u",&ctx->physics2.fieldname[0]);CHKERRQ(ierr);
113*84ff8c8bSHong Zhang   user->a = 1;
114*84ff8c8bSHong Zhang   ierr = PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");CHKERRQ(ierr);
115*84ff8c8bSHong Zhang     ierr = PetscOptionsReal("-physics_advect_a","Speed","",user->a,&user->a,NULL);CHKERRQ(ierr);
116*84ff8c8bSHong Zhang   ierr = PetscOptionsEnd();CHKERRQ(ierr);
117*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
118*84ff8c8bSHong Zhang }
119*84ff8c8bSHong Zhang /* --------------------------------- Shallow Water ----------------------------------- */
120*84ff8c8bSHong Zhang 
121*84ff8c8bSHong Zhang typedef struct {
122*84ff8c8bSHong Zhang   PetscReal gravity;
123*84ff8c8bSHong Zhang } ShallowCtx;
124*84ff8c8bSHong Zhang 
125*84ff8c8bSHong Zhang PETSC_STATIC_INLINE void ShallowFlux(ShallowCtx *phys,const PetscScalar *u,PetscScalar *f)
126*84ff8c8bSHong Zhang {
127*84ff8c8bSHong Zhang   f[0] = u[1];
128*84ff8c8bSHong Zhang   f[1] = PetscSqr(u[1])/u[0] + 0.5*phys->gravity*PetscSqr(u[0]);
129*84ff8c8bSHong Zhang }
130*84ff8c8bSHong Zhang 
131*84ff8c8bSHong Zhang static PetscErrorCode PhysicsRiemann_Shallow_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
132*84ff8c8bSHong Zhang {
133*84ff8c8bSHong Zhang   ShallowCtx                *phys = (ShallowCtx*)vctx;
134*84ff8c8bSHong Zhang   PetscScalar               g    = phys->gravity,ustar[2],cL,cR,c,cstar;
135*84ff8c8bSHong Zhang   struct {PetscScalar h,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]},star;
136*84ff8c8bSHong Zhang   PetscInt                  i;
137*84ff8c8bSHong Zhang 
138*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
139*84ff8c8bSHong Zhang   if (!(L.h > 0 && R.h > 0)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative");
140*84ff8c8bSHong Zhang   cL = PetscSqrtScalar(g*L.h);
141*84ff8c8bSHong Zhang   cR = PetscSqrtScalar(g*R.h);
142*84ff8c8bSHong Zhang   c  = PetscMax(cL,cR);
143*84ff8c8bSHong Zhang   {
144*84ff8c8bSHong Zhang     /* Solve for star state */
145*84ff8c8bSHong Zhang     const PetscInt maxits = 50;
146*84ff8c8bSHong Zhang     PetscScalar tmp,res,res0=0,h0,h = 0.5*(L.h + R.h); /* initial guess */
147*84ff8c8bSHong Zhang     h0 = h;
148*84ff8c8bSHong Zhang     for (i=0; i<maxits; i++) {
149*84ff8c8bSHong Zhang       PetscScalar fr,fl,dfr,dfl;
150*84ff8c8bSHong Zhang       fl = (L.h < h)
151*84ff8c8bSHong Zhang         ? PetscSqrtScalar(0.5*g*(h*h - L.h*L.h)*(1/L.h - 1/h)) /* shock */
152*84ff8c8bSHong Zhang         : 2*PetscSqrtScalar(g*h) - 2*PetscSqrtScalar(g*L.h);   /* rarefaction */
153*84ff8c8bSHong Zhang       fr = (R.h < h)
154*84ff8c8bSHong Zhang         ? PetscSqrtScalar(0.5*g*(h*h - R.h*R.h)*(1/R.h - 1/h)) /* shock */
155*84ff8c8bSHong Zhang         : 2*PetscSqrtScalar(g*h) - 2*PetscSqrtScalar(g*R.h);   /* rarefaction */
156*84ff8c8bSHong Zhang       res = R.u - L.u + fr + fl;
157*84ff8c8bSHong Zhang       if (PetscIsInfOrNanScalar(res)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinity or Not-a-Number generated in computation");
158*84ff8c8bSHong Zhang       if (PetscAbsScalar(res) < PETSC_SQRT_MACHINE_EPSILON || (i > 0 && PetscAbsScalar(h-h0) < PETSC_SQRT_MACHINE_EPSILON)) {
159*84ff8c8bSHong Zhang         star.h = h;
160*84ff8c8bSHong Zhang         star.u = L.u - fl;
161*84ff8c8bSHong Zhang         goto converged;
162*84ff8c8bSHong Zhang       } else if (i > 0 && PetscAbsScalar(res) >= PetscAbsScalar(res0)) {        /* Line search */
163*84ff8c8bSHong Zhang         h = 0.8*h0 + 0.2*h;
164*84ff8c8bSHong Zhang         continue;
165*84ff8c8bSHong Zhang       }
166*84ff8c8bSHong Zhang       /* Accept the last step and take another */
167*84ff8c8bSHong Zhang       res0 = res;
168*84ff8c8bSHong Zhang       h0 = h;
169*84ff8c8bSHong Zhang       dfl = (L.h < h) ? 0.5/fl*0.5*g*(-L.h*L.h/(h*h) - 1 + 2*h/L.h) : PetscSqrtScalar(g/h);
170*84ff8c8bSHong Zhang       dfr = (R.h < h) ? 0.5/fr*0.5*g*(-R.h*R.h/(h*h) - 1 + 2*h/R.h) : PetscSqrtScalar(g/h);
171*84ff8c8bSHong Zhang       tmp = h - res/(dfr+dfl);
172*84ff8c8bSHong Zhang       if (tmp <= 0) h /= 2;   /* Guard against Newton shooting off to a negative thickness */
173*84ff8c8bSHong Zhang       else h = tmp;
174*84ff8c8bSHong Zhang       if (!((h > 0) && PetscIsNormalScalar(h))) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FP,"non-normal iterate h=%g",(double)h);
175*84ff8c8bSHong Zhang     }
176*84ff8c8bSHong Zhang     SETERRQ1(PETSC_COMM_SELF,1,"Newton iteration for star.h diverged after %D iterations",i);
177*84ff8c8bSHong Zhang   }
178*84ff8c8bSHong Zhang converged:
179*84ff8c8bSHong Zhang   cstar = PetscSqrtScalar(g*star.h);
180*84ff8c8bSHong Zhang   if (L.u-cL < 0 && 0 < star.u-cstar) { /* 1-wave is sonic rarefaction */
181*84ff8c8bSHong Zhang     PetscScalar ufan[2];
182*84ff8c8bSHong Zhang     ufan[0] = 1/g*PetscSqr(L.u/3 + 2./3*cL);
183*84ff8c8bSHong Zhang     ufan[1] = PetscSqrtScalar(g*ufan[0])*ufan[0];
184*84ff8c8bSHong Zhang     ShallowFlux(phys,ufan,flux);
185*84ff8c8bSHong Zhang   } else if (star.u+cstar < 0 && 0 < R.u+cR) { /* 2-wave is sonic rarefaction */
186*84ff8c8bSHong Zhang     PetscScalar ufan[2];
187*84ff8c8bSHong Zhang     ufan[0] = 1/g*PetscSqr(R.u/3 - 2./3*cR);
188*84ff8c8bSHong Zhang     ufan[1] = -PetscSqrtScalar(g*ufan[0])*ufan[0];
189*84ff8c8bSHong Zhang     ShallowFlux(phys,ufan,flux);
190*84ff8c8bSHong Zhang   } else if ((L.h >= star.h && L.u-c >= 0) || (L.h<star.h && (star.h*star.u-L.h*L.u)/(star.h-L.h) > 0)) {
191*84ff8c8bSHong Zhang     /* 1-wave is right-travelling shock (supersonic) */
192*84ff8c8bSHong Zhang     ShallowFlux(phys,uL,flux);
193*84ff8c8bSHong Zhang   } else if ((star.h <= R.h && R.u+c <= 0) || (star.h>R.h && (R.h*R.u-star.h*star.h)/(R.h-star.h) < 0)) {
194*84ff8c8bSHong Zhang     /* 2-wave is left-travelling shock (supersonic) */
195*84ff8c8bSHong Zhang     ShallowFlux(phys,uR,flux);
196*84ff8c8bSHong Zhang   } else {
197*84ff8c8bSHong Zhang     ustar[0] = star.h;
198*84ff8c8bSHong Zhang     ustar[1] = star.h*star.u;
199*84ff8c8bSHong Zhang     ShallowFlux(phys,ustar,flux);
200*84ff8c8bSHong Zhang   }
201*84ff8c8bSHong Zhang   *maxspeed = MaxAbs(MaxAbs(star.u-cstar,star.u+cstar),MaxAbs(L.u-cL,R.u+cR));
202*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
203*84ff8c8bSHong Zhang }
204*84ff8c8bSHong Zhang 
205*84ff8c8bSHong Zhang static PetscErrorCode PhysicsRiemann_Shallow_Rusanov(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
206*84ff8c8bSHong Zhang {
207*84ff8c8bSHong Zhang   ShallowCtx                *phys = (ShallowCtx*)vctx;
208*84ff8c8bSHong Zhang   PetscScalar               g = phys->gravity,fL[2],fR[2],s;
209*84ff8c8bSHong Zhang   struct {PetscScalar h,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]};
210*84ff8c8bSHong Zhang   PetscReal                 tol = 1e-6;
211*84ff8c8bSHong Zhang 
212*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
213*84ff8c8bSHong Zhang   /* Positivity preserving modification*/
214*84ff8c8bSHong Zhang   if (L.h < tol) L.u = 0.0;
215*84ff8c8bSHong Zhang   if (R.h < tol) R.u = 0.0;
216*84ff8c8bSHong Zhang 
217*84ff8c8bSHong Zhang   /*simple pos preserve limiter*/
218*84ff8c8bSHong Zhang   if (L.h < 0) L.h = 0;
219*84ff8c8bSHong Zhang   if (R.h < 0) R.h = 0;
220*84ff8c8bSHong Zhang 
221*84ff8c8bSHong Zhang   ShallowFlux(phys,uL,fL);
222*84ff8c8bSHong Zhang   ShallowFlux(phys,uR,fR);
223*84ff8c8bSHong Zhang 
224*84ff8c8bSHong Zhang   s         = PetscMax(PetscAbs(L.u)+PetscSqrtScalar(g*L.h),PetscAbs(R.u)+PetscSqrtScalar(g*R.h));
225*84ff8c8bSHong Zhang   flux[0]   = 0.5*(fL[0] + fR[0]) + 0.5*s*(L.h - R.h);
226*84ff8c8bSHong Zhang   flux[1]   = 0.5*(fL[1] + fR[1]) + 0.5*s*(uL[1] - uR[1]);
227*84ff8c8bSHong Zhang   *maxspeed = s;
228*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
229*84ff8c8bSHong Zhang }
230*84ff8c8bSHong Zhang 
231*84ff8c8bSHong Zhang static PetscErrorCode PhysicsCharacteristic_Conservative(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
232*84ff8c8bSHong Zhang {
233*84ff8c8bSHong Zhang   PetscInt i,j;
234*84ff8c8bSHong Zhang 
235*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
236*84ff8c8bSHong Zhang   for (i=0; i<m; i++) {
237*84ff8c8bSHong Zhang     for (j=0; j<m; j++) Xi[i*m+j] = X[i*m+j] = (PetscScalar)(i==j);
238*84ff8c8bSHong Zhang     speeds[i] = PETSC_MAX_REAL; /* Indicates invalid */
239*84ff8c8bSHong Zhang   }
240*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
241*84ff8c8bSHong Zhang }
242*84ff8c8bSHong Zhang 
243*84ff8c8bSHong Zhang static PetscErrorCode PhysicsCharacteristic_Shallow(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
244*84ff8c8bSHong Zhang {
245*84ff8c8bSHong Zhang   ShallowCtx     *phys = (ShallowCtx*)vctx;
246*84ff8c8bSHong Zhang   PetscReal      c;
247*84ff8c8bSHong Zhang   PetscErrorCode ierr;
248*84ff8c8bSHong Zhang   PetscReal      tol = 1e-6;
249*84ff8c8bSHong Zhang 
250*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
251*84ff8c8bSHong Zhang   c           = PetscSqrtScalar(u[0]*phys->gravity);
252*84ff8c8bSHong Zhang 
253*84ff8c8bSHong Zhang   if (u[0] < tol) { /*Use conservative variables*/
254*84ff8c8bSHong Zhang     X[0*2+0]  = 1;
255*84ff8c8bSHong Zhang     X[0*2+1]  = 0;
256*84ff8c8bSHong Zhang     X[1*2+0]  = 0;
257*84ff8c8bSHong Zhang     X[1*2+1]  = 1;
258*84ff8c8bSHong Zhang   } else {
259*84ff8c8bSHong Zhang     speeds[0] = u[1]/u[0] - c;
260*84ff8c8bSHong Zhang     speeds[1] = u[1]/u[0] + c;
261*84ff8c8bSHong Zhang     X[0*2+0]  = 1;
262*84ff8c8bSHong Zhang     X[0*2+1]  = speeds[0];
263*84ff8c8bSHong Zhang     X[1*2+0]  = 1;
264*84ff8c8bSHong Zhang     X[1*2+1]  = speeds[1];
265*84ff8c8bSHong Zhang   }
266*84ff8c8bSHong Zhang 
267*84ff8c8bSHong Zhang   ierr = PetscArraycpy(Xi,X,4);CHKERRQ(ierr);
268*84ff8c8bSHong Zhang   ierr = PetscKernel_A_gets_inverse_A_2(Xi,0,PETSC_FALSE,NULL);CHKERRQ(ierr);
269*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
270*84ff8c8bSHong Zhang }
271*84ff8c8bSHong Zhang 
272*84ff8c8bSHong Zhang static PetscErrorCode PhysicsSample_Shallow(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
273*84ff8c8bSHong Zhang {
274*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
275*84ff8c8bSHong Zhang   if (t > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Exact solutions not implemented for t > 0");
276*84ff8c8bSHong Zhang   switch (initial) {
277*84ff8c8bSHong Zhang     case 0:
278*84ff8c8bSHong Zhang       u[0] = (x < 0) ? 2 : 1; /* Standard Dam Break Problem */
279*84ff8c8bSHong Zhang       u[1] = (x < 0) ? 0 : 0;
280*84ff8c8bSHong Zhang       break;
281*84ff8c8bSHong Zhang     case 1:
282*84ff8c8bSHong Zhang       u[0] = (x < 10) ?   1 : 0.1; /*The Next 5 problems are standard Riemann problem tests */
283*84ff8c8bSHong Zhang       u[1] = (x < 10) ? 2.5 : 0;
284*84ff8c8bSHong Zhang       break;
285*84ff8c8bSHong Zhang     case 2:
286*84ff8c8bSHong Zhang       u[0] = (x < 25) ?  1 : 1;
287*84ff8c8bSHong Zhang       u[1] = (x < 25) ? -5 : 5;
288*84ff8c8bSHong Zhang       break;
289*84ff8c8bSHong Zhang     case 3:
290*84ff8c8bSHong Zhang       u[0] = (x < 20) ?  1 : 1e-12;
291*84ff8c8bSHong Zhang       u[1] = (x < 20) ?  0 : 0;
292*84ff8c8bSHong Zhang       break;
293*84ff8c8bSHong Zhang     case 4:
294*84ff8c8bSHong Zhang       u[0] = (x < 30) ? 1e-12 : 1;
295*84ff8c8bSHong Zhang       u[1] = (x < 30) ? 0 : 0;
296*84ff8c8bSHong Zhang       break;
297*84ff8c8bSHong Zhang     case 5:
298*84ff8c8bSHong Zhang       u[0] = (x < 25) ?  0.1 : 0.1;
299*84ff8c8bSHong Zhang       u[1] = (x < 25) ? -0.3 : 0.3;
300*84ff8c8bSHong Zhang       break;
301*84ff8c8bSHong Zhang     case 6:
302*84ff8c8bSHong Zhang       u[0] = 1+0.5*PetscSinReal(2*PETSC_PI*x);
303*84ff8c8bSHong Zhang       u[1] = 1*u[0];
304*84ff8c8bSHong Zhang       break;
305*84ff8c8bSHong Zhang     case 7:
306*84ff8c8bSHong Zhang       if (x < -0.1) {
307*84ff8c8bSHong Zhang        u[0] = 1e-9;
308*84ff8c8bSHong Zhang        u[1] = 0.0;
309*84ff8c8bSHong Zhang       } else if (x < 0.1) {
310*84ff8c8bSHong Zhang        u[0] = 1.0;
311*84ff8c8bSHong Zhang        u[1] = 0.0;
312*84ff8c8bSHong Zhang       } else {
313*84ff8c8bSHong Zhang        u[0] = 1e-9;
314*84ff8c8bSHong Zhang        u[1] = 0.0;
315*84ff8c8bSHong Zhang       }
316*84ff8c8bSHong Zhang       break;
317*84ff8c8bSHong Zhang     case 8:
318*84ff8c8bSHong Zhang      if (x < -0.1) {
319*84ff8c8bSHong Zhang        u[0] = 2;
320*84ff8c8bSHong Zhang        u[1] = 0.0;
321*84ff8c8bSHong Zhang       } else if (x < 0.1) {
322*84ff8c8bSHong Zhang        u[0] = 3.0;
323*84ff8c8bSHong Zhang        u[1] = 0.0;
324*84ff8c8bSHong Zhang       } else {
325*84ff8c8bSHong Zhang        u[0] = 2;
326*84ff8c8bSHong Zhang        u[1] = 0.0;
327*84ff8c8bSHong Zhang       }
328*84ff8c8bSHong Zhang       break;
329*84ff8c8bSHong Zhang     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition");
330*84ff8c8bSHong Zhang   }
331*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
332*84ff8c8bSHong Zhang }
333*84ff8c8bSHong Zhang 
334*84ff8c8bSHong Zhang /* Implements inflow conditions for the given initial conditions. Which conditions are actually enforced depends on
335*84ff8c8bSHong Zhang    on the results of PhysicsSetInflowType_Shallow. */
336*84ff8c8bSHong Zhang static PetscErrorCode PhysicsInflow_Shallow(void *vctx,PetscReal t,PetscReal x,PetscReal *u)
337*84ff8c8bSHong Zhang {
338*84ff8c8bSHong Zhang   FVCtx          *ctx = (FVCtx*)vctx;
339*84ff8c8bSHong Zhang 
340*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
341*84ff8c8bSHong Zhang   if (ctx->bctype == FVBC_INFLOW) {
342*84ff8c8bSHong Zhang     switch (ctx->initial) {
343*84ff8c8bSHong Zhang       case 0:
344*84ff8c8bSHong Zhang       case 1:
345*84ff8c8bSHong Zhang       case 2:
346*84ff8c8bSHong Zhang       case 3:
347*84ff8c8bSHong Zhang       case 4:
348*84ff8c8bSHong Zhang       case 5:
349*84ff8c8bSHong Zhang         u[0] = 0; u[1] = 0.0; /* Left boundary conditions */
350*84ff8c8bSHong Zhang         u[2] = 0; u[3] = 0.0; /* Right boundary conditions */
351*84ff8c8bSHong Zhang         break;
352*84ff8c8bSHong Zhang       case 6:
353*84ff8c8bSHong Zhang         u[0] = 0; u[1] = 0.0; /* Left boundary conditions */
354*84ff8c8bSHong Zhang         u[2] = 0; u[3] = 0.0; /* Right boundary conditions */
355*84ff8c8bSHong Zhang         break;
356*84ff8c8bSHong Zhang       case 7:
357*84ff8c8bSHong Zhang         u[0] = 0; u[1] = 0.0; /* Left boundary conditions */
358*84ff8c8bSHong Zhang         u[2] = 0; u[3] = 0.0; /* Right boundary conditions */
359*84ff8c8bSHong Zhang         break;
360*84ff8c8bSHong Zhang       case 8:
361*84ff8c8bSHong Zhang         u[0] = 0; u[1] = 1.0; /* Left boundary conditions */
362*84ff8c8bSHong Zhang         u[2] = 0; u[3] = -1.0;/* Right boundary conditions */
363*84ff8c8bSHong Zhang         break;
364*84ff8c8bSHong Zhang       default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition");
365*84ff8c8bSHong Zhang     }
366*84ff8c8bSHong Zhang   }
367*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
368*84ff8c8bSHong Zhang }
369*84ff8c8bSHong Zhang 
370*84ff8c8bSHong Zhang /* Selects which boundary conditions are marked as inflow and which as outflow when FVBC_INFLOW is selected. */
371*84ff8c8bSHong Zhang static PetscErrorCode PhysicsSetInflowType_Shallow(FVCtx *ctx)
372*84ff8c8bSHong Zhang {
373*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
374*84ff8c8bSHong Zhang   switch (ctx->initial) {
375*84ff8c8bSHong Zhang     case 0:
376*84ff8c8bSHong Zhang     case 1:
377*84ff8c8bSHong Zhang     case 2:
378*84ff8c8bSHong Zhang     case 3:
379*84ff8c8bSHong Zhang     case 4:
380*84ff8c8bSHong Zhang     case 5:
381*84ff8c8bSHong Zhang     case 6:
382*84ff8c8bSHong Zhang     case 7:
383*84ff8c8bSHong Zhang     case 8: /* Fix left and right momentum, height is outflow */
384*84ff8c8bSHong Zhang       ctx->physics2.bcinflowindex[0] = PETSC_FALSE;
385*84ff8c8bSHong Zhang       ctx->physics2.bcinflowindex[1] = PETSC_TRUE;
386*84ff8c8bSHong Zhang       ctx->physics2.bcinflowindex[2] = PETSC_FALSE;
387*84ff8c8bSHong Zhang       ctx->physics2.bcinflowindex[3] = PETSC_TRUE;
388*84ff8c8bSHong Zhang       break;
389*84ff8c8bSHong Zhang     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition");
390*84ff8c8bSHong Zhang   }
391*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
392*84ff8c8bSHong Zhang }
393*84ff8c8bSHong Zhang 
394*84ff8c8bSHong Zhang static PetscErrorCode PhysicsCreate_Shallow(FVCtx *ctx)
395*84ff8c8bSHong Zhang {
396*84ff8c8bSHong Zhang   PetscErrorCode    ierr;
397*84ff8c8bSHong Zhang   ShallowCtx        *user;
398*84ff8c8bSHong Zhang   PetscFunctionList rlist = 0,rclist = 0;
399*84ff8c8bSHong Zhang   char              rname[256] = "rusanov",rcname[256] = "characteristic";
400*84ff8c8bSHong Zhang 
401*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
402*84ff8c8bSHong Zhang   ierr = PetscNew(&user);CHKERRQ(ierr);
403*84ff8c8bSHong Zhang   ctx->physics2.sample2         = PhysicsSample_Shallow;
404*84ff8c8bSHong Zhang   ctx->physics2.inflow          = PhysicsInflow_Shallow;
405*84ff8c8bSHong Zhang   ctx->physics2.destroy         = PhysicsDestroy_SimpleFree;
406*84ff8c8bSHong Zhang   ctx->physics2.riemann2        = PhysicsRiemann_Shallow_Rusanov;
407*84ff8c8bSHong Zhang   ctx->physics2.characteristic2 = PhysicsCharacteristic_Shallow;
408*84ff8c8bSHong Zhang   ctx->physics2.user            = user;
409*84ff8c8bSHong Zhang   ctx->physics2.dof             = 2;
410*84ff8c8bSHong Zhang 
411*84ff8c8bSHong Zhang   PetscMalloc1(2*(ctx->physics2.dof),&ctx->physics2.bcinflowindex);
412*84ff8c8bSHong Zhang   PhysicsSetInflowType_Shallow(ctx);
413*84ff8c8bSHong Zhang 
414*84ff8c8bSHong Zhang   ierr = PetscStrallocpy("height",&ctx->physics2.fieldname[0]);CHKERRQ(ierr);
415*84ff8c8bSHong Zhang   ierr = PetscStrallocpy("momentum",&ctx->physics2.fieldname[1]);CHKERRQ(ierr);
416*84ff8c8bSHong Zhang 
417*84ff8c8bSHong Zhang   user->gravity = 9.81;
418*84ff8c8bSHong Zhang 
419*84ff8c8bSHong Zhang   ierr = RiemannListAdd_2WaySplit(&rlist,"exact",  PhysicsRiemann_Shallow_Exact);CHKERRQ(ierr);
420*84ff8c8bSHong Zhang   ierr = RiemannListAdd_2WaySplit(&rlist,"rusanov",PhysicsRiemann_Shallow_Rusanov);CHKERRQ(ierr);
421*84ff8c8bSHong Zhang   ierr = ReconstructListAdd_2WaySplit(&rclist,"characteristic",PhysicsCharacteristic_Shallow);CHKERRQ(ierr);
422*84ff8c8bSHong Zhang   ierr = ReconstructListAdd_2WaySplit(&rclist,"conservative",PhysicsCharacteristic_Conservative);CHKERRQ(ierr);
423*84ff8c8bSHong Zhang   ierr = PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for Shallow","");CHKERRQ(ierr);
424*84ff8c8bSHong Zhang     ierr = PetscOptionsReal("-physics_shallow_gravity","Gravity","",user->gravity,&user->gravity,NULL);CHKERRQ(ierr);
425*84ff8c8bSHong Zhang     ierr = PetscOptionsFList("-physics_shallow_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL);CHKERRQ(ierr);
426*84ff8c8bSHong Zhang     ierr = PetscOptionsFList("-physics_shallow_reconstruct","Reconstruction","",rclist,rcname,rcname,sizeof(rcname),NULL);CHKERRQ(ierr);
427*84ff8c8bSHong Zhang   ierr = PetscOptionsEnd();CHKERRQ(ierr);
428*84ff8c8bSHong Zhang   ierr = RiemannListFind_2WaySplit(rlist,rname,&ctx->physics2.riemann2);CHKERRQ(ierr);
429*84ff8c8bSHong Zhang   ierr = ReconstructListFind_2WaySplit(rclist,rcname,&ctx->physics2.characteristic2);CHKERRQ(ierr);
430*84ff8c8bSHong Zhang   ierr = PetscFunctionListDestroy(&rlist);CHKERRQ(ierr);
431*84ff8c8bSHong Zhang   ierr = PetscFunctionListDestroy(&rclist);CHKERRQ(ierr);
432*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
433*84ff8c8bSHong Zhang }
434*84ff8c8bSHong Zhang 
435*84ff8c8bSHong Zhang PetscErrorCode FVSample_2WaySplit(FVCtx *ctx,DM da,PetscReal time,Vec U)
436*84ff8c8bSHong Zhang {
437*84ff8c8bSHong Zhang   PetscErrorCode  ierr;
438*84ff8c8bSHong Zhang   PetscScalar     *u,*uj,xj,xi;
439*84ff8c8bSHong Zhang   PetscInt        i,j,k,dof,xs,xm,Mx;
440*84ff8c8bSHong Zhang   const PetscInt  N = 200;
441*84ff8c8bSHong Zhang   PetscReal       hs,hf;
442*84ff8c8bSHong Zhang 
443*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
444*84ff8c8bSHong Zhang   if (!ctx->physics2.sample2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Physics has not provided a sampling function");
445*84ff8c8bSHong Zhang   ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr);
446*84ff8c8bSHong Zhang   ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr);
447*84ff8c8bSHong Zhang   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
448*84ff8c8bSHong Zhang   ierr = PetscMalloc1(dof,&uj);CHKERRQ(ierr);
449*84ff8c8bSHong Zhang   hs   = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
450*84ff8c8bSHong Zhang   hf   = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
451*84ff8c8bSHong Zhang   for (i=xs; i<xs+xm; i++) {
452*84ff8c8bSHong Zhang     if (i < ctx->sf) {
453*84ff8c8bSHong Zhang       xi = ctx->xmin+0.5*hs+i*hs;
454*84ff8c8bSHong Zhang       /* Integrate over cell i using trapezoid rule with N points. */
455*84ff8c8bSHong Zhang       for (k=0; k<dof; k++) u[i*dof+k] = 0;
456*84ff8c8bSHong Zhang       for (j=0; j<N+1; j++) {
457*84ff8c8bSHong Zhang         xj = xi+hs*(j-N/2)/(PetscReal)N;
458*84ff8c8bSHong Zhang         ierr = (*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr);
459*84ff8c8bSHong Zhang         for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N;
460*84ff8c8bSHong Zhang       }
461*84ff8c8bSHong Zhang     } else if (i < ctx->fs) {
462*84ff8c8bSHong Zhang       xi = ctx->xmin+ctx->sf*hs+0.5*hf+(i-ctx->sf)*hf;
463*84ff8c8bSHong Zhang       /* Integrate over cell i using trapezoid rule with N points. */
464*84ff8c8bSHong Zhang       for (k=0; k<dof; k++) u[i*dof+k] = 0;
465*84ff8c8bSHong Zhang       for (j=0; j<N+1; j++) {
466*84ff8c8bSHong Zhang         xj = xi+hf*(j-N/2)/(PetscReal)N;
467*84ff8c8bSHong Zhang         ierr = (*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr);
468*84ff8c8bSHong Zhang         for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N;
469*84ff8c8bSHong Zhang       }
470*84ff8c8bSHong Zhang     } else {
471*84ff8c8bSHong Zhang       xi = ctx->xmin+ctx->sf*hs+(ctx->fs-ctx->sf)*hf+0.5*hs+(i-ctx->fs)*hs;
472*84ff8c8bSHong Zhang       /* Integrate over cell i using trapezoid rule with N points. */
473*84ff8c8bSHong Zhang       for (k=0; k<dof; k++) u[i*dof+k] = 0;
474*84ff8c8bSHong Zhang       for (j=0; j<N+1; j++) {
475*84ff8c8bSHong Zhang         xj = xi+hs*(j-N/2)/(PetscReal)N;
476*84ff8c8bSHong Zhang         ierr = (*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr);
477*84ff8c8bSHong Zhang         for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N;
478*84ff8c8bSHong Zhang       }
479*84ff8c8bSHong Zhang     }
480*84ff8c8bSHong Zhang   }
481*84ff8c8bSHong Zhang   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
482*84ff8c8bSHong Zhang   ierr = PetscFree(uj);CHKERRQ(ierr);
483*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
484*84ff8c8bSHong Zhang }
485*84ff8c8bSHong Zhang 
486*84ff8c8bSHong Zhang static PetscErrorCode SolutionErrorNorms_2WaySplit(FVCtx *ctx,DM da,PetscReal t,Vec X,PetscReal *nrm1)
487*84ff8c8bSHong Zhang {
488*84ff8c8bSHong Zhang   PetscErrorCode    ierr;
489*84ff8c8bSHong Zhang   Vec               Y;
490*84ff8c8bSHong Zhang   PetscInt          i,Mx;
491*84ff8c8bSHong Zhang   const PetscScalar *ptr_X,*ptr_Y;
492*84ff8c8bSHong Zhang   PetscReal         hs,hf;
493*84ff8c8bSHong Zhang 
494*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
495*84ff8c8bSHong Zhang   ierr = VecGetSize(X,&Mx);CHKERRQ(ierr);
496*84ff8c8bSHong Zhang   ierr = VecDuplicate(X,&Y);CHKERRQ(ierr);
497*84ff8c8bSHong Zhang   ierr = FVSample_2WaySplit(ctx,da,t,Y);CHKERRQ(ierr);
498*84ff8c8bSHong Zhang   hs   = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
499*84ff8c8bSHong Zhang   hf   = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
500*84ff8c8bSHong Zhang   ierr = VecGetArrayRead(X,&ptr_X);CHKERRQ(ierr);
501*84ff8c8bSHong Zhang   ierr = VecGetArrayRead(Y,&ptr_Y);CHKERRQ(ierr);
502*84ff8c8bSHong Zhang   for (i=0; i<Mx; i++) {
503*84ff8c8bSHong Zhang     if (i < ctx->sf || i > ctx->fs -1) *nrm1 +=  hs*PetscAbs(ptr_X[i]-ptr_Y[i]);
504*84ff8c8bSHong Zhang     else *nrm1 += hf*PetscAbs(ptr_X[i]-ptr_Y[i]);
505*84ff8c8bSHong Zhang   }
506*84ff8c8bSHong Zhang   ierr = VecRestoreArrayRead(X,&ptr_X);CHKERRQ(ierr);
507*84ff8c8bSHong Zhang   ierr = VecRestoreArrayRead(Y,&ptr_Y);CHKERRQ(ierr);
508*84ff8c8bSHong Zhang   ierr = VecDestroy(&Y);CHKERRQ(ierr);
509*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
510*84ff8c8bSHong Zhang }
511*84ff8c8bSHong Zhang 
512*84ff8c8bSHong Zhang PetscErrorCode FVRHSFunction_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
513*84ff8c8bSHong Zhang {
514*84ff8c8bSHong Zhang   FVCtx          *ctx = (FVCtx*)vctx;
515*84ff8c8bSHong Zhang   PetscErrorCode ierr;
516*84ff8c8bSHong Zhang   PetscInt       i,j,k,Mx,dof,xs,xm,sf = ctx->sf,fs = ctx->fs;
517*84ff8c8bSHong Zhang   PetscReal      hxf,hxs,cfl_idt = 0;
518*84ff8c8bSHong Zhang   PetscScalar    *x,*f,*slope;
519*84ff8c8bSHong Zhang   Vec            Xloc;
520*84ff8c8bSHong Zhang   DM             da;
521*84ff8c8bSHong Zhang 
522*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
523*84ff8c8bSHong Zhang   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
524*84ff8c8bSHong Zhang   ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr);                          /* Xloc contains ghost points                                     */
525*84ff8c8bSHong Zhang   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                              */
526*84ff8c8bSHong Zhang   hxs  = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
527*84ff8c8bSHong Zhang   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
528*84ff8c8bSHong Zhang   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);       /* X is solution vector which does not contain ghost points       */
529*84ff8c8bSHong Zhang   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
530*84ff8c8bSHong Zhang 
531*84ff8c8bSHong Zhang   ierr = VecZeroEntries(F);CHKERRQ(ierr);                                   /* F is the right hand side function corresponds to center points */
532*84ff8c8bSHong Zhang 
533*84ff8c8bSHong Zhang   ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr);
534*84ff8c8bSHong Zhang   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
535*84ff8c8bSHong Zhang   ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);                  /* contains ghost points                                           */
536*84ff8c8bSHong Zhang   ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr);
537*84ff8c8bSHong Zhang 
538*84ff8c8bSHong Zhang   if (ctx->bctype == FVBC_OUTFLOW) {
539*84ff8c8bSHong Zhang     for (i=xs-2; i<0; i++) {
540*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
541*84ff8c8bSHong Zhang     }
542*84ff8c8bSHong Zhang     for (i=Mx; i<xs+xm+2; i++) {
543*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
544*84ff8c8bSHong Zhang     }
545*84ff8c8bSHong Zhang   }
546*84ff8c8bSHong Zhang 
547*84ff8c8bSHong Zhang   if (ctx->bctype == FVBC_INFLOW) {
548*84ff8c8bSHong Zhang     /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253
549*84ff8c8bSHong Zhang     pages 137-138 for the scheme. */
550*84ff8c8bSHong Zhang     if (xs==0) { /* Left Boundary */
551*84ff8c8bSHong Zhang       ctx->physics2.inflow(ctx,time,ctx->xmin,ctx->ub);
552*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
553*84ff8c8bSHong Zhang         if (ctx->physics2.bcinflowindex[j]) {
554*84ff8c8bSHong Zhang           for (i=-2; i<0; i++) x[i*dof+j] = 2.0*ctx->ub[j]-x[-(i+1)*dof+j];
555*84ff8c8bSHong Zhang         } else {
556*84ff8c8bSHong Zhang           for (i=-2; i<0; i++) x[i*dof+j] = x[j]; /* Outflow */
557*84ff8c8bSHong Zhang         }
558*84ff8c8bSHong Zhang       }
559*84ff8c8bSHong Zhang     }
560*84ff8c8bSHong Zhang     if (xs+xm==Mx) { /* Right Boundary */
561*84ff8c8bSHong Zhang       ctx->physics2.inflow(ctx,time,ctx->xmax,ctx->ub);
562*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
563*84ff8c8bSHong Zhang         if (ctx->physics2.bcinflowindex[dof+j]) {
564*84ff8c8bSHong Zhang           for (i=Mx; i<Mx+2; i++) x[i*dof+j] = 2.0*ctx->ub[dof+j]-x[(2*Mx-(i+1))*dof+j];
565*84ff8c8bSHong Zhang         } else {
566*84ff8c8bSHong Zhang           for (i=Mx; i<Mx+2; i++) x[i*dof+j] = x[(Mx-1)*dof+j]; /* Outflow */
567*84ff8c8bSHong Zhang         }
568*84ff8c8bSHong Zhang       }
569*84ff8c8bSHong Zhang     }
570*84ff8c8bSHong Zhang   }
571*84ff8c8bSHong Zhang 
572*84ff8c8bSHong Zhang   for (i=xs-1; i<xs+xm+1; i++) {
573*84ff8c8bSHong Zhang     struct _LimitInfo info;
574*84ff8c8bSHong Zhang     PetscScalar       *cjmpL,*cjmpR;
575*84ff8c8bSHong Zhang     /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
576*84ff8c8bSHong Zhang     ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr);
577*84ff8c8bSHong Zhang     /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
578*84ff8c8bSHong Zhang     ierr  = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr);
579*84ff8c8bSHong Zhang     cjmpL = &ctx->cjmpLR[0];
580*84ff8c8bSHong Zhang     cjmpR = &ctx->cjmpLR[dof];
581*84ff8c8bSHong Zhang     for (j=0; j<dof; j++) {
582*84ff8c8bSHong Zhang       PetscScalar jmpL,jmpR;
583*84ff8c8bSHong Zhang       jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
584*84ff8c8bSHong Zhang       jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
585*84ff8c8bSHong Zhang       for (k=0; k<dof; k++) {
586*84ff8c8bSHong Zhang         cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
587*84ff8c8bSHong Zhang         cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
588*84ff8c8bSHong Zhang       }
589*84ff8c8bSHong Zhang     }
590*84ff8c8bSHong Zhang     /* Apply limiter to the left and right characteristic jumps */
591*84ff8c8bSHong Zhang     info.m  = dof;
592*84ff8c8bSHong Zhang     info.hxs = hxs;
593*84ff8c8bSHong Zhang     info.hxf = hxf;
594*84ff8c8bSHong Zhang     (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope);
595*84ff8c8bSHong Zhang     for (j=0; j<dof; j++) {
596*84ff8c8bSHong Zhang       PetscScalar tmp = 0;
597*84ff8c8bSHong Zhang       for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
598*84ff8c8bSHong Zhang       slope[i*dof+j] = tmp;
599*84ff8c8bSHong Zhang     }
600*84ff8c8bSHong Zhang   }
601*84ff8c8bSHong Zhang 
602*84ff8c8bSHong Zhang   for (i=xs; i<xs+xm+1; i++) {
603*84ff8c8bSHong Zhang     PetscReal   maxspeed;
604*84ff8c8bSHong Zhang     PetscScalar *uL,*uR;
605*84ff8c8bSHong Zhang     uL = &ctx->uLR[0];
606*84ff8c8bSHong Zhang     uR = &ctx->uLR[dof];
607*84ff8c8bSHong Zhang     if (i < sf) { /* slow region */
608*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
609*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
610*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
611*84ff8c8bSHong Zhang       }
612*84ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
613*84ff8c8bSHong Zhang       if (i > xs) {
614*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs;
615*84ff8c8bSHong Zhang       }
616*84ff8c8bSHong Zhang       if (i < xs+xm) {
617*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs;
618*84ff8c8bSHong Zhang       }
619*84ff8c8bSHong Zhang     } else if (i == sf) { /* interface between the slow region and the fast region */
620*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
621*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
622*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
623*84ff8c8bSHong Zhang       }
624*84ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
625*84ff8c8bSHong Zhang       if (i > xs) {
626*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs;
627*84ff8c8bSHong Zhang       }
628*84ff8c8bSHong Zhang       if (i < xs+xm) {
629*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxf;
630*84ff8c8bSHong Zhang       }
631*84ff8c8bSHong Zhang     } else if (i > sf && i < fs) { /* fast region */
632*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
633*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
634*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
635*84ff8c8bSHong Zhang       }
636*84ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
637*84ff8c8bSHong Zhang       if (i > xs) {
638*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf;
639*84ff8c8bSHong Zhang       }
640*84ff8c8bSHong Zhang       if (i < xs+xm) {
641*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxf;
642*84ff8c8bSHong Zhang       }
643*84ff8c8bSHong Zhang     } else if (i == fs) { /* interface between the fast region and the slow region */
644*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
645*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
646*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
647*84ff8c8bSHong Zhang       }
648*84ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
649*84ff8c8bSHong Zhang       if (i > xs) {
650*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf;
651*84ff8c8bSHong Zhang       }
652*84ff8c8bSHong Zhang       if (i < xs+xm) {
653*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs;
654*84ff8c8bSHong Zhang       }
655*84ff8c8bSHong Zhang     } else { /* slow region */
656*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
657*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
658*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
659*84ff8c8bSHong Zhang       }
660*84ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
661*84ff8c8bSHong Zhang       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */
662*84ff8c8bSHong Zhang       if (i > xs) {
663*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs;
664*84ff8c8bSHong Zhang       }
665*84ff8c8bSHong Zhang       if (i < xs+xm) {
666*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs;
667*84ff8c8bSHong Zhang       }
668*84ff8c8bSHong Zhang     }
669*84ff8c8bSHong Zhang   }
670*84ff8c8bSHong Zhang   ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr);
671*84ff8c8bSHong Zhang   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
672*84ff8c8bSHong Zhang   ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
673*84ff8c8bSHong Zhang   ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr);
674*84ff8c8bSHong Zhang   ierr = MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr);
675*84ff8c8bSHong Zhang   if (0) {
676*84ff8c8bSHong Zhang     /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */
677*84ff8c8bSHong Zhang     PetscReal dt,tnow;
678*84ff8c8bSHong Zhang     ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
679*84ff8c8bSHong Zhang     ierr = TSGetTime(ts,&tnow);CHKERRQ(ierr);
680*84ff8c8bSHong Zhang     if (dt > 0.5/ctx->cfl_idt) {
681*84ff8c8bSHong Zhang       if (1) {
682*84ff8c8bSHong Zhang         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);
683*84ff8c8bSHong Zhang       } else SETERRQ2(PETSC_COMM_SELF,1,"Stability constraint exceeded, %g > %g",(double)dt,(double)(ctx->cfl/ctx->cfl_idt));
684*84ff8c8bSHong Zhang     }
685*84ff8c8bSHong Zhang   }
686*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
687*84ff8c8bSHong Zhang }
688*84ff8c8bSHong Zhang 
689*84ff8c8bSHong Zhang /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */
690*84ff8c8bSHong Zhang PetscErrorCode FVRHSFunctionslow_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
691*84ff8c8bSHong Zhang {
692*84ff8c8bSHong Zhang   FVCtx          *ctx = (FVCtx*)vctx;
693*84ff8c8bSHong Zhang   PetscErrorCode ierr;
694*84ff8c8bSHong Zhang   PetscInt       i,j,k,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth;
695*84ff8c8bSHong Zhang   PetscReal      hxs,hxf,cfl_idt = 0;
696*84ff8c8bSHong Zhang   PetscScalar    *x,*f,*slope;
697*84ff8c8bSHong Zhang   Vec            Xloc;
698*84ff8c8bSHong Zhang   DM             da;
699*84ff8c8bSHong Zhang 
700*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
701*84ff8c8bSHong Zhang   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
702*84ff8c8bSHong Zhang   ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr);
703*84ff8c8bSHong Zhang   ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr);
704*84ff8c8bSHong Zhang   hxs  = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
705*84ff8c8bSHong Zhang   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
706*84ff8c8bSHong Zhang   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
707*84ff8c8bSHong Zhang   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
708*84ff8c8bSHong Zhang   ierr = VecZeroEntries(F);CHKERRQ(ierr);
709*84ff8c8bSHong Zhang   ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr);
710*84ff8c8bSHong Zhang   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
711*84ff8c8bSHong Zhang   ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
712*84ff8c8bSHong Zhang   ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr);
713*84ff8c8bSHong Zhang 
714*84ff8c8bSHong Zhang   if (ctx->bctype == FVBC_OUTFLOW) {
715*84ff8c8bSHong Zhang     for (i=xs-2; i<0; i++) {
716*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
717*84ff8c8bSHong Zhang     }
718*84ff8c8bSHong Zhang     for (i=Mx; i<xs+xm+2; i++) {
719*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
720*84ff8c8bSHong Zhang     }
721*84ff8c8bSHong Zhang   }
722*84ff8c8bSHong Zhang   if (ctx->bctype == FVBC_INFLOW) {
723*84ff8c8bSHong Zhang     /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253
724*84ff8c8bSHong Zhang     pages 137-138 for the scheme. */
725*84ff8c8bSHong Zhang     if (xs==0) { /* Left Boundary */
726*84ff8c8bSHong Zhang       ctx->physics2.inflow(ctx,time,ctx->xmin,ctx->ub);
727*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
728*84ff8c8bSHong Zhang         if (ctx->physics2.bcinflowindex[j]==PETSC_TRUE) {
729*84ff8c8bSHong Zhang           for (i=-2; i<0; i++) x[i*dof+j] = 2.0*ctx->ub[j]-x[-(i+1)*dof+j];
730*84ff8c8bSHong Zhang         } else {
731*84ff8c8bSHong Zhang           for (i=-2; i<0; i++) x[i*dof+j] = x[j]; /* Outflow */
732*84ff8c8bSHong Zhang         }
733*84ff8c8bSHong Zhang       }
734*84ff8c8bSHong Zhang     }
735*84ff8c8bSHong Zhang     if (xs+xm==Mx) { /* Right Boundary */
736*84ff8c8bSHong Zhang       ctx->physics2.inflow(ctx,time,ctx->xmax,ctx->ub);
737*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
738*84ff8c8bSHong Zhang         if (ctx->physics2.bcinflowindex[dof+j]==PETSC_TRUE) {
739*84ff8c8bSHong Zhang           for (i=Mx; i<Mx+2; i++) x[i*dof+j] = 2.0*ctx->ub[dof+j]-x[(2*Mx-(i+1))*dof+j];
740*84ff8c8bSHong Zhang         } else {
741*84ff8c8bSHong Zhang           for (i=Mx; i<Mx+2; i++) x[i*dof+j] = x[(Mx-1)*dof+j]; /* Outflow */
742*84ff8c8bSHong Zhang         }
743*84ff8c8bSHong Zhang       }
744*84ff8c8bSHong Zhang     }
745*84ff8c8bSHong Zhang   }
746*84ff8c8bSHong Zhang   for (i=xs-1; i<xs+xm+1; i++) {
747*84ff8c8bSHong Zhang     struct _LimitInfo info;
748*84ff8c8bSHong Zhang     PetscScalar       *cjmpL,*cjmpR;
749*84ff8c8bSHong Zhang     if (i < sf-lsbwidth+1 || i > fs+rsbwidth-2) { /* slow components and the first and last fast components */
750*84ff8c8bSHong Zhang       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
751*84ff8c8bSHong Zhang       ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr);
752*84ff8c8bSHong Zhang       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
753*84ff8c8bSHong Zhang       ierr  = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr);
754*84ff8c8bSHong Zhang       cjmpL = &ctx->cjmpLR[0];
755*84ff8c8bSHong Zhang       cjmpR = &ctx->cjmpLR[dof];
756*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
757*84ff8c8bSHong Zhang         PetscScalar jmpL,jmpR;
758*84ff8c8bSHong Zhang         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
759*84ff8c8bSHong Zhang         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
760*84ff8c8bSHong Zhang         for (k=0; k<dof; k++) {
761*84ff8c8bSHong Zhang           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
762*84ff8c8bSHong Zhang           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
763*84ff8c8bSHong Zhang         }
764*84ff8c8bSHong Zhang       }
765*84ff8c8bSHong Zhang       /* Apply limiter to the left and right characteristic jumps */
766*84ff8c8bSHong Zhang       info.m  = dof;
767*84ff8c8bSHong Zhang       info.hxs = hxs;
768*84ff8c8bSHong Zhang       info.hxf = hxf;
769*84ff8c8bSHong Zhang       (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope);
770*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
771*84ff8c8bSHong Zhang         PetscScalar tmp = 0;
772*84ff8c8bSHong Zhang         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
773*84ff8c8bSHong Zhang           slope[i*dof+j] = tmp;
774*84ff8c8bSHong Zhang       }
775*84ff8c8bSHong Zhang     }
776*84ff8c8bSHong Zhang   }
777*84ff8c8bSHong Zhang 
778*84ff8c8bSHong Zhang   for (i=xs; i<xs+xm+1; i++) {
779*84ff8c8bSHong Zhang     PetscReal   maxspeed;
780*84ff8c8bSHong Zhang     PetscScalar *uL,*uR;
781*84ff8c8bSHong Zhang     uL = &ctx->uLR[0];
782*84ff8c8bSHong Zhang     uR = &ctx->uLR[dof];
783*84ff8c8bSHong Zhang     if (i < sf-lsbwidth) { /* slow region */
784*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
785*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
786*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
787*84ff8c8bSHong Zhang       }
788*84ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
789*84ff8c8bSHong Zhang       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */
790*84ff8c8bSHong Zhang       if (i > xs) {
791*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
792*84ff8c8bSHong Zhang       }
793*84ff8c8bSHong Zhang       if (i < xs+xm) {
794*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
795*84ff8c8bSHong Zhang         islow++;
796*84ff8c8bSHong Zhang       }
797*84ff8c8bSHong Zhang     }
798*84ff8c8bSHong Zhang     if (i == sf-lsbwidth) { /* interface between the slow region and the fast region */
799*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
800*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
801*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
802*84ff8c8bSHong Zhang       }
803*84ff8c8bSHong Zhang       ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
804*84ff8c8bSHong Zhang       if (i > xs) {
805*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
806*84ff8c8bSHong Zhang       }
807*84ff8c8bSHong Zhang     }
808*84ff8c8bSHong Zhang     if (i == fs+rsbwidth) { /* slow region */
809*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
810*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
811*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
812*84ff8c8bSHong Zhang       }
813*84ff8c8bSHong Zhang       ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
814*84ff8c8bSHong Zhang       if (i < xs+xm) {
815*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
816*84ff8c8bSHong Zhang         islow++;
817*84ff8c8bSHong Zhang       }
818*84ff8c8bSHong Zhang     }
819*84ff8c8bSHong Zhang     if (i > fs+rsbwidth) { /* slow region */
820*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
821*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
822*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
823*84ff8c8bSHong Zhang       }
824*84ff8c8bSHong Zhang       ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
825*84ff8c8bSHong Zhang       if (i > xs) {
826*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
827*84ff8c8bSHong Zhang       }
828*84ff8c8bSHong Zhang       if (i < xs+xm) {
829*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
830*84ff8c8bSHong Zhang         islow++;
831*84ff8c8bSHong Zhang       }
832*84ff8c8bSHong Zhang     }
833*84ff8c8bSHong Zhang   }
834*84ff8c8bSHong Zhang   ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr);
835*84ff8c8bSHong Zhang   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
836*84ff8c8bSHong Zhang   ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
837*84ff8c8bSHong Zhang   ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr);
838*84ff8c8bSHong Zhang   ierr = MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr);
839*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
840*84ff8c8bSHong Zhang }
841*84ff8c8bSHong Zhang 
842*84ff8c8bSHong Zhang PetscErrorCode FVRHSFunctionslowbuffer_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
843*84ff8c8bSHong Zhang {
844*84ff8c8bSHong Zhang   FVCtx          *ctx = (FVCtx*)vctx;
845*84ff8c8bSHong Zhang   PetscErrorCode ierr;
846*84ff8c8bSHong Zhang   PetscInt       i,j,k,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth;
847*84ff8c8bSHong Zhang   PetscReal      hxs,hxf;
848*84ff8c8bSHong Zhang   PetscScalar    *x,*f,*slope;
849*84ff8c8bSHong Zhang   Vec            Xloc;
850*84ff8c8bSHong Zhang   DM             da;
851*84ff8c8bSHong Zhang 
852*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
853*84ff8c8bSHong Zhang   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
854*84ff8c8bSHong Zhang   ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr);
855*84ff8c8bSHong Zhang   ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr);
856*84ff8c8bSHong Zhang   hxs  = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
857*84ff8c8bSHong Zhang   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
858*84ff8c8bSHong Zhang   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
859*84ff8c8bSHong Zhang   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
860*84ff8c8bSHong Zhang   ierr = VecZeroEntries(F);CHKERRQ(ierr);
861*84ff8c8bSHong Zhang   ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr);
862*84ff8c8bSHong Zhang   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
863*84ff8c8bSHong Zhang   ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
864*84ff8c8bSHong Zhang   ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr);
865*84ff8c8bSHong Zhang 
866*84ff8c8bSHong Zhang   if (ctx->bctype == FVBC_OUTFLOW) {
867*84ff8c8bSHong Zhang     for (i=xs-2; i<0; i++) {
868*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
869*84ff8c8bSHong Zhang     }
870*84ff8c8bSHong Zhang     for (i=Mx; i<xs+xm+2; i++) {
871*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
872*84ff8c8bSHong Zhang     }
873*84ff8c8bSHong Zhang   }
874*84ff8c8bSHong Zhang   if (ctx->bctype == FVBC_INFLOW) {
875*84ff8c8bSHong Zhang     /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253
876*84ff8c8bSHong Zhang     pages 137-138 for the scheme. */
877*84ff8c8bSHong Zhang     if (xs==0) { /* Left Boundary */
878*84ff8c8bSHong Zhang       ctx->physics2.inflow(ctx,time,ctx->xmin,ctx->ub);
879*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
880*84ff8c8bSHong Zhang         if (ctx->physics2.bcinflowindex[j]==PETSC_TRUE) {
881*84ff8c8bSHong Zhang           for (i=-2; i<0; i++) x[i*dof+j] = 2.0*ctx->ub[j]-x[-(i+1)*dof+j];
882*84ff8c8bSHong Zhang         } else {
883*84ff8c8bSHong Zhang           for (i=-2; i<0; i++) x[i*dof+j] = x[j]; /* Outflow */
884*84ff8c8bSHong Zhang         }
885*84ff8c8bSHong Zhang       }
886*84ff8c8bSHong Zhang     }
887*84ff8c8bSHong Zhang     if (xs+xm==Mx) { /* Right Boundary */
888*84ff8c8bSHong Zhang       ctx->physics2.inflow(ctx,time,ctx->xmax,ctx->ub);
889*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
890*84ff8c8bSHong Zhang         if (ctx->physics2.bcinflowindex[dof+j]==PETSC_TRUE) {
891*84ff8c8bSHong Zhang           for (i=Mx; i<Mx+2; i++) x[i*dof+j] = 2.0*ctx->ub[dof+j]-x[(2*Mx-(i+1))*dof+j];
892*84ff8c8bSHong Zhang         } else {
893*84ff8c8bSHong Zhang           for (i=Mx; i<Mx+2; i++) x[i*dof+j] = x[(Mx-1)*dof+j]; /* Outflow */
894*84ff8c8bSHong Zhang         }
895*84ff8c8bSHong Zhang       }
896*84ff8c8bSHong Zhang     }
897*84ff8c8bSHong Zhang   }
898*84ff8c8bSHong Zhang   for (i=xs-1; i<xs+xm+1; i++) {
899*84ff8c8bSHong Zhang     struct _LimitInfo info;
900*84ff8c8bSHong Zhang     PetscScalar       *cjmpL,*cjmpR;
901*84ff8c8bSHong Zhang     if ((i > sf-lsbwidth-2 && i < sf+1) || (i > fs-2 && i < fs+rsbwidth+1)) {
902*84ff8c8bSHong Zhang       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
903*84ff8c8bSHong Zhang       ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr);
904*84ff8c8bSHong Zhang       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
905*84ff8c8bSHong Zhang       ierr  = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr);
906*84ff8c8bSHong Zhang       cjmpL = &ctx->cjmpLR[0];
907*84ff8c8bSHong Zhang       cjmpR = &ctx->cjmpLR[dof];
908*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
909*84ff8c8bSHong Zhang         PetscScalar jmpL,jmpR;
910*84ff8c8bSHong Zhang         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
911*84ff8c8bSHong Zhang         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
912*84ff8c8bSHong Zhang         for (k=0; k<dof; k++) {
913*84ff8c8bSHong Zhang           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
914*84ff8c8bSHong Zhang           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
915*84ff8c8bSHong Zhang         }
916*84ff8c8bSHong Zhang       }
917*84ff8c8bSHong Zhang       /* Apply limiter to the left and right characteristic jumps */
918*84ff8c8bSHong Zhang       info.m  = dof;
919*84ff8c8bSHong Zhang       info.hxs = hxs;
920*84ff8c8bSHong Zhang       info.hxf = hxf;
921*84ff8c8bSHong Zhang       (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope);
922*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
923*84ff8c8bSHong Zhang         PetscScalar tmp = 0;
924*84ff8c8bSHong Zhang         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
925*84ff8c8bSHong Zhang           slope[i*dof+j] = tmp;
926*84ff8c8bSHong Zhang       }
927*84ff8c8bSHong Zhang     }
928*84ff8c8bSHong Zhang   }
929*84ff8c8bSHong Zhang 
930*84ff8c8bSHong Zhang   for (i=xs; i<xs+xm+1; i++) {
931*84ff8c8bSHong Zhang     PetscReal   maxspeed;
932*84ff8c8bSHong Zhang     PetscScalar *uL,*uR;
933*84ff8c8bSHong Zhang     uL = &ctx->uLR[0];
934*84ff8c8bSHong Zhang     uR = &ctx->uLR[dof];
935*84ff8c8bSHong Zhang     if (i == sf-lsbwidth) {
936*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
937*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
938*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
939*84ff8c8bSHong Zhang       }
940*84ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
941*84ff8c8bSHong Zhang       if (i < xs+xm) {
942*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
943*84ff8c8bSHong Zhang         islow++;
944*84ff8c8bSHong Zhang       }
945*84ff8c8bSHong Zhang     }
946*84ff8c8bSHong Zhang     if (i > sf-lsbwidth && i < sf) {
947*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
948*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
949*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
950*84ff8c8bSHong Zhang       }
951*84ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
952*84ff8c8bSHong Zhang       if (i > xs) {
953*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
954*84ff8c8bSHong Zhang       }
955*84ff8c8bSHong Zhang       if (i < xs+xm) {
956*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
957*84ff8c8bSHong Zhang         islow++;
958*84ff8c8bSHong Zhang       }
959*84ff8c8bSHong Zhang     }
960*84ff8c8bSHong Zhang     if (i == sf) { /* interface between the slow region and the fast region */
961*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
962*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
963*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
964*84ff8c8bSHong Zhang       }
965*84ff8c8bSHong Zhang       ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
966*84ff8c8bSHong Zhang       if (i > xs) {
967*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
968*84ff8c8bSHong Zhang       }
969*84ff8c8bSHong Zhang     }
970*84ff8c8bSHong Zhang     if (i == fs) { /* interface between the fast region and the slow region */
971*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
972*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
973*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
974*84ff8c8bSHong Zhang       }
975*84ff8c8bSHong Zhang       ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
976*84ff8c8bSHong Zhang       if (i < xs+xm) {
977*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
978*84ff8c8bSHong Zhang         islow++;
979*84ff8c8bSHong Zhang       }
980*84ff8c8bSHong Zhang     }
981*84ff8c8bSHong Zhang     if (i > fs && i < fs+rsbwidth) {
982*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
983*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
984*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
985*84ff8c8bSHong Zhang       }
986*84ff8c8bSHong Zhang       ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
987*84ff8c8bSHong Zhang       if (i > xs) {
988*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
989*84ff8c8bSHong Zhang       }
990*84ff8c8bSHong Zhang       if (i < xs+xm) {
991*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
992*84ff8c8bSHong Zhang         islow++;
993*84ff8c8bSHong Zhang       }
994*84ff8c8bSHong Zhang     }
995*84ff8c8bSHong Zhang     if (i == fs+rsbwidth) {
996*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
997*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
998*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
999*84ff8c8bSHong Zhang       }
1000*84ff8c8bSHong Zhang       ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
1001*84ff8c8bSHong Zhang       if (i > xs) {
1002*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
1003*84ff8c8bSHong Zhang       }
1004*84ff8c8bSHong Zhang     }
1005*84ff8c8bSHong Zhang   }
1006*84ff8c8bSHong Zhang   ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr);
1007*84ff8c8bSHong Zhang   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
1008*84ff8c8bSHong Zhang   ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
1009*84ff8c8bSHong Zhang   ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr);
1010*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
1011*84ff8c8bSHong Zhang }
1012*84ff8c8bSHong Zhang 
1013*84ff8c8bSHong Zhang /* --------------------------------- Finite Volume Solver for fast  parts ----------------------------------- */
1014*84ff8c8bSHong Zhang PetscErrorCode FVRHSFunctionfast_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
1015*84ff8c8bSHong Zhang {
1016*84ff8c8bSHong Zhang   FVCtx          *ctx = (FVCtx*)vctx;
1017*84ff8c8bSHong Zhang   PetscErrorCode ierr;
1018*84ff8c8bSHong Zhang   PetscInt       i,j,k,Mx,dof,xs,xm,ifast = 0,sf = ctx->sf,fs = ctx->fs;
1019*84ff8c8bSHong Zhang   PetscReal      hxs,hxf;
1020*84ff8c8bSHong Zhang   PetscScalar    *x,*f,*slope;
1021*84ff8c8bSHong Zhang   Vec            Xloc;
1022*84ff8c8bSHong Zhang   DM             da;
1023*84ff8c8bSHong Zhang 
1024*84ff8c8bSHong Zhang   PetscFunctionBeginUser;
1025*84ff8c8bSHong Zhang   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
1026*84ff8c8bSHong Zhang   ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr);
1027*84ff8c8bSHong Zhang   ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr);
1028*84ff8c8bSHong Zhang   hxs  = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
1029*84ff8c8bSHong Zhang   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
1030*84ff8c8bSHong Zhang   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
1031*84ff8c8bSHong Zhang   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
1032*84ff8c8bSHong Zhang   ierr = VecZeroEntries(F);CHKERRQ(ierr);
1033*84ff8c8bSHong Zhang   ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr);
1034*84ff8c8bSHong Zhang   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
1035*84ff8c8bSHong Zhang   ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
1036*84ff8c8bSHong Zhang   ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr);
1037*84ff8c8bSHong Zhang 
1038*84ff8c8bSHong Zhang   if (ctx->bctype == FVBC_OUTFLOW) {
1039*84ff8c8bSHong Zhang     for (i=xs-2; i<0; i++) {
1040*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
1041*84ff8c8bSHong Zhang     }
1042*84ff8c8bSHong Zhang     for (i=Mx; i<xs+xm+2; i++) {
1043*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
1044*84ff8c8bSHong Zhang     }
1045*84ff8c8bSHong Zhang   }
1046*84ff8c8bSHong Zhang   if (ctx->bctype == FVBC_INFLOW) {
1047*84ff8c8bSHong Zhang     /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253
1048*84ff8c8bSHong Zhang     pages 137-138 for the scheme. */
1049*84ff8c8bSHong Zhang     if (xs==0) { /* Left Boundary */
1050*84ff8c8bSHong Zhang       ctx->physics2.inflow(ctx,time,ctx->xmin,ctx->ub);
1051*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
1052*84ff8c8bSHong Zhang         if (ctx->physics2.bcinflowindex[j]==PETSC_TRUE) {
1053*84ff8c8bSHong Zhang           for (i=-2; i<0; i++) x[i*dof+j] = 2.0*ctx->ub[j]-x[-(i+1)*dof+j];
1054*84ff8c8bSHong Zhang         } else {
1055*84ff8c8bSHong Zhang           for (i=-2; i<0; i++) x[i*dof+j] = x[j]; /* Outflow */
1056*84ff8c8bSHong Zhang         }
1057*84ff8c8bSHong Zhang       }
1058*84ff8c8bSHong Zhang     }
1059*84ff8c8bSHong Zhang     if (xs+xm==Mx) { /* Right Boundary */
1060*84ff8c8bSHong Zhang       ctx->physics2.inflow(ctx,time,ctx->xmax,ctx->ub);
1061*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
1062*84ff8c8bSHong Zhang         if (ctx->physics2.bcinflowindex[dof+j]==PETSC_TRUE) {
1063*84ff8c8bSHong Zhang           for (i=Mx; i<Mx+2; i++) x[i*dof+j] = 2.0*ctx->ub[dof+j]-x[(2*Mx-(i+1))*dof+j];
1064*84ff8c8bSHong Zhang         } else {
1065*84ff8c8bSHong Zhang           for (i=Mx; i<Mx+2; i++) x[i*dof+j] = x[(Mx-1)*dof+j]; /* Outflow */
1066*84ff8c8bSHong Zhang         }
1067*84ff8c8bSHong Zhang       }
1068*84ff8c8bSHong Zhang     }
1069*84ff8c8bSHong Zhang   }
1070*84ff8c8bSHong Zhang   for (i=xs-1; i<xs+xm+1; i++) {
1071*84ff8c8bSHong Zhang     struct _LimitInfo info;
1072*84ff8c8bSHong Zhang     PetscScalar       *cjmpL,*cjmpR;
1073*84ff8c8bSHong Zhang     if (i > sf-2 && i < fs+1) {
1074*84ff8c8bSHong Zhang       ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr);
1075*84ff8c8bSHong Zhang       ierr  = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr);
1076*84ff8c8bSHong Zhang       cjmpL = &ctx->cjmpLR[0];
1077*84ff8c8bSHong Zhang       cjmpR = &ctx->cjmpLR[dof];
1078*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
1079*84ff8c8bSHong Zhang         PetscScalar jmpL,jmpR;
1080*84ff8c8bSHong Zhang         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
1081*84ff8c8bSHong Zhang         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
1082*84ff8c8bSHong Zhang         for (k=0; k<dof; k++) {
1083*84ff8c8bSHong Zhang           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
1084*84ff8c8bSHong Zhang           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
1085*84ff8c8bSHong Zhang         }
1086*84ff8c8bSHong Zhang       }
1087*84ff8c8bSHong Zhang       /* Apply limiter to the left and right characteristic jumps */
1088*84ff8c8bSHong Zhang       info.m  = dof;
1089*84ff8c8bSHong Zhang       info.hxs = hxs;
1090*84ff8c8bSHong Zhang       info.hxf = hxf;
1091*84ff8c8bSHong Zhang       (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope);
1092*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
1093*84ff8c8bSHong Zhang       PetscScalar tmp = 0;
1094*84ff8c8bSHong Zhang       for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
1095*84ff8c8bSHong Zhang         slope[i*dof+j] = tmp;
1096*84ff8c8bSHong Zhang       }
1097*84ff8c8bSHong Zhang     }
1098*84ff8c8bSHong Zhang   }
1099*84ff8c8bSHong Zhang 
1100*84ff8c8bSHong Zhang   for (i=xs; i<xs+xm+1; i++) {
1101*84ff8c8bSHong Zhang     PetscReal   maxspeed;
1102*84ff8c8bSHong Zhang     PetscScalar *uL,*uR;
1103*84ff8c8bSHong Zhang     uL = &ctx->uLR[0];
1104*84ff8c8bSHong Zhang     uR = &ctx->uLR[dof];
1105*84ff8c8bSHong Zhang     if (i == sf) { /* interface between the slow region and the fast region */
1106*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
1107*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
1108*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
1109*84ff8c8bSHong Zhang       }
1110*84ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
1111*84ff8c8bSHong Zhang       if (i < xs+xm) {
1112*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf;
1113*84ff8c8bSHong Zhang         ifast++;
1114*84ff8c8bSHong Zhang       }
1115*84ff8c8bSHong Zhang     }
1116*84ff8c8bSHong Zhang     if (i > sf && i < fs) { /* fast region */
1117*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
1118*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
1119*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
1120*84ff8c8bSHong Zhang       }
1121*84ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
1122*84ff8c8bSHong Zhang       if (i > xs) {
1123*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf;
1124*84ff8c8bSHong Zhang       }
1125*84ff8c8bSHong Zhang       if (i < xs+xm) {
1126*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf;
1127*84ff8c8bSHong Zhang         ifast++;
1128*84ff8c8bSHong Zhang       }
1129*84ff8c8bSHong Zhang     }
1130*84ff8c8bSHong Zhang     if (i == fs) { /* interface between the fast region and the slow region */
1131*84ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
1132*84ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
1133*84ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
1134*84ff8c8bSHong Zhang       }
1135*84ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
1136*84ff8c8bSHong Zhang       if (i > xs) {
1137*84ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf;
1138*84ff8c8bSHong Zhang       }
1139*84ff8c8bSHong Zhang     }
1140*84ff8c8bSHong Zhang   }
1141*84ff8c8bSHong Zhang   ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr);
1142*84ff8c8bSHong Zhang   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
1143*84ff8c8bSHong Zhang   ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
1144*84ff8c8bSHong Zhang   ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr);
1145*84ff8c8bSHong Zhang   PetscFunctionReturn(0);
1146*84ff8c8bSHong Zhang }
1147*84ff8c8bSHong Zhang 
1148*84ff8c8bSHong Zhang int main(int argc,char *argv[])
1149*84ff8c8bSHong Zhang {
1150*84ff8c8bSHong Zhang   char              lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m";
1151*84ff8c8bSHong Zhang   PetscFunctionList limiters   = 0,physics = 0;
1152*84ff8c8bSHong Zhang   MPI_Comm          comm;
1153*84ff8c8bSHong Zhang   TS                ts;
1154*84ff8c8bSHong Zhang   DM                da;
1155*84ff8c8bSHong Zhang   Vec               X,X0,R;
1156*84ff8c8bSHong Zhang   FVCtx             ctx;
1157*84ff8c8bSHong Zhang   PetscInt          i,k,dof,xs,xm,Mx,draw = 0,count_slow,count_fast,islow = 0,ifast =0,islowbuffer = 0,*index_slow,*index_fast,*index_slowbuffer;
1158*84ff8c8bSHong Zhang   PetscBool         view_final = PETSC_FALSE;
1159*84ff8c8bSHong Zhang   PetscReal         ptime,maxtime;
1160*84ff8c8bSHong Zhang   PetscErrorCode    ierr;
1161*84ff8c8bSHong Zhang 
1162*84ff8c8bSHong Zhang   ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
1163*84ff8c8bSHong Zhang   comm = PETSC_COMM_WORLD;
1164*84ff8c8bSHong Zhang   ierr = PetscMemzero(&ctx,sizeof(ctx));CHKERRQ(ierr);
1165*84ff8c8bSHong Zhang 
1166*84ff8c8bSHong Zhang   /* Register limiters to be available on the command line */
1167*84ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&limiters,"upwind"              ,Limit2_Upwind);CHKERRQ(ierr);
1168*84ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&limiters,"lax-wendroff"        ,Limit2_LaxWendroff);CHKERRQ(ierr);
1169*84ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&limiters,"beam-warming"        ,Limit2_BeamWarming);CHKERRQ(ierr);
1170*84ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&limiters,"fromm"               ,Limit2_Fromm);CHKERRQ(ierr);
1171*84ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&limiters,"minmod"              ,Limit2_Minmod);CHKERRQ(ierr);
1172*84ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&limiters,"superbee"            ,Limit2_Superbee);CHKERRQ(ierr);
1173*84ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&limiters,"mc"                  ,Limit2_MC);CHKERRQ(ierr);
1174*84ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&limiters,"koren3"              ,Limit2_Koren3);CHKERRQ(ierr);
1175*84ff8c8bSHong Zhang 
1176*84ff8c8bSHong Zhang   /* Register physical models to be available on the command line */
1177*84ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&physics,"shallow"         ,PhysicsCreate_Shallow);CHKERRQ(ierr);
1178*84ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&physics,"advect"          ,PhysicsCreate_Advect);CHKERRQ(ierr);
1179*84ff8c8bSHong Zhang 
1180*84ff8c8bSHong Zhang   ctx.comm    = comm;
1181*84ff8c8bSHong Zhang   ctx.cfl     = 0.9;
1182*84ff8c8bSHong Zhang   ctx.bctype  = FVBC_PERIODIC;
1183*84ff8c8bSHong Zhang   ctx.xmin    = -1.0;
1184*84ff8c8bSHong Zhang   ctx.xmax    = 1.0;
1185*84ff8c8bSHong Zhang   ctx.initial = 1;
1186*84ff8c8bSHong Zhang   ctx.hratio  = 2;
1187*84ff8c8bSHong Zhang   maxtime     = 10.0;
1188*84ff8c8bSHong Zhang   ctx.simulation = PETSC_FALSE;
1189*84ff8c8bSHong Zhang   ierr = PetscOptionsBegin(comm,NULL,"Finite Volume solver options","");CHKERRQ(ierr);
1190*84ff8c8bSHong Zhang   ierr = PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL);CHKERRQ(ierr);
1191*84ff8c8bSHong Zhang   ierr = PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL);CHKERRQ(ierr);
1192*84ff8c8bSHong Zhang   ierr = PetscOptionsFList("-limit","Name of flux imiter to use","",limiters,lname,lname,sizeof(lname),NULL);CHKERRQ(ierr);
1193*84ff8c8bSHong Zhang   ierr = PetscOptionsFList("-physics","Name of physics model to use","",physics,physname,physname,sizeof(physname),NULL);CHKERRQ(ierr);
1194*84ff8c8bSHong Zhang   ierr = PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL);CHKERRQ(ierr);
1195*84ff8c8bSHong Zhang   ierr = PetscOptionsString("-view_final","Write final solution in ASCII MATLAB format to given file name","",final_fname,final_fname,sizeof(final_fname),&view_final);CHKERRQ(ierr);
1196*84ff8c8bSHong Zhang   ierr = PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL);CHKERRQ(ierr);
1197*84ff8c8bSHong Zhang   ierr = PetscOptionsBool("-exact","Compare errors with exact solution","",ctx.exact,&ctx.exact,NULL);CHKERRQ(ierr);
1198*84ff8c8bSHong Zhang   ierr = PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL);CHKERRQ(ierr);
1199*84ff8c8bSHong Zhang   ierr = PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL);CHKERRQ(ierr);
1200*84ff8c8bSHong Zhang   ierr = PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL);CHKERRQ(ierr);
1201*84ff8c8bSHong Zhang   ierr = PetscOptionsInt("-hratio","Spacing ratio","",ctx.hratio,&ctx.hratio,NULL);CHKERRQ(ierr);
1202*84ff8c8bSHong Zhang   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1203*84ff8c8bSHong Zhang 
1204*84ff8c8bSHong Zhang   /* Choose the limiter from the list of registered limiters */
1205*84ff8c8bSHong Zhang   ierr = PetscFunctionListFind(limiters,lname,&ctx.limit2);CHKERRQ(ierr);
1206*84ff8c8bSHong Zhang   if (!ctx.limit2) SETERRQ1(PETSC_COMM_SELF,1,"Limiter '%s' not found",lname);
1207*84ff8c8bSHong Zhang 
1208*84ff8c8bSHong Zhang   /* Choose the physics from the list of registered models */
1209*84ff8c8bSHong Zhang   {
1210*84ff8c8bSHong Zhang     PetscErrorCode (*r)(FVCtx*);
1211*84ff8c8bSHong Zhang     ierr = PetscFunctionListFind(physics,physname,&r);CHKERRQ(ierr);
1212*84ff8c8bSHong Zhang     if (!r) SETERRQ1(PETSC_COMM_SELF,1,"Physics '%s' not found",physname);
1213*84ff8c8bSHong Zhang     /* Create the physics, will set the number of fields and their names */
1214*84ff8c8bSHong Zhang     ierr = (*r)(&ctx);CHKERRQ(ierr);
1215*84ff8c8bSHong Zhang   }
1216*84ff8c8bSHong Zhang 
1217*84ff8c8bSHong Zhang   /* Create a DMDA to manage the parallel grid */
1218*84ff8c8bSHong Zhang   ierr = DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics2.dof,2,NULL,&da);CHKERRQ(ierr);
1219*84ff8c8bSHong Zhang   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
1220*84ff8c8bSHong Zhang   ierr = DMSetUp(da);CHKERRQ(ierr);
1221*84ff8c8bSHong Zhang   /* Inform the DMDA of the field names provided by the physics. */
1222*84ff8c8bSHong Zhang   /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
1223*84ff8c8bSHong Zhang   for (i=0; i<ctx.physics2.dof; i++) {
1224*84ff8c8bSHong Zhang     ierr = DMDASetFieldName(da,i,ctx.physics2.fieldname[i]);CHKERRQ(ierr);
1225*84ff8c8bSHong Zhang   }
1226*84ff8c8bSHong Zhang   ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr);
1227*84ff8c8bSHong Zhang   ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr);
1228*84ff8c8bSHong Zhang 
1229*84ff8c8bSHong Zhang   /* Set coordinates of cell centers */
1230*84ff8c8bSHong Zhang   ierr = DMDASetUniformCoordinates(da,ctx.xmin+0.5*(ctx.xmax-ctx.xmin)/Mx,ctx.xmax+0.5*(ctx.xmax-ctx.xmin)/Mx,0,0,0,0);CHKERRQ(ierr);
1231*84ff8c8bSHong Zhang 
1232*84ff8c8bSHong Zhang   /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
1233*84ff8c8bSHong Zhang   ierr = PetscMalloc4(dof*dof,&ctx.R,dof*dof,&ctx.Rinv,2*dof,&ctx.cjmpLR,1*dof,&ctx.cslope);CHKERRQ(ierr);
1234*84ff8c8bSHong Zhang   ierr = PetscMalloc3(2*dof,&ctx.uLR,dof,&ctx.flux,dof,&ctx.speeds);CHKERRQ(ierr);
1235*84ff8c8bSHong Zhang   ierr = PetscMalloc1(2*dof,&ctx.ub);CHKERRQ(ierr);
1236*84ff8c8bSHong Zhang 
1237*84ff8c8bSHong Zhang   /* Create a vector to store the solution and to save the initial state */
1238*84ff8c8bSHong Zhang   ierr = DMCreateGlobalVector(da,&X);CHKERRQ(ierr);
1239*84ff8c8bSHong Zhang   ierr = VecDuplicate(X,&X0);CHKERRQ(ierr);
1240*84ff8c8bSHong Zhang   ierr = VecDuplicate(X,&R);CHKERRQ(ierr);
1241*84ff8c8bSHong Zhang 
1242*84ff8c8bSHong Zhang   /* create index for slow parts and fast parts,
1243*84ff8c8bSHong Zhang      count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */
1244*84ff8c8bSHong Zhang   count_slow = Mx/(1.0+ctx.hratio/3.0);
1245*84ff8c8bSHong Zhang   if (count_slow%2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Please adjust grid size Mx (-da_grid_x) and hratio (-hratio) so that Mx/(1+hratio/3) is even");
1246*84ff8c8bSHong Zhang   count_fast = Mx-count_slow;
1247*84ff8c8bSHong Zhang   ctx.sf = count_slow/2;
1248*84ff8c8bSHong Zhang   ctx.fs = ctx.sf+count_fast;
1249*84ff8c8bSHong Zhang   ierr = PetscMalloc1(xm*dof,&index_slow);CHKERRQ(ierr);
1250*84ff8c8bSHong Zhang   ierr = PetscMalloc1(xm*dof,&index_fast);CHKERRQ(ierr);
1251*84ff8c8bSHong Zhang   ierr = PetscMalloc1(8*dof,&index_slowbuffer);CHKERRQ(ierr);
1252*84ff8c8bSHong Zhang   ctx.lsbwidth = 4;
1253*84ff8c8bSHong Zhang   ctx.rsbwidth = 4;
1254*84ff8c8bSHong Zhang 
1255*84ff8c8bSHong Zhang   for (i=xs; i<xs+xm; i++) {
1256*84ff8c8bSHong Zhang     if (i < ctx.sf-ctx.lsbwidth || i > ctx.fs+ctx.rsbwidth-1)
1257*84ff8c8bSHong Zhang       for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k;
1258*84ff8c8bSHong Zhang     else if ((i >= ctx.sf-ctx.lsbwidth && i < ctx.sf) || (i > ctx.fs-1 && i <= ctx.fs+ctx.rsbwidth-1))
1259*84ff8c8bSHong Zhang       for (k=0; k<dof; k++) index_slowbuffer[islowbuffer++] = i*dof+k;
1260*84ff8c8bSHong Zhang     else
1261*84ff8c8bSHong Zhang       for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k;
1262*84ff8c8bSHong Zhang   }
1263*84ff8c8bSHong Zhang   ierr = ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss);CHKERRQ(ierr);
1264*84ff8c8bSHong Zhang   ierr = ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf);CHKERRQ(ierr);
1265*84ff8c8bSHong Zhang   ierr = ISCreateGeneral(PETSC_COMM_WORLD,islowbuffer,index_slowbuffer,PETSC_COPY_VALUES,&ctx.issb);CHKERRQ(ierr);
1266*84ff8c8bSHong Zhang 
1267*84ff8c8bSHong Zhang   /* Create a time-stepping object */
1268*84ff8c8bSHong Zhang   ierr = TSCreate(comm,&ts);CHKERRQ(ierr);
1269*84ff8c8bSHong Zhang   ierr = TSSetDM(ts,da);CHKERRQ(ierr);
1270*84ff8c8bSHong Zhang   ierr = TSSetRHSFunction(ts,R,FVRHSFunction_2WaySplit,&ctx);CHKERRQ(ierr);
1271*84ff8c8bSHong Zhang   ierr = TSRHSSplitSetIS(ts,"slow",ctx.iss);CHKERRQ(ierr);
1272*84ff8c8bSHong Zhang   ierr = TSRHSSplitSetIS(ts,"slowbuffer",ctx.issb);CHKERRQ(ierr);
1273*84ff8c8bSHong Zhang   ierr = TSRHSSplitSetIS(ts,"fast",ctx.isf);CHKERRQ(ierr);
1274*84ff8c8bSHong Zhang   ierr = TSRHSSplitSetRHSFunction(ts,"slow",NULL,FVRHSFunctionslow_2WaySplit,&ctx);CHKERRQ(ierr);
1275*84ff8c8bSHong Zhang   ierr = TSRHSSplitSetRHSFunction(ts,"fast",NULL,FVRHSFunctionfast_2WaySplit,&ctx);CHKERRQ(ierr);
1276*84ff8c8bSHong Zhang   ierr = TSRHSSplitSetRHSFunction(ts,"slowbuffer",NULL,FVRHSFunctionslowbuffer_2WaySplit,&ctx);CHKERRQ(ierr);
1277*84ff8c8bSHong Zhang 
1278*84ff8c8bSHong Zhang   ierr = TSSetType(ts,TSMPRK);CHKERRQ(ierr);
1279*84ff8c8bSHong Zhang   ierr = TSSetMaxTime(ts,maxtime);CHKERRQ(ierr);
1280*84ff8c8bSHong Zhang   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
1281*84ff8c8bSHong Zhang 
1282*84ff8c8bSHong Zhang   /* Compute initial conditions and starting time step */
1283*84ff8c8bSHong Zhang   ierr = FVSample_2WaySplit(&ctx,da,0,X0);CHKERRQ(ierr);
1284*84ff8c8bSHong Zhang   ierr = FVRHSFunction_2WaySplit(ts,0,X0,X,(void*)&ctx);CHKERRQ(ierr); /* Initial function evaluation, only used to determine max speed */
1285*84ff8c8bSHong Zhang   ierr = VecCopy(X0,X);CHKERRQ(ierr);                        /* The function value was not used so we set X=X0 again */
1286*84ff8c8bSHong Zhang   ierr = TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt);CHKERRQ(ierr);
1287*84ff8c8bSHong Zhang   ierr = TSSetFromOptions(ts);CHKERRQ(ierr); /* Take runtime options */
1288*84ff8c8bSHong Zhang   ierr = SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1289*84ff8c8bSHong Zhang   {
1290*84ff8c8bSHong Zhang     PetscInt          steps;
1291*84ff8c8bSHong Zhang     PetscScalar       mass_initial,mass_final,mass_difference,mass_differenceg;
1292*84ff8c8bSHong Zhang     const PetscScalar *ptr_X,*ptr_X0;
1293*84ff8c8bSHong Zhang     const PetscReal   hs = (ctx.xmax-ctx.xmin)*3.0/4.0/count_slow;
1294*84ff8c8bSHong Zhang     const PetscReal   hf = (ctx.xmax-ctx.xmin)/4.0/count_fast;
1295*84ff8c8bSHong Zhang 
1296*84ff8c8bSHong Zhang     ierr = TSSolve(ts,X);CHKERRQ(ierr);
1297*84ff8c8bSHong Zhang     ierr = TSGetSolveTime(ts,&ptime);CHKERRQ(ierr);
1298*84ff8c8bSHong Zhang     ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
1299*84ff8c8bSHong Zhang     /* calculate the total mass at initial time and final time */
1300*84ff8c8bSHong Zhang     mass_initial = 0.0;
1301*84ff8c8bSHong Zhang     mass_final   = 0.0;
1302*84ff8c8bSHong Zhang     ierr = DMDAVecGetArrayRead(da,X0,(void*)&ptr_X0);CHKERRQ(ierr);
1303*84ff8c8bSHong Zhang     ierr = DMDAVecGetArrayRead(da,X,(void*)&ptr_X);CHKERRQ(ierr);
1304*84ff8c8bSHong Zhang     for (i=xs;i<xs+xm;i++) {
1305*84ff8c8bSHong Zhang       if (i < ctx.sf || i > ctx.fs-1) {
1306*84ff8c8bSHong Zhang         for (k=0; k<dof; k++) {
1307*84ff8c8bSHong Zhang           mass_initial = mass_initial + hs*ptr_X0[i*dof+k];
1308*84ff8c8bSHong Zhang           mass_final = mass_final + hs*ptr_X[i*dof+k];
1309*84ff8c8bSHong Zhang         }
1310*84ff8c8bSHong Zhang       } else {
1311*84ff8c8bSHong Zhang         for (k=0; k<dof; k++) {
1312*84ff8c8bSHong Zhang           mass_initial = mass_initial + hf*ptr_X0[i*dof+k];
1313*84ff8c8bSHong Zhang           mass_final = mass_final + hf*ptr_X[i*dof+k];
1314*84ff8c8bSHong Zhang         }
1315*84ff8c8bSHong Zhang       }
1316*84ff8c8bSHong Zhang     }
1317*84ff8c8bSHong Zhang     ierr = DMDAVecRestoreArrayRead(da,X0,(void*)&ptr_X0);CHKERRQ(ierr);
1318*84ff8c8bSHong Zhang     ierr = DMDAVecRestoreArrayRead(da,X,(void*)&ptr_X);CHKERRQ(ierr);
1319*84ff8c8bSHong Zhang     mass_difference = mass_final - mass_initial;
1320*84ff8c8bSHong Zhang     ierr = MPI_Allreduce(&mass_difference,&mass_differenceg,1,MPIU_SCALAR,MPIU_SUM,comm);CHKERRMPI(ierr);
1321*84ff8c8bSHong Zhang     ierr = PetscPrintf(comm,"Mass difference %g\n",(double)mass_differenceg);CHKERRQ(ierr);
1322*84ff8c8bSHong Zhang     ierr = PetscPrintf(comm,"Final time %g, steps %D\n",(double)ptime,steps);CHKERRQ(ierr);
1323*84ff8c8bSHong Zhang     ierr = PetscPrintf(comm,"Maximum allowable stepsize according to CFL %g\n",(double)(1.0/ctx.cfl_idt));CHKERRQ(ierr);
1324*84ff8c8bSHong Zhang     if (ctx.exact) {
1325*84ff8c8bSHong Zhang       PetscReal nrm1=0;
1326*84ff8c8bSHong Zhang       ierr = SolutionErrorNorms_2WaySplit(&ctx,da,ptime,X,&nrm1);CHKERRQ(ierr);
1327*84ff8c8bSHong Zhang       ierr = PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1);CHKERRQ(ierr);
1328*84ff8c8bSHong Zhang     }
1329*84ff8c8bSHong Zhang     if (ctx.simulation) {
1330*84ff8c8bSHong Zhang       PetscReal    nrm1=0;
1331*84ff8c8bSHong Zhang       PetscViewer  fd;
1332*84ff8c8bSHong Zhang       char         filename[PETSC_MAX_PATH_LEN] = "binaryoutput";
1333*84ff8c8bSHong Zhang       Vec          XR;
1334*84ff8c8bSHong Zhang       PetscBool    flg;
1335*84ff8c8bSHong Zhang       const PetscScalar  *ptr_XR;
1336*84ff8c8bSHong Zhang       ierr = PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg);CHKERRQ(ierr);
1337*84ff8c8bSHong Zhang       if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
1338*84ff8c8bSHong Zhang       ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&fd);CHKERRQ(ierr);
1339*84ff8c8bSHong Zhang       ierr = VecDuplicate(X0,&XR);CHKERRQ(ierr);
1340*84ff8c8bSHong Zhang       ierr = VecLoad(XR,fd);CHKERRQ(ierr);
1341*84ff8c8bSHong Zhang       ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
1342*84ff8c8bSHong Zhang       ierr = VecGetArrayRead(X,&ptr_X);CHKERRQ(ierr);
1343*84ff8c8bSHong Zhang       ierr = VecGetArrayRead(XR,&ptr_XR);CHKERRQ(ierr);
1344*84ff8c8bSHong Zhang       for (i=xs;i<xs+xm;i++) {
1345*84ff8c8bSHong Zhang         if (i < ctx.sf || i > ctx.fs-1)
1346*84ff8c8bSHong Zhang           for (k=0; k<dof; k++) nrm1 = nrm1 + hs*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]);
1347*84ff8c8bSHong Zhang         else
1348*84ff8c8bSHong Zhang           for (k=0; k<dof; k++) nrm1 = nrm1 + hf*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]);
1349*84ff8c8bSHong Zhang       }
1350*84ff8c8bSHong Zhang       ierr = VecRestoreArrayRead(X,&ptr_X);CHKERRQ(ierr);
1351*84ff8c8bSHong Zhang       ierr = VecRestoreArrayRead(XR,&ptr_XR);CHKERRQ(ierr);
1352*84ff8c8bSHong Zhang       ierr = PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1);CHKERRQ(ierr);
1353*84ff8c8bSHong Zhang       ierr = VecDestroy(&XR);CHKERRQ(ierr);
1354*84ff8c8bSHong Zhang     }
1355*84ff8c8bSHong Zhang   }
1356*84ff8c8bSHong Zhang 
1357*84ff8c8bSHong Zhang   ierr = SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1358*84ff8c8bSHong Zhang   if (draw & 0x1) {ierr = VecView(X0,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);}
1359*84ff8c8bSHong Zhang   if (draw & 0x2) {ierr = VecView(X,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);}
1360*84ff8c8bSHong Zhang   if (draw & 0x4) {
1361*84ff8c8bSHong Zhang     Vec Y;
1362*84ff8c8bSHong Zhang     ierr = VecDuplicate(X,&Y);CHKERRQ(ierr);
1363*84ff8c8bSHong Zhang     ierr = FVSample_2WaySplit(&ctx,da,ptime,Y);CHKERRQ(ierr);
1364*84ff8c8bSHong Zhang     ierr = VecAYPX(Y,-1,X);CHKERRQ(ierr);
1365*84ff8c8bSHong Zhang     ierr = VecView(Y,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);
1366*84ff8c8bSHong Zhang     ierr = VecDestroy(&Y);CHKERRQ(ierr);
1367*84ff8c8bSHong Zhang   }
1368*84ff8c8bSHong Zhang 
1369*84ff8c8bSHong Zhang   if (view_final) {
1370*84ff8c8bSHong Zhang     PetscViewer viewer;
1371*84ff8c8bSHong Zhang     ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer);CHKERRQ(ierr);
1372*84ff8c8bSHong Zhang     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
1373*84ff8c8bSHong Zhang     ierr = VecView(X,viewer);CHKERRQ(ierr);
1374*84ff8c8bSHong Zhang     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1375*84ff8c8bSHong Zhang     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1376*84ff8c8bSHong Zhang   }
1377*84ff8c8bSHong Zhang 
1378*84ff8c8bSHong Zhang   /* Clean up */
1379*84ff8c8bSHong Zhang   ierr = (*ctx.physics2.destroy)(ctx.physics2.user);CHKERRQ(ierr);
1380*84ff8c8bSHong Zhang   for (i=0; i<ctx.physics2.dof; i++) {ierr = PetscFree(ctx.physics2.fieldname[i]);CHKERRQ(ierr);}
1381*84ff8c8bSHong Zhang   ierr = PetscFree(ctx.physics2.bcinflowindex);CHKERRQ(ierr);
1382*84ff8c8bSHong Zhang   ierr = PetscFree(ctx.ub);
1383*84ff8c8bSHong Zhang   ierr = PetscFree4(ctx.R,ctx.Rinv,ctx.cjmpLR,ctx.cslope);CHKERRQ(ierr);
1384*84ff8c8bSHong Zhang   ierr = PetscFree3(ctx.uLR,ctx.flux,ctx.speeds);CHKERRQ(ierr);
1385*84ff8c8bSHong Zhang   ierr = VecDestroy(&X);CHKERRQ(ierr);
1386*84ff8c8bSHong Zhang   ierr = VecDestroy(&X0);CHKERRQ(ierr);
1387*84ff8c8bSHong Zhang   ierr = VecDestroy(&R);CHKERRQ(ierr);
1388*84ff8c8bSHong Zhang   ierr = DMDestroy(&da);CHKERRQ(ierr);
1389*84ff8c8bSHong Zhang   ierr = TSDestroy(&ts);CHKERRQ(ierr);
1390*84ff8c8bSHong Zhang   ierr = ISDestroy(&ctx.iss);CHKERRQ(ierr);
1391*84ff8c8bSHong Zhang   ierr = ISDestroy(&ctx.isf);CHKERRQ(ierr);
1392*84ff8c8bSHong Zhang   ierr = ISDestroy(&ctx.issb);CHKERRQ(ierr);
1393*84ff8c8bSHong Zhang   ierr = PetscFree(index_slow);CHKERRQ(ierr);
1394*84ff8c8bSHong Zhang   ierr = PetscFree(index_fast);CHKERRQ(ierr);
1395*84ff8c8bSHong Zhang   ierr = PetscFree(index_slowbuffer);CHKERRQ(ierr);
1396*84ff8c8bSHong Zhang   ierr = PetscFunctionListDestroy(&limiters);CHKERRQ(ierr);
1397*84ff8c8bSHong Zhang   ierr = PetscFunctionListDestroy(&physics);CHKERRQ(ierr);
1398*84ff8c8bSHong Zhang   ierr = PetscFinalize();
1399*84ff8c8bSHong Zhang   return ierr;
1400*84ff8c8bSHong Zhang }
1401*84ff8c8bSHong Zhang 
1402*84ff8c8bSHong Zhang /*TEST
1403*84ff8c8bSHong Zhang 
1404*84ff8c8bSHong Zhang     build:
1405*84ff8c8bSHong Zhang       requires: !complex !single
1406*84ff8c8bSHong Zhang       depends: finitevolume1d.c
1407*84ff8c8bSHong Zhang 
1408*84ff8c8bSHong Zhang     test:
1409*84ff8c8bSHong Zhang       suffix: 1
1410*84ff8c8bSHong Zhang       args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22
1411*84ff8c8bSHong Zhang       output_file: output/ex4_1.out
1412*84ff8c8bSHong Zhang 
1413*84ff8c8bSHong Zhang     test:
1414*84ff8c8bSHong Zhang       suffix: 2
1415*84ff8c8bSHong Zhang       args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 1
1416*84ff8c8bSHong Zhang       output_file: output/ex4_1.out
1417*84ff8c8bSHong Zhang 
1418*84ff8c8bSHong Zhang     test:
1419*84ff8c8bSHong Zhang       suffix: 3
1420*84ff8c8bSHong Zhang       args: -da_grid_x 40 -initial 1 -hratio 2 -limit mc -ts_dt 0.1 -ts_max_steps 24 -ts_max_time 7.0 -ts_type mprk -ts_mprk_type 2a22 -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0
1421*84ff8c8bSHong Zhang       output_file: output/ex4_3.out
1422*84ff8c8bSHong Zhang 
1423*84ff8c8bSHong Zhang     test:
1424*84ff8c8bSHong Zhang       suffix: 4
1425*84ff8c8bSHong Zhang       nsize: 2
1426*84ff8c8bSHong Zhang       args: -da_grid_x 40 -initial 1 -hratio 2 -limit mc -ts_dt 0.1 -ts_max_steps 24 -ts_max_time 7.0 -ts_type mprk -ts_mprk_type 2a22 -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 1
1427*84ff8c8bSHong Zhang       output_file: output/ex4_3.out
1428*84ff8c8bSHong Zhang 
1429*84ff8c8bSHong Zhang     test:
1430*84ff8c8bSHong Zhang       suffix: 5
1431*84ff8c8bSHong Zhang       nsize: 4
1432*84ff8c8bSHong Zhang       args: args: -da_grid_x 40 -initial 1 -hratio 2 -limit mc -ts_dt 0.1 -ts_max_steps 24 -ts_max_time 7.0 -ts_type mprk -ts_mprk_type 2a22 -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 1
1433*84ff8c8bSHong Zhang       output_file: output/ex4_3.out
1434*84ff8c8bSHong Zhang TEST*/
1435