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