xref: /petsc/src/ts/tutorials/multirate/ex4.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
184ff8c8bSHong Zhang /*
284ff8c8bSHong Zhang   Note:
384ff8c8bSHong Zhang     -hratio is the ratio between mesh size of corse grids and fine grids
484ff8c8bSHong Zhang */
584ff8c8bSHong Zhang 
684ff8c8bSHong Zhang static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n"
784ff8c8bSHong Zhang   "  advect      - Constant coefficient scalar advection\n"
884ff8c8bSHong Zhang   "                u_t       + (a*u)_x               = 0\n"
984ff8c8bSHong Zhang   "  shallow     - 1D Shallow water equations (Saint Venant System)\n"
1084ff8c8bSHong Zhang   "                h_t + (q)_x = 0 \n"
1184ff8c8bSHong Zhang   "                q_t + (\frac{q^2}{h} + g/2*h^2)_x = 0  \n"
1284ff8c8bSHong Zhang   "                where, h(x,t) denotes the height of the water and q(x,t) the momentum.\n"
1384ff8c8bSHong Zhang   "  for this toy problem, we choose different meshsizes for different sub-domains (slow-fast-slow), say\n"
1484ff8c8bSHong Zhang   "                hxs  = hratio*hxf \n"
1584ff8c8bSHong Zhang   "  where hxs and hxf are the grid spacings for coarse and fine grids respectively.\n"
1684ff8c8bSHong Zhang   "  exact       - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n"
1784ff8c8bSHong Zhang   "                the states across shocks and rarefactions\n"
1884ff8c8bSHong Zhang   "  simulation  - use reference solution which is generated by smaller time step size to be true solution,\n"
1984ff8c8bSHong Zhang   "                also the reference solution should be generated by user and stored in a binary file.\n"
2084ff8c8bSHong Zhang   "  characteristic - Limit the characteristic variables, this is usually preferred (default)\n"
2184ff8c8bSHong Zhang   "  bc_type     - Boundary condition for the problem, options are: periodic, outflow, inflow "
2284ff8c8bSHong Zhang   "Several problem descriptions (initial data, physics specific features, boundary data) can be chosen with -initial N\n\n"
2384ff8c8bSHong Zhang   "The problem size should be set with -da_grid_x M\n\n";
2484ff8c8bSHong Zhang 
2584ff8c8bSHong Zhang /*
2684ff8c8bSHong Zhang   Example:
2784ff8c8bSHong 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
2884ff8c8bSHong 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
2984ff8c8bSHong 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
3084ff8c8bSHong 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
3184ff8c8bSHong 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
3284ff8c8bSHong Zhang 
3384ff8c8bSHong Zhang   Contributed by: Aidan Hamilton <aidan@udel.edu>
3484ff8c8bSHong Zhang */
3584ff8c8bSHong Zhang 
3684ff8c8bSHong Zhang #include <petscts.h>
3784ff8c8bSHong Zhang #include <petscdm.h>
3884ff8c8bSHong Zhang #include <petscdmda.h>
3984ff8c8bSHong Zhang #include <petscdraw.h>
4084ff8c8bSHong Zhang #include "finitevolume1d.h"
4184ff8c8bSHong Zhang #include <petsc/private/kernels/blockinvert.h>
4284ff8c8bSHong Zhang 
439fbee547SJacob Faibussowitsch static inline PetscReal RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax) { PetscReal range = xmax-xmin; return xmin +PetscFmodReal(range+PetscFmodReal(a,range),range); }
449fbee547SJacob Faibussowitsch static inline PetscReal MaxAbs(PetscReal a,PetscReal b) { return (PetscAbs(a) > PetscAbs(b)) ? a : b; }
4584ff8c8bSHong Zhang /* --------------------------------- Advection ----------------------------------- */
4684ff8c8bSHong Zhang typedef struct {
4784ff8c8bSHong Zhang   PetscReal a;                  /* advective velocity */
4884ff8c8bSHong Zhang } AdvectCtx;
4984ff8c8bSHong Zhang 
5084ff8c8bSHong Zhang static PetscErrorCode PhysicsRiemann_Advect(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
5184ff8c8bSHong Zhang {
5284ff8c8bSHong Zhang   AdvectCtx *ctx = (AdvectCtx*)vctx;
5384ff8c8bSHong Zhang   PetscReal speed;
5484ff8c8bSHong Zhang 
5584ff8c8bSHong Zhang   PetscFunctionBeginUser;
5684ff8c8bSHong Zhang   speed     = ctx->a;
5784ff8c8bSHong Zhang   flux[0]   = PetscMax(0,speed)*uL[0] + PetscMin(0,speed)*uR[0];
5884ff8c8bSHong Zhang   *maxspeed = speed;
5984ff8c8bSHong Zhang   PetscFunctionReturn(0);
6084ff8c8bSHong Zhang }
6184ff8c8bSHong Zhang 
6284ff8c8bSHong Zhang static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
6384ff8c8bSHong Zhang {
6484ff8c8bSHong Zhang   AdvectCtx *ctx = (AdvectCtx*)vctx;
6584ff8c8bSHong Zhang 
6684ff8c8bSHong Zhang   PetscFunctionBeginUser;
6784ff8c8bSHong Zhang   X[0]      = 1.;
6884ff8c8bSHong Zhang   Xi[0]     = 1.;
6984ff8c8bSHong Zhang   speeds[0] = ctx->a;
7084ff8c8bSHong Zhang   PetscFunctionReturn(0);
7184ff8c8bSHong Zhang }
7284ff8c8bSHong Zhang 
7384ff8c8bSHong Zhang static PetscErrorCode PhysicsSample_Advect(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
7484ff8c8bSHong Zhang {
7584ff8c8bSHong Zhang   AdvectCtx *ctx = (AdvectCtx*)vctx;
7684ff8c8bSHong Zhang   PetscReal a    = ctx->a,x0;
7784ff8c8bSHong Zhang 
7884ff8c8bSHong Zhang   PetscFunctionBeginUser;
7984ff8c8bSHong Zhang   switch (bctype) {
8084ff8c8bSHong Zhang     case FVBC_OUTFLOW:   x0 = x-a*t; break;
8184ff8c8bSHong Zhang     case FVBC_PERIODIC: x0 = RangeMod(x-a*t,xmin,xmax); break;
8284ff8c8bSHong Zhang     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown BCType");
8384ff8c8bSHong Zhang   }
8484ff8c8bSHong Zhang   switch (initial) {
8584ff8c8bSHong Zhang     case 0: u[0] = (x0 < 0) ? 1 : -1; break;
8684ff8c8bSHong Zhang     case 1: u[0] = (x0 < 0) ? -1 : 1; break;
8784ff8c8bSHong Zhang     case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break;
8884ff8c8bSHong Zhang     case 3: u[0] = PetscSinReal(2*PETSC_PI*x0); break;
8984ff8c8bSHong Zhang     case 4: u[0] = PetscAbs(x0); break;
9084ff8c8bSHong Zhang     case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2*PETSC_PI*x0)); break;
9184ff8c8bSHong Zhang     case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2-x0 : 0)); break;
9284ff8c8bSHong Zhang     case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI*x0),10.0);break;
9384ff8c8bSHong Zhang     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition");
9484ff8c8bSHong Zhang   }
9584ff8c8bSHong Zhang   PetscFunctionReturn(0);
9684ff8c8bSHong Zhang }
9784ff8c8bSHong Zhang 
9884ff8c8bSHong Zhang static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx)
9984ff8c8bSHong Zhang {
10084ff8c8bSHong Zhang   PetscErrorCode ierr;
10184ff8c8bSHong Zhang   AdvectCtx      *user;
10284ff8c8bSHong Zhang 
10384ff8c8bSHong Zhang   PetscFunctionBeginUser;
10484ff8c8bSHong Zhang   ierr = PetscNew(&user);CHKERRQ(ierr);
10584ff8c8bSHong Zhang   ctx->physics2.sample2         = PhysicsSample_Advect;
10684ff8c8bSHong Zhang   ctx->physics2.riemann2        = PhysicsRiemann_Advect;
10784ff8c8bSHong Zhang   ctx->physics2.characteristic2 = PhysicsCharacteristic_Advect;
10884ff8c8bSHong Zhang   ctx->physics2.destroy         = PhysicsDestroy_SimpleFree;
10984ff8c8bSHong Zhang   ctx->physics2.user            = user;
11084ff8c8bSHong Zhang   ctx->physics2.dof             = 1;
11184ff8c8bSHong Zhang 
11284ff8c8bSHong Zhang   ierr = PetscStrallocpy("u",&ctx->physics2.fieldname[0]);CHKERRQ(ierr);
11384ff8c8bSHong Zhang   user->a = 1;
11484ff8c8bSHong Zhang   ierr = PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");CHKERRQ(ierr);
11584ff8c8bSHong Zhang     ierr = PetscOptionsReal("-physics_advect_a","Speed","",user->a,&user->a,NULL);CHKERRQ(ierr);
11684ff8c8bSHong Zhang   ierr = PetscOptionsEnd();CHKERRQ(ierr);
11784ff8c8bSHong Zhang   PetscFunctionReturn(0);
11884ff8c8bSHong Zhang }
11984ff8c8bSHong Zhang /* --------------------------------- Shallow Water ----------------------------------- */
12084ff8c8bSHong Zhang 
12184ff8c8bSHong Zhang typedef struct {
12284ff8c8bSHong Zhang   PetscReal gravity;
12384ff8c8bSHong Zhang } ShallowCtx;
12484ff8c8bSHong Zhang 
1259fbee547SJacob Faibussowitsch static inline void ShallowFlux(ShallowCtx *phys,const PetscScalar *u,PetscScalar *f)
12684ff8c8bSHong Zhang {
12784ff8c8bSHong Zhang   f[0] = u[1];
12884ff8c8bSHong Zhang   f[1] = PetscSqr(u[1])/u[0] + 0.5*phys->gravity*PetscSqr(u[0]);
12984ff8c8bSHong Zhang }
13084ff8c8bSHong Zhang 
13184ff8c8bSHong Zhang static PetscErrorCode PhysicsRiemann_Shallow_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
13284ff8c8bSHong Zhang {
13384ff8c8bSHong Zhang   ShallowCtx                *phys = (ShallowCtx*)vctx;
13484ff8c8bSHong Zhang   PetscScalar               g    = phys->gravity,ustar[2],cL,cR,c,cstar;
13584ff8c8bSHong Zhang   struct {PetscScalar h,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]},star;
13684ff8c8bSHong Zhang   PetscInt                  i;
13784ff8c8bSHong Zhang 
13884ff8c8bSHong Zhang   PetscFunctionBeginUser;
139*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!(L.h > 0 && R.h > 0),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative");
14084ff8c8bSHong Zhang   cL = PetscSqrtScalar(g*L.h);
14184ff8c8bSHong Zhang   cR = PetscSqrtScalar(g*R.h);
14284ff8c8bSHong Zhang   c  = PetscMax(cL,cR);
14384ff8c8bSHong Zhang   {
14484ff8c8bSHong Zhang     /* Solve for star state */
14584ff8c8bSHong Zhang     const PetscInt maxits = 50;
14684ff8c8bSHong Zhang     PetscScalar tmp,res,res0=0,h0,h = 0.5*(L.h + R.h); /* initial guess */
14784ff8c8bSHong Zhang     h0 = h;
14884ff8c8bSHong Zhang     for (i=0; i<maxits; i++) {
14984ff8c8bSHong Zhang       PetscScalar fr,fl,dfr,dfl;
15084ff8c8bSHong Zhang       fl = (L.h < h)
15184ff8c8bSHong Zhang         ? PetscSqrtScalar(0.5*g*(h*h - L.h*L.h)*(1/L.h - 1/h)) /* shock */
15284ff8c8bSHong Zhang         : 2*PetscSqrtScalar(g*h) - 2*PetscSqrtScalar(g*L.h);   /* rarefaction */
15384ff8c8bSHong Zhang       fr = (R.h < h)
15484ff8c8bSHong Zhang         ? PetscSqrtScalar(0.5*g*(h*h - R.h*R.h)*(1/R.h - 1/h)) /* shock */
15584ff8c8bSHong Zhang         : 2*PetscSqrtScalar(g*h) - 2*PetscSqrtScalar(g*R.h);   /* rarefaction */
15684ff8c8bSHong Zhang       res = R.u - L.u + fr + fl;
157*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(PetscIsInfOrNanScalar(res),PETSC_COMM_SELF,PETSC_ERR_FP,"Infinity or Not-a-Number generated in computation");
15884ff8c8bSHong Zhang       if (PetscAbsScalar(res) < PETSC_SQRT_MACHINE_EPSILON || (i > 0 && PetscAbsScalar(h-h0) < PETSC_SQRT_MACHINE_EPSILON)) {
15984ff8c8bSHong Zhang         star.h = h;
16084ff8c8bSHong Zhang         star.u = L.u - fl;
16184ff8c8bSHong Zhang         goto converged;
16284ff8c8bSHong Zhang       } else if (i > 0 && PetscAbsScalar(res) >= PetscAbsScalar(res0)) {        /* Line search */
16384ff8c8bSHong Zhang         h = 0.8*h0 + 0.2*h;
16484ff8c8bSHong Zhang         continue;
16584ff8c8bSHong Zhang       }
16684ff8c8bSHong Zhang       /* Accept the last step and take another */
16784ff8c8bSHong Zhang       res0 = res;
16884ff8c8bSHong Zhang       h0 = h;
16984ff8c8bSHong 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);
17084ff8c8bSHong 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);
17184ff8c8bSHong Zhang       tmp = h - res/(dfr+dfl);
17284ff8c8bSHong Zhang       if (tmp <= 0) h /= 2;   /* Guard against Newton shooting off to a negative thickness */
17384ff8c8bSHong Zhang       else h = tmp;
174*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!((h > 0) && PetscIsNormalScalar(h)),PETSC_COMM_SELF,PETSC_ERR_FP,"non-normal iterate h=%g",(double)h);
17584ff8c8bSHong Zhang     }
17698921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF,1,"Newton iteration for star.h diverged after %D iterations",i);
17784ff8c8bSHong Zhang   }
17884ff8c8bSHong Zhang converged:
17984ff8c8bSHong Zhang   cstar = PetscSqrtScalar(g*star.h);
18084ff8c8bSHong Zhang   if (L.u-cL < 0 && 0 < star.u-cstar) { /* 1-wave is sonic rarefaction */
18184ff8c8bSHong Zhang     PetscScalar ufan[2];
18284ff8c8bSHong Zhang     ufan[0] = 1/g*PetscSqr(L.u/3 + 2./3*cL);
18384ff8c8bSHong Zhang     ufan[1] = PetscSqrtScalar(g*ufan[0])*ufan[0];
18484ff8c8bSHong Zhang     ShallowFlux(phys,ufan,flux);
18584ff8c8bSHong Zhang   } else if (star.u+cstar < 0 && 0 < R.u+cR) { /* 2-wave is sonic rarefaction */
18684ff8c8bSHong Zhang     PetscScalar ufan[2];
18784ff8c8bSHong Zhang     ufan[0] = 1/g*PetscSqr(R.u/3 - 2./3*cR);
18884ff8c8bSHong Zhang     ufan[1] = -PetscSqrtScalar(g*ufan[0])*ufan[0];
18984ff8c8bSHong Zhang     ShallowFlux(phys,ufan,flux);
19084ff8c8bSHong 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)) {
19184ff8c8bSHong Zhang     /* 1-wave is right-travelling shock (supersonic) */
19284ff8c8bSHong Zhang     ShallowFlux(phys,uL,flux);
19384ff8c8bSHong 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)) {
19484ff8c8bSHong Zhang     /* 2-wave is left-travelling shock (supersonic) */
19584ff8c8bSHong Zhang     ShallowFlux(phys,uR,flux);
19684ff8c8bSHong Zhang   } else {
19784ff8c8bSHong Zhang     ustar[0] = star.h;
19884ff8c8bSHong Zhang     ustar[1] = star.h*star.u;
19984ff8c8bSHong Zhang     ShallowFlux(phys,ustar,flux);
20084ff8c8bSHong Zhang   }
20184ff8c8bSHong Zhang   *maxspeed = MaxAbs(MaxAbs(star.u-cstar,star.u+cstar),MaxAbs(L.u-cL,R.u+cR));
20284ff8c8bSHong Zhang   PetscFunctionReturn(0);
20384ff8c8bSHong Zhang }
20484ff8c8bSHong Zhang 
20584ff8c8bSHong Zhang static PetscErrorCode PhysicsRiemann_Shallow_Rusanov(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
20684ff8c8bSHong Zhang {
20784ff8c8bSHong Zhang   ShallowCtx                *phys = (ShallowCtx*)vctx;
20884ff8c8bSHong Zhang   PetscScalar               g = phys->gravity,fL[2],fR[2],s;
20984ff8c8bSHong Zhang   struct {PetscScalar h,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]};
21084ff8c8bSHong Zhang   PetscReal                 tol = 1e-6;
21184ff8c8bSHong Zhang 
21284ff8c8bSHong Zhang   PetscFunctionBeginUser;
21384ff8c8bSHong Zhang   /* Positivity preserving modification*/
21484ff8c8bSHong Zhang   if (L.h < tol) L.u = 0.0;
21584ff8c8bSHong Zhang   if (R.h < tol) R.u = 0.0;
21684ff8c8bSHong Zhang 
21784ff8c8bSHong Zhang   /*simple pos preserve limiter*/
21884ff8c8bSHong Zhang   if (L.h < 0) L.h = 0;
21984ff8c8bSHong Zhang   if (R.h < 0) R.h = 0;
22084ff8c8bSHong Zhang 
22184ff8c8bSHong Zhang   ShallowFlux(phys,uL,fL);
22284ff8c8bSHong Zhang   ShallowFlux(phys,uR,fR);
22384ff8c8bSHong Zhang 
22484ff8c8bSHong Zhang   s         = PetscMax(PetscAbs(L.u)+PetscSqrtScalar(g*L.h),PetscAbs(R.u)+PetscSqrtScalar(g*R.h));
22584ff8c8bSHong Zhang   flux[0]   = 0.5*(fL[0] + fR[0]) + 0.5*s*(L.h - R.h);
22684ff8c8bSHong Zhang   flux[1]   = 0.5*(fL[1] + fR[1]) + 0.5*s*(uL[1] - uR[1]);
22784ff8c8bSHong Zhang   *maxspeed = s;
22884ff8c8bSHong Zhang   PetscFunctionReturn(0);
22984ff8c8bSHong Zhang }
23084ff8c8bSHong Zhang 
23184ff8c8bSHong Zhang static PetscErrorCode PhysicsCharacteristic_Conservative(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
23284ff8c8bSHong Zhang {
23384ff8c8bSHong Zhang   PetscInt i,j;
23484ff8c8bSHong Zhang 
23584ff8c8bSHong Zhang   PetscFunctionBeginUser;
23684ff8c8bSHong Zhang   for (i=0; i<m; i++) {
23784ff8c8bSHong Zhang     for (j=0; j<m; j++) Xi[i*m+j] = X[i*m+j] = (PetscScalar)(i==j);
23884ff8c8bSHong Zhang     speeds[i] = PETSC_MAX_REAL; /* Indicates invalid */
23984ff8c8bSHong Zhang   }
24084ff8c8bSHong Zhang   PetscFunctionReturn(0);
24184ff8c8bSHong Zhang }
24284ff8c8bSHong Zhang 
24384ff8c8bSHong Zhang static PetscErrorCode PhysicsCharacteristic_Shallow(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
24484ff8c8bSHong Zhang {
24584ff8c8bSHong Zhang   ShallowCtx     *phys = (ShallowCtx*)vctx;
24684ff8c8bSHong Zhang   PetscReal      c;
24784ff8c8bSHong Zhang   PetscErrorCode ierr;
24884ff8c8bSHong Zhang   PetscReal      tol = 1e-6;
24984ff8c8bSHong Zhang 
25084ff8c8bSHong Zhang   PetscFunctionBeginUser;
25184ff8c8bSHong Zhang   c           = PetscSqrtScalar(u[0]*phys->gravity);
25284ff8c8bSHong Zhang 
25384ff8c8bSHong Zhang   if (u[0] < tol) { /*Use conservative variables*/
25484ff8c8bSHong Zhang     X[0*2+0]  = 1;
25584ff8c8bSHong Zhang     X[0*2+1]  = 0;
25684ff8c8bSHong Zhang     X[1*2+0]  = 0;
25784ff8c8bSHong Zhang     X[1*2+1]  = 1;
25884ff8c8bSHong Zhang   } else {
25984ff8c8bSHong Zhang     speeds[0] = u[1]/u[0] - c;
26084ff8c8bSHong Zhang     speeds[1] = u[1]/u[0] + c;
26184ff8c8bSHong Zhang     X[0*2+0]  = 1;
26284ff8c8bSHong Zhang     X[0*2+1]  = speeds[0];
26384ff8c8bSHong Zhang     X[1*2+0]  = 1;
26484ff8c8bSHong Zhang     X[1*2+1]  = speeds[1];
26584ff8c8bSHong Zhang   }
26684ff8c8bSHong Zhang 
26784ff8c8bSHong Zhang   ierr = PetscArraycpy(Xi,X,4);CHKERRQ(ierr);
26884ff8c8bSHong Zhang   ierr = PetscKernel_A_gets_inverse_A_2(Xi,0,PETSC_FALSE,NULL);CHKERRQ(ierr);
26984ff8c8bSHong Zhang   PetscFunctionReturn(0);
27084ff8c8bSHong Zhang }
27184ff8c8bSHong Zhang 
27284ff8c8bSHong Zhang static PetscErrorCode PhysicsSample_Shallow(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
27384ff8c8bSHong Zhang {
27484ff8c8bSHong Zhang   PetscFunctionBeginUser;
275*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(t > 0,PETSC_COMM_SELF,PETSC_ERR_SUP,"Exact solutions not implemented for t > 0");
27684ff8c8bSHong Zhang   switch (initial) {
27784ff8c8bSHong Zhang     case 0:
27884ff8c8bSHong Zhang       u[0] = (x < 0) ? 2 : 1; /* Standard Dam Break Problem */
27984ff8c8bSHong Zhang       u[1] = (x < 0) ? 0 : 0;
28084ff8c8bSHong Zhang       break;
28184ff8c8bSHong Zhang     case 1:
28284ff8c8bSHong Zhang       u[0] = (x < 10) ?   1 : 0.1; /*The Next 5 problems are standard Riemann problem tests */
28384ff8c8bSHong Zhang       u[1] = (x < 10) ? 2.5 : 0;
28484ff8c8bSHong Zhang       break;
28584ff8c8bSHong Zhang     case 2:
28684ff8c8bSHong Zhang       u[0] = (x < 25) ?  1 : 1;
28784ff8c8bSHong Zhang       u[1] = (x < 25) ? -5 : 5;
28884ff8c8bSHong Zhang       break;
28984ff8c8bSHong Zhang     case 3:
29084ff8c8bSHong Zhang       u[0] = (x < 20) ?  1 : 1e-12;
29184ff8c8bSHong Zhang       u[1] = (x < 20) ?  0 : 0;
29284ff8c8bSHong Zhang       break;
29384ff8c8bSHong Zhang     case 4:
29484ff8c8bSHong Zhang       u[0] = (x < 30) ? 1e-12 : 1;
29584ff8c8bSHong Zhang       u[1] = (x < 30) ? 0 : 0;
29684ff8c8bSHong Zhang       break;
29784ff8c8bSHong Zhang     case 5:
29884ff8c8bSHong Zhang       u[0] = (x < 25) ?  0.1 : 0.1;
29984ff8c8bSHong Zhang       u[1] = (x < 25) ? -0.3 : 0.3;
30084ff8c8bSHong Zhang       break;
30184ff8c8bSHong Zhang     case 6:
30284ff8c8bSHong Zhang       u[0] = 1+0.5*PetscSinReal(2*PETSC_PI*x);
30384ff8c8bSHong Zhang       u[1] = 1*u[0];
30484ff8c8bSHong Zhang       break;
30584ff8c8bSHong Zhang     case 7:
30684ff8c8bSHong Zhang       if (x < -0.1) {
30784ff8c8bSHong Zhang        u[0] = 1e-9;
30884ff8c8bSHong Zhang        u[1] = 0.0;
30984ff8c8bSHong Zhang       } else if (x < 0.1) {
31084ff8c8bSHong Zhang        u[0] = 1.0;
31184ff8c8bSHong Zhang        u[1] = 0.0;
31284ff8c8bSHong Zhang       } else {
31384ff8c8bSHong Zhang        u[0] = 1e-9;
31484ff8c8bSHong Zhang        u[1] = 0.0;
31584ff8c8bSHong Zhang       }
31684ff8c8bSHong Zhang       break;
31784ff8c8bSHong Zhang     case 8:
31884ff8c8bSHong Zhang      if (x < -0.1) {
31984ff8c8bSHong Zhang        u[0] = 2;
32084ff8c8bSHong Zhang        u[1] = 0.0;
32184ff8c8bSHong Zhang       } else if (x < 0.1) {
32284ff8c8bSHong Zhang        u[0] = 3.0;
32384ff8c8bSHong Zhang        u[1] = 0.0;
32484ff8c8bSHong Zhang       } else {
32584ff8c8bSHong Zhang        u[0] = 2;
32684ff8c8bSHong Zhang        u[1] = 0.0;
32784ff8c8bSHong Zhang       }
32884ff8c8bSHong Zhang       break;
32984ff8c8bSHong Zhang     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition");
33084ff8c8bSHong Zhang   }
33184ff8c8bSHong Zhang   PetscFunctionReturn(0);
33284ff8c8bSHong Zhang }
33384ff8c8bSHong Zhang 
33484ff8c8bSHong Zhang /* Implements inflow conditions for the given initial conditions. Which conditions are actually enforced depends on
33584ff8c8bSHong Zhang    on the results of PhysicsSetInflowType_Shallow. */
33684ff8c8bSHong Zhang static PetscErrorCode PhysicsInflow_Shallow(void *vctx,PetscReal t,PetscReal x,PetscReal *u)
33784ff8c8bSHong Zhang {
33884ff8c8bSHong Zhang   FVCtx          *ctx = (FVCtx*)vctx;
33984ff8c8bSHong Zhang 
34084ff8c8bSHong Zhang   PetscFunctionBeginUser;
34184ff8c8bSHong Zhang   if (ctx->bctype == FVBC_INFLOW) {
34284ff8c8bSHong Zhang     switch (ctx->initial) {
34384ff8c8bSHong Zhang       case 0:
34484ff8c8bSHong Zhang       case 1:
34584ff8c8bSHong Zhang       case 2:
34684ff8c8bSHong Zhang       case 3:
34784ff8c8bSHong Zhang       case 4:
34884ff8c8bSHong Zhang       case 5:
34984ff8c8bSHong Zhang         u[0] = 0; u[1] = 0.0; /* Left boundary conditions */
35084ff8c8bSHong Zhang         u[2] = 0; u[3] = 0.0; /* Right boundary conditions */
35184ff8c8bSHong Zhang         break;
35284ff8c8bSHong Zhang       case 6:
35384ff8c8bSHong Zhang         u[0] = 0; u[1] = 0.0; /* Left boundary conditions */
35484ff8c8bSHong Zhang         u[2] = 0; u[3] = 0.0; /* Right boundary conditions */
35584ff8c8bSHong Zhang         break;
35684ff8c8bSHong Zhang       case 7:
35784ff8c8bSHong Zhang         u[0] = 0; u[1] = 0.0; /* Left boundary conditions */
35884ff8c8bSHong Zhang         u[2] = 0; u[3] = 0.0; /* Right boundary conditions */
35984ff8c8bSHong Zhang         break;
36084ff8c8bSHong Zhang       case 8:
36184ff8c8bSHong Zhang         u[0] = 0; u[1] = 1.0; /* Left boundary conditions */
36284ff8c8bSHong Zhang         u[2] = 0; u[3] = -1.0;/* Right boundary conditions */
36384ff8c8bSHong Zhang         break;
36484ff8c8bSHong Zhang       default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition");
36584ff8c8bSHong Zhang     }
36684ff8c8bSHong Zhang   }
36784ff8c8bSHong Zhang   PetscFunctionReturn(0);
36884ff8c8bSHong Zhang }
36984ff8c8bSHong Zhang 
37084ff8c8bSHong Zhang /* Selects which boundary conditions are marked as inflow and which as outflow when FVBC_INFLOW is selected. */
37184ff8c8bSHong Zhang static PetscErrorCode PhysicsSetInflowType_Shallow(FVCtx *ctx)
37284ff8c8bSHong Zhang {
37384ff8c8bSHong Zhang   PetscFunctionBeginUser;
37484ff8c8bSHong Zhang   switch (ctx->initial) {
37584ff8c8bSHong Zhang     case 0:
37684ff8c8bSHong Zhang     case 1:
37784ff8c8bSHong Zhang     case 2:
37884ff8c8bSHong Zhang     case 3:
37984ff8c8bSHong Zhang     case 4:
38084ff8c8bSHong Zhang     case 5:
38184ff8c8bSHong Zhang     case 6:
38284ff8c8bSHong Zhang     case 7:
38384ff8c8bSHong Zhang     case 8: /* Fix left and right momentum, height is outflow */
38484ff8c8bSHong Zhang       ctx->physics2.bcinflowindex[0] = PETSC_FALSE;
38584ff8c8bSHong Zhang       ctx->physics2.bcinflowindex[1] = PETSC_TRUE;
38684ff8c8bSHong Zhang       ctx->physics2.bcinflowindex[2] = PETSC_FALSE;
38784ff8c8bSHong Zhang       ctx->physics2.bcinflowindex[3] = PETSC_TRUE;
38884ff8c8bSHong Zhang       break;
38984ff8c8bSHong Zhang     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition");
39084ff8c8bSHong Zhang   }
39184ff8c8bSHong Zhang   PetscFunctionReturn(0);
39284ff8c8bSHong Zhang }
39384ff8c8bSHong Zhang 
39484ff8c8bSHong Zhang static PetscErrorCode PhysicsCreate_Shallow(FVCtx *ctx)
39584ff8c8bSHong Zhang {
39684ff8c8bSHong Zhang   PetscErrorCode    ierr;
39784ff8c8bSHong Zhang   ShallowCtx        *user;
39884ff8c8bSHong Zhang   PetscFunctionList rlist = 0,rclist = 0;
39984ff8c8bSHong Zhang   char              rname[256] = "rusanov",rcname[256] = "characteristic";
40084ff8c8bSHong Zhang 
40184ff8c8bSHong Zhang   PetscFunctionBeginUser;
40284ff8c8bSHong Zhang   ierr = PetscNew(&user);CHKERRQ(ierr);
40384ff8c8bSHong Zhang   ctx->physics2.sample2         = PhysicsSample_Shallow;
40484ff8c8bSHong Zhang   ctx->physics2.inflow          = PhysicsInflow_Shallow;
40584ff8c8bSHong Zhang   ctx->physics2.destroy         = PhysicsDestroy_SimpleFree;
40684ff8c8bSHong Zhang   ctx->physics2.riemann2        = PhysicsRiemann_Shallow_Rusanov;
40784ff8c8bSHong Zhang   ctx->physics2.characteristic2 = PhysicsCharacteristic_Shallow;
40884ff8c8bSHong Zhang   ctx->physics2.user            = user;
40984ff8c8bSHong Zhang   ctx->physics2.dof             = 2;
41084ff8c8bSHong Zhang 
41184ff8c8bSHong Zhang   PetscMalloc1(2*(ctx->physics2.dof),&ctx->physics2.bcinflowindex);
41284ff8c8bSHong Zhang   PhysicsSetInflowType_Shallow(ctx);
41384ff8c8bSHong Zhang 
41484ff8c8bSHong Zhang   ierr = PetscStrallocpy("height",&ctx->physics2.fieldname[0]);CHKERRQ(ierr);
41584ff8c8bSHong Zhang   ierr = PetscStrallocpy("momentum",&ctx->physics2.fieldname[1]);CHKERRQ(ierr);
41684ff8c8bSHong Zhang 
41784ff8c8bSHong Zhang   user->gravity = 9.81;
41884ff8c8bSHong Zhang 
41984ff8c8bSHong Zhang   ierr = RiemannListAdd_2WaySplit(&rlist,"exact",  PhysicsRiemann_Shallow_Exact);CHKERRQ(ierr);
42084ff8c8bSHong Zhang   ierr = RiemannListAdd_2WaySplit(&rlist,"rusanov",PhysicsRiemann_Shallow_Rusanov);CHKERRQ(ierr);
42184ff8c8bSHong Zhang   ierr = ReconstructListAdd_2WaySplit(&rclist,"characteristic",PhysicsCharacteristic_Shallow);CHKERRQ(ierr);
42284ff8c8bSHong Zhang   ierr = ReconstructListAdd_2WaySplit(&rclist,"conservative",PhysicsCharacteristic_Conservative);CHKERRQ(ierr);
42384ff8c8bSHong Zhang   ierr = PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for Shallow","");CHKERRQ(ierr);
42484ff8c8bSHong Zhang     ierr = PetscOptionsReal("-physics_shallow_gravity","Gravity","",user->gravity,&user->gravity,NULL);CHKERRQ(ierr);
42584ff8c8bSHong Zhang     ierr = PetscOptionsFList("-physics_shallow_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL);CHKERRQ(ierr);
42684ff8c8bSHong Zhang     ierr = PetscOptionsFList("-physics_shallow_reconstruct","Reconstruction","",rclist,rcname,rcname,sizeof(rcname),NULL);CHKERRQ(ierr);
42784ff8c8bSHong Zhang   ierr = PetscOptionsEnd();CHKERRQ(ierr);
42884ff8c8bSHong Zhang   ierr = RiemannListFind_2WaySplit(rlist,rname,&ctx->physics2.riemann2);CHKERRQ(ierr);
42984ff8c8bSHong Zhang   ierr = ReconstructListFind_2WaySplit(rclist,rcname,&ctx->physics2.characteristic2);CHKERRQ(ierr);
43084ff8c8bSHong Zhang   ierr = PetscFunctionListDestroy(&rlist);CHKERRQ(ierr);
43184ff8c8bSHong Zhang   ierr = PetscFunctionListDestroy(&rclist);CHKERRQ(ierr);
43284ff8c8bSHong Zhang   PetscFunctionReturn(0);
43384ff8c8bSHong Zhang }
43484ff8c8bSHong Zhang 
43584ff8c8bSHong Zhang PetscErrorCode FVSample_2WaySplit(FVCtx *ctx,DM da,PetscReal time,Vec U)
43684ff8c8bSHong Zhang {
43784ff8c8bSHong Zhang   PetscErrorCode  ierr;
43884ff8c8bSHong Zhang   PetscScalar     *u,*uj,xj,xi;
43984ff8c8bSHong Zhang   PetscInt        i,j,k,dof,xs,xm,Mx;
44084ff8c8bSHong Zhang   const PetscInt  N = 200;
44184ff8c8bSHong Zhang   PetscReal       hs,hf;
44284ff8c8bSHong Zhang 
44384ff8c8bSHong Zhang   PetscFunctionBeginUser;
444*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!ctx->physics2.sample2,PETSC_COMM_SELF,PETSC_ERR_SUP,"Physics has not provided a sampling function");
44584ff8c8bSHong Zhang   ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr);
44684ff8c8bSHong Zhang   ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr);
44784ff8c8bSHong Zhang   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
44884ff8c8bSHong Zhang   ierr = PetscMalloc1(dof,&uj);CHKERRQ(ierr);
44984ff8c8bSHong Zhang   hs   = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
45084ff8c8bSHong Zhang   hf   = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
45184ff8c8bSHong Zhang   for (i=xs; i<xs+xm; i++) {
45284ff8c8bSHong Zhang     if (i < ctx->sf) {
45384ff8c8bSHong Zhang       xi = ctx->xmin+0.5*hs+i*hs;
45484ff8c8bSHong Zhang       /* Integrate over cell i using trapezoid rule with N points. */
45584ff8c8bSHong Zhang       for (k=0; k<dof; k++) u[i*dof+k] = 0;
45684ff8c8bSHong Zhang       for (j=0; j<N+1; j++) {
45784ff8c8bSHong Zhang         xj = xi+hs*(j-N/2)/(PetscReal)N;
45884ff8c8bSHong Zhang         ierr = (*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr);
45984ff8c8bSHong Zhang         for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N;
46084ff8c8bSHong Zhang       }
46184ff8c8bSHong Zhang     } else if (i < ctx->fs) {
46284ff8c8bSHong Zhang       xi = ctx->xmin+ctx->sf*hs+0.5*hf+(i-ctx->sf)*hf;
46384ff8c8bSHong Zhang       /* Integrate over cell i using trapezoid rule with N points. */
46484ff8c8bSHong Zhang       for (k=0; k<dof; k++) u[i*dof+k] = 0;
46584ff8c8bSHong Zhang       for (j=0; j<N+1; j++) {
46684ff8c8bSHong Zhang         xj = xi+hf*(j-N/2)/(PetscReal)N;
46784ff8c8bSHong Zhang         ierr = (*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr);
46884ff8c8bSHong Zhang         for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N;
46984ff8c8bSHong Zhang       }
47084ff8c8bSHong Zhang     } else {
47184ff8c8bSHong Zhang       xi = ctx->xmin+ctx->sf*hs+(ctx->fs-ctx->sf)*hf+0.5*hs+(i-ctx->fs)*hs;
47284ff8c8bSHong Zhang       /* Integrate over cell i using trapezoid rule with N points. */
47384ff8c8bSHong Zhang       for (k=0; k<dof; k++) u[i*dof+k] = 0;
47484ff8c8bSHong Zhang       for (j=0; j<N+1; j++) {
47584ff8c8bSHong Zhang         xj = xi+hs*(j-N/2)/(PetscReal)N;
47684ff8c8bSHong Zhang         ierr = (*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr);
47784ff8c8bSHong Zhang         for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N;
47884ff8c8bSHong Zhang       }
47984ff8c8bSHong Zhang     }
48084ff8c8bSHong Zhang   }
48184ff8c8bSHong Zhang   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
48284ff8c8bSHong Zhang   ierr = PetscFree(uj);CHKERRQ(ierr);
48384ff8c8bSHong Zhang   PetscFunctionReturn(0);
48484ff8c8bSHong Zhang }
48584ff8c8bSHong Zhang 
48684ff8c8bSHong Zhang static PetscErrorCode SolutionErrorNorms_2WaySplit(FVCtx *ctx,DM da,PetscReal t,Vec X,PetscReal *nrm1)
48784ff8c8bSHong Zhang {
48884ff8c8bSHong Zhang   PetscErrorCode    ierr;
48984ff8c8bSHong Zhang   Vec               Y;
49084ff8c8bSHong Zhang   PetscInt          i,Mx;
49184ff8c8bSHong Zhang   const PetscScalar *ptr_X,*ptr_Y;
49284ff8c8bSHong Zhang   PetscReal         hs,hf;
49384ff8c8bSHong Zhang 
49484ff8c8bSHong Zhang   PetscFunctionBeginUser;
49584ff8c8bSHong Zhang   ierr = VecGetSize(X,&Mx);CHKERRQ(ierr);
49684ff8c8bSHong Zhang   ierr = VecDuplicate(X,&Y);CHKERRQ(ierr);
49784ff8c8bSHong Zhang   ierr = FVSample_2WaySplit(ctx,da,t,Y);CHKERRQ(ierr);
49884ff8c8bSHong Zhang   hs   = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
49984ff8c8bSHong Zhang   hf   = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
50084ff8c8bSHong Zhang   ierr = VecGetArrayRead(X,&ptr_X);CHKERRQ(ierr);
50184ff8c8bSHong Zhang   ierr = VecGetArrayRead(Y,&ptr_Y);CHKERRQ(ierr);
50284ff8c8bSHong Zhang   for (i=0; i<Mx; i++) {
50384ff8c8bSHong Zhang     if (i < ctx->sf || i > ctx->fs -1) *nrm1 +=  hs*PetscAbs(ptr_X[i]-ptr_Y[i]);
50484ff8c8bSHong Zhang     else *nrm1 += hf*PetscAbs(ptr_X[i]-ptr_Y[i]);
50584ff8c8bSHong Zhang   }
50684ff8c8bSHong Zhang   ierr = VecRestoreArrayRead(X,&ptr_X);CHKERRQ(ierr);
50784ff8c8bSHong Zhang   ierr = VecRestoreArrayRead(Y,&ptr_Y);CHKERRQ(ierr);
50884ff8c8bSHong Zhang   ierr = VecDestroy(&Y);CHKERRQ(ierr);
50984ff8c8bSHong Zhang   PetscFunctionReturn(0);
51084ff8c8bSHong Zhang }
51184ff8c8bSHong Zhang 
51284ff8c8bSHong Zhang PetscErrorCode FVRHSFunction_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
51384ff8c8bSHong Zhang {
51484ff8c8bSHong Zhang   FVCtx          *ctx = (FVCtx*)vctx;
51584ff8c8bSHong Zhang   PetscErrorCode ierr;
51684ff8c8bSHong Zhang   PetscInt       i,j,k,Mx,dof,xs,xm,sf = ctx->sf,fs = ctx->fs;
51784ff8c8bSHong Zhang   PetscReal      hxf,hxs,cfl_idt = 0;
51884ff8c8bSHong Zhang   PetscScalar    *x,*f,*slope;
51984ff8c8bSHong Zhang   Vec            Xloc;
52084ff8c8bSHong Zhang   DM             da;
52184ff8c8bSHong Zhang 
52284ff8c8bSHong Zhang   PetscFunctionBeginUser;
52384ff8c8bSHong Zhang   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
52484ff8c8bSHong Zhang   ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr);                          /* Xloc contains ghost points                                     */
52584ff8c8bSHong 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                              */
52684ff8c8bSHong Zhang   hxs  = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
52784ff8c8bSHong Zhang   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
52884ff8c8bSHong Zhang   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);       /* X is solution vector which does not contain ghost points       */
52984ff8c8bSHong Zhang   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
53084ff8c8bSHong Zhang 
53184ff8c8bSHong Zhang   ierr = VecZeroEntries(F);CHKERRQ(ierr);                                   /* F is the right hand side function corresponds to center points */
53284ff8c8bSHong Zhang 
53384ff8c8bSHong Zhang   ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr);
53484ff8c8bSHong Zhang   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
53584ff8c8bSHong Zhang   ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);                  /* contains ghost points                                           */
53684ff8c8bSHong Zhang   ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr);
53784ff8c8bSHong Zhang 
53884ff8c8bSHong Zhang   if (ctx->bctype == FVBC_OUTFLOW) {
53984ff8c8bSHong Zhang     for (i=xs-2; i<0; i++) {
54084ff8c8bSHong Zhang       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
54184ff8c8bSHong Zhang     }
54284ff8c8bSHong Zhang     for (i=Mx; i<xs+xm+2; i++) {
54384ff8c8bSHong Zhang       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
54484ff8c8bSHong Zhang     }
54584ff8c8bSHong Zhang   }
54684ff8c8bSHong Zhang 
54784ff8c8bSHong Zhang   if (ctx->bctype == FVBC_INFLOW) {
54884ff8c8bSHong Zhang     /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253
54984ff8c8bSHong Zhang     pages 137-138 for the scheme. */
55084ff8c8bSHong Zhang     if (xs==0) { /* Left Boundary */
55184ff8c8bSHong Zhang       ctx->physics2.inflow(ctx,time,ctx->xmin,ctx->ub);
55284ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
55384ff8c8bSHong Zhang         if (ctx->physics2.bcinflowindex[j]) {
55484ff8c8bSHong Zhang           for (i=-2; i<0; i++) x[i*dof+j] = 2.0*ctx->ub[j]-x[-(i+1)*dof+j];
55584ff8c8bSHong Zhang         } else {
55684ff8c8bSHong Zhang           for (i=-2; i<0; i++) x[i*dof+j] = x[j]; /* Outflow */
55784ff8c8bSHong Zhang         }
55884ff8c8bSHong Zhang       }
55984ff8c8bSHong Zhang     }
56084ff8c8bSHong Zhang     if (xs+xm==Mx) { /* Right Boundary */
56184ff8c8bSHong Zhang       ctx->physics2.inflow(ctx,time,ctx->xmax,ctx->ub);
56284ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
56384ff8c8bSHong Zhang         if (ctx->physics2.bcinflowindex[dof+j]) {
56484ff8c8bSHong 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];
56584ff8c8bSHong Zhang         } else {
56684ff8c8bSHong Zhang           for (i=Mx; i<Mx+2; i++) x[i*dof+j] = x[(Mx-1)*dof+j]; /* Outflow */
56784ff8c8bSHong Zhang         }
56884ff8c8bSHong Zhang       }
56984ff8c8bSHong Zhang     }
57084ff8c8bSHong Zhang   }
57184ff8c8bSHong Zhang 
57284ff8c8bSHong Zhang   for (i=xs-1; i<xs+xm+1; i++) {
57384ff8c8bSHong Zhang     struct _LimitInfo info;
57484ff8c8bSHong Zhang     PetscScalar       *cjmpL,*cjmpR;
57584ff8c8bSHong Zhang     /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
57684ff8c8bSHong Zhang     ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr);
57784ff8c8bSHong Zhang     /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
57884ff8c8bSHong Zhang     ierr  = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr);
57984ff8c8bSHong Zhang     cjmpL = &ctx->cjmpLR[0];
58084ff8c8bSHong Zhang     cjmpR = &ctx->cjmpLR[dof];
58184ff8c8bSHong Zhang     for (j=0; j<dof; j++) {
58284ff8c8bSHong Zhang       PetscScalar jmpL,jmpR;
58384ff8c8bSHong Zhang       jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
58484ff8c8bSHong Zhang       jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
58584ff8c8bSHong Zhang       for (k=0; k<dof; k++) {
58684ff8c8bSHong Zhang         cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
58784ff8c8bSHong Zhang         cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
58884ff8c8bSHong Zhang       }
58984ff8c8bSHong Zhang     }
59084ff8c8bSHong Zhang     /* Apply limiter to the left and right characteristic jumps */
59184ff8c8bSHong Zhang     info.m  = dof;
59284ff8c8bSHong Zhang     info.hxs = hxs;
59384ff8c8bSHong Zhang     info.hxf = hxf;
59484ff8c8bSHong Zhang     (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope);
59584ff8c8bSHong Zhang     for (j=0; j<dof; j++) {
59684ff8c8bSHong Zhang       PetscScalar tmp = 0;
59784ff8c8bSHong Zhang       for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
59884ff8c8bSHong Zhang       slope[i*dof+j] = tmp;
59984ff8c8bSHong Zhang     }
60084ff8c8bSHong Zhang   }
60184ff8c8bSHong Zhang 
60284ff8c8bSHong Zhang   for (i=xs; i<xs+xm+1; i++) {
60384ff8c8bSHong Zhang     PetscReal   maxspeed;
60484ff8c8bSHong Zhang     PetscScalar *uL,*uR;
60584ff8c8bSHong Zhang     uL = &ctx->uLR[0];
60684ff8c8bSHong Zhang     uR = &ctx->uLR[dof];
60784ff8c8bSHong Zhang     if (i < sf) { /* slow region */
60884ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
60984ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
61084ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
61184ff8c8bSHong Zhang       }
61284ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
61384ff8c8bSHong Zhang       if (i > xs) {
61484ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs;
61584ff8c8bSHong Zhang       }
61684ff8c8bSHong Zhang       if (i < xs+xm) {
61784ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs;
61884ff8c8bSHong Zhang       }
61984ff8c8bSHong Zhang     } else if (i == sf) { /* interface between the slow region and the fast region */
62084ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
62184ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
62284ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
62384ff8c8bSHong Zhang       }
62484ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
62584ff8c8bSHong Zhang       if (i > xs) {
62684ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs;
62784ff8c8bSHong Zhang       }
62884ff8c8bSHong Zhang       if (i < xs+xm) {
62984ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxf;
63084ff8c8bSHong Zhang       }
63184ff8c8bSHong Zhang     } else if (i > sf && i < fs) { /* fast region */
63284ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
63384ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
63484ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
63584ff8c8bSHong Zhang       }
63684ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
63784ff8c8bSHong Zhang       if (i > xs) {
63884ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf;
63984ff8c8bSHong Zhang       }
64084ff8c8bSHong Zhang       if (i < xs+xm) {
64184ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxf;
64284ff8c8bSHong Zhang       }
64384ff8c8bSHong Zhang     } else if (i == fs) { /* interface between the fast region and the slow region */
64484ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
64584ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
64684ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
64784ff8c8bSHong Zhang       }
64884ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
64984ff8c8bSHong Zhang       if (i > xs) {
65084ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf;
65184ff8c8bSHong Zhang       }
65284ff8c8bSHong Zhang       if (i < xs+xm) {
65384ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs;
65484ff8c8bSHong Zhang       }
65584ff8c8bSHong Zhang     } else { /* slow region */
65684ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
65784ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
65884ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
65984ff8c8bSHong Zhang       }
66084ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
66184ff8c8bSHong Zhang       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */
66284ff8c8bSHong Zhang       if (i > xs) {
66384ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs;
66484ff8c8bSHong Zhang       }
66584ff8c8bSHong Zhang       if (i < xs+xm) {
66684ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs;
66784ff8c8bSHong Zhang       }
66884ff8c8bSHong Zhang     }
66984ff8c8bSHong Zhang   }
67084ff8c8bSHong Zhang   ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr);
67184ff8c8bSHong Zhang   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
67284ff8c8bSHong Zhang   ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
67384ff8c8bSHong Zhang   ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr);
67484ff8c8bSHong Zhang   ierr = MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr);
67584ff8c8bSHong Zhang   if (0) {
67684ff8c8bSHong Zhang     /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */
67784ff8c8bSHong Zhang     PetscReal dt,tnow;
67884ff8c8bSHong Zhang     ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
67984ff8c8bSHong Zhang     ierr = TSGetTime(ts,&tnow);CHKERRQ(ierr);
68084ff8c8bSHong Zhang     if (dt > 0.5/ctx->cfl_idt) {
68184ff8c8bSHong Zhang       if (1) {
68284ff8c8bSHong 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);
68398921bdaSJacob Faibussowitsch       } else SETERRQ(PETSC_COMM_SELF,1,"Stability constraint exceeded, %g > %g",(double)dt,(double)(ctx->cfl/ctx->cfl_idt));
68484ff8c8bSHong Zhang     }
68584ff8c8bSHong Zhang   }
68684ff8c8bSHong Zhang   PetscFunctionReturn(0);
68784ff8c8bSHong Zhang }
68884ff8c8bSHong Zhang 
68984ff8c8bSHong Zhang /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */
69084ff8c8bSHong Zhang PetscErrorCode FVRHSFunctionslow_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
69184ff8c8bSHong Zhang {
69284ff8c8bSHong Zhang   FVCtx          *ctx = (FVCtx*)vctx;
69384ff8c8bSHong Zhang   PetscErrorCode ierr;
69484ff8c8bSHong Zhang   PetscInt       i,j,k,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth;
69584ff8c8bSHong Zhang   PetscReal      hxs,hxf,cfl_idt = 0;
69684ff8c8bSHong Zhang   PetscScalar    *x,*f,*slope;
69784ff8c8bSHong Zhang   Vec            Xloc;
69884ff8c8bSHong Zhang   DM             da;
69984ff8c8bSHong Zhang 
70084ff8c8bSHong Zhang   PetscFunctionBeginUser;
70184ff8c8bSHong Zhang   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
70284ff8c8bSHong Zhang   ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr);
70384ff8c8bSHong Zhang   ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr);
70484ff8c8bSHong Zhang   hxs  = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
70584ff8c8bSHong Zhang   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
70684ff8c8bSHong Zhang   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
70784ff8c8bSHong Zhang   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
70884ff8c8bSHong Zhang   ierr = VecZeroEntries(F);CHKERRQ(ierr);
70984ff8c8bSHong Zhang   ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr);
71084ff8c8bSHong Zhang   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
71184ff8c8bSHong Zhang   ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
71284ff8c8bSHong Zhang   ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr);
71384ff8c8bSHong Zhang 
71484ff8c8bSHong Zhang   if (ctx->bctype == FVBC_OUTFLOW) {
71584ff8c8bSHong Zhang     for (i=xs-2; i<0; i++) {
71684ff8c8bSHong Zhang       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
71784ff8c8bSHong Zhang     }
71884ff8c8bSHong Zhang     for (i=Mx; i<xs+xm+2; i++) {
71984ff8c8bSHong Zhang       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
72084ff8c8bSHong Zhang     }
72184ff8c8bSHong Zhang   }
72284ff8c8bSHong Zhang   if (ctx->bctype == FVBC_INFLOW) {
72384ff8c8bSHong Zhang     /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253
72484ff8c8bSHong Zhang     pages 137-138 for the scheme. */
72584ff8c8bSHong Zhang     if (xs==0) { /* Left Boundary */
72684ff8c8bSHong Zhang       ctx->physics2.inflow(ctx,time,ctx->xmin,ctx->ub);
72784ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
72884ff8c8bSHong Zhang         if (ctx->physics2.bcinflowindex[j]==PETSC_TRUE) {
72984ff8c8bSHong Zhang           for (i=-2; i<0; i++) x[i*dof+j] = 2.0*ctx->ub[j]-x[-(i+1)*dof+j];
73084ff8c8bSHong Zhang         } else {
73184ff8c8bSHong Zhang           for (i=-2; i<0; i++) x[i*dof+j] = x[j]; /* Outflow */
73284ff8c8bSHong Zhang         }
73384ff8c8bSHong Zhang       }
73484ff8c8bSHong Zhang     }
73584ff8c8bSHong Zhang     if (xs+xm==Mx) { /* Right Boundary */
73684ff8c8bSHong Zhang       ctx->physics2.inflow(ctx,time,ctx->xmax,ctx->ub);
73784ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
73884ff8c8bSHong Zhang         if (ctx->physics2.bcinflowindex[dof+j]==PETSC_TRUE) {
73984ff8c8bSHong 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];
74084ff8c8bSHong Zhang         } else {
74184ff8c8bSHong Zhang           for (i=Mx; i<Mx+2; i++) x[i*dof+j] = x[(Mx-1)*dof+j]; /* Outflow */
74284ff8c8bSHong Zhang         }
74384ff8c8bSHong Zhang       }
74484ff8c8bSHong Zhang     }
74584ff8c8bSHong Zhang   }
74684ff8c8bSHong Zhang   for (i=xs-1; i<xs+xm+1; i++) {
74784ff8c8bSHong Zhang     struct _LimitInfo info;
74884ff8c8bSHong Zhang     PetscScalar       *cjmpL,*cjmpR;
74984ff8c8bSHong Zhang     if (i < sf-lsbwidth+1 || i > fs+rsbwidth-2) { /* slow components and the first and last fast components */
75084ff8c8bSHong Zhang       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
75184ff8c8bSHong Zhang       ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr);
75284ff8c8bSHong Zhang       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
75384ff8c8bSHong Zhang       ierr  = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr);
75484ff8c8bSHong Zhang       cjmpL = &ctx->cjmpLR[0];
75584ff8c8bSHong Zhang       cjmpR = &ctx->cjmpLR[dof];
75684ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
75784ff8c8bSHong Zhang         PetscScalar jmpL,jmpR;
75884ff8c8bSHong Zhang         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
75984ff8c8bSHong Zhang         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
76084ff8c8bSHong Zhang         for (k=0; k<dof; k++) {
76184ff8c8bSHong Zhang           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
76284ff8c8bSHong Zhang           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
76384ff8c8bSHong Zhang         }
76484ff8c8bSHong Zhang       }
76584ff8c8bSHong Zhang       /* Apply limiter to the left and right characteristic jumps */
76684ff8c8bSHong Zhang       info.m  = dof;
76784ff8c8bSHong Zhang       info.hxs = hxs;
76884ff8c8bSHong Zhang       info.hxf = hxf;
76984ff8c8bSHong Zhang       (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope);
77084ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
77184ff8c8bSHong Zhang         PetscScalar tmp = 0;
77284ff8c8bSHong Zhang         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
77384ff8c8bSHong Zhang           slope[i*dof+j] = tmp;
77484ff8c8bSHong Zhang       }
77584ff8c8bSHong Zhang     }
77684ff8c8bSHong Zhang   }
77784ff8c8bSHong Zhang 
77884ff8c8bSHong Zhang   for (i=xs; i<xs+xm+1; i++) {
77984ff8c8bSHong Zhang     PetscReal   maxspeed;
78084ff8c8bSHong Zhang     PetscScalar *uL,*uR;
78184ff8c8bSHong Zhang     uL = &ctx->uLR[0];
78284ff8c8bSHong Zhang     uR = &ctx->uLR[dof];
78384ff8c8bSHong Zhang     if (i < sf-lsbwidth) { /* slow region */
78484ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
78584ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
78684ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
78784ff8c8bSHong Zhang       }
78884ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
78984ff8c8bSHong Zhang       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */
79084ff8c8bSHong Zhang       if (i > xs) {
79184ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
79284ff8c8bSHong Zhang       }
79384ff8c8bSHong Zhang       if (i < xs+xm) {
79484ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
79584ff8c8bSHong Zhang         islow++;
79684ff8c8bSHong Zhang       }
79784ff8c8bSHong Zhang     }
79884ff8c8bSHong Zhang     if (i == sf-lsbwidth) { /* interface between the slow region and the fast region */
79984ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
80084ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
80184ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
80284ff8c8bSHong Zhang       }
80384ff8c8bSHong Zhang       ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
80484ff8c8bSHong Zhang       if (i > xs) {
80584ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
80684ff8c8bSHong Zhang       }
80784ff8c8bSHong Zhang     }
80884ff8c8bSHong Zhang     if (i == fs+rsbwidth) { /* slow region */
80984ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
81084ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
81184ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
81284ff8c8bSHong Zhang       }
81384ff8c8bSHong Zhang       ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
81484ff8c8bSHong Zhang       if (i < xs+xm) {
81584ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
81684ff8c8bSHong Zhang         islow++;
81784ff8c8bSHong Zhang       }
81884ff8c8bSHong Zhang     }
81984ff8c8bSHong Zhang     if (i > fs+rsbwidth) { /* slow region */
82084ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
82184ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
82284ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
82384ff8c8bSHong Zhang       }
82484ff8c8bSHong Zhang       ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
82584ff8c8bSHong Zhang       if (i > xs) {
82684ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
82784ff8c8bSHong Zhang       }
82884ff8c8bSHong Zhang       if (i < xs+xm) {
82984ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
83084ff8c8bSHong Zhang         islow++;
83184ff8c8bSHong Zhang       }
83284ff8c8bSHong Zhang     }
83384ff8c8bSHong Zhang   }
83484ff8c8bSHong Zhang   ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr);
83584ff8c8bSHong Zhang   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
83684ff8c8bSHong Zhang   ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
83784ff8c8bSHong Zhang   ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr);
83884ff8c8bSHong Zhang   ierr = MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr);
83984ff8c8bSHong Zhang   PetscFunctionReturn(0);
84084ff8c8bSHong Zhang }
84184ff8c8bSHong Zhang 
84284ff8c8bSHong Zhang PetscErrorCode FVRHSFunctionslowbuffer_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
84384ff8c8bSHong Zhang {
84484ff8c8bSHong Zhang   FVCtx          *ctx = (FVCtx*)vctx;
84584ff8c8bSHong Zhang   PetscErrorCode ierr;
84684ff8c8bSHong Zhang   PetscInt       i,j,k,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth;
84784ff8c8bSHong Zhang   PetscReal      hxs,hxf;
84884ff8c8bSHong Zhang   PetscScalar    *x,*f,*slope;
84984ff8c8bSHong Zhang   Vec            Xloc;
85084ff8c8bSHong Zhang   DM             da;
85184ff8c8bSHong Zhang 
85284ff8c8bSHong Zhang   PetscFunctionBeginUser;
85384ff8c8bSHong Zhang   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
85484ff8c8bSHong Zhang   ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr);
85584ff8c8bSHong Zhang   ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr);
85684ff8c8bSHong Zhang   hxs  = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
85784ff8c8bSHong Zhang   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
85884ff8c8bSHong Zhang   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
85984ff8c8bSHong Zhang   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
86084ff8c8bSHong Zhang   ierr = VecZeroEntries(F);CHKERRQ(ierr);
86184ff8c8bSHong Zhang   ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr);
86284ff8c8bSHong Zhang   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
86384ff8c8bSHong Zhang   ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
86484ff8c8bSHong Zhang   ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr);
86584ff8c8bSHong Zhang 
86684ff8c8bSHong Zhang   if (ctx->bctype == FVBC_OUTFLOW) {
86784ff8c8bSHong Zhang     for (i=xs-2; i<0; i++) {
86884ff8c8bSHong Zhang       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
86984ff8c8bSHong Zhang     }
87084ff8c8bSHong Zhang     for (i=Mx; i<xs+xm+2; i++) {
87184ff8c8bSHong Zhang       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
87284ff8c8bSHong Zhang     }
87384ff8c8bSHong Zhang   }
87484ff8c8bSHong Zhang   if (ctx->bctype == FVBC_INFLOW) {
87584ff8c8bSHong Zhang     /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253
87684ff8c8bSHong Zhang     pages 137-138 for the scheme. */
87784ff8c8bSHong Zhang     if (xs==0) { /* Left Boundary */
87884ff8c8bSHong Zhang       ctx->physics2.inflow(ctx,time,ctx->xmin,ctx->ub);
87984ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
88084ff8c8bSHong Zhang         if (ctx->physics2.bcinflowindex[j]==PETSC_TRUE) {
88184ff8c8bSHong Zhang           for (i=-2; i<0; i++) x[i*dof+j] = 2.0*ctx->ub[j]-x[-(i+1)*dof+j];
88284ff8c8bSHong Zhang         } else {
88384ff8c8bSHong Zhang           for (i=-2; i<0; i++) x[i*dof+j] = x[j]; /* Outflow */
88484ff8c8bSHong Zhang         }
88584ff8c8bSHong Zhang       }
88684ff8c8bSHong Zhang     }
88784ff8c8bSHong Zhang     if (xs+xm==Mx) { /* Right Boundary */
88884ff8c8bSHong Zhang       ctx->physics2.inflow(ctx,time,ctx->xmax,ctx->ub);
88984ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
89084ff8c8bSHong Zhang         if (ctx->physics2.bcinflowindex[dof+j]==PETSC_TRUE) {
89184ff8c8bSHong 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];
89284ff8c8bSHong Zhang         } else {
89384ff8c8bSHong Zhang           for (i=Mx; i<Mx+2; i++) x[i*dof+j] = x[(Mx-1)*dof+j]; /* Outflow */
89484ff8c8bSHong Zhang         }
89584ff8c8bSHong Zhang       }
89684ff8c8bSHong Zhang     }
89784ff8c8bSHong Zhang   }
89884ff8c8bSHong Zhang   for (i=xs-1; i<xs+xm+1; i++) {
89984ff8c8bSHong Zhang     struct _LimitInfo info;
90084ff8c8bSHong Zhang     PetscScalar       *cjmpL,*cjmpR;
90184ff8c8bSHong Zhang     if ((i > sf-lsbwidth-2 && i < sf+1) || (i > fs-2 && i < fs+rsbwidth+1)) {
90284ff8c8bSHong Zhang       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
90384ff8c8bSHong Zhang       ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr);
90484ff8c8bSHong Zhang       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
90584ff8c8bSHong Zhang       ierr  = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr);
90684ff8c8bSHong Zhang       cjmpL = &ctx->cjmpLR[0];
90784ff8c8bSHong Zhang       cjmpR = &ctx->cjmpLR[dof];
90884ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
90984ff8c8bSHong Zhang         PetscScalar jmpL,jmpR;
91084ff8c8bSHong Zhang         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
91184ff8c8bSHong Zhang         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
91284ff8c8bSHong Zhang         for (k=0; k<dof; k++) {
91384ff8c8bSHong Zhang           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
91484ff8c8bSHong Zhang           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
91584ff8c8bSHong Zhang         }
91684ff8c8bSHong Zhang       }
91784ff8c8bSHong Zhang       /* Apply limiter to the left and right characteristic jumps */
91884ff8c8bSHong Zhang       info.m  = dof;
91984ff8c8bSHong Zhang       info.hxs = hxs;
92084ff8c8bSHong Zhang       info.hxf = hxf;
92184ff8c8bSHong Zhang       (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope);
92284ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
92384ff8c8bSHong Zhang         PetscScalar tmp = 0;
92484ff8c8bSHong Zhang         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
92584ff8c8bSHong Zhang           slope[i*dof+j] = tmp;
92684ff8c8bSHong Zhang       }
92784ff8c8bSHong Zhang     }
92884ff8c8bSHong Zhang   }
92984ff8c8bSHong Zhang 
93084ff8c8bSHong Zhang   for (i=xs; i<xs+xm+1; i++) {
93184ff8c8bSHong Zhang     PetscReal   maxspeed;
93284ff8c8bSHong Zhang     PetscScalar *uL,*uR;
93384ff8c8bSHong Zhang     uL = &ctx->uLR[0];
93484ff8c8bSHong Zhang     uR = &ctx->uLR[dof];
93584ff8c8bSHong Zhang     if (i == sf-lsbwidth) {
93684ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
93784ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
93884ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
93984ff8c8bSHong Zhang       }
94084ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
94184ff8c8bSHong Zhang       if (i < xs+xm) {
94284ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
94384ff8c8bSHong Zhang         islow++;
94484ff8c8bSHong Zhang       }
94584ff8c8bSHong Zhang     }
94684ff8c8bSHong Zhang     if (i > sf-lsbwidth && i < sf) {
94784ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
94884ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
94984ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
95084ff8c8bSHong Zhang       }
95184ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
95284ff8c8bSHong Zhang       if (i > xs) {
95384ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
95484ff8c8bSHong Zhang       }
95584ff8c8bSHong Zhang       if (i < xs+xm) {
95684ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
95784ff8c8bSHong Zhang         islow++;
95884ff8c8bSHong Zhang       }
95984ff8c8bSHong Zhang     }
96084ff8c8bSHong Zhang     if (i == sf) { /* interface between the slow region and the fast region */
96184ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
96284ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
96384ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
96484ff8c8bSHong Zhang       }
96584ff8c8bSHong Zhang       ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
96684ff8c8bSHong Zhang       if (i > xs) {
96784ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
96884ff8c8bSHong Zhang       }
96984ff8c8bSHong Zhang     }
97084ff8c8bSHong Zhang     if (i == fs) { /* interface between the fast region and the slow region */
97184ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
97284ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
97384ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
97484ff8c8bSHong Zhang       }
97584ff8c8bSHong Zhang       ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
97684ff8c8bSHong Zhang       if (i < xs+xm) {
97784ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
97884ff8c8bSHong Zhang         islow++;
97984ff8c8bSHong Zhang       }
98084ff8c8bSHong Zhang     }
98184ff8c8bSHong Zhang     if (i > fs && i < fs+rsbwidth) {
98284ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
98384ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
98484ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
98584ff8c8bSHong Zhang       }
98684ff8c8bSHong Zhang       ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
98784ff8c8bSHong Zhang       if (i > xs) {
98884ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
98984ff8c8bSHong Zhang       }
99084ff8c8bSHong Zhang       if (i < xs+xm) {
99184ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
99284ff8c8bSHong Zhang         islow++;
99384ff8c8bSHong Zhang       }
99484ff8c8bSHong Zhang     }
99584ff8c8bSHong Zhang     if (i == fs+rsbwidth) {
99684ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
99784ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
99884ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
99984ff8c8bSHong Zhang       }
100084ff8c8bSHong Zhang       ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
100184ff8c8bSHong Zhang       if (i > xs) {
100284ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
100384ff8c8bSHong Zhang       }
100484ff8c8bSHong Zhang     }
100584ff8c8bSHong Zhang   }
100684ff8c8bSHong Zhang   ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr);
100784ff8c8bSHong Zhang   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
100884ff8c8bSHong Zhang   ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
100984ff8c8bSHong Zhang   ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr);
101084ff8c8bSHong Zhang   PetscFunctionReturn(0);
101184ff8c8bSHong Zhang }
101284ff8c8bSHong Zhang 
101384ff8c8bSHong Zhang /* --------------------------------- Finite Volume Solver for fast  parts ----------------------------------- */
101484ff8c8bSHong Zhang PetscErrorCode FVRHSFunctionfast_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
101584ff8c8bSHong Zhang {
101684ff8c8bSHong Zhang   FVCtx          *ctx = (FVCtx*)vctx;
101784ff8c8bSHong Zhang   PetscErrorCode ierr;
101884ff8c8bSHong Zhang   PetscInt       i,j,k,Mx,dof,xs,xm,ifast = 0,sf = ctx->sf,fs = ctx->fs;
101984ff8c8bSHong Zhang   PetscReal      hxs,hxf;
102084ff8c8bSHong Zhang   PetscScalar    *x,*f,*slope;
102184ff8c8bSHong Zhang   Vec            Xloc;
102284ff8c8bSHong Zhang   DM             da;
102384ff8c8bSHong Zhang 
102484ff8c8bSHong Zhang   PetscFunctionBeginUser;
102584ff8c8bSHong Zhang   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
102684ff8c8bSHong Zhang   ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr);
102784ff8c8bSHong Zhang   ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr);
102884ff8c8bSHong Zhang   hxs  = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
102984ff8c8bSHong Zhang   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
103084ff8c8bSHong Zhang   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
103184ff8c8bSHong Zhang   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
103284ff8c8bSHong Zhang   ierr = VecZeroEntries(F);CHKERRQ(ierr);
103384ff8c8bSHong Zhang   ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr);
103484ff8c8bSHong Zhang   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
103584ff8c8bSHong Zhang   ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
103684ff8c8bSHong Zhang   ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr);
103784ff8c8bSHong Zhang 
103884ff8c8bSHong Zhang   if (ctx->bctype == FVBC_OUTFLOW) {
103984ff8c8bSHong Zhang     for (i=xs-2; i<0; i++) {
104084ff8c8bSHong Zhang       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
104184ff8c8bSHong Zhang     }
104284ff8c8bSHong Zhang     for (i=Mx; i<xs+xm+2; i++) {
104384ff8c8bSHong Zhang       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
104484ff8c8bSHong Zhang     }
104584ff8c8bSHong Zhang   }
104684ff8c8bSHong Zhang   if (ctx->bctype == FVBC_INFLOW) {
104784ff8c8bSHong Zhang     /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253
104884ff8c8bSHong Zhang     pages 137-138 for the scheme. */
104984ff8c8bSHong Zhang     if (xs==0) { /* Left Boundary */
105084ff8c8bSHong Zhang       ctx->physics2.inflow(ctx,time,ctx->xmin,ctx->ub);
105184ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
105284ff8c8bSHong Zhang         if (ctx->physics2.bcinflowindex[j]==PETSC_TRUE) {
105384ff8c8bSHong Zhang           for (i=-2; i<0; i++) x[i*dof+j] = 2.0*ctx->ub[j]-x[-(i+1)*dof+j];
105484ff8c8bSHong Zhang         } else {
105584ff8c8bSHong Zhang           for (i=-2; i<0; i++) x[i*dof+j] = x[j]; /* Outflow */
105684ff8c8bSHong Zhang         }
105784ff8c8bSHong Zhang       }
105884ff8c8bSHong Zhang     }
105984ff8c8bSHong Zhang     if (xs+xm==Mx) { /* Right Boundary */
106084ff8c8bSHong Zhang       ctx->physics2.inflow(ctx,time,ctx->xmax,ctx->ub);
106184ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
106284ff8c8bSHong Zhang         if (ctx->physics2.bcinflowindex[dof+j]==PETSC_TRUE) {
106384ff8c8bSHong 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];
106484ff8c8bSHong Zhang         } else {
106584ff8c8bSHong Zhang           for (i=Mx; i<Mx+2; i++) x[i*dof+j] = x[(Mx-1)*dof+j]; /* Outflow */
106684ff8c8bSHong Zhang         }
106784ff8c8bSHong Zhang       }
106884ff8c8bSHong Zhang     }
106984ff8c8bSHong Zhang   }
107084ff8c8bSHong Zhang   for (i=xs-1; i<xs+xm+1; i++) {
107184ff8c8bSHong Zhang     struct _LimitInfo info;
107284ff8c8bSHong Zhang     PetscScalar       *cjmpL,*cjmpR;
107384ff8c8bSHong Zhang     if (i > sf-2 && i < fs+1) {
107484ff8c8bSHong Zhang       ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr);
107584ff8c8bSHong Zhang       ierr  = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr);
107684ff8c8bSHong Zhang       cjmpL = &ctx->cjmpLR[0];
107784ff8c8bSHong Zhang       cjmpR = &ctx->cjmpLR[dof];
107884ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
107984ff8c8bSHong Zhang         PetscScalar jmpL,jmpR;
108084ff8c8bSHong Zhang         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
108184ff8c8bSHong Zhang         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
108284ff8c8bSHong Zhang         for (k=0; k<dof; k++) {
108384ff8c8bSHong Zhang           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
108484ff8c8bSHong Zhang           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
108584ff8c8bSHong Zhang         }
108684ff8c8bSHong Zhang       }
108784ff8c8bSHong Zhang       /* Apply limiter to the left and right characteristic jumps */
108884ff8c8bSHong Zhang       info.m  = dof;
108984ff8c8bSHong Zhang       info.hxs = hxs;
109084ff8c8bSHong Zhang       info.hxf = hxf;
109184ff8c8bSHong Zhang       (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope);
109284ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
109384ff8c8bSHong Zhang       PetscScalar tmp = 0;
109484ff8c8bSHong Zhang       for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
109584ff8c8bSHong Zhang         slope[i*dof+j] = tmp;
109684ff8c8bSHong Zhang       }
109784ff8c8bSHong Zhang     }
109884ff8c8bSHong Zhang   }
109984ff8c8bSHong Zhang 
110084ff8c8bSHong Zhang   for (i=xs; i<xs+xm+1; i++) {
110184ff8c8bSHong Zhang     PetscReal   maxspeed;
110284ff8c8bSHong Zhang     PetscScalar *uL,*uR;
110384ff8c8bSHong Zhang     uL = &ctx->uLR[0];
110484ff8c8bSHong Zhang     uR = &ctx->uLR[dof];
110584ff8c8bSHong Zhang     if (i == sf) { /* interface between the slow region and the fast region */
110684ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
110784ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
110884ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
110984ff8c8bSHong Zhang       }
111084ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
111184ff8c8bSHong Zhang       if (i < xs+xm) {
111284ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf;
111384ff8c8bSHong Zhang         ifast++;
111484ff8c8bSHong Zhang       }
111584ff8c8bSHong Zhang     }
111684ff8c8bSHong Zhang     if (i > sf && i < fs) { /* fast region */
111784ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
111884ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
111984ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
112084ff8c8bSHong Zhang       }
112184ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
112284ff8c8bSHong Zhang       if (i > xs) {
112384ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf;
112484ff8c8bSHong Zhang       }
112584ff8c8bSHong Zhang       if (i < xs+xm) {
112684ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf;
112784ff8c8bSHong Zhang         ifast++;
112884ff8c8bSHong Zhang       }
112984ff8c8bSHong Zhang     }
113084ff8c8bSHong Zhang     if (i == fs) { /* interface between the fast region and the slow region */
113184ff8c8bSHong Zhang       for (j=0; j<dof; j++) {
113284ff8c8bSHong Zhang         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
113384ff8c8bSHong Zhang         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
113484ff8c8bSHong Zhang       }
113584ff8c8bSHong Zhang       ierr    = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr);
113684ff8c8bSHong Zhang       if (i > xs) {
113784ff8c8bSHong Zhang         for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf;
113884ff8c8bSHong Zhang       }
113984ff8c8bSHong Zhang     }
114084ff8c8bSHong Zhang   }
114184ff8c8bSHong Zhang   ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr);
114284ff8c8bSHong Zhang   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
114384ff8c8bSHong Zhang   ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
114484ff8c8bSHong Zhang   ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr);
114584ff8c8bSHong Zhang   PetscFunctionReturn(0);
114684ff8c8bSHong Zhang }
114784ff8c8bSHong Zhang 
114884ff8c8bSHong Zhang int main(int argc,char *argv[])
114984ff8c8bSHong Zhang {
115084ff8c8bSHong Zhang   char              lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m";
115184ff8c8bSHong Zhang   PetscFunctionList limiters   = 0,physics = 0;
115284ff8c8bSHong Zhang   MPI_Comm          comm;
115384ff8c8bSHong Zhang   TS                ts;
115484ff8c8bSHong Zhang   DM                da;
115584ff8c8bSHong Zhang   Vec               X,X0,R;
115684ff8c8bSHong Zhang   FVCtx             ctx;
115784ff8c8bSHong 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;
115884ff8c8bSHong Zhang   PetscBool         view_final = PETSC_FALSE;
115984ff8c8bSHong Zhang   PetscReal         ptime,maxtime;
116084ff8c8bSHong Zhang   PetscErrorCode    ierr;
116184ff8c8bSHong Zhang 
116284ff8c8bSHong Zhang   ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
116384ff8c8bSHong Zhang   comm = PETSC_COMM_WORLD;
116484ff8c8bSHong Zhang   ierr = PetscMemzero(&ctx,sizeof(ctx));CHKERRQ(ierr);
116584ff8c8bSHong Zhang 
116684ff8c8bSHong Zhang   /* Register limiters to be available on the command line */
116784ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&limiters,"upwind"              ,Limit2_Upwind);CHKERRQ(ierr);
116884ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&limiters,"lax-wendroff"        ,Limit2_LaxWendroff);CHKERRQ(ierr);
116984ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&limiters,"beam-warming"        ,Limit2_BeamWarming);CHKERRQ(ierr);
117084ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&limiters,"fromm"               ,Limit2_Fromm);CHKERRQ(ierr);
117184ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&limiters,"minmod"              ,Limit2_Minmod);CHKERRQ(ierr);
117284ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&limiters,"superbee"            ,Limit2_Superbee);CHKERRQ(ierr);
117384ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&limiters,"mc"                  ,Limit2_MC);CHKERRQ(ierr);
117484ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&limiters,"koren3"              ,Limit2_Koren3);CHKERRQ(ierr);
117584ff8c8bSHong Zhang 
117684ff8c8bSHong Zhang   /* Register physical models to be available on the command line */
117784ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&physics,"shallow"         ,PhysicsCreate_Shallow);CHKERRQ(ierr);
117884ff8c8bSHong Zhang   ierr = PetscFunctionListAdd(&physics,"advect"          ,PhysicsCreate_Advect);CHKERRQ(ierr);
117984ff8c8bSHong Zhang 
118084ff8c8bSHong Zhang   ctx.comm    = comm;
118184ff8c8bSHong Zhang   ctx.cfl     = 0.9;
118284ff8c8bSHong Zhang   ctx.bctype  = FVBC_PERIODIC;
118384ff8c8bSHong Zhang   ctx.xmin    = -1.0;
118484ff8c8bSHong Zhang   ctx.xmax    = 1.0;
118584ff8c8bSHong Zhang   ctx.initial = 1;
118684ff8c8bSHong Zhang   ctx.hratio  = 2;
118784ff8c8bSHong Zhang   maxtime     = 10.0;
118884ff8c8bSHong Zhang   ctx.simulation = PETSC_FALSE;
118984ff8c8bSHong Zhang   ierr = PetscOptionsBegin(comm,NULL,"Finite Volume solver options","");CHKERRQ(ierr);
119084ff8c8bSHong Zhang   ierr = PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL);CHKERRQ(ierr);
119184ff8c8bSHong Zhang   ierr = PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL);CHKERRQ(ierr);
119284ff8c8bSHong Zhang   ierr = PetscOptionsFList("-limit","Name of flux imiter to use","",limiters,lname,lname,sizeof(lname),NULL);CHKERRQ(ierr);
119384ff8c8bSHong Zhang   ierr = PetscOptionsFList("-physics","Name of physics model to use","",physics,physname,physname,sizeof(physname),NULL);CHKERRQ(ierr);
119484ff8c8bSHong Zhang   ierr = PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL);CHKERRQ(ierr);
119584ff8c8bSHong 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);
119684ff8c8bSHong Zhang   ierr = PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL);CHKERRQ(ierr);
119784ff8c8bSHong Zhang   ierr = PetscOptionsBool("-exact","Compare errors with exact solution","",ctx.exact,&ctx.exact,NULL);CHKERRQ(ierr);
119884ff8c8bSHong Zhang   ierr = PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL);CHKERRQ(ierr);
119984ff8c8bSHong Zhang   ierr = PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL);CHKERRQ(ierr);
120084ff8c8bSHong Zhang   ierr = PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL);CHKERRQ(ierr);
120184ff8c8bSHong Zhang   ierr = PetscOptionsInt("-hratio","Spacing ratio","",ctx.hratio,&ctx.hratio,NULL);CHKERRQ(ierr);
120284ff8c8bSHong Zhang   ierr = PetscOptionsEnd();CHKERRQ(ierr);
120384ff8c8bSHong Zhang 
120484ff8c8bSHong Zhang   /* Choose the limiter from the list of registered limiters */
120584ff8c8bSHong Zhang   ierr = PetscFunctionListFind(limiters,lname,&ctx.limit2);CHKERRQ(ierr);
1206*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!ctx.limit2,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Limiter '%s' not found",lname);
120784ff8c8bSHong Zhang 
120884ff8c8bSHong Zhang   /* Choose the physics from the list of registered models */
120984ff8c8bSHong Zhang   {
121084ff8c8bSHong Zhang     PetscErrorCode (*r)(FVCtx*);
121184ff8c8bSHong Zhang     ierr = PetscFunctionListFind(physics,physname,&r);CHKERRQ(ierr);
1212*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Physics '%s' not found",physname);
121384ff8c8bSHong Zhang     /* Create the physics, will set the number of fields and their names */
121484ff8c8bSHong Zhang     ierr = (*r)(&ctx);CHKERRQ(ierr);
121584ff8c8bSHong Zhang   }
121684ff8c8bSHong Zhang 
121784ff8c8bSHong Zhang   /* Create a DMDA to manage the parallel grid */
121884ff8c8bSHong Zhang   ierr = DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics2.dof,2,NULL,&da);CHKERRQ(ierr);
121984ff8c8bSHong Zhang   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
122084ff8c8bSHong Zhang   ierr = DMSetUp(da);CHKERRQ(ierr);
122184ff8c8bSHong Zhang   /* Inform the DMDA of the field names provided by the physics. */
122284ff8c8bSHong Zhang   /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
122384ff8c8bSHong Zhang   for (i=0; i<ctx.physics2.dof; i++) {
122484ff8c8bSHong Zhang     ierr = DMDASetFieldName(da,i,ctx.physics2.fieldname[i]);CHKERRQ(ierr);
122584ff8c8bSHong Zhang   }
122684ff8c8bSHong Zhang   ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr);
122784ff8c8bSHong Zhang   ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr);
122884ff8c8bSHong Zhang 
122984ff8c8bSHong Zhang   /* Set coordinates of cell centers */
123084ff8c8bSHong 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);
123184ff8c8bSHong Zhang 
123284ff8c8bSHong Zhang   /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
123384ff8c8bSHong Zhang   ierr = PetscMalloc4(dof*dof,&ctx.R,dof*dof,&ctx.Rinv,2*dof,&ctx.cjmpLR,1*dof,&ctx.cslope);CHKERRQ(ierr);
123484ff8c8bSHong Zhang   ierr = PetscMalloc3(2*dof,&ctx.uLR,dof,&ctx.flux,dof,&ctx.speeds);CHKERRQ(ierr);
123584ff8c8bSHong Zhang   ierr = PetscMalloc1(2*dof,&ctx.ub);CHKERRQ(ierr);
123684ff8c8bSHong Zhang 
123784ff8c8bSHong Zhang   /* Create a vector to store the solution and to save the initial state */
123884ff8c8bSHong Zhang   ierr = DMCreateGlobalVector(da,&X);CHKERRQ(ierr);
123984ff8c8bSHong Zhang   ierr = VecDuplicate(X,&X0);CHKERRQ(ierr);
124084ff8c8bSHong Zhang   ierr = VecDuplicate(X,&R);CHKERRQ(ierr);
124184ff8c8bSHong Zhang 
124284ff8c8bSHong Zhang   /* create index for slow parts and fast parts,
124384ff8c8bSHong Zhang      count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */
124484ff8c8bSHong Zhang   count_slow = Mx/(1.0+ctx.hratio/3.0);
1245*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(count_slow%2,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");
124684ff8c8bSHong Zhang   count_fast = Mx-count_slow;
124784ff8c8bSHong Zhang   ctx.sf = count_slow/2;
124884ff8c8bSHong Zhang   ctx.fs = ctx.sf+count_fast;
124984ff8c8bSHong Zhang   ierr = PetscMalloc1(xm*dof,&index_slow);CHKERRQ(ierr);
125084ff8c8bSHong Zhang   ierr = PetscMalloc1(xm*dof,&index_fast);CHKERRQ(ierr);
125184ff8c8bSHong Zhang   ierr = PetscMalloc1(8*dof,&index_slowbuffer);CHKERRQ(ierr);
125284ff8c8bSHong Zhang   ctx.lsbwidth = 4;
125384ff8c8bSHong Zhang   ctx.rsbwidth = 4;
125484ff8c8bSHong Zhang 
125584ff8c8bSHong Zhang   for (i=xs; i<xs+xm; i++) {
125684ff8c8bSHong Zhang     if (i < ctx.sf-ctx.lsbwidth || i > ctx.fs+ctx.rsbwidth-1)
125784ff8c8bSHong Zhang       for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k;
125884ff8c8bSHong Zhang     else if ((i >= ctx.sf-ctx.lsbwidth && i < ctx.sf) || (i > ctx.fs-1 && i <= ctx.fs+ctx.rsbwidth-1))
125984ff8c8bSHong Zhang       for (k=0; k<dof; k++) index_slowbuffer[islowbuffer++] = i*dof+k;
126084ff8c8bSHong Zhang     else
126184ff8c8bSHong Zhang       for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k;
126284ff8c8bSHong Zhang   }
126384ff8c8bSHong Zhang   ierr = ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss);CHKERRQ(ierr);
126484ff8c8bSHong Zhang   ierr = ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf);CHKERRQ(ierr);
126584ff8c8bSHong Zhang   ierr = ISCreateGeneral(PETSC_COMM_WORLD,islowbuffer,index_slowbuffer,PETSC_COPY_VALUES,&ctx.issb);CHKERRQ(ierr);
126684ff8c8bSHong Zhang 
126784ff8c8bSHong Zhang   /* Create a time-stepping object */
126884ff8c8bSHong Zhang   ierr = TSCreate(comm,&ts);CHKERRQ(ierr);
126984ff8c8bSHong Zhang   ierr = TSSetDM(ts,da);CHKERRQ(ierr);
127084ff8c8bSHong Zhang   ierr = TSSetRHSFunction(ts,R,FVRHSFunction_2WaySplit,&ctx);CHKERRQ(ierr);
127184ff8c8bSHong Zhang   ierr = TSRHSSplitSetIS(ts,"slow",ctx.iss);CHKERRQ(ierr);
127284ff8c8bSHong Zhang   ierr = TSRHSSplitSetIS(ts,"slowbuffer",ctx.issb);CHKERRQ(ierr);
127384ff8c8bSHong Zhang   ierr = TSRHSSplitSetIS(ts,"fast",ctx.isf);CHKERRQ(ierr);
127484ff8c8bSHong Zhang   ierr = TSRHSSplitSetRHSFunction(ts,"slow",NULL,FVRHSFunctionslow_2WaySplit,&ctx);CHKERRQ(ierr);
127584ff8c8bSHong Zhang   ierr = TSRHSSplitSetRHSFunction(ts,"fast",NULL,FVRHSFunctionfast_2WaySplit,&ctx);CHKERRQ(ierr);
127684ff8c8bSHong Zhang   ierr = TSRHSSplitSetRHSFunction(ts,"slowbuffer",NULL,FVRHSFunctionslowbuffer_2WaySplit,&ctx);CHKERRQ(ierr);
127784ff8c8bSHong Zhang 
127884ff8c8bSHong Zhang   ierr = TSSetType(ts,TSMPRK);CHKERRQ(ierr);
127984ff8c8bSHong Zhang   ierr = TSSetMaxTime(ts,maxtime);CHKERRQ(ierr);
128084ff8c8bSHong Zhang   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
128184ff8c8bSHong Zhang 
128284ff8c8bSHong Zhang   /* Compute initial conditions and starting time step */
128384ff8c8bSHong Zhang   ierr = FVSample_2WaySplit(&ctx,da,0,X0);CHKERRQ(ierr);
128484ff8c8bSHong Zhang   ierr = FVRHSFunction_2WaySplit(ts,0,X0,X,(void*)&ctx);CHKERRQ(ierr); /* Initial function evaluation, only used to determine max speed */
128584ff8c8bSHong Zhang   ierr = VecCopy(X0,X);CHKERRQ(ierr);                        /* The function value was not used so we set X=X0 again */
128684ff8c8bSHong Zhang   ierr = TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt);CHKERRQ(ierr);
128784ff8c8bSHong Zhang   ierr = TSSetFromOptions(ts);CHKERRQ(ierr); /* Take runtime options */
128884ff8c8bSHong Zhang   ierr = SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
128984ff8c8bSHong Zhang   {
129084ff8c8bSHong Zhang     PetscInt          steps;
129184ff8c8bSHong Zhang     PetscScalar       mass_initial,mass_final,mass_difference,mass_differenceg;
129284ff8c8bSHong Zhang     const PetscScalar *ptr_X,*ptr_X0;
129384ff8c8bSHong Zhang     const PetscReal   hs = (ctx.xmax-ctx.xmin)*3.0/4.0/count_slow;
129484ff8c8bSHong Zhang     const PetscReal   hf = (ctx.xmax-ctx.xmin)/4.0/count_fast;
129584ff8c8bSHong Zhang 
129684ff8c8bSHong Zhang     ierr = TSSolve(ts,X);CHKERRQ(ierr);
129784ff8c8bSHong Zhang     ierr = TSGetSolveTime(ts,&ptime);CHKERRQ(ierr);
129884ff8c8bSHong Zhang     ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
129984ff8c8bSHong Zhang     /* calculate the total mass at initial time and final time */
130084ff8c8bSHong Zhang     mass_initial = 0.0;
130184ff8c8bSHong Zhang     mass_final   = 0.0;
130284ff8c8bSHong Zhang     ierr = DMDAVecGetArrayRead(da,X0,(void*)&ptr_X0);CHKERRQ(ierr);
130384ff8c8bSHong Zhang     ierr = DMDAVecGetArrayRead(da,X,(void*)&ptr_X);CHKERRQ(ierr);
130484ff8c8bSHong Zhang     for (i=xs;i<xs+xm;i++) {
130584ff8c8bSHong Zhang       if (i < ctx.sf || i > ctx.fs-1) {
130684ff8c8bSHong Zhang         for (k=0; k<dof; k++) {
130784ff8c8bSHong Zhang           mass_initial = mass_initial + hs*ptr_X0[i*dof+k];
130884ff8c8bSHong Zhang           mass_final = mass_final + hs*ptr_X[i*dof+k];
130984ff8c8bSHong Zhang         }
131084ff8c8bSHong Zhang       } else {
131184ff8c8bSHong Zhang         for (k=0; k<dof; k++) {
131284ff8c8bSHong Zhang           mass_initial = mass_initial + hf*ptr_X0[i*dof+k];
131384ff8c8bSHong Zhang           mass_final = mass_final + hf*ptr_X[i*dof+k];
131484ff8c8bSHong Zhang         }
131584ff8c8bSHong Zhang       }
131684ff8c8bSHong Zhang     }
131784ff8c8bSHong Zhang     ierr = DMDAVecRestoreArrayRead(da,X0,(void*)&ptr_X0);CHKERRQ(ierr);
131884ff8c8bSHong Zhang     ierr = DMDAVecRestoreArrayRead(da,X,(void*)&ptr_X);CHKERRQ(ierr);
131984ff8c8bSHong Zhang     mass_difference = mass_final - mass_initial;
132084ff8c8bSHong Zhang     ierr = MPI_Allreduce(&mass_difference,&mass_differenceg,1,MPIU_SCALAR,MPIU_SUM,comm);CHKERRMPI(ierr);
132184ff8c8bSHong Zhang     ierr = PetscPrintf(comm,"Mass difference %g\n",(double)mass_differenceg);CHKERRQ(ierr);
132284ff8c8bSHong Zhang     ierr = PetscPrintf(comm,"Final time %g, steps %D\n",(double)ptime,steps);CHKERRQ(ierr);
132384ff8c8bSHong Zhang     ierr = PetscPrintf(comm,"Maximum allowable stepsize according to CFL %g\n",(double)(1.0/ctx.cfl_idt));CHKERRQ(ierr);
132484ff8c8bSHong Zhang     if (ctx.exact) {
132584ff8c8bSHong Zhang       PetscReal nrm1=0;
132684ff8c8bSHong Zhang       ierr = SolutionErrorNorms_2WaySplit(&ctx,da,ptime,X,&nrm1);CHKERRQ(ierr);
132784ff8c8bSHong Zhang       ierr = PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1);CHKERRQ(ierr);
132884ff8c8bSHong Zhang     }
132984ff8c8bSHong Zhang     if (ctx.simulation) {
133084ff8c8bSHong Zhang       PetscReal    nrm1=0;
133184ff8c8bSHong Zhang       PetscViewer  fd;
133284ff8c8bSHong Zhang       char         filename[PETSC_MAX_PATH_LEN] = "binaryoutput";
133384ff8c8bSHong Zhang       Vec          XR;
133484ff8c8bSHong Zhang       PetscBool    flg;
133584ff8c8bSHong Zhang       const PetscScalar  *ptr_XR;
133684ff8c8bSHong Zhang       ierr = PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg);CHKERRQ(ierr);
1337*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
133884ff8c8bSHong Zhang       ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&fd);CHKERRQ(ierr);
133984ff8c8bSHong Zhang       ierr = VecDuplicate(X0,&XR);CHKERRQ(ierr);
134084ff8c8bSHong Zhang       ierr = VecLoad(XR,fd);CHKERRQ(ierr);
134184ff8c8bSHong Zhang       ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
134284ff8c8bSHong Zhang       ierr = VecGetArrayRead(X,&ptr_X);CHKERRQ(ierr);
134384ff8c8bSHong Zhang       ierr = VecGetArrayRead(XR,&ptr_XR);CHKERRQ(ierr);
134484ff8c8bSHong Zhang       for (i=xs;i<xs+xm;i++) {
134584ff8c8bSHong Zhang         if (i < ctx.sf || i > ctx.fs-1)
134684ff8c8bSHong Zhang           for (k=0; k<dof; k++) nrm1 = nrm1 + hs*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]);
134784ff8c8bSHong Zhang         else
134884ff8c8bSHong Zhang           for (k=0; k<dof; k++) nrm1 = nrm1 + hf*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]);
134984ff8c8bSHong Zhang       }
135084ff8c8bSHong Zhang       ierr = VecRestoreArrayRead(X,&ptr_X);CHKERRQ(ierr);
135184ff8c8bSHong Zhang       ierr = VecRestoreArrayRead(XR,&ptr_XR);CHKERRQ(ierr);
135284ff8c8bSHong Zhang       ierr = PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1);CHKERRQ(ierr);
135384ff8c8bSHong Zhang       ierr = VecDestroy(&XR);CHKERRQ(ierr);
135484ff8c8bSHong Zhang     }
135584ff8c8bSHong Zhang   }
135684ff8c8bSHong Zhang 
135784ff8c8bSHong Zhang   ierr = SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
135884ff8c8bSHong Zhang   if (draw & 0x1) {ierr = VecView(X0,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);}
135984ff8c8bSHong Zhang   if (draw & 0x2) {ierr = VecView(X,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);}
136084ff8c8bSHong Zhang   if (draw & 0x4) {
136184ff8c8bSHong Zhang     Vec Y;
136284ff8c8bSHong Zhang     ierr = VecDuplicate(X,&Y);CHKERRQ(ierr);
136384ff8c8bSHong Zhang     ierr = FVSample_2WaySplit(&ctx,da,ptime,Y);CHKERRQ(ierr);
136484ff8c8bSHong Zhang     ierr = VecAYPX(Y,-1,X);CHKERRQ(ierr);
136584ff8c8bSHong Zhang     ierr = VecView(Y,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);
136684ff8c8bSHong Zhang     ierr = VecDestroy(&Y);CHKERRQ(ierr);
136784ff8c8bSHong Zhang   }
136884ff8c8bSHong Zhang 
136984ff8c8bSHong Zhang   if (view_final) {
137084ff8c8bSHong Zhang     PetscViewer viewer;
137184ff8c8bSHong Zhang     ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer);CHKERRQ(ierr);
137284ff8c8bSHong Zhang     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
137384ff8c8bSHong Zhang     ierr = VecView(X,viewer);CHKERRQ(ierr);
137484ff8c8bSHong Zhang     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
137584ff8c8bSHong Zhang     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
137684ff8c8bSHong Zhang   }
137784ff8c8bSHong Zhang 
137884ff8c8bSHong Zhang   /* Clean up */
137984ff8c8bSHong Zhang   ierr = (*ctx.physics2.destroy)(ctx.physics2.user);CHKERRQ(ierr);
138084ff8c8bSHong Zhang   for (i=0; i<ctx.physics2.dof; i++) {ierr = PetscFree(ctx.physics2.fieldname[i]);CHKERRQ(ierr);}
138184ff8c8bSHong Zhang   ierr = PetscFree(ctx.physics2.bcinflowindex);CHKERRQ(ierr);
138284ff8c8bSHong Zhang   ierr = PetscFree(ctx.ub);
138384ff8c8bSHong Zhang   ierr = PetscFree4(ctx.R,ctx.Rinv,ctx.cjmpLR,ctx.cslope);CHKERRQ(ierr);
138484ff8c8bSHong Zhang   ierr = PetscFree3(ctx.uLR,ctx.flux,ctx.speeds);CHKERRQ(ierr);
138584ff8c8bSHong Zhang   ierr = VecDestroy(&X);CHKERRQ(ierr);
138684ff8c8bSHong Zhang   ierr = VecDestroy(&X0);CHKERRQ(ierr);
138784ff8c8bSHong Zhang   ierr = VecDestroy(&R);CHKERRQ(ierr);
138884ff8c8bSHong Zhang   ierr = DMDestroy(&da);CHKERRQ(ierr);
138984ff8c8bSHong Zhang   ierr = TSDestroy(&ts);CHKERRQ(ierr);
139084ff8c8bSHong Zhang   ierr = ISDestroy(&ctx.iss);CHKERRQ(ierr);
139184ff8c8bSHong Zhang   ierr = ISDestroy(&ctx.isf);CHKERRQ(ierr);
139284ff8c8bSHong Zhang   ierr = ISDestroy(&ctx.issb);CHKERRQ(ierr);
139384ff8c8bSHong Zhang   ierr = PetscFree(index_slow);CHKERRQ(ierr);
139484ff8c8bSHong Zhang   ierr = PetscFree(index_fast);CHKERRQ(ierr);
139584ff8c8bSHong Zhang   ierr = PetscFree(index_slowbuffer);CHKERRQ(ierr);
139684ff8c8bSHong Zhang   ierr = PetscFunctionListDestroy(&limiters);CHKERRQ(ierr);
139784ff8c8bSHong Zhang   ierr = PetscFunctionListDestroy(&physics);CHKERRQ(ierr);
139884ff8c8bSHong Zhang   ierr = PetscFinalize();
139984ff8c8bSHong Zhang   return ierr;
140084ff8c8bSHong Zhang }
140184ff8c8bSHong Zhang 
140284ff8c8bSHong Zhang /*TEST
140384ff8c8bSHong Zhang 
140484ff8c8bSHong Zhang     build:
140584ff8c8bSHong Zhang       requires: !complex !single
140684ff8c8bSHong Zhang       depends: finitevolume1d.c
140784ff8c8bSHong Zhang 
140884ff8c8bSHong Zhang     test:
140984ff8c8bSHong Zhang       suffix: 1
141084ff8c8bSHong 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
141184ff8c8bSHong Zhang       output_file: output/ex4_1.out
141284ff8c8bSHong Zhang 
141384ff8c8bSHong Zhang     test:
141484ff8c8bSHong Zhang       suffix: 2
141584ff8c8bSHong 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
141684ff8c8bSHong Zhang       output_file: output/ex4_1.out
141784ff8c8bSHong Zhang 
141884ff8c8bSHong Zhang     test:
141984ff8c8bSHong Zhang       suffix: 3
142084ff8c8bSHong 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
142184ff8c8bSHong Zhang       output_file: output/ex4_3.out
142284ff8c8bSHong Zhang 
142384ff8c8bSHong Zhang     test:
142484ff8c8bSHong Zhang       suffix: 4
142584ff8c8bSHong Zhang       nsize: 2
142684ff8c8bSHong 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
142784ff8c8bSHong Zhang       output_file: output/ex4_3.out
142884ff8c8bSHong Zhang 
142984ff8c8bSHong Zhang     test:
143084ff8c8bSHong Zhang       suffix: 5
143184ff8c8bSHong Zhang       nsize: 4
143284ff8c8bSHong 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
143384ff8c8bSHong Zhang       output_file: output/ex4_3.out
143484ff8c8bSHong Zhang TEST*/
1435