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