xref: /petsc/src/ts/tutorials/multirate/ex6.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown /*
2c4762a1bSJed Brown   Note:
3c4762a1bSJed Brown     -hratio is the ratio between mesh size of corse grids and fine grids
4c4762a1bSJed Brown     -ts_rk_dtratio is the ratio between the large stepsize and the small stepsize
5c4762a1bSJed Brown */
6c4762a1bSJed Brown 
7c4762a1bSJed Brown static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n"
8c4762a1bSJed Brown   "  advection   - Constant coefficient scalar advection\n"
9c4762a1bSJed Brown   "                u_t       + (a*u)_x               = 0\n"
10c4762a1bSJed Brown   "  for this toy problem, we choose different meshsizes for different sub-domains (slow-fast-slow), say\n"
11c4762a1bSJed Brown   "                hxs  = hratio*hxf \n"
12c4762a1bSJed Brown   "  where hxs and hxf are the grid spacings for coarse and fine grids respectively.\n"
13c4762a1bSJed Brown   "  exact       - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n"
14c4762a1bSJed Brown   "                the states across shocks and rarefactions\n"
15c4762a1bSJed Brown   "  simulation  - use reference solution which is generated by smaller time step size to be true solution,\n"
16c4762a1bSJed Brown   "                also the reference solution should be generated by user and stored in a binary file.\n"
17c4762a1bSJed Brown   "  characteristic - Limit the characteristic variables, this is usually preferred (default)\n"
18c4762a1bSJed Brown   "Several initial conditions can be chosen with -initial N\n\n"
19c4762a1bSJed Brown   "The problem size should be set with -da_grid_x M\n\n";
20c4762a1bSJed Brown 
21c4762a1bSJed Brown #include <petscts.h>
22c4762a1bSJed Brown #include <petscdm.h>
23c4762a1bSJed Brown #include <petscdmda.h>
24c4762a1bSJed Brown #include <petscdraw.h>
25c4762a1bSJed Brown #include "finitevolume1d.h"
26c4762a1bSJed Brown 
279fbee547SJacob Faibussowitsch static inline PetscReal RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax) { PetscReal range = xmax-xmin; return xmin +PetscFmodReal(range+PetscFmodReal(a,range),range); }
28c4762a1bSJed Brown 
29c4762a1bSJed Brown /* --------------------------------- Advection ----------------------------------- */
30c4762a1bSJed Brown typedef struct {
31c4762a1bSJed Brown   PetscReal a;                  /* advective velocity */
32c4762a1bSJed Brown } AdvectCtx;
33c4762a1bSJed Brown 
34c4762a1bSJed Brown static PetscErrorCode PhysicsRiemann_Advect(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
35c4762a1bSJed Brown {
36c4762a1bSJed Brown   AdvectCtx *ctx = (AdvectCtx*)vctx;
37c4762a1bSJed Brown   PetscReal speed;
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   PetscFunctionBeginUser;
40c4762a1bSJed Brown   speed     = ctx->a;
41c4762a1bSJed Brown   flux[0]   = PetscMax(0,speed)*uL[0] + PetscMin(0,speed)*uR[0];
42c4762a1bSJed Brown   *maxspeed = speed;
43c4762a1bSJed Brown   PetscFunctionReturn(0);
44c4762a1bSJed Brown }
45c4762a1bSJed Brown 
46c4762a1bSJed Brown static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
47c4762a1bSJed Brown {
48c4762a1bSJed Brown   AdvectCtx *ctx = (AdvectCtx*)vctx;
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   PetscFunctionBeginUser;
51c4762a1bSJed Brown   X[0]      = 1.;
52c4762a1bSJed Brown   Xi[0]     = 1.;
53c4762a1bSJed Brown   speeds[0] = ctx->a;
54c4762a1bSJed Brown   PetscFunctionReturn(0);
55c4762a1bSJed Brown }
56c4762a1bSJed Brown 
57c4762a1bSJed Brown static PetscErrorCode PhysicsSample_Advect(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
58c4762a1bSJed Brown {
59c4762a1bSJed Brown   AdvectCtx *ctx = (AdvectCtx*)vctx;
60c4762a1bSJed Brown   PetscReal a    = ctx->a,x0;
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   PetscFunctionBeginUser;
63c4762a1bSJed Brown   switch (bctype) {
64c4762a1bSJed Brown     case FVBC_OUTFLOW:   x0 = x-a*t; break;
65c4762a1bSJed Brown     case FVBC_PERIODIC: x0 = RangeMod(x-a*t,xmin,xmax); break;
66c4762a1bSJed Brown     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown BCType");
67c4762a1bSJed Brown   }
68c4762a1bSJed Brown   switch (initial) {
69c4762a1bSJed Brown     case 0: u[0] = (x0 < 0) ? 1 : -1; break;
70c4762a1bSJed Brown     case 1: u[0] = (x0 < 0) ? -1 : 1; break;
71c4762a1bSJed Brown     case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break;
72c4762a1bSJed Brown     case 3: u[0] = PetscSinReal(2*PETSC_PI*x0); break;
73c4762a1bSJed Brown     case 4: u[0] = PetscAbs(x0); break;
74c4762a1bSJed Brown     case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2*PETSC_PI*x0)); break;
75c4762a1bSJed Brown     case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2-x0 : 0)); break;
76c4762a1bSJed Brown     case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI*x0),10.0);break;
77c4762a1bSJed Brown     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition");
78c4762a1bSJed Brown   }
79c4762a1bSJed Brown   PetscFunctionReturn(0);
80c4762a1bSJed Brown }
81c4762a1bSJed Brown 
82c4762a1bSJed Brown static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx)
83c4762a1bSJed Brown {
84c4762a1bSJed Brown   PetscErrorCode ierr;
85c4762a1bSJed Brown   AdvectCtx      *user;
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   PetscFunctionBeginUser;
885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNew(&user));
89c4762a1bSJed Brown   ctx->physics2.sample2         = PhysicsSample_Advect;
90c4762a1bSJed Brown   ctx->physics2.riemann2        = PhysicsRiemann_Advect;
91c4762a1bSJed Brown   ctx->physics2.characteristic2 = PhysicsCharacteristic_Advect;
92c4762a1bSJed Brown   ctx->physics2.destroy         = PhysicsDestroy_SimpleFree;
93c4762a1bSJed Brown   ctx->physics2.user            = user;
94c4762a1bSJed Brown   ctx->physics2.dof             = 1;
955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrallocpy("u",&ctx->physics2.fieldname[0]));
96c4762a1bSJed Brown   user->a = 1;
97c4762a1bSJed Brown   ierr = PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");CHKERRQ(ierr);
98c4762a1bSJed Brown   {
995f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-physics_advect_a","Speed","",user->a,&user->a,NULL));
100c4762a1bSJed Brown   }
101c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
102c4762a1bSJed Brown   PetscFunctionReturn(0);
103c4762a1bSJed Brown }
104c4762a1bSJed Brown 
105c4762a1bSJed Brown PetscErrorCode FVSample_2WaySplit(FVCtx *ctx,DM da,PetscReal time,Vec U)
106c4762a1bSJed Brown {
107c4762a1bSJed Brown   PetscScalar     *u,*uj,xj,xi;
108c4762a1bSJed Brown   PetscInt        i,j,k,dof,xs,xm,Mx;
109c4762a1bSJed Brown   const PetscInt  N = 200;
110c4762a1bSJed Brown   PetscReal       hs,hf;
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   PetscFunctionBeginUser;
1133c633725SBarry Smith   PetscCheck(ctx->physics2.sample2,PETSC_COMM_SELF,PETSC_ERR_SUP,"Physics has not provided a sampling function");
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,U,&u));
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(dof,&uj));
118c4762a1bSJed Brown   hs   = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
119c4762a1bSJed Brown   hf   = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
120c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
121c4762a1bSJed Brown     if (i < ctx->sf) {
122c4762a1bSJed Brown       xi = ctx->xmin+0.5*hs+i*hs;
123c4762a1bSJed Brown       /* Integrate over cell i using trapezoid rule with N points. */
124c4762a1bSJed Brown       for (k=0; k<dof; k++) u[i*dof+k] = 0;
125c4762a1bSJed Brown       for (j=0; j<N+1; j++) {
126c4762a1bSJed Brown         xj = xi+hs*(j-N/2)/(PetscReal)N;
1275f80ce2aSJacob Faibussowitsch         CHKERRQ((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj));
128c4762a1bSJed Brown         for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N;
129c4762a1bSJed Brown       }
130c4762a1bSJed Brown     } else if (i < ctx->fs) {
131c4762a1bSJed Brown       xi = ctx->xmin+ctx->sf*hs+0.5*hf+(i-ctx->sf)*hf;
132c4762a1bSJed Brown       /* Integrate over cell i using trapezoid rule with N points. */
133c4762a1bSJed Brown       for (k=0; k<dof; k++) u[i*dof+k] = 0;
134c4762a1bSJed Brown       for (j=0; j<N+1; j++) {
135c4762a1bSJed Brown         xj = xi+hf*(j-N/2)/(PetscReal)N;
1365f80ce2aSJacob Faibussowitsch         CHKERRQ((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj));
137c4762a1bSJed Brown         for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N;
138c4762a1bSJed Brown       }
139c4762a1bSJed Brown     } else {
140c4762a1bSJed Brown       xi = ctx->xmin+ctx->sf*hs+(ctx->fs-ctx->sf)*hf+0.5*hs+(i-ctx->fs)*hs;
141c4762a1bSJed Brown       /* Integrate over cell i using trapezoid rule with N points. */
142c4762a1bSJed Brown       for (k=0; k<dof; k++) u[i*dof+k] = 0;
143c4762a1bSJed Brown       for (j=0; j<N+1; j++) {
144c4762a1bSJed Brown         xj = xi+hs*(j-N/2)/(PetscReal)N;
1455f80ce2aSJacob Faibussowitsch         CHKERRQ((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj));
146c4762a1bSJed Brown         for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N;
147c4762a1bSJed Brown       }
148c4762a1bSJed Brown     }
149c4762a1bSJed Brown   }
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,U,&u));
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(uj));
152c4762a1bSJed Brown   PetscFunctionReturn(0);
153c4762a1bSJed Brown }
154c4762a1bSJed Brown 
155c4762a1bSJed Brown static PetscErrorCode SolutionErrorNorms_2WaySplit(FVCtx *ctx,DM da,PetscReal t,Vec X,PetscReal *nrm1)
156c4762a1bSJed Brown {
157c4762a1bSJed Brown   Vec               Y;
158c4762a1bSJed Brown   PetscInt          i,Mx;
159c4762a1bSJed Brown   const PetscScalar *ptr_X,*ptr_Y;
160c4762a1bSJed Brown   PetscReal         hs,hf;
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   PetscFunctionBeginUser;
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(X,&Mx));
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(X,&Y));
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(FVSample_2WaySplit(ctx,da,t,Y));
166c4762a1bSJed Brown   hs   = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
167c4762a1bSJed Brown   hf   = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
1685f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&ptr_X));
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Y,&ptr_Y));
170c4762a1bSJed Brown   for (i=0; i<Mx; i++) {
171c4762a1bSJed Brown     if (i < ctx->sf || i > ctx->fs -1) *nrm1 +=  hs*PetscAbs(ptr_X[i]-ptr_Y[i]);
172c4762a1bSJed Brown     else *nrm1 += hf*PetscAbs(ptr_X[i]-ptr_Y[i]);
173c4762a1bSJed Brown   }
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&ptr_X));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Y,&ptr_Y));
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Y));
177c4762a1bSJed Brown   PetscFunctionReturn(0);
178c4762a1bSJed Brown }
179c4762a1bSJed Brown 
180c4762a1bSJed Brown PetscErrorCode FVRHSFunction_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
181c4762a1bSJed Brown {
182c4762a1bSJed Brown   FVCtx          *ctx = (FVCtx*)vctx;
183c4762a1bSJed Brown   PetscInt       i,j,k,Mx,dof,xs,xm,sf = ctx->sf,fs = ctx->fs;
184c4762a1bSJed Brown   PetscReal      hxf,hxs,cfl_idt = 0;
185c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
186c4762a1bSJed Brown   Vec            Xloc;
187c4762a1bSJed Brown   DM             da;
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   PetscFunctionBeginUser;
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
1915f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&Xloc));                          /* Xloc contains ghost points                                     */
1925f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));   /* Mx is the number of center points                              */
193c4762a1bSJed Brown   hxs  = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
194c4762a1bSJed Brown   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
1955f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));       /* X is solution vector which does not contain ghost points       */
1965f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc));
197c4762a1bSJed Brown 
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(F));                                   /* F is the right hand side function corresponds to center points */
199c4762a1bSJed Brown 
2005f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,Xloc,&x));
2015f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,F,&f));
2025f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope));                  /* contains ghost points                                           */
203c4762a1bSJed Brown 
2045f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
207c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
208c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
209c4762a1bSJed Brown     }
210c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
211c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
212c4762a1bSJed Brown     }
213c4762a1bSJed Brown   }
214c4762a1bSJed Brown   for (i=xs-1; i<xs+xm+1; i++) {
215c4762a1bSJed Brown     struct _LimitInfo info;
216c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
217c4762a1bSJed Brown     /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
2185f80ce2aSJacob Faibussowitsch     CHKERRQ((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds));
219c4762a1bSJed Brown     /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
2205f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof));
221c4762a1bSJed Brown     cjmpL = &ctx->cjmpLR[0];
222c4762a1bSJed Brown     cjmpR = &ctx->cjmpLR[dof];
223c4762a1bSJed Brown     for (j=0; j<dof; j++) {
224c4762a1bSJed Brown       PetscScalar jmpL,jmpR;
225c4762a1bSJed Brown       jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
226c4762a1bSJed Brown       jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
227c4762a1bSJed Brown       for (k=0; k<dof; k++) {
228c4762a1bSJed Brown         cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
229c4762a1bSJed Brown         cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
230c4762a1bSJed Brown       }
231c4762a1bSJed Brown     }
232c4762a1bSJed Brown     /* Apply limiter to the left and right characteristic jumps */
233c4762a1bSJed Brown     info.m  = dof;
234c4762a1bSJed Brown     info.hxs = hxs;
235c4762a1bSJed Brown     info.hxf = hxf;
236c4762a1bSJed Brown     (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope);
237c4762a1bSJed Brown     for (j=0; j<dof; j++) {
238c4762a1bSJed Brown       PetscScalar tmp = 0;
239c4762a1bSJed Brown       for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
240c4762a1bSJed Brown       slope[i*dof+j] = tmp;
241c4762a1bSJed Brown     }
242c4762a1bSJed Brown   }
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
245c4762a1bSJed Brown     PetscReal   maxspeed;
246c4762a1bSJed Brown     PetscScalar *uL,*uR;
247c4762a1bSJed Brown     uL = &ctx->uLR[0];
248c4762a1bSJed Brown     uR = &ctx->uLR[dof];
249c4762a1bSJed Brown     if (i < sf) { /* slow region */
250c4762a1bSJed Brown       for (j=0; j<dof; j++) {
251c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
252c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
253c4762a1bSJed Brown       }
2545f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
255c4762a1bSJed Brown       if (i > xs) {
256c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs;
257c4762a1bSJed Brown       }
258c4762a1bSJed Brown       if (i < xs+xm) {
259c4762a1bSJed Brown         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs;
260c4762a1bSJed Brown       }
261c4762a1bSJed Brown     } else if (i == sf) { /* interface between the slow region and the fast region */
262c4762a1bSJed Brown       for (j=0; j<dof; j++) {
263c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
264c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
265c4762a1bSJed Brown       }
2665f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
267c4762a1bSJed Brown       if (i > xs) {
268c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs;
269c4762a1bSJed Brown       }
270c4762a1bSJed Brown       if (i < xs+xm) {
271c4762a1bSJed Brown         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxf;
272c4762a1bSJed Brown       }
273c4762a1bSJed Brown     } else if (i > sf && i < fs) { /* fast region */
274c4762a1bSJed Brown       for (j=0; j<dof; j++) {
275c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
276c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
277c4762a1bSJed Brown       }
2785f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
279c4762a1bSJed Brown       if (i > xs) {
280c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf;
281c4762a1bSJed Brown       }
282c4762a1bSJed Brown       if (i < xs+xm) {
283c4762a1bSJed Brown         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxf;
284c4762a1bSJed Brown       }
285c4762a1bSJed Brown     } else if (i == fs) { /* interface between the fast region and the slow region */
286c4762a1bSJed Brown       for (j=0; j<dof; j++) {
287c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
288c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
289c4762a1bSJed Brown       }
2905f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
291c4762a1bSJed Brown       if (i > xs) {
292c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf;
293c4762a1bSJed Brown       }
294c4762a1bSJed Brown       if (i < xs+xm) {
295c4762a1bSJed Brown         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs;
296c4762a1bSJed Brown       }
297c4762a1bSJed Brown     } else { /* slow region */
298c4762a1bSJed Brown       for (j=0; j<dof; j++) {
299c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
300c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
301c4762a1bSJed Brown       }
3025f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
303c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */
304c4762a1bSJed Brown       if (i > xs) {
305c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs;
306c4762a1bSJed Brown       }
307c4762a1bSJed Brown       if (i < xs+xm) {
308c4762a1bSJed Brown         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs;
309c4762a1bSJed Brown       }
310c4762a1bSJed Brown     }
311c4762a1bSJed Brown   }
3125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x));
3135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,F,&f));
3145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope));
3155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&Xloc));
3165f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da)));
317c4762a1bSJed Brown   if (0) {
318c4762a1bSJed Brown     /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */
319c4762a1bSJed Brown     PetscReal dt,tnow;
3205f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetTimeStep(ts,&dt));
3215f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetTime(ts,&tnow));
322c4762a1bSJed Brown     if (dt > 0.5/ctx->cfl_idt) {
323c4762a1bSJed Brown       if (1) {
3245f80ce2aSJacob 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)));
32598921bdaSJacob Faibussowitsch       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Stability constraint exceeded, %g > %g",(double)dt,(double)(ctx->cfl/ctx->cfl_idt));
326c4762a1bSJed Brown     }
327c4762a1bSJed Brown   }
328c4762a1bSJed Brown   PetscFunctionReturn(0);
329c4762a1bSJed Brown }
330c4762a1bSJed Brown 
331c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */
332c4762a1bSJed Brown PetscErrorCode FVRHSFunctionslow_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
333c4762a1bSJed Brown {
334c4762a1bSJed Brown   FVCtx          *ctx = (FVCtx*)vctx;
335c4762a1bSJed Brown   PetscInt       i,j,k,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth;
336c4762a1bSJed Brown   PetscReal      hxs,hxf,cfl_idt = 0;
337c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
338c4762a1bSJed Brown   Vec            Xloc;
339c4762a1bSJed Brown   DM             da;
340c4762a1bSJed Brown 
341c4762a1bSJed Brown   PetscFunctionBeginUser;
3425f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
3435f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&Xloc));
3445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
345c4762a1bSJed Brown   hxs  = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
346c4762a1bSJed Brown   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
3475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));
3485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc));
3495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(F));
3505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,Xloc,&x));
3515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
3525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope));
3535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
354c4762a1bSJed Brown 
355c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
356c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
357c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
358c4762a1bSJed Brown     }
359c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
360c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
361c4762a1bSJed Brown     }
362c4762a1bSJed Brown   }
363c4762a1bSJed Brown   for (i=xs-1; i<xs+xm+1; i++) {
364c4762a1bSJed Brown     struct _LimitInfo info;
365c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
366c4762a1bSJed Brown     if (i < sf-lsbwidth+1 || i > fs+rsbwidth-2) { /* slow components and the first and last fast components */
367c4762a1bSJed Brown       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
3685f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds));
369c4762a1bSJed Brown       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
3705f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof));
371c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
372c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
373c4762a1bSJed Brown       for (j=0; j<dof; j++) {
374c4762a1bSJed Brown         PetscScalar jmpL,jmpR;
375c4762a1bSJed Brown         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
376c4762a1bSJed Brown         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
377c4762a1bSJed Brown         for (k=0; k<dof; k++) {
378c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
379c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
380c4762a1bSJed Brown         }
381c4762a1bSJed Brown       }
382c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
383c4762a1bSJed Brown       info.m  = dof;
384c4762a1bSJed Brown       info.hxs = hxs;
385c4762a1bSJed Brown       info.hxf = hxf;
386c4762a1bSJed Brown       (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope);
387c4762a1bSJed Brown       for (j=0; j<dof; j++) {
388c4762a1bSJed Brown         PetscScalar tmp = 0;
389c4762a1bSJed Brown         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
390c4762a1bSJed Brown           slope[i*dof+j] = tmp;
391c4762a1bSJed Brown       }
392c4762a1bSJed Brown     }
393c4762a1bSJed Brown   }
394c4762a1bSJed Brown 
395c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
396c4762a1bSJed Brown     PetscReal   maxspeed;
397c4762a1bSJed Brown     PetscScalar *uL,*uR;
398c4762a1bSJed Brown     uL = &ctx->uLR[0];
399c4762a1bSJed Brown     uR = &ctx->uLR[dof];
400c4762a1bSJed Brown     if (i < sf-lsbwidth) { /* slow region */
401c4762a1bSJed Brown       for (j=0; j<dof; j++) {
402c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
403c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
404c4762a1bSJed Brown       }
4055f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
406c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */
407c4762a1bSJed Brown       if (i > xs) {
408c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
409c4762a1bSJed Brown       }
410c4762a1bSJed Brown       if (i < xs+xm) {
411c4762a1bSJed Brown         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
412c4762a1bSJed Brown         islow++;
413c4762a1bSJed Brown       }
414c4762a1bSJed Brown     }
415c4762a1bSJed Brown     if (i == sf-lsbwidth) { /* interface between the slow region and the fast region */
416c4762a1bSJed Brown       for (j=0; j<dof; j++) {
417c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
418c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
419c4762a1bSJed Brown       }
4205f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
421c4762a1bSJed Brown       if (i > xs) {
422c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
423c4762a1bSJed Brown       }
424c4762a1bSJed Brown     }
425c4762a1bSJed Brown     if (i == fs+rsbwidth) { /* slow region */
426c4762a1bSJed Brown       for (j=0; j<dof; j++) {
427c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
428c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
429c4762a1bSJed Brown       }
4305f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
431c4762a1bSJed Brown       if (i < xs+xm) {
432c4762a1bSJed Brown         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
433c4762a1bSJed Brown         islow++;
434c4762a1bSJed Brown       }
435c4762a1bSJed Brown     }
436c4762a1bSJed Brown     if (i > fs+rsbwidth) { /* slow region */
437c4762a1bSJed Brown       for (j=0; j<dof; j++) {
438c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
439c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
440c4762a1bSJed Brown       }
4415f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
442c4762a1bSJed Brown       if (i > xs) {
443c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
444c4762a1bSJed Brown       }
445c4762a1bSJed Brown       if (i < xs+xm) {
446c4762a1bSJed Brown         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
447c4762a1bSJed Brown         islow++;
448c4762a1bSJed Brown       }
449c4762a1bSJed Brown     }
450c4762a1bSJed Brown   }
4515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x));
4525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
4535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope));
4545f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&Xloc));
4555f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da)));
456c4762a1bSJed Brown   PetscFunctionReturn(0);
457c4762a1bSJed Brown }
458c4762a1bSJed Brown 
459c4762a1bSJed Brown PetscErrorCode FVRHSFunctionslowbuffer_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
460c4762a1bSJed Brown {
461c4762a1bSJed Brown   FVCtx          *ctx = (FVCtx*)vctx;
462c4762a1bSJed Brown   PetscInt       i,j,k,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth;
463c4762a1bSJed Brown   PetscReal      hxs,hxf;
464c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
465c4762a1bSJed Brown   Vec            Xloc;
466c4762a1bSJed Brown   DM             da;
467c4762a1bSJed Brown 
468c4762a1bSJed Brown   PetscFunctionBeginUser;
4695f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
4705f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&Xloc));
4715f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
472c4762a1bSJed Brown   hxs  = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
473c4762a1bSJed Brown   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
4745f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));
4755f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc));
4765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(F));
4775f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,Xloc,&x));
4785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
4795f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope));
4805f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
481c4762a1bSJed Brown 
482c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
483c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
484c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
485c4762a1bSJed Brown     }
486c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
487c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
488c4762a1bSJed Brown     }
489c4762a1bSJed Brown   }
490c4762a1bSJed Brown   for (i=xs-1; i<xs+xm+1; i++) {
491c4762a1bSJed Brown     struct _LimitInfo info;
492c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
493c4762a1bSJed Brown     if ((i > sf-lsbwidth-2 && i < sf+1) || (i > fs-2 && i < fs+rsbwidth+1)) {
494c4762a1bSJed Brown       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
4955f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds));
496c4762a1bSJed Brown       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
4975f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof));
498c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
499c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
500c4762a1bSJed Brown       for (j=0; j<dof; j++) {
501c4762a1bSJed Brown         PetscScalar jmpL,jmpR;
502c4762a1bSJed Brown         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
503c4762a1bSJed Brown         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
504c4762a1bSJed Brown         for (k=0; k<dof; k++) {
505c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
506c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
507c4762a1bSJed Brown         }
508c4762a1bSJed Brown       }
509c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
510c4762a1bSJed Brown       info.m  = dof;
511c4762a1bSJed Brown       info.hxs = hxs;
512c4762a1bSJed Brown       info.hxf = hxf;
513c4762a1bSJed Brown       (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope);
514c4762a1bSJed Brown       for (j=0; j<dof; j++) {
515c4762a1bSJed Brown         PetscScalar tmp = 0;
516c4762a1bSJed Brown         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
517c4762a1bSJed Brown           slope[i*dof+j] = tmp;
518c4762a1bSJed Brown       }
519c4762a1bSJed Brown     }
520c4762a1bSJed Brown   }
521c4762a1bSJed Brown 
522c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
523c4762a1bSJed Brown     PetscReal   maxspeed;
524c4762a1bSJed Brown     PetscScalar *uL,*uR;
525c4762a1bSJed Brown     uL = &ctx->uLR[0];
526c4762a1bSJed Brown     uR = &ctx->uLR[dof];
527c4762a1bSJed Brown     if (i == sf-lsbwidth) {
528c4762a1bSJed Brown       for (j=0; j<dof; j++) {
529c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
530c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
531c4762a1bSJed Brown       }
5325f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
533c4762a1bSJed Brown       if (i < xs+xm) {
534c4762a1bSJed Brown         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
535c4762a1bSJed Brown         islow++;
536c4762a1bSJed Brown       }
537c4762a1bSJed Brown     }
538c4762a1bSJed Brown     if (i > sf-lsbwidth && i < sf) {
539c4762a1bSJed Brown       for (j=0; j<dof; j++) {
540c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
541c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
542c4762a1bSJed Brown       }
5435f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
544c4762a1bSJed Brown       if (i > xs) {
545c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
546c4762a1bSJed Brown       }
547c4762a1bSJed Brown       if (i < xs+xm) {
548c4762a1bSJed Brown         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
549c4762a1bSJed Brown         islow++;
550c4762a1bSJed Brown       }
551c4762a1bSJed Brown     }
552c4762a1bSJed Brown     if (i == sf) { /* interface between the slow region and the fast region */
553c4762a1bSJed Brown       for (j=0; j<dof; j++) {
554c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
555c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
556c4762a1bSJed Brown       }
5575f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
558c4762a1bSJed Brown       if (i > xs) {
559c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
560c4762a1bSJed Brown       }
561c4762a1bSJed Brown     }
562c4762a1bSJed Brown     if (i == fs) { /* interface between the fast region and the slow region */
563c4762a1bSJed Brown       for (j=0; j<dof; j++) {
564c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
565c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
566c4762a1bSJed Brown       }
5675f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
568c4762a1bSJed Brown       if (i < xs+xm) {
569c4762a1bSJed Brown         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
570c4762a1bSJed Brown         islow++;
571c4762a1bSJed Brown       }
572c4762a1bSJed Brown     }
573c4762a1bSJed Brown     if (i > fs && i < fs+rsbwidth) {
574c4762a1bSJed Brown       for (j=0; j<dof; j++) {
575c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
576c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
577c4762a1bSJed Brown       }
5785f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
579c4762a1bSJed Brown       if (i > xs) {
580c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
581c4762a1bSJed Brown       }
582c4762a1bSJed Brown       if (i < xs+xm) {
583c4762a1bSJed Brown         for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs;
584c4762a1bSJed Brown         islow++;
585c4762a1bSJed Brown       }
586c4762a1bSJed Brown     }
587c4762a1bSJed Brown     if (i == fs+rsbwidth) {
588c4762a1bSJed Brown       for (j=0; j<dof; j++) {
589c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
590c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
591c4762a1bSJed Brown       }
5925f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
593c4762a1bSJed Brown       if (i > xs) {
594c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs;
595c4762a1bSJed Brown       }
596c4762a1bSJed Brown     }
597c4762a1bSJed Brown   }
5985f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x));
5995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
6005f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope));
6015f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&Xloc));
602c4762a1bSJed Brown   PetscFunctionReturn(0);
603c4762a1bSJed Brown }
604c4762a1bSJed Brown 
605c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for fast  parts ----------------------------------- */
606c4762a1bSJed Brown PetscErrorCode FVRHSFunctionfast_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
607c4762a1bSJed Brown {
608c4762a1bSJed Brown   FVCtx          *ctx = (FVCtx*)vctx;
609c4762a1bSJed Brown   PetscInt       i,j,k,Mx,dof,xs,xm,ifast = 0,sf = ctx->sf,fs = ctx->fs;
610c4762a1bSJed Brown   PetscReal      hxs,hxf;
611c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
612c4762a1bSJed Brown   Vec            Xloc;
613c4762a1bSJed Brown   DM             da;
614c4762a1bSJed Brown 
615c4762a1bSJed Brown   PetscFunctionBeginUser;
6165f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
6175f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&Xloc));
6185f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
619c4762a1bSJed Brown   hxs  = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf;
620c4762a1bSJed Brown   hxf  = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf);
6215f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));
6225f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc));
6235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(F));
6245f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,Xloc,&x));
6255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
6265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope));
6275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
628c4762a1bSJed Brown 
629c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
630c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
631c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
632c4762a1bSJed Brown     }
633c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
634c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
635c4762a1bSJed Brown     }
636c4762a1bSJed Brown   }
637c4762a1bSJed Brown   for (i=xs-1; i<xs+xm+1; i++) {
638c4762a1bSJed Brown     struct _LimitInfo info;
639c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
640c4762a1bSJed Brown     if (i > sf-2 && i < fs+1) {
6415f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds));
6425f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof));
643c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
644c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
645c4762a1bSJed Brown       for (j=0; j<dof; j++) {
646c4762a1bSJed Brown         PetscScalar jmpL,jmpR;
647c4762a1bSJed Brown         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
648c4762a1bSJed Brown         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
649c4762a1bSJed Brown         for (k=0; k<dof; k++) {
650c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
651c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
652c4762a1bSJed Brown         }
653c4762a1bSJed Brown       }
654c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
655c4762a1bSJed Brown       info.m  = dof;
656c4762a1bSJed Brown       info.hxs = hxs;
657c4762a1bSJed Brown       info.hxf = hxf;
658c4762a1bSJed Brown       (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope);
659c4762a1bSJed Brown       for (j=0; j<dof; j++) {
660c4762a1bSJed Brown       PetscScalar tmp = 0;
661c4762a1bSJed Brown       for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
662c4762a1bSJed Brown         slope[i*dof+j] = tmp;
663c4762a1bSJed Brown       }
664c4762a1bSJed Brown     }
665c4762a1bSJed Brown   }
666c4762a1bSJed Brown 
667c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
668c4762a1bSJed Brown     PetscReal   maxspeed;
669c4762a1bSJed Brown     PetscScalar *uL,*uR;
670c4762a1bSJed Brown     uL = &ctx->uLR[0];
671c4762a1bSJed Brown     uR = &ctx->uLR[dof];
672c4762a1bSJed Brown     if (i == sf) { /* interface between the slow region and the fast region */
673c4762a1bSJed Brown       for (j=0; j<dof; j++) {
674c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2;
675c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
676c4762a1bSJed Brown       }
6775f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
678c4762a1bSJed Brown       if (i < xs+xm) {
679c4762a1bSJed Brown         for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf;
680c4762a1bSJed Brown         ifast++;
681c4762a1bSJed Brown       }
682c4762a1bSJed Brown     }
683c4762a1bSJed Brown     if (i > sf && i < fs) { /* fast region */
684c4762a1bSJed Brown       for (j=0; j<dof; j++) {
685c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
686c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2;
687c4762a1bSJed Brown       }
6885f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
689c4762a1bSJed Brown       if (i > xs) {
690c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf;
691c4762a1bSJed Brown       }
692c4762a1bSJed Brown       if (i < xs+xm) {
693c4762a1bSJed Brown         for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf;
694c4762a1bSJed Brown         ifast++;
695c4762a1bSJed Brown       }
696c4762a1bSJed Brown     }
697c4762a1bSJed Brown     if (i == fs) { /* interface between the fast region and the slow region */
698c4762a1bSJed Brown       for (j=0; j<dof; j++) {
699c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2;
700c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2;
701c4762a1bSJed Brown       }
7025f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed));
703c4762a1bSJed Brown       if (i > xs) {
704c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf;
705c4762a1bSJed Brown       }
706c4762a1bSJed Brown     }
707c4762a1bSJed Brown   }
7085f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x));
7095f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
7105f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope));
7115f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&Xloc));
712c4762a1bSJed Brown   PetscFunctionReturn(0);
713c4762a1bSJed Brown }
714c4762a1bSJed Brown 
715c4762a1bSJed Brown int main(int argc,char *argv[])
716c4762a1bSJed Brown {
717c4762a1bSJed Brown   char              lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m";
718c4762a1bSJed Brown   PetscFunctionList limiters   = 0,physics = 0;
719c4762a1bSJed Brown   MPI_Comm          comm;
720c4762a1bSJed Brown   TS                ts;
721c4762a1bSJed Brown   DM                da;
722c4762a1bSJed Brown   Vec               X,X0,R;
723c4762a1bSJed Brown   FVCtx             ctx;
724c4762a1bSJed Brown   PetscInt          i,k,dof,xs,xm,Mx,draw = 0,count_slow,count_fast,islow = 0,ifast =0,islowbuffer = 0,*index_slow,*index_fast,*index_slowbuffer;
725c4762a1bSJed Brown   PetscBool         view_final = PETSC_FALSE;
726c4762a1bSJed Brown   PetscReal         ptime;
727c4762a1bSJed Brown   PetscErrorCode    ierr;
728c4762a1bSJed Brown 
729*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,0,help));
730c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
7315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMemzero(&ctx,sizeof(ctx)));
732c4762a1bSJed Brown 
733c4762a1bSJed Brown   /* Register limiters to be available on the command line */
7345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"upwind"              ,Limit2_Upwind));
7355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"lax-wendroff"        ,Limit2_LaxWendroff));
7365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"beam-warming"        ,Limit2_BeamWarming));
7375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"fromm"               ,Limit2_Fromm));
7385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"minmod"              ,Limit2_Minmod));
7395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"superbee"            ,Limit2_Superbee));
7405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"mc"                  ,Limit2_MC));
7415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"koren3"              ,Limit2_Koren3));
742c4762a1bSJed Brown 
743c4762a1bSJed Brown   /* Register physical models to be available on the command line */
7445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&physics,"advect"          ,PhysicsCreate_Advect));
745c4762a1bSJed Brown 
746c4762a1bSJed Brown   ctx.comm = comm;
747c4762a1bSJed Brown   ctx.cfl  = 0.9;
748c4762a1bSJed Brown   ctx.bctype = FVBC_PERIODIC;
749c4762a1bSJed Brown   ctx.xmin = -1.0;
750c4762a1bSJed Brown   ctx.xmax = 1.0;
751c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm,NULL,"Finite Volume solver options","");CHKERRQ(ierr);
7525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL));
7535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL));
7545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsFList("-limit","Name of flux imiter to use","",limiters,lname,lname,sizeof(lname),NULL));
7555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL));
7565f80ce2aSJacob 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));
7575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL));
7585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-exact","Compare errors with exact solution","",ctx.exact,&ctx.exact,NULL));
7595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL));
7605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL));
7615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL));
7625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-hratio","Spacing ratio","",ctx.hratio,&ctx.hratio,NULL));
763c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
764c4762a1bSJed Brown 
765c4762a1bSJed Brown   /* Choose the limiter from the list of registered limiters */
7665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListFind(limiters,lname,&ctx.limit2));
7673c633725SBarry Smith   PetscCheck(ctx.limit2,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Limiter '%s' not found",lname);
768c4762a1bSJed Brown 
769c4762a1bSJed Brown   /* Choose the physics from the list of registered models */
770c4762a1bSJed Brown   {
771c4762a1bSJed Brown     PetscErrorCode (*r)(FVCtx*);
7725f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFunctionListFind(physics,physname,&r));
7733c633725SBarry Smith     PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Physics '%s' not found",physname);
774c4762a1bSJed Brown     /* Create the physics, will set the number of fields and their names */
7755f80ce2aSJacob Faibussowitsch     CHKERRQ((*r)(&ctx));
776c4762a1bSJed Brown   }
777c4762a1bSJed Brown 
778c4762a1bSJed Brown   /* Create a DMDA to manage the parallel grid */
7795f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics2.dof,2,NULL,&da));
7805f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
7815f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
782c4762a1bSJed Brown   /* Inform the DMDA of the field names provided by the physics. */
783c4762a1bSJed Brown   /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
784c4762a1bSJed Brown   for (i=0; i<ctx.physics2.dof; i++) {
7855f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDASetFieldName(da,i,ctx.physics2.fieldname[i]));
786c4762a1bSJed Brown   }
7875f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
7885f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
789c4762a1bSJed Brown 
790c4762a1bSJed Brown   /* Set coordinates of cell centers */
7915f80ce2aSJacob 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));
792c4762a1bSJed Brown 
793c4762a1bSJed Brown   /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
7945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc4(dof*dof,&ctx.R,dof*dof,&ctx.Rinv,2*dof,&ctx.cjmpLR,1*dof,&ctx.cslope));
7955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3(2*dof,&ctx.uLR,dof,&ctx.flux,dof,&ctx.speeds));
796c4762a1bSJed Brown 
797c4762a1bSJed Brown   /* Create a vector to store the solution and to save the initial state */
7985f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&X));
7995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(X,&X0));
8005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(X,&R));
801c4762a1bSJed Brown 
802c4762a1bSJed Brown   /* create index for slow parts and fast parts,
803c4762a1bSJed Brown      count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */
804c4762a1bSJed Brown   count_slow = Mx/(1.0+ctx.hratio/3.0);
8052c71b3e2SJacob Faibussowitsch   PetscCheckFalse(count_slow%2,PETSC_COMM_WORLD,PETSC_ERR_USER,"Please adjust grid size Mx (-da_grid_x) and hratio (-hratio) so that Mx/(1+hartio/3) is even");
806c4762a1bSJed Brown   count_fast = Mx-count_slow;
807c4762a1bSJed Brown   ctx.sf = count_slow/2;
808c4762a1bSJed Brown   ctx.fs = ctx.sf+count_fast;
8095f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(xm*dof,&index_slow));
8105f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(xm*dof,&index_fast));
8115f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(6*dof,&index_slowbuffer));
812c4762a1bSJed Brown   if (((AdvectCtx*)ctx.physics2.user)->a > 0) {
813c4762a1bSJed Brown     ctx.lsbwidth = 2;
814c4762a1bSJed Brown     ctx.rsbwidth = 4;
815c4762a1bSJed Brown   } else {
816c4762a1bSJed Brown     ctx.lsbwidth = 4;
817c4762a1bSJed Brown     ctx.rsbwidth = 2;
818c4762a1bSJed Brown   }
819c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
820c4762a1bSJed Brown     if (i < ctx.sf-ctx.lsbwidth || i > ctx.fs+ctx.rsbwidth-1)
821c4762a1bSJed Brown       for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k;
822c4762a1bSJed Brown     else if ((i >= ctx.sf-ctx.lsbwidth && i < ctx.sf) || (i > ctx.fs-1 && i <= ctx.fs+ctx.rsbwidth-1))
823c4762a1bSJed Brown       for (k=0; k<dof; k++) index_slowbuffer[islowbuffer++] = i*dof+k;
824c4762a1bSJed Brown     else
825c4762a1bSJed Brown       for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k;
826c4762a1bSJed Brown   }
8275f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss));
8285f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf));
8295f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,islowbuffer,index_slowbuffer,PETSC_COPY_VALUES,&ctx.issb));
830c4762a1bSJed Brown 
831c4762a1bSJed Brown   /* Create a time-stepping object */
8325f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(comm,&ts));
8335f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts,da));
8345f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,R,FVRHSFunction_2WaySplit,&ctx));
8355f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetIS(ts,"slow",ctx.iss));
8365f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetIS(ts,"slowbuffer",ctx.issb));
8375f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetIS(ts,"fast",ctx.isf));
8385f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetRHSFunction(ts,"slow",NULL,FVRHSFunctionslow_2WaySplit,&ctx));
8395f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetRHSFunction(ts,"fast",NULL,FVRHSFunctionfast_2WaySplit,&ctx));
8405f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetRHSFunction(ts,"slowbuffer",NULL,FVRHSFunctionslowbuffer_2WaySplit,&ctx));
841c4762a1bSJed Brown 
8425f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSSSP));
8435f80ce2aSJacob Faibussowitsch   /*CHKERRQ(TSSetType(ts,TSMPRK));*/
8445f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,10));
8455f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
846c4762a1bSJed Brown 
847c4762a1bSJed Brown   /* Compute initial conditions and starting time step */
8485f80ce2aSJacob Faibussowitsch   CHKERRQ(FVSample_2WaySplit(&ctx,da,0,X0));
8495f80ce2aSJacob Faibussowitsch   CHKERRQ(FVRHSFunction_2WaySplit(ts,0,X0,X,(void*)&ctx)); /* Initial function evaluation, only used to determine max speed */
8505f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(X0,X));                        /* The function value was not used so we set X=X0 again */
8515f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt));
8525f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts)); /* Take runtime options */
8535f80ce2aSJacob Faibussowitsch   CHKERRQ(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD));
854c4762a1bSJed Brown   {
855c4762a1bSJed Brown     PetscInt          steps;
856c4762a1bSJed Brown     PetscScalar       mass_initial,mass_final,mass_difference,mass_differenceg;
857c4762a1bSJed Brown     const PetscScalar *ptr_X,*ptr_X0;
858c4762a1bSJed Brown     const PetscReal   hs = (ctx.xmax-ctx.xmin)*3.0/4.0/count_slow;
859c4762a1bSJed Brown     const PetscReal   hf = (ctx.xmax-ctx.xmin)/4.0/count_fast;
860c4762a1bSJed Brown 
8615f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSolve(ts,X));
8625f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetSolveTime(ts,&ptime));
8635f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetStepNumber(ts,&steps));
864c4762a1bSJed Brown     /* calculate the total mass at initial time and final time */
865c4762a1bSJed Brown     mass_initial = 0.0;
866c4762a1bSJed Brown     mass_final   = 0.0;
8675f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArrayRead(da,X0,(void*)&ptr_X0));
8685f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArrayRead(da,X,(void*)&ptr_X));
869c4762a1bSJed Brown     for (i=xs;i<xs+xm;i++) {
870c4762a1bSJed Brown       if (i < ctx.sf || i > ctx.fs-1) {
871c4762a1bSJed Brown         for (k=0; k<dof; k++) {
872c4762a1bSJed Brown           mass_initial = mass_initial + hs*ptr_X0[i*dof+k];
873c4762a1bSJed Brown           mass_final = mass_final + hs*ptr_X[i*dof+k];
874c4762a1bSJed Brown         }
875c4762a1bSJed Brown       } else {
876c4762a1bSJed Brown         for (k=0; k<dof; k++) {
877c4762a1bSJed Brown           mass_initial = mass_initial + hf*ptr_X0[i*dof+k];
878c4762a1bSJed Brown           mass_final = mass_final + hf*ptr_X[i*dof+k];
879c4762a1bSJed Brown         }
880c4762a1bSJed Brown       }
881c4762a1bSJed Brown     }
8825f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArrayRead(da,X0,(void*)&ptr_X0));
8835f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArrayRead(da,X,(void*)&ptr_X));
884c4762a1bSJed Brown     mass_difference = mass_final - mass_initial;
8855f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Allreduce(&mass_difference,&mass_differenceg,1,MPIU_SCALAR,MPIU_SUM,comm));
8865f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm,"Mass difference %g\n",(double)mass_differenceg));
8875f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm,"Final time %g, steps %D\n",(double)ptime,steps));
8885f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm,"Maximum allowable stepsize according to CFL %g\n",(double)(1.0/ctx.cfl_idt)));
889c4762a1bSJed Brown     if (ctx.exact) {
890c4762a1bSJed Brown       PetscReal nrm1=0;
8915f80ce2aSJacob Faibussowitsch       CHKERRQ(SolutionErrorNorms_2WaySplit(&ctx,da,ptime,X,&nrm1));
8925f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1));
893c4762a1bSJed Brown     }
894c4762a1bSJed Brown     if (ctx.simulation) {
895c4762a1bSJed Brown       PetscReal    nrm1=0;
896c4762a1bSJed Brown       PetscViewer  fd;
897c4762a1bSJed Brown       char         filename[PETSC_MAX_PATH_LEN] = "binaryoutput";
898c4762a1bSJed Brown       Vec          XR;
899c4762a1bSJed Brown       PetscBool    flg;
900c4762a1bSJed Brown       const PetscScalar  *ptr_XR;
9015f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg));
9023c633725SBarry Smith       PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
9035f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&fd));
9045f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDuplicate(X0,&XR));
9055f80ce2aSJacob Faibussowitsch       CHKERRQ(VecLoad(XR,fd));
9065f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerDestroy(&fd));
9075f80ce2aSJacob Faibussowitsch       CHKERRQ(VecGetArrayRead(X,&ptr_X));
9085f80ce2aSJacob Faibussowitsch       CHKERRQ(VecGetArrayRead(XR,&ptr_XR));
909c4762a1bSJed Brown       for (i=xs;i<xs+xm;i++) {
910c4762a1bSJed Brown         if (i < ctx.sf || i > ctx.fs-1)
911c4762a1bSJed Brown           for (k=0; k<dof; k++) nrm1 = nrm1 + hs*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]);
912c4762a1bSJed Brown         else
913c4762a1bSJed Brown           for (k=0; k<dof; k++) nrm1 = nrm1 + hf*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]);
914c4762a1bSJed Brown       }
9155f80ce2aSJacob Faibussowitsch       CHKERRQ(VecRestoreArrayRead(X,&ptr_X));
9165f80ce2aSJacob Faibussowitsch       CHKERRQ(VecRestoreArrayRead(XR,&ptr_XR));
9175f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1));
9185f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&XR));
919c4762a1bSJed Brown     }
920c4762a1bSJed Brown   }
921c4762a1bSJed Brown 
9225f80ce2aSJacob Faibussowitsch   CHKERRQ(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD));
9235f80ce2aSJacob Faibussowitsch   if (draw & 0x1) CHKERRQ(VecView(X0,PETSC_VIEWER_DRAW_WORLD));
9245f80ce2aSJacob Faibussowitsch   if (draw & 0x2) CHKERRQ(VecView(X,PETSC_VIEWER_DRAW_WORLD));
925c4762a1bSJed Brown   if (draw & 0x4) {
926c4762a1bSJed Brown     Vec Y;
9275f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(X,&Y));
9285f80ce2aSJacob Faibussowitsch     CHKERRQ(FVSample_2WaySplit(&ctx,da,ptime,Y));
9295f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAYPX(Y,-1,X));
9305f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(Y,PETSC_VIEWER_DRAW_WORLD));
9315f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&Y));
932c4762a1bSJed Brown   }
933c4762a1bSJed Brown 
934c4762a1bSJed Brown   if (view_final) {
935c4762a1bSJed Brown     PetscViewer viewer;
9365f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer));
9375f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
9385f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(X,viewer));
9395f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerPopFormat(viewer));
9405f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&viewer));
941c4762a1bSJed Brown   }
942c4762a1bSJed Brown 
943c4762a1bSJed Brown   /* Clean up */
9445f80ce2aSJacob Faibussowitsch   CHKERRQ((*ctx.physics2.destroy)(ctx.physics2.user));
9455f80ce2aSJacob Faibussowitsch   for (i=0; i<ctx.physics2.dof; i++) CHKERRQ(PetscFree(ctx.physics2.fieldname[i]));
9465f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree4(ctx.R,ctx.Rinv,ctx.cjmpLR,ctx.cslope));
9475f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(ctx.uLR,ctx.flux,ctx.speeds));
9485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X));
9495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X0));
9505f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&R));
9515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
9525f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
9535f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&ctx.iss));
9545f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&ctx.isf));
9555f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&ctx.issb));
9565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(index_slow));
9575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(index_fast));
9585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(index_slowbuffer));
9595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListDestroy(&limiters));
9605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListDestroy(&physics));
961*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
962*b122ec5aSJacob Faibussowitsch   return 0;
963c4762a1bSJed Brown }
964c4762a1bSJed Brown 
965c4762a1bSJed Brown /*TEST
966c4762a1bSJed Brown 
967c4762a1bSJed Brown     build:
968f56ea12dSJed Brown       requires: !complex
969c4762a1bSJed Brown       depends: finitevolume1d.c
970c4762a1bSJed Brown 
971c4762a1bSJed Brown     test:
972c4762a1bSJed Brown       suffix: 1
973c4762a1bSJed Brown       args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22
974c4762a1bSJed Brown 
975c4762a1bSJed Brown     test:
976c4762a1bSJed Brown       suffix: 2
977c4762a1bSJed Brown       args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 0
978c4762a1bSJed Brown       output_file: output/ex6_1.out
979c4762a1bSJed Brown 
980c4762a1bSJed Brown TEST*/
981