xref: /petsc/src/ts/tutorials/multirate/ex8.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 
209fbee547SJacob Faibussowitsch 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;
815f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNew(&user));
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;
885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrallocpy("u",&ctx->physics2.fieldname[0]));
89c4762a1bSJed Brown   user->a = 1;
90c4762a1bSJed Brown   ierr = PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");CHKERRQ(ierr);
91c4762a1bSJed Brown   {
925f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-physics_advect_a","Speed","",user->a,&user->a,NULL));
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   PetscScalar     *u,*uj,xj,xi;
101c4762a1bSJed Brown   PetscInt        i,j,k,dof,xs,xm,Mx;
102c4762a1bSJed Brown   const PetscInt  N = 200;
103c4762a1bSJed Brown   PetscReal       hs,hm,hf;
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   PetscFunctionBeginUser;
1063c633725SBarry Smith   PetscCheck(ctx->physics2.sample2,PETSC_COMM_SELF,PETSC_ERR_SUP,"Physics has not provided a sampling function");
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,U,&u));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(dof,&uj));
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   hs  = (ctx->xmax-ctx->xmin)/8.0/ctx->sm;
113c4762a1bSJed Brown   hm  = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm);
114c4762a1bSJed Brown   hf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf);
115c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
116c4762a1bSJed Brown     if (i < ctx->sm) {
117c4762a1bSJed Brown       xi = ctx->xmin+0.5*hs+i*hs;
118c4762a1bSJed Brown       /* Integrate over cell i using trapezoid rule with N points. */
119c4762a1bSJed Brown       for (k=0; k<dof; k++) u[i*dof+k] = 0;
120c4762a1bSJed Brown       for (j=0; j<N+1; j++) {
121c4762a1bSJed Brown         xj = xi+hs*(j-N/2)/(PetscReal)N;
1225f80ce2aSJacob Faibussowitsch         CHKERRQ((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj));
123c4762a1bSJed Brown         for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N;
124c4762a1bSJed Brown       }
125c4762a1bSJed Brown     } else if (i < ctx->mf) {
126c4762a1bSJed Brown       xi = ctx->xmin+ctx->sm*hs+0.5*hm+(i-ctx->sm)*hm;
127c4762a1bSJed Brown       /* Integrate over cell i using trapezoid rule with N points. */
128c4762a1bSJed Brown       for (k=0; k<dof; k++) u[i*dof+k] = 0;
129c4762a1bSJed Brown       for (j=0; j<N+1; j++) {
130c4762a1bSJed Brown         xj = xi+hm*(j-N/2)/(PetscReal)N;
1315f80ce2aSJacob Faibussowitsch         CHKERRQ((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj));
132c4762a1bSJed Brown         for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N;
133c4762a1bSJed Brown       }
134c4762a1bSJed Brown     } else if (i < ctx->fm) {
135c4762a1bSJed Brown       xi = ctx->xmin+ctx->sm*hs+(ctx->mf-ctx->sm)*hm+0.5*hf+(i-ctx->mf)*hf;
136c4762a1bSJed Brown       /* Integrate over cell i using trapezoid rule with N points. */
137c4762a1bSJed Brown       for (k=0; k<dof; k++) u[i*dof+k] = 0;
138c4762a1bSJed Brown       for (j=0; j<N+1; j++) {
139c4762a1bSJed Brown         xj = xi+hf*(j-N/2)/(PetscReal)N;
1405f80ce2aSJacob Faibussowitsch         CHKERRQ((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj));
141c4762a1bSJed Brown         for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N;
142c4762a1bSJed Brown       }
143c4762a1bSJed Brown     } else if (i < ctx->ms) {
144c4762a1bSJed Brown       xi = ctx->xmin+ctx->sm*hs+(ctx->mf-ctx->sm)*hm+(ctx->fm-ctx->mf)*hf+0.5*hm+(i-ctx->fm)*hm;
145c4762a1bSJed Brown       /* Integrate over cell i using trapezoid rule with N points. */
146c4762a1bSJed Brown       for (k=0; k<dof; k++) u[i*dof+k] = 0;
147c4762a1bSJed Brown       for (j=0; j<N+1; j++) {
148c4762a1bSJed Brown         xj = xi+hm*(j-N/2)/(PetscReal)N;
1495f80ce2aSJacob Faibussowitsch         CHKERRQ((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj));
150c4762a1bSJed Brown         for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N;
151c4762a1bSJed Brown       }
152c4762a1bSJed Brown     } else {
153c4762a1bSJed 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;
154c4762a1bSJed Brown       /* Integrate over cell i using trapezoid rule with N points. */
155c4762a1bSJed Brown       for (k=0; k<dof; k++) u[i*dof+k] = 0;
156c4762a1bSJed Brown       for (j=0; j<N+1; j++) {
157c4762a1bSJed Brown         xj = xi+hs*(j-N/2)/(PetscReal)N;
1585f80ce2aSJacob Faibussowitsch         CHKERRQ((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj));
159c4762a1bSJed Brown         for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N;
160c4762a1bSJed Brown       }
161c4762a1bSJed Brown     }
162c4762a1bSJed Brown   }
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,U,&u));
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(uj));
165c4762a1bSJed Brown   PetscFunctionReturn(0);
166c4762a1bSJed Brown }
167c4762a1bSJed Brown 
168c4762a1bSJed Brown static PetscErrorCode SolutionErrorNorms_3WaySplit(FVCtx *ctx,DM da,PetscReal t,Vec X,PetscReal *nrm1)
169c4762a1bSJed Brown {
170c4762a1bSJed Brown   Vec               Y;
171c4762a1bSJed Brown   PetscInt          i,Mx;
172c4762a1bSJed Brown   const PetscScalar *ptr_X,*ptr_Y;
173c4762a1bSJed Brown   PetscReal         hs,hm,hf;
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   PetscFunctionBeginUser;
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(X,&Mx));
177c4762a1bSJed Brown   hs   = (ctx->xmax-ctx->xmin)/8.0/ctx->sm;
178c4762a1bSJed Brown   hm   = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm);
179c4762a1bSJed Brown   hf   = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf);
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(X,&Y));
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(FVSample_3WaySplit(ctx,da,t,Y));
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&ptr_X));
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Y,&ptr_Y));
184c4762a1bSJed Brown   for (i=0;i<Mx;i++) {
185c4762a1bSJed Brown     if (i < ctx->sm  || i > ctx->ms - 1) *nrm1 +=  hs*PetscAbs(ptr_X[i]-ptr_Y[i]);
186c4762a1bSJed Brown     else if (i < ctx->mf  || i > ctx->fm - 1) *nrm1 +=  hm*PetscAbs(ptr_X[i]-ptr_Y[i]);
187c4762a1bSJed Brown     else *nrm1 += hf*PetscAbs(ptr_X[i]-ptr_Y[i]);
188c4762a1bSJed Brown   }
1895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&ptr_X));
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Y,&ptr_Y));
1915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Y));
192c4762a1bSJed Brown   PetscFunctionReturn(0);
193c4762a1bSJed Brown }
194c4762a1bSJed Brown 
195c4762a1bSJed Brown PetscErrorCode FVRHSFunction_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
196c4762a1bSJed Brown {
197c4762a1bSJed Brown   FVCtx          *ctx = (FVCtx*)vctx;
198c4762a1bSJed Brown   PetscInt       i,j,k,Mx,dof,xs,xm,sm = ctx->sm,mf = ctx->mf,fm = ctx->fm,ms = ctx->ms;
199c4762a1bSJed Brown   PetscReal      hxf,hxm,hxs,cfl_idt = 0;
200c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
201c4762a1bSJed Brown   Vec            Xloc;
202c4762a1bSJed Brown   DM             da;
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   PetscFunctionBeginUser;
2055f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
2065f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&Xloc));                          /* Xloc contains ghost points                                     */
2075f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));   /* Mx is the number of center points                              */
208c4762a1bSJed Brown   hxs  = (ctx->xmax-ctx->xmin)/8.0/ctx->sm;
209c4762a1bSJed Brown   hxm  = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm);
210c4762a1bSJed Brown   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf);
2115f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));       /* X is solution vector which does not contain ghost points       */
2125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc));
213c4762a1bSJed Brown 
2145f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(F));                                   /* F is the right hand side function corresponds to center points */
215c4762a1bSJed Brown 
2165f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,Xloc,&x));
2175f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,F,&f));
2185f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope));                  /* contains ghost points                                           */
219c4762a1bSJed Brown 
2205f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
223c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
224c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
225c4762a1bSJed Brown     }
226c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
227c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
228c4762a1bSJed Brown     }
229c4762a1bSJed Brown   }
230c4762a1bSJed Brown   for (i=xs-1; i<xs+xm+1; i++) {
231c4762a1bSJed Brown     struct _LimitInfo info;
232c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
233c4762a1bSJed Brown     /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
2345f80ce2aSJacob Faibussowitsch     CHKERRQ((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds));
235c4762a1bSJed Brown     /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
2365f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof));
237c4762a1bSJed Brown     cjmpL = &ctx->cjmpLR[0];
238c4762a1bSJed Brown     cjmpR = &ctx->cjmpLR[dof];
239c4762a1bSJed Brown     for (j=0; j<dof; j++) {
240c4762a1bSJed Brown       PetscScalar jmpL,jmpR;
241c4762a1bSJed Brown       jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
242c4762a1bSJed Brown       jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
243c4762a1bSJed Brown       for (k=0; k<dof; k++) {
244c4762a1bSJed Brown         cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
245c4762a1bSJed Brown         cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
246c4762a1bSJed Brown       }
247c4762a1bSJed Brown     }
248c4762a1bSJed Brown     /* Apply limiter to the left and right characteristic jumps */
249c4762a1bSJed Brown     info.m  = dof;
250c4762a1bSJed Brown     info.hxs = hxs;
251c4762a1bSJed Brown     info.hxm = hxm;
252c4762a1bSJed Brown     info.hxf = hxf;
253c4762a1bSJed Brown     (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope);
254c4762a1bSJed Brown     for (j=0; j<dof; j++) {
255c4762a1bSJed Brown       PetscScalar tmp = 0;
256c4762a1bSJed Brown       for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
257c4762a1bSJed Brown       slope[i*dof+j] = tmp;
258c4762a1bSJed Brown     }
259c4762a1bSJed Brown   }
260c4762a1bSJed Brown 
261c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
262c4762a1bSJed Brown     PetscReal   maxspeed;
263c4762a1bSJed Brown     PetscScalar *uL,*uR;
264c4762a1bSJed Brown     uL = &ctx->uLR[0];
265c4762a1bSJed Brown     uR = &ctx->uLR[dof];
266c4762a1bSJed Brown     if (i < sm || i > ms) { /* slow region */
267c4762a1bSJed Brown       for (j=0; j<dof; j++) {
268c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
269c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
270c4762a1bSJed Brown       }
2715f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
272c4762a1bSJed Brown       if (i > xs) {
273c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs;
274c4762a1bSJed Brown       }
275c4762a1bSJed Brown       if (i < xs+xm) {
276c4762a1bSJed Brown         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs;
277c4762a1bSJed Brown       }
278c4762a1bSJed Brown     } else if (i == sm) { /* interface between slow and medium component */
279c4762a1bSJed Brown       for (j=0; j<dof; j++) {
280c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
281c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2;
282c4762a1bSJed Brown       }
2835f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
284c4762a1bSJed Brown       if (i > xs) {
285c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs;
286c4762a1bSJed Brown       }
287c4762a1bSJed Brown       if (i < xs+xm) {
288c4762a1bSJed Brown         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxm;
289c4762a1bSJed Brown       }
290c4762a1bSJed Brown     } else if (i == ms) { /* interface between medium and slow regions */
291c4762a1bSJed Brown       for (j=0; j<dof; j++) {
292c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2;
293c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
294c4762a1bSJed Brown       }
2955f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
296c4762a1bSJed Brown       if (i > xs) {
297c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxm;
298c4762a1bSJed Brown       }
299c4762a1bSJed Brown       if (i < xs+xm) {
300c4762a1bSJed Brown         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs;
301c4762a1bSJed Brown       }
302c4762a1bSJed Brown     } else if (i < mf || i > fm) { /* medium region */
303c4762a1bSJed Brown       for (j=0; j<dof; j++) {
304c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2;
305c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2;
306c4762a1bSJed Brown       }
3075f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
308c4762a1bSJed Brown       if (i > xs) {
309c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxm;
310c4762a1bSJed Brown       }
311c4762a1bSJed Brown       if (i < xs+xm) {
312c4762a1bSJed Brown         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxm;
313c4762a1bSJed Brown       }
314c4762a1bSJed Brown     } else if (i == mf) { /* interface between medium and fast regions */
315c4762a1bSJed Brown       for (j=0; j<dof; j++) {
316c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2;
317c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
318c4762a1bSJed Brown       }
3195f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
320c4762a1bSJed Brown       if (i > xs) {
321c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxm;
322c4762a1bSJed Brown       }
323c4762a1bSJed Brown       if (i < xs+xm) {
324c4762a1bSJed Brown         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxf;
325c4762a1bSJed Brown       }
326c4762a1bSJed Brown     } else if (i == fm) { /* interface between fast and medium regions */
327c4762a1bSJed Brown       for (j=0; j<dof; j++) {
328c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
329c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2;
330c4762a1bSJed Brown       }
3315f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
332c4762a1bSJed Brown       if (i > xs) {
333c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf;
334c4762a1bSJed Brown       }
335c4762a1bSJed Brown       if (i < xs+xm) {
336c4762a1bSJed Brown         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxm;
337c4762a1bSJed Brown       }
338c4762a1bSJed Brown     } else { /* fast region */
339c4762a1bSJed Brown       for (j=0; j<dof; j++) {
340c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
341c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
342c4762a1bSJed Brown       }
3435f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
344c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */
345c4762a1bSJed Brown       if (i > xs) {
346c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf;
347c4762a1bSJed Brown       }
348c4762a1bSJed Brown       if (i < xs+xm) {
349c4762a1bSJed Brown         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxf;
350c4762a1bSJed Brown       }
351c4762a1bSJed Brown     }
352c4762a1bSJed Brown   }
3535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x));
3545f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,F,&f));
3555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope));
3565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&Xloc));
3575f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da)));
358c4762a1bSJed Brown   if (0) {
359c4762a1bSJed Brown     /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */
360c4762a1bSJed Brown     PetscReal dt,tnow;
3615f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetTimeStep(ts,&dt));
3625f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetTime(ts,&tnow));
363c4762a1bSJed Brown     if (dt > 0.5/ctx->cfl_idt) {
3645f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(ctx->comm,"Stability constraint exceeded at t=%g, dt %g > %g\n",(double)tnow,(double)dt,(double)(0.5/ctx->cfl_idt)));
365c4762a1bSJed Brown     }
366c4762a1bSJed Brown   }
367c4762a1bSJed Brown   PetscFunctionReturn(0);
368c4762a1bSJed Brown }
369c4762a1bSJed Brown 
370c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */
371c4762a1bSJed Brown PetscErrorCode FVRHSFunctionslow_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
372c4762a1bSJed Brown {
373c4762a1bSJed Brown   FVCtx          *ctx = (FVCtx*)vctx;
374c4762a1bSJed Brown   PetscInt       i,j,k,Mx,dof,xs,xm,islow = 0,sm = ctx->sm,ms = ctx->ms,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth;
375c4762a1bSJed Brown   PetscReal      hxs,hxm,hxf,cfl_idt = 0;
376c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
377c4762a1bSJed Brown   Vec            Xloc;
378c4762a1bSJed Brown   DM             da;
379c4762a1bSJed Brown 
380c4762a1bSJed Brown   PetscFunctionBeginUser;
3815f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
3825f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&Xloc));
3835f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
384c4762a1bSJed Brown   hxs  = (ctx->xmax-ctx->xmin)/8.0/ctx->sm;
385c4762a1bSJed Brown   hxm  = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm);
386c4762a1bSJed Brown   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf);
3875f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));
3885f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc));
3895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(F));
3905f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,Xloc,&x));
3915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
3925f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope));
3935f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
394c4762a1bSJed Brown 
395c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
396c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
397c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
398c4762a1bSJed Brown     }
399c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
400c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
401c4762a1bSJed Brown     }
402c4762a1bSJed Brown   }
403c4762a1bSJed Brown   for (i=xs-1; i<xs+xm+1; i++) {
404c4762a1bSJed Brown     struct _LimitInfo info;
405c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
406c4762a1bSJed Brown     if (i < sm-lsbwidth+1 || i > ms+rsbwidth-2) { /* slow components and the first and last fast components */
407c4762a1bSJed Brown       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
4085f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds));
409c4762a1bSJed Brown       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
4105f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof));
411c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
412c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
413c4762a1bSJed Brown       for (j=0; j<dof; j++) {
414c4762a1bSJed Brown         PetscScalar jmpL,jmpR;
415c4762a1bSJed Brown         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
416c4762a1bSJed Brown         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
417c4762a1bSJed Brown         for (k=0; k<dof; k++) {
418c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
419c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
420c4762a1bSJed Brown         }
421c4762a1bSJed Brown       }
422c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
423c4762a1bSJed Brown       info.m  = dof;
424c4762a1bSJed Brown       info.hxs = hxs;
425c4762a1bSJed Brown       info.hxm = hxm;
426c4762a1bSJed Brown       info.hxf = hxf;
427c4762a1bSJed Brown       (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope);
428c4762a1bSJed Brown       for (j=0; j<dof; j++) {
429c4762a1bSJed Brown         PetscScalar tmp = 0;
430c4762a1bSJed Brown         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
431c4762a1bSJed Brown           slope[i*dof+j] = tmp;
432c4762a1bSJed Brown       }
433c4762a1bSJed Brown     }
434c4762a1bSJed Brown   }
435c4762a1bSJed Brown 
436c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
437c4762a1bSJed Brown     PetscReal   maxspeed;
438c4762a1bSJed Brown     PetscScalar *uL,*uR;
439c4762a1bSJed Brown     uL = &ctx->uLR[0];
440c4762a1bSJed Brown     uR = &ctx->uLR[dof];
441c4762a1bSJed Brown     if (i < sm-lsbwidth) { /* slow region */
442c4762a1bSJed Brown       for (j=0; j<dof; j++) {
443c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
444c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
445c4762a1bSJed Brown       }
4465f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
447c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */
448c4762a1bSJed Brown       if (i > xs) {
449c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
450c4762a1bSJed Brown       }
451c4762a1bSJed Brown       if (i < xs+xm) {
452c4762a1bSJed Brown         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
453c4762a1bSJed Brown         islow++;
454c4762a1bSJed Brown       }
455c4762a1bSJed Brown     }
456c4762a1bSJed Brown     if (i == sm-lsbwidth) { /* interface between slow and medium regions */
457c4762a1bSJed Brown       for (j=0; j<dof; j++) {
458c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
459c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
460c4762a1bSJed Brown       }
4615f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
462c4762a1bSJed Brown       if (i > xs) {
463c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
464c4762a1bSJed Brown       }
465c4762a1bSJed Brown     }
466c4762a1bSJed Brown     if (i == ms+rsbwidth) { /* interface between medium and slow regions */
467c4762a1bSJed Brown       for (j=0; j<dof; j++) {
468c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
469c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
470c4762a1bSJed Brown       }
4715f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
472c4762a1bSJed Brown       if (i < xs+xm) {
473c4762a1bSJed Brown         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
474c4762a1bSJed Brown         islow++;
475c4762a1bSJed Brown       }
476c4762a1bSJed Brown     }
477c4762a1bSJed Brown     if (i > ms+rsbwidth) { /* slow region */
478c4762a1bSJed Brown       for (j=0; j<dof; j++) {
479c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
480c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
481c4762a1bSJed Brown       }
4825f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
483c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */
484c4762a1bSJed Brown       if (i > xs) {
485c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
486c4762a1bSJed Brown       }
487c4762a1bSJed Brown       if (i < xs+xm) {
488c4762a1bSJed Brown         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
489c4762a1bSJed Brown         islow++;
490c4762a1bSJed Brown       }
491c4762a1bSJed Brown     }
492c4762a1bSJed Brown   }
4935f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x));
4945f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
4955f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope));
4965f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&Xloc));
4975f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da)));
498c4762a1bSJed Brown   PetscFunctionReturn(0);
499c4762a1bSJed Brown }
500c4762a1bSJed Brown 
501c4762a1bSJed Brown PetscErrorCode FVRHSFunctionslowbuffer_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
502c4762a1bSJed Brown {
503c4762a1bSJed Brown   FVCtx          *ctx = (FVCtx*)vctx;
504c4762a1bSJed Brown   PetscInt       i,j,k,Mx,dof,xs,xm,islowbuffer = 0,sm = ctx->sm,ms = ctx->ms,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth;
505c4762a1bSJed Brown   PetscReal      hxs,hxm,hxf;
506c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
507c4762a1bSJed Brown   Vec            Xloc;
508c4762a1bSJed Brown   DM             da;
509c4762a1bSJed Brown 
510c4762a1bSJed Brown   PetscFunctionBeginUser;
5115f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
5125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&Xloc));
5135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
514c4762a1bSJed Brown   hxs  = (ctx->xmax-ctx->xmin)/8.0/ctx->sm;
515c4762a1bSJed Brown   hxm  = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm);
516c4762a1bSJed Brown   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf);
5175f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));
5185f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc));
5195f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(F));
5205f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,Xloc,&x));
5215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
5225f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope));
5235f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
524c4762a1bSJed Brown 
525c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
526c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
527c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
528c4762a1bSJed Brown     }
529c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
530c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
531c4762a1bSJed Brown     }
532c4762a1bSJed Brown   }
533c4762a1bSJed Brown   for (i=xs-1; i<xs+xm+1; i++) {
534c4762a1bSJed Brown     struct _LimitInfo info;
535c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
536c4762a1bSJed Brown     if ((i > sm-lsbwidth-2 && i < sm+1) || (i > ms-2 && i < ms+rsbwidth+1)) {
537c4762a1bSJed Brown       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
5385f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds));
539c4762a1bSJed Brown       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
5405f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof));
541c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
542c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
543c4762a1bSJed Brown       for (j=0; j<dof; j++) {
544c4762a1bSJed Brown         PetscScalar jmpL,jmpR;
545c4762a1bSJed Brown         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
546c4762a1bSJed Brown         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
547c4762a1bSJed Brown         for (k=0; k<dof; k++) {
548c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
549c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
550c4762a1bSJed Brown         }
551c4762a1bSJed Brown       }
552c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
553c4762a1bSJed Brown       info.m  = dof;
554c4762a1bSJed Brown       info.hxs = hxs;
555c4762a1bSJed Brown       info.hxm = hxm;
556c4762a1bSJed Brown       info.hxf = hxf;
557c4762a1bSJed Brown       (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope);
558c4762a1bSJed Brown       for (j=0; j<dof; j++) {
559c4762a1bSJed Brown         PetscScalar tmp = 0;
560c4762a1bSJed Brown         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
561c4762a1bSJed Brown           slope[i*dof+j] = tmp;
562c4762a1bSJed Brown       }
563c4762a1bSJed Brown     }
564c4762a1bSJed Brown   }
565c4762a1bSJed Brown 
566c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
567c4762a1bSJed Brown     PetscReal   maxspeed;
568c4762a1bSJed Brown     PetscScalar *uL,*uR;
569c4762a1bSJed Brown     uL = &ctx->uLR[0];
570c4762a1bSJed Brown     uR = &ctx->uLR[dof];
571c4762a1bSJed Brown     if (i == sm-lsbwidth) {
572c4762a1bSJed Brown       for (j=0; j<dof; j++) {
573c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
574c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
575c4762a1bSJed Brown       }
5765f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
577c4762a1bSJed Brown       if (i < xs+xm) {
578c4762a1bSJed Brown         for (j=0; j<dof; j++) f[islowbuffer*dof+j] += ctx->flux[j]/hxs;
579c4762a1bSJed Brown         islowbuffer++;
580c4762a1bSJed Brown       }
581c4762a1bSJed Brown     }
582c4762a1bSJed Brown     if (i > sm-lsbwidth && i < sm) {
583c4762a1bSJed Brown       for (j=0; j<dof; j++) {
584c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
585c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
586c4762a1bSJed Brown       }
5875f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
588c4762a1bSJed Brown       if (i > xs) {
589c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(islowbuffer-1)*dof+j] -= ctx->flux[j]/hxs;
590c4762a1bSJed Brown       }
591c4762a1bSJed Brown       if (i < xs+xm) {
592c4762a1bSJed Brown         for (j=0; j<dof; j++) f[islowbuffer*dof+j] += ctx->flux[j]/hxs;
593c4762a1bSJed Brown         islowbuffer++;
594c4762a1bSJed Brown       }
595c4762a1bSJed Brown     }
596c4762a1bSJed Brown     if (i == sm) { /* interface between the slow region and the medium region */
597c4762a1bSJed Brown       for (j=0; j<dof; j++) {
598c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
599c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2;
600c4762a1bSJed Brown       }
6015f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
602c4762a1bSJed Brown       if (i > xs) {
603c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(islowbuffer-1)*dof+j] -= ctx->flux[j]/hxs;
604c4762a1bSJed Brown       }
605c4762a1bSJed Brown     }
606c4762a1bSJed Brown     if (i == ms) { /* interface between the medium region and the slow region */
607c4762a1bSJed Brown       for (j=0; j<dof; j++) {
608c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2;
609c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
610c4762a1bSJed Brown       }
6115f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
612c4762a1bSJed Brown       if (i < xs+xm) {
613c4762a1bSJed Brown         for (j=0; j<dof; j++) f[islowbuffer*dof+j] += ctx->flux[j]/hxs;
614c4762a1bSJed Brown         islowbuffer++;
615c4762a1bSJed Brown       }
616c4762a1bSJed Brown     }
617c4762a1bSJed Brown     if (i > ms && i < ms+rsbwidth) {
618c4762a1bSJed Brown       for (j=0; j<dof; j++) {
619c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
620c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
621c4762a1bSJed Brown       }
6225f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
623c4762a1bSJed Brown       if (i > xs) {
624c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(islowbuffer-1)*dof+j] -= ctx->flux[j]/hxs;
625c4762a1bSJed Brown       }
626c4762a1bSJed Brown       if (i < xs+xm) {
627c4762a1bSJed Brown         for (j=0; j<dof; j++) f[islowbuffer*dof+j] += ctx->flux[j]/hxs;
628c4762a1bSJed Brown         islowbuffer++;
629c4762a1bSJed Brown       }
630c4762a1bSJed Brown     }
631c4762a1bSJed Brown     if (i == ms+rsbwidth) {
632c4762a1bSJed Brown       for (j=0; j<dof; j++) {
633c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
634c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
635c4762a1bSJed Brown       }
6365f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
637c4762a1bSJed Brown       if (i > xs) {
638c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(islowbuffer-1)*dof+j] -= ctx->flux[j]/hxs;
639c4762a1bSJed Brown       }
640c4762a1bSJed Brown     }
641c4762a1bSJed Brown   }
6425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x));
6435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
6445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope));
6455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&Xloc));
646c4762a1bSJed Brown   PetscFunctionReturn(0);
647c4762a1bSJed Brown }
648c4762a1bSJed Brown 
649c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for medium components ----------------------------------- */
650c4762a1bSJed Brown PetscErrorCode FVRHSFunctionmedium_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
651c4762a1bSJed Brown {
652c4762a1bSJed Brown   FVCtx          *ctx = (FVCtx*)vctx;
653c4762a1bSJed 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;
654c4762a1bSJed Brown   PetscReal      hxs,hxm,hxf;
655c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
656c4762a1bSJed Brown   Vec            Xloc;
657c4762a1bSJed Brown   DM             da;
658c4762a1bSJed Brown 
659c4762a1bSJed Brown   PetscFunctionBeginUser;
6605f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
6615f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&Xloc));
6625f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
663c4762a1bSJed Brown   hxs  = (ctx->xmax-ctx->xmin)/8.0/ctx->sm;
664c4762a1bSJed Brown   hxm  = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm);
665c4762a1bSJed Brown   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf);
6665f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));
6675f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc));
6685f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(F));
6695f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,Xloc,&x));
6705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
6715f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope));
6725f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
673c4762a1bSJed Brown 
674c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
675c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
676c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
677c4762a1bSJed Brown     }
678c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
679c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
680c4762a1bSJed Brown     }
681c4762a1bSJed Brown   }
682c4762a1bSJed Brown   for (i=xs-1; i<xs+xm+1; i++) {
683c4762a1bSJed Brown     struct _LimitInfo info;
684c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
685c4762a1bSJed 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 */
686c4762a1bSJed Brown       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
6875f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds));
688c4762a1bSJed Brown       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
6895f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof));
690c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
691c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
692c4762a1bSJed Brown       for (j=0; j<dof; j++) {
693c4762a1bSJed Brown         PetscScalar jmpL,jmpR;
694c4762a1bSJed Brown         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
695c4762a1bSJed Brown         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
696c4762a1bSJed Brown         for (k=0; k<dof; k++) {
697c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
698c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
699c4762a1bSJed Brown         }
700c4762a1bSJed Brown       }
701c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
702c4762a1bSJed Brown       info.m  = dof;
703c4762a1bSJed Brown       info.hxs = hxs;
704c4762a1bSJed Brown       info.hxm = hxm;
705c4762a1bSJed Brown       info.hxf = hxf;
706c4762a1bSJed Brown       (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope);
707c4762a1bSJed Brown       for (j=0; j<dof; j++) {
708c4762a1bSJed Brown         PetscScalar tmp = 0;
709c4762a1bSJed Brown         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
710c4762a1bSJed Brown           slope[i*dof+j] = tmp;
711c4762a1bSJed Brown       }
712c4762a1bSJed Brown     }
713c4762a1bSJed Brown   }
714c4762a1bSJed Brown 
715c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
716c4762a1bSJed Brown     PetscReal   maxspeed;
717c4762a1bSJed Brown     PetscScalar *uL,*uR;
718c4762a1bSJed Brown     uL = &ctx->uLR[0];
719c4762a1bSJed Brown     uR = &ctx->uLR[dof];
720c4762a1bSJed Brown     if (i == sm) { /* interface between slow and medium regions */
721c4762a1bSJed Brown       for (j=0; j<dof; j++) {
722c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
723c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2;
724c4762a1bSJed Brown       }
7255f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
726c4762a1bSJed Brown       if (i < xs+xm) {
727c4762a1bSJed Brown         for (j=0; j<dof; j++) f[imedium*dof+j] += ctx->flux[j]/hxm;
728c4762a1bSJed Brown         imedium++;
729c4762a1bSJed Brown       }
730c4762a1bSJed Brown     }
731c4762a1bSJed Brown     if (i > sm && i < mf-lmbwidth) { /* medium region */
732c4762a1bSJed Brown       for (j=0; j<dof; j++) {
733c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2;
734c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2;
735c4762a1bSJed Brown       }
7365f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
737c4762a1bSJed Brown       if (i > xs) {
738c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(imedium-1)*dof+j] -= ctx->flux[j]/hxm;
739c4762a1bSJed Brown       }
740c4762a1bSJed Brown       if (i < xs+xm) {
741c4762a1bSJed Brown         for (j=0; j<dof; j++) f[imedium*dof+j] += ctx->flux[j]/hxm;
742c4762a1bSJed Brown         imedium++;
743c4762a1bSJed Brown       }
744c4762a1bSJed Brown     }
745c4762a1bSJed Brown     if (i == mf-lmbwidth) { /* interface between medium and fast regions */
746c4762a1bSJed Brown       for (j=0; j<dof; j++) {
747c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2;
748c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2;
749c4762a1bSJed Brown       }
7505f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
751c4762a1bSJed Brown       if (i > xs) {
752c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(imedium-1)*dof+j] -= ctx->flux[j]/hxm;
753c4762a1bSJed Brown       }
754c4762a1bSJed Brown     }
755c4762a1bSJed Brown     if (i == fm+rmbwidth) { /* interface between fast and medium regions */
756c4762a1bSJed Brown       for (j=0; j<dof; j++) {
757c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2;
758c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2;
759c4762a1bSJed Brown       }
7605f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
761c4762a1bSJed Brown       if (i < xs+xm) {
762c4762a1bSJed Brown         for (j=0; j<dof; j++) f[imedium*dof+j] += ctx->flux[j]/hxm;
763c4762a1bSJed Brown         imedium++;
764c4762a1bSJed Brown       }
765c4762a1bSJed Brown     }
766c4762a1bSJed Brown     if (i > fm+rmbwidth && i < ms) { /* medium region */
767c4762a1bSJed Brown       for (j=0; j<dof; j++) {
768c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2;
769c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2;
770c4762a1bSJed Brown       }
7715f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
772c4762a1bSJed Brown       if (i > xs) {
773c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(imedium-1)*dof+j] -= ctx->flux[j]/hxm;
774c4762a1bSJed Brown       }
775c4762a1bSJed Brown       if (i < xs+xm) {
776c4762a1bSJed Brown         for (j=0; j<dof; j++) f[imedium*dof+j] += ctx->flux[j]/hxm;
777c4762a1bSJed Brown         imedium++;
778c4762a1bSJed Brown       }
779c4762a1bSJed Brown     }
780c4762a1bSJed Brown     if (i == ms) { /* interface between medium and slow regions */
781c4762a1bSJed Brown       for (j=0; j<dof; j++) {
782c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2;
783c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
784c4762a1bSJed Brown       }
7855f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
786c4762a1bSJed Brown       if (i > xs) {
787c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(imedium-1)*dof+j] -= ctx->flux[j]/hxm;
788c4762a1bSJed Brown       }
789c4762a1bSJed Brown     }
790c4762a1bSJed Brown   }
7915f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x));
7925f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
7935f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope));
7945f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&Xloc));
795c4762a1bSJed Brown   PetscFunctionReturn(0);
796c4762a1bSJed Brown }
797c4762a1bSJed Brown 
798c4762a1bSJed Brown PetscErrorCode FVRHSFunctionmediumbuffer_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
799c4762a1bSJed Brown {
800c4762a1bSJed Brown   FVCtx          *ctx = (FVCtx*)vctx;
801c4762a1bSJed Brown   PetscInt       i,j,k,Mx,dof,xs,xm,imediumbuffer = 0,mf = ctx->mf,fm = ctx->fm,lmbwidth = ctx->lmbwidth,rmbwidth = ctx->rmbwidth;
802c4762a1bSJed Brown   PetscReal      hxs,hxm,hxf;
803c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
804c4762a1bSJed Brown   Vec            Xloc;
805c4762a1bSJed Brown   DM             da;
806c4762a1bSJed Brown 
807c4762a1bSJed Brown   PetscFunctionBeginUser;
8085f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
8095f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&Xloc));
8105f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
811c4762a1bSJed Brown   hxs  = (ctx->xmax-ctx->xmin)/8.0/ctx->sm;
812c4762a1bSJed Brown   hxm  = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm);
813c4762a1bSJed Brown   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf);
8145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));
8155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc));
8165f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(F));
8175f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,Xloc,&x));
8185f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
8195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope));
8205f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
821c4762a1bSJed Brown 
822c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
823c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
824c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
825c4762a1bSJed Brown     }
826c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
827c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
828c4762a1bSJed Brown     }
829c4762a1bSJed Brown   }
830c4762a1bSJed Brown   for (i=xs-1; i<xs+xm+1; i++) {
831c4762a1bSJed Brown     struct _LimitInfo info;
832c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
833c4762a1bSJed Brown     if ((i > mf-lmbwidth-2 && i < mf+1) || (i > fm-2 && i < fm+rmbwidth+1)) {
834c4762a1bSJed Brown       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
8355f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds));
836c4762a1bSJed Brown       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
8375f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof));
838c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
839c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
840c4762a1bSJed Brown       for (j=0; j<dof; j++) {
841c4762a1bSJed Brown         PetscScalar jmpL,jmpR;
842c4762a1bSJed Brown         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
843c4762a1bSJed Brown         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
844c4762a1bSJed Brown         for (k=0; k<dof; k++) {
845c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
846c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
847c4762a1bSJed Brown         }
848c4762a1bSJed Brown       }
849c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
850c4762a1bSJed Brown       info.m  = dof;
851c4762a1bSJed Brown       info.hxs = hxs;
852c4762a1bSJed Brown       info.hxm = hxm;
853c4762a1bSJed Brown       info.hxf = hxf;
854c4762a1bSJed Brown       (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope);
855c4762a1bSJed Brown       for (j=0; j<dof; j++) {
856c4762a1bSJed Brown         PetscScalar tmp = 0;
857c4762a1bSJed Brown         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
858c4762a1bSJed Brown           slope[i*dof+j] = tmp;
859c4762a1bSJed Brown       }
860c4762a1bSJed Brown     }
861c4762a1bSJed Brown   }
862c4762a1bSJed Brown 
863c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
864c4762a1bSJed Brown     PetscReal   maxspeed;
865c4762a1bSJed Brown     PetscScalar *uL,*uR;
866c4762a1bSJed Brown     uL = &ctx->uLR[0];
867c4762a1bSJed Brown     uR = &ctx->uLR[dof];
868c4762a1bSJed Brown     if (i == mf-lmbwidth) {
869c4762a1bSJed Brown       for (j=0; j<dof; j++) {
870c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2;
871c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2;
872c4762a1bSJed Brown       }
8735f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
874c4762a1bSJed Brown       if (i < xs+xm) {
875c4762a1bSJed Brown         for (j=0; j<dof; j++) f[imediumbuffer*dof+j] += ctx->flux[j]/hxm;
876c4762a1bSJed Brown         imediumbuffer++;
877c4762a1bSJed Brown       }
878c4762a1bSJed Brown     }
879c4762a1bSJed Brown     if (i > mf-lmbwidth && i < mf) {
880c4762a1bSJed Brown       for (j=0; j<dof; j++) {
881c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2;
882c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2;
883c4762a1bSJed Brown       }
8845f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
885c4762a1bSJed Brown       if (i > xs) {
886c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(imediumbuffer-1)*dof+j] -= ctx->flux[j]/hxm;
887c4762a1bSJed Brown       }
888c4762a1bSJed Brown       if (i < xs+xm) {
889c4762a1bSJed Brown         for (j=0; j<dof; j++) f[imediumbuffer*dof+j] += ctx->flux[j]/hxm;
890c4762a1bSJed Brown         imediumbuffer++;
891c4762a1bSJed Brown       }
892c4762a1bSJed Brown     }
893c4762a1bSJed Brown     if (i == mf) { /* interface between the medium region and the fast region */
894c4762a1bSJed Brown       for (j=0; j<dof; j++) {
895c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2;
896c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
897c4762a1bSJed Brown       }
8985f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
899c4762a1bSJed Brown       if (i > xs) {
900c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(imediumbuffer-1)*dof+j] -= ctx->flux[j]/hxm;
901c4762a1bSJed Brown       }
902c4762a1bSJed Brown     }
903c4762a1bSJed Brown     if (i == fm) { /* interface between the fast region and the medium region */
904c4762a1bSJed Brown       for (j=0; j<dof; j++) {
905c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
906c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2;
907c4762a1bSJed Brown       }
9085f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
909c4762a1bSJed Brown       if (i < xs+xm) {
910c4762a1bSJed Brown         for (j=0; j<dof; j++) f[imediumbuffer*dof+j] += ctx->flux[j]/hxm;
911c4762a1bSJed Brown         imediumbuffer++;
912c4762a1bSJed Brown       }
913c4762a1bSJed Brown     }
914c4762a1bSJed Brown     if (i > fm && i < fm+rmbwidth) {
915c4762a1bSJed Brown       for (j=0; j<dof; j++) {
916c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2;
917c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2;
918c4762a1bSJed Brown       }
9195f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
920c4762a1bSJed Brown       if (i > xs) {
921c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(imediumbuffer-1)*dof+j] -= ctx->flux[j]/hxm;
922c4762a1bSJed Brown       }
923c4762a1bSJed Brown       if (i < xs+xm) {
924c4762a1bSJed Brown         for (j=0; j<dof; j++) f[imediumbuffer*dof+j] += ctx->flux[j]/hxm;
925c4762a1bSJed Brown         imediumbuffer++;
926c4762a1bSJed Brown       }
927c4762a1bSJed Brown     }
928c4762a1bSJed Brown     if (i == fm+rmbwidth) {
929c4762a1bSJed Brown       for (j=0; j<dof; j++) {
930c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2;
931c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2;
932c4762a1bSJed Brown       }
9335f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
934c4762a1bSJed Brown       if (i > xs) {
935c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(imediumbuffer-1)*dof+j] -= ctx->flux[j]/hxm;
936c4762a1bSJed Brown       }
937c4762a1bSJed Brown     }
938c4762a1bSJed Brown   }
9395f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x));
9405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
9415f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope));
9425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&Xloc));
943c4762a1bSJed Brown   PetscFunctionReturn(0);
944c4762a1bSJed Brown }
945c4762a1bSJed Brown 
946c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for fast  parts ----------------------------------- */
947c4762a1bSJed Brown PetscErrorCode FVRHSFunctionfast_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
948c4762a1bSJed Brown {
949c4762a1bSJed Brown   FVCtx          *ctx = (FVCtx*)vctx;
950c4762a1bSJed Brown   PetscInt       i,j,k,Mx,dof,xs,xm,ifast = 0,mf = ctx->mf,fm = ctx->fm;
951c4762a1bSJed Brown   PetscReal      hxs,hxm,hxf;
952c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
953c4762a1bSJed Brown   Vec            Xloc;
954c4762a1bSJed Brown   DM             da;
955c4762a1bSJed Brown 
956c4762a1bSJed Brown   PetscFunctionBeginUser;
9575f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
9585f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&Xloc));
9595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
960c4762a1bSJed Brown   hxs  = (ctx->xmax-ctx->xmin)/8.0/ctx->sm;
961c4762a1bSJed Brown   hxm  = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm);
962c4762a1bSJed Brown   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf);
9635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));
9645f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc));
9655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(F));
9665f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,Xloc,&x));
9675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
9685f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope));
9695f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
970c4762a1bSJed Brown 
971c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
972c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
973c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
974c4762a1bSJed Brown     }
975c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
976c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
977c4762a1bSJed Brown     }
978c4762a1bSJed Brown   }
979c4762a1bSJed 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 */
980c4762a1bSJed Brown     struct _LimitInfo info;
981c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
982c4762a1bSJed Brown     if (i > mf-2 && i < fm+1) {
9835f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds));
9845f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof));
985c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
986c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
987c4762a1bSJed Brown       for (j=0; j<dof; j++) {
988c4762a1bSJed Brown         PetscScalar jmpL,jmpR;
989c4762a1bSJed Brown         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
990c4762a1bSJed Brown         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
991c4762a1bSJed Brown         for (k=0; k<dof; k++) {
992c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
993c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
994c4762a1bSJed Brown         }
995c4762a1bSJed Brown       }
996c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
997c4762a1bSJed Brown       info.m  = dof;
998c4762a1bSJed Brown       info.hxs = hxs;
999c4762a1bSJed Brown       info.hxm = hxm;
1000c4762a1bSJed Brown       info.hxf = hxf;
1001c4762a1bSJed Brown       (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope);
1002c4762a1bSJed Brown       for (j=0; j<dof; j++) {
1003c4762a1bSJed Brown       PetscScalar tmp = 0;
1004c4762a1bSJed Brown       for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
1005c4762a1bSJed Brown         slope[i*dof+j] = tmp;
1006c4762a1bSJed Brown       }
1007c4762a1bSJed Brown     }
1008c4762a1bSJed Brown   }
1009c4762a1bSJed Brown 
1010c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
1011c4762a1bSJed Brown     PetscReal   maxspeed;
1012c4762a1bSJed Brown     PetscScalar *uL,*uR;
1013c4762a1bSJed Brown     uL = &ctx->uLR[0];
1014c4762a1bSJed Brown     uR = &ctx->uLR[dof];
1015c4762a1bSJed Brown     if (i == mf) { /* interface between medium and fast regions */
1016c4762a1bSJed Brown       for (j=0; j<dof; j++) {
1017c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2;
1018c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
1019c4762a1bSJed Brown       }
10205f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
1021c4762a1bSJed Brown       if (i < xs+xm) {
1022c4762a1bSJed Brown         for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf;
1023c4762a1bSJed Brown         ifast++;
1024c4762a1bSJed Brown       }
1025c4762a1bSJed Brown     }
1026c4762a1bSJed Brown     if (i > mf && i < fm) { /* fast region */
1027c4762a1bSJed Brown       for (j=0; j<dof; j++) {
1028c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
1029c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
1030c4762a1bSJed Brown       }
10315f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
1032c4762a1bSJed Brown       if (i > xs) {
1033c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf;
1034c4762a1bSJed Brown       }
1035c4762a1bSJed Brown       if (i < xs+xm) {
1036c4762a1bSJed Brown         for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf;
1037c4762a1bSJed Brown         ifast++;
1038c4762a1bSJed Brown       }
1039c4762a1bSJed Brown     }
1040c4762a1bSJed Brown     if (i == fm) { /* interface between fast and medium regions */
1041c4762a1bSJed Brown       for (j=0; j<dof; j++) {
1042c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
1043c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2;
1044c4762a1bSJed Brown       }
10455f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
1046c4762a1bSJed Brown       if (i > xs) {
1047c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf;
1048c4762a1bSJed Brown       }
1049c4762a1bSJed Brown     }
1050c4762a1bSJed Brown   }
10515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x));
10525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
10535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope));
10545f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&Xloc));
1055c4762a1bSJed Brown   PetscFunctionReturn(0);
1056c4762a1bSJed Brown }
1057c4762a1bSJed Brown 
1058c4762a1bSJed Brown int main(int argc,char *argv[])
1059c4762a1bSJed Brown {
1060c4762a1bSJed Brown   char              lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m";
1061c4762a1bSJed Brown   PetscFunctionList limiters   = 0,physics = 0;
1062c4762a1bSJed Brown   MPI_Comm          comm;
1063c4762a1bSJed Brown   TS                ts;
1064c4762a1bSJed Brown   DM                da;
1065c4762a1bSJed Brown   Vec               X,X0,R;
1066c4762a1bSJed Brown   FVCtx             ctx;
1067c4762a1bSJed 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;
1068c4762a1bSJed Brown   PetscBool         view_final = PETSC_FALSE;
1069c4762a1bSJed Brown   PetscReal         ptime;
1070c4762a1bSJed Brown   PetscErrorCode    ierr;
1071c4762a1bSJed Brown 
1072*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,0,help));
1073c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
10745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMemzero(&ctx,sizeof(ctx)));
1075c4762a1bSJed Brown 
1076c4762a1bSJed Brown   /* Register limiters to be available on the command line */
10775f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"upwind"              ,Limit3_Upwind));
10785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"lax-wendroff"        ,Limit3_LaxWendroff));
10795f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"beam-warming"        ,Limit3_BeamWarming));
10805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"fromm"               ,Limit3_Fromm));
10815f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"minmod"              ,Limit3_Minmod));
10825f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"superbee"            ,Limit3_Superbee));
10835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"mc"                  ,Limit3_MC));
10845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"koren3"              ,Limit3_Koren3));
1085c4762a1bSJed Brown 
1086c4762a1bSJed Brown   /* Register physical models to be available on the command line */
10875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&physics,"advect"          ,PhysicsCreate_Advect));
1088c4762a1bSJed Brown 
1089c4762a1bSJed Brown   ctx.comm = comm;
1090c4762a1bSJed Brown   ctx.cfl  = 0.9;
1091c4762a1bSJed Brown   ctx.bctype = FVBC_PERIODIC;
1092c4762a1bSJed Brown   ctx.xmin = -1.0;
1093c4762a1bSJed Brown   ctx.xmax = 1.0;
1094c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm,NULL,"Finite Volume solver options","");CHKERRQ(ierr);
10955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL));
10965f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL));
10975f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsFList("-limit","Name of flux imiter to use","",limiters,lname,lname,sizeof(lname),NULL));
10985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL));
10995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsString("-view_final","Write final solution in ASCII MATLAB format to given file name","",final_fname,final_fname,sizeof(final_fname),&view_final));
11005f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL));
11015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-exact","Compare errors with exact solution","",ctx.exact,&ctx.exact,NULL));
11025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL));
11035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL));
11045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL));
11055f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-hratio","Spacing ratio","",ctx.hratio,&ctx.hratio,NULL));
1106c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1107c4762a1bSJed Brown 
1108c4762a1bSJed Brown   /* Choose the limiter from the list of registered limiters */
11095f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListFind(limiters,lname,&ctx.limit3));
11103c633725SBarry Smith   PetscCheck(ctx.limit3,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Limiter '%s' not found",lname);
1111c4762a1bSJed Brown 
1112c4762a1bSJed Brown   /* Choose the physics from the list of registered models */
1113c4762a1bSJed Brown   {
1114c4762a1bSJed Brown     PetscErrorCode (*r)(FVCtx*);
11155f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFunctionListFind(physics,physname,&r));
11163c633725SBarry Smith     PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Physics '%s' not found",physname);
1117c4762a1bSJed Brown     /* Create the physics, will set the number of fields and their names */
11185f80ce2aSJacob Faibussowitsch     CHKERRQ((*r)(&ctx));
1119c4762a1bSJed Brown   }
1120c4762a1bSJed Brown 
1121c4762a1bSJed Brown   /* Create a DMDA to manage the parallel grid */
11225f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics2.dof,2,NULL,&da));
11235f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
11245f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
1125c4762a1bSJed Brown   /* Inform the DMDA of the field names provided by the physics. */
1126c4762a1bSJed Brown   /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
1127c4762a1bSJed Brown   for (i=0; i<ctx.physics2.dof; i++) {
11285f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDASetFieldName(da,i,ctx.physics2.fieldname[i]));
1129c4762a1bSJed Brown   }
11305f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
11315f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
1132c4762a1bSJed Brown 
1133c4762a1bSJed Brown   /* Set coordinates of cell centers */
11345f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(da,ctx.xmin+0.5*(ctx.xmax-ctx.xmin)/Mx,ctx.xmax+0.5*(ctx.xmax-ctx.xmin)/Mx,0,0,0,0));
1135c4762a1bSJed Brown 
1136c4762a1bSJed Brown   /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
11375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc4(dof*dof,&ctx.R,dof*dof,&ctx.Rinv,2*dof,&ctx.cjmpLR,1*dof,&ctx.cslope));
11385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3(2*dof,&ctx.uLR,dof,&ctx.flux,dof,&ctx.speeds));
1139c4762a1bSJed Brown 
1140c4762a1bSJed Brown   /* Create a vector to store the solution and to save the initial state */
11415f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&X));
11425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(X,&X0));
11435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(X,&R));
1144c4762a1bSJed Brown 
1145c4762a1bSJed Brown   /* create index for slow parts and fast parts,
1146c4762a1bSJed Brown      count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */
1147c4762a1bSJed Brown   count_slow = Mx/(1+ctx.hratio)/(1+ctx.hratio);
1148c4762a1bSJed Brown   count_medium = 2*ctx.hratio*count_slow;
11492c71b3e2SJacob Faibussowitsch   PetscCheckFalse(count_slow%2 || count_medium%2,PETSC_COMM_WORLD,PETSC_ERR_USER,"Please adjust grid size Mx (-da_grid_x) and hratio (-hratio) so that Mx/(1+hartio)^2 and Mx*2*hratio/(1+hratio)^2 is even");
1150c4762a1bSJed Brown   count_fast = ctx.hratio*ctx.hratio*count_slow;
1151c4762a1bSJed Brown   ctx.sm = count_slow/2;
1152c4762a1bSJed Brown   ctx.mf = ctx.sm + count_medium/2;
1153c4762a1bSJed Brown   ctx.fm = ctx.mf + count_fast;
1154c4762a1bSJed Brown   ctx.ms = ctx.fm + count_medium/2;
11555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(xm*dof,&index_slow));
11565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(xm*dof,&index_medium));
11575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(xm*dof,&index_fast));
11585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(6*dof,&index_slowbuffer));
11595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(6*dof,&index_mediumbuffer));
1160c4762a1bSJed Brown   if (((AdvectCtx*)ctx.physics2.user)->a > 0) {
1161c4762a1bSJed Brown     ctx.lsbwidth = 2;
1162c4762a1bSJed Brown     ctx.rsbwidth = 4;
1163c4762a1bSJed Brown     ctx.lmbwidth = 2;
1164c4762a1bSJed Brown     ctx.rmbwidth = 4;
1165c4762a1bSJed Brown   } else {
1166c4762a1bSJed Brown     ctx.lsbwidth = 4;
1167c4762a1bSJed Brown     ctx.rsbwidth = 2;
1168c4762a1bSJed Brown     ctx.lmbwidth = 4;
1169c4762a1bSJed Brown     ctx.rmbwidth = 2;
1170c4762a1bSJed Brown   }
1171c4762a1bSJed Brown 
1172c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
1173c4762a1bSJed Brown     if (i < ctx.sm-ctx.lsbwidth || i > ctx.ms+ctx.rsbwidth-1)
1174c4762a1bSJed Brown       for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k;
1175c4762a1bSJed Brown     else if ((i >= ctx.sm-ctx.lsbwidth && i < ctx.sm) || (i > ctx.ms-1 && i <= ctx.ms+ctx.rsbwidth-1))
1176c4762a1bSJed Brown       for (k=0; k<dof; k++) index_slowbuffer[islowbuffer++] = i*dof+k;
1177c4762a1bSJed Brown     else if (i < ctx.mf-ctx.lmbwidth || i > ctx.fm+ctx.rmbwidth-1)
1178c4762a1bSJed Brown       for (k=0; k<dof; k++) index_medium[imedium++] = i*dof+k;
1179c4762a1bSJed Brown     else if ((i >= ctx.mf-ctx.lmbwidth && i < ctx.mf) || (i > ctx.fm-1 && i <= ctx.fm+ctx.rmbwidth-1))
1180c4762a1bSJed Brown       for (k=0; k<dof; k++) index_mediumbuffer[imediumbuffer++] = i*dof+k;
1181c4762a1bSJed Brown     else
1182c4762a1bSJed Brown       for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k;
1183c4762a1bSJed Brown   }
11845f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss));
11855f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,imedium,index_medium,PETSC_COPY_VALUES,&ctx.ism));
11865f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf));
11875f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,islowbuffer,index_slowbuffer,PETSC_COPY_VALUES,&ctx.issb));
11885f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,imediumbuffer,index_mediumbuffer,PETSC_COPY_VALUES,&ctx.ismb));
1189c4762a1bSJed Brown 
1190c4762a1bSJed Brown   /* Create a time-stepping object */
11915f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(comm,&ts));
11925f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts,da));
11935f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,R,FVRHSFunction_3WaySplit,&ctx));
11945f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetIS(ts,"slow",ctx.iss));
11955f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetIS(ts,"medium",ctx.ism));
11965f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetIS(ts,"fast",ctx.isf));
11975f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetIS(ts,"slowbuffer",ctx.issb));
11985f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetIS(ts,"mediumbuffer",ctx.ismb));
11995f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetRHSFunction(ts,"slow",NULL,FVRHSFunctionslow_3WaySplit,&ctx));
12005f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetRHSFunction(ts,"medium",NULL,FVRHSFunctionmedium_3WaySplit,&ctx));
12015f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetRHSFunction(ts,"fast",NULL,FVRHSFunctionfast_3WaySplit,&ctx));
12025f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetRHSFunction(ts,"slowbuffer",NULL,FVRHSFunctionslowbuffer_3WaySplit,&ctx));
12035f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetRHSFunction(ts,"mediumbuffer",NULL,FVRHSFunctionmediumbuffer_3WaySplit,&ctx));
1204c4762a1bSJed Brown 
12055f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSSSP));
12065f80ce2aSJacob Faibussowitsch   /*CHKERRQ(TSSetType(ts,TSMPRK));*/
12075f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,10));
12085f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
1209c4762a1bSJed Brown 
1210c4762a1bSJed Brown   /* Compute initial conditions and starting time step */
12115f80ce2aSJacob Faibussowitsch   CHKERRQ(FVSample_3WaySplit(&ctx,da,0,X0));
12125f80ce2aSJacob Faibussowitsch   CHKERRQ(FVRHSFunction_3WaySplit(ts,0,X0,X,(void*)&ctx)); /* Initial function evaluation, only used to determine max speed */
12135f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(X0,X));                        /* The function value was not used so we set X=X0 again */
12145f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt));
12155f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts)); /* Take runtime options */
12165f80ce2aSJacob Faibussowitsch   CHKERRQ(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD));
1217c4762a1bSJed Brown   {
1218c4762a1bSJed Brown     PetscInt          steps;
1219c4762a1bSJed Brown     PetscScalar       mass_initial,mass_final,mass_difference,mass_differenceg;
1220c4762a1bSJed Brown     const PetscScalar *ptr_X,*ptr_X0;
1221c4762a1bSJed Brown     const PetscReal   hs = (ctx.xmax-ctx.xmin)/4.0/count_slow;
1222c4762a1bSJed Brown     const PetscReal   hm = (ctx.xmax-ctx.xmin)/2.0/count_medium;
1223c4762a1bSJed Brown     const PetscReal   hf = (ctx.xmax-ctx.xmin)/4.0/count_fast;
1224c4762a1bSJed Brown 
12255f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSolve(ts,X));
12265f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetSolveTime(ts,&ptime));
12275f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetStepNumber(ts,&steps));
1228c4762a1bSJed Brown     /* calculate the total mass at initial time and final time */
1229c4762a1bSJed Brown     mass_initial = 0.0;
1230c4762a1bSJed Brown     mass_final   = 0.0;
12315f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArrayRead(da,X0,(void*)&ptr_X0));
12325f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArrayRead(da,X,(void*)&ptr_X));
1233c4762a1bSJed Brown     for (i=xs;i<xs+xm;i++) {
1234c4762a1bSJed Brown       if (i < ctx.sm || i > ctx.ms-1)
1235c4762a1bSJed Brown         for (k=0; k<dof; k++) {
1236c4762a1bSJed Brown           mass_initial = mass_initial + hs*ptr_X0[i*dof+k];
1237c4762a1bSJed Brown           mass_final = mass_final + hs*ptr_X[i*dof+k];
1238c4762a1bSJed Brown         }
1239c4762a1bSJed Brown       else if (i < ctx.mf || i > ctx.fm-1)
1240c4762a1bSJed Brown         for (k=0; k<dof; k++) {
1241c4762a1bSJed Brown           mass_initial = mass_initial + hm*ptr_X0[i*dof+k];
1242c4762a1bSJed Brown           mass_final = mass_final + hm*ptr_X[i*dof+k];
1243c4762a1bSJed Brown         }
1244c4762a1bSJed Brown       else {
1245c4762a1bSJed Brown         for (k=0; k<dof; k++) {
1246c4762a1bSJed Brown           mass_initial = mass_initial + hf*ptr_X0[i*dof+k];
1247c4762a1bSJed Brown           mass_final = mass_final + hf*ptr_X[i*dof+k];
1248c4762a1bSJed Brown         }
1249c4762a1bSJed Brown       }
1250c4762a1bSJed Brown     }
12515f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArrayRead(da,X0,(void*)&ptr_X0));
12525f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArrayRead(da,X,(void*)&ptr_X));
1253c4762a1bSJed Brown     mass_difference = mass_final - mass_initial;
12545f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Allreduce(&mass_difference,&mass_differenceg,1,MPIU_SCALAR,MPIU_SUM,comm));
12555f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm,"Mass difference %g\n",(double)mass_differenceg));
12565f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm,"Final time %g, steps %D\n",(double)ptime,steps));
12575f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm,"Maximum allowable stepsize according to CFL %g\n",(double)(1.0/ctx.cfl_idt)));
1258c4762a1bSJed Brown     if (ctx.exact) {
1259c4762a1bSJed Brown       PetscReal nrm1=0;
12605f80ce2aSJacob Faibussowitsch       CHKERRQ(SolutionErrorNorms_3WaySplit(&ctx,da,ptime,X,&nrm1));
12615f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1));
1262c4762a1bSJed Brown     }
1263c4762a1bSJed Brown     if (ctx.simulation) {
1264c4762a1bSJed Brown       PetscReal    nrm1=0;
1265c4762a1bSJed Brown       PetscViewer  fd;
1266c4762a1bSJed Brown       char         filename[PETSC_MAX_PATH_LEN] = "binaryoutput";
1267c4762a1bSJed Brown       Vec          XR;
1268c4762a1bSJed Brown       PetscBool    flg;
1269c4762a1bSJed Brown       const PetscScalar  *ptr_XR;
12705f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg));
12713c633725SBarry Smith       PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
12725f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&fd));
12735f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDuplicate(X0,&XR));
12745f80ce2aSJacob Faibussowitsch       CHKERRQ(VecLoad(XR,fd));
12755f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerDestroy(&fd));
12765f80ce2aSJacob Faibussowitsch       CHKERRQ(VecGetArrayRead(X,&ptr_X));
12775f80ce2aSJacob Faibussowitsch       CHKERRQ(VecGetArrayRead(XR,&ptr_XR));
1278c4762a1bSJed Brown       for (i=xs;i<xs+xm;i++) {
1279c4762a1bSJed Brown         if (i < ctx.sm || i > ctx.ms-1)
1280c4762a1bSJed Brown           for (k=0; k<dof; k++) nrm1 = nrm1 + hs*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]);
1281c4762a1bSJed Brown         else if (i < ctx.mf || i < ctx.fm-1)
1282c4762a1bSJed Brown           for (k=0; k<dof; k++) nrm1 = nrm1 + hm*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]);
1283c4762a1bSJed Brown         else
1284c4762a1bSJed Brown           for (k=0; k<dof; k++) nrm1 = nrm1 + hf*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]);
1285c4762a1bSJed Brown       }
12865f80ce2aSJacob Faibussowitsch       CHKERRQ(VecRestoreArrayRead(X,&ptr_X));
12875f80ce2aSJacob Faibussowitsch       CHKERRQ(VecRestoreArrayRead(XR,&ptr_XR));
12885f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1));
12895f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&XR));
1290c4762a1bSJed Brown     }
1291c4762a1bSJed Brown   }
1292c4762a1bSJed Brown 
12935f80ce2aSJacob Faibussowitsch   CHKERRQ(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD));
12945f80ce2aSJacob Faibussowitsch   if (draw & 0x1) CHKERRQ(VecView(X0,PETSC_VIEWER_DRAW_WORLD));
12955f80ce2aSJacob Faibussowitsch   if (draw & 0x2) CHKERRQ(VecView(X,PETSC_VIEWER_DRAW_WORLD));
1296c4762a1bSJed Brown   if (draw & 0x4) {
1297c4762a1bSJed Brown     Vec Y;
12985f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(X,&Y));
12995f80ce2aSJacob Faibussowitsch     CHKERRQ(FVSample_3WaySplit(&ctx,da,ptime,Y));
13005f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAYPX(Y,-1,X));
13015f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(Y,PETSC_VIEWER_DRAW_WORLD));
13025f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&Y));
1303c4762a1bSJed Brown   }
1304c4762a1bSJed Brown 
1305c4762a1bSJed Brown   if (view_final) {
1306c4762a1bSJed Brown     PetscViewer viewer;
13075f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer));
13085f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
13095f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(X,viewer));
13105f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerPopFormat(viewer));
13115f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&viewer));
1312c4762a1bSJed Brown   }
1313c4762a1bSJed Brown 
1314c4762a1bSJed Brown   /* Clean up */
13155f80ce2aSJacob Faibussowitsch   CHKERRQ((*ctx.physics2.destroy)(ctx.physics2.user));
13165f80ce2aSJacob Faibussowitsch   for (i=0; i<ctx.physics2.dof; i++) CHKERRQ(PetscFree(ctx.physics2.fieldname[i]));
13175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree4(ctx.R,ctx.Rinv,ctx.cjmpLR,ctx.cslope));
13185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(ctx.uLR,ctx.flux,ctx.speeds));
13195f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X));
13205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X0));
13215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&R));
13225f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
13235f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
13245f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&ctx.iss));
13255f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&ctx.ism));
13265f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&ctx.isf));
13275f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&ctx.issb));
13285f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&ctx.ismb));
13295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(index_slow));
13305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(index_medium));
13315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(index_fast));
13325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(index_slowbuffer));
13335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(index_mediumbuffer));
13345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListDestroy(&limiters));
13355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListDestroy(&physics));
1336*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
1337*b122ec5aSJacob Faibussowitsch   return 0;
1338c4762a1bSJed Brown }
1339c4762a1bSJed Brown 
1340c4762a1bSJed Brown /*TEST
1341c4762a1bSJed Brown 
1342c4762a1bSJed Brown     build:
1343f56ea12dSJed Brown       requires: !complex
1344c4762a1bSJed Brown       depends: finitevolume1d.c
1345c4762a1bSJed Brown 
1346c4762a1bSJed Brown     test:
1347c4762a1bSJed Brown       suffix: 1
1348c4762a1bSJed 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
1349c4762a1bSJed Brown 
1350c4762a1bSJed Brown     test:
1351c4762a1bSJed Brown       suffix: 2
1352c4762a1bSJed 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
1353c4762a1bSJed Brown 
1354c4762a1bSJed Brown TEST*/
1355