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