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