xref: /petsc/src/ts/tutorials/multirate/ex5.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown /*
2c4762a1bSJed Brown   Note:
3c4762a1bSJed Brown   To performance a convergence study, a reference solution can be generated with a small stepsize and stored into a binary file, e.g. -ts_type ssp -ts_dt 1e-5 -ts_max_steps 30000 -ts_view_solution binary:reference.bin
4c4762a1bSJed Brown   Errors can be computed in the following runs with -simulation -f reference.bin
5c4762a1bSJed Brown 
6c4762a1bSJed Brown   Multirate RK options:
7c4762a1bSJed Brown   -ts_rk_dtratio is the ratio between larger time step size and small time step size
8c4762a1bSJed Brown   -ts_rk_multirate_type has three choices: none (for single step RK)
9c4762a1bSJed Brown                                            combined (for multirate method and user just need to provide the same RHS with the single step RK)
10c4762a1bSJed Brown                                            partitioned (for multiraet method and user need to provide two RHS, one is for fast components and the orther is for slow components
11c4762a1bSJed Brown */
12c4762a1bSJed Brown 
13c4762a1bSJed Brown static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n"
14c4762a1bSJed Brown   " advection   - Variable coefficient scalar advection\n"
15c4762a1bSJed Brown   "                u_t       + (a*u)_x               = 0\n"
16c4762a1bSJed Brown   " a is a piecewise function with two values say a0 and a1 (both a0 and a1 are constant).\n"
17c4762a1bSJed Brown   " in this toy problem we assume a=a0 when x<0 and a=a1 when x>0 with a0<a1 which has the same discontinuous point with initial condition.\n"
18c4762a1bSJed Brown   " we don't provide an exact solution, so you need to generate reference solution to do convergnece staudy,\n"
19c4762a1bSJed Brown   " more precisely, you need firstly generate a reference to a binary file say file.bin, then on commend line,\n"
20c4762a1bSJed Brown   " you should type -simulation -f file.bin.\n"
21c4762a1bSJed Brown   " you can choose the number of grids by -da_grid_x.\n"
22c4762a1bSJed Brown   " you can choose the value of a by -physics_advect_a1 and -physics_advect_a2.\n";
23c4762a1bSJed Brown 
24c4762a1bSJed Brown #include <petscts.h>
25c4762a1bSJed Brown #include <petscdm.h>
26c4762a1bSJed Brown #include <petscdmda.h>
27c4762a1bSJed Brown #include <petscdraw.h>
28c4762a1bSJed Brown #include <petsc/private/tsimpl.h>
29c4762a1bSJed Brown 
30c4762a1bSJed Brown #include <petsc/private/kernels/blockinvert.h> /* For the Kernel_*_gets_* stuff for BAIJ */
31c4762a1bSJed Brown 
32c4762a1bSJed Brown #include "finitevolume1d.h"
33c4762a1bSJed Brown 
349fbee547SJacob Faibussowitsch static inline PetscReal RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax) { PetscReal range = xmax-xmin; return xmin +PetscFmodReal(range+PetscFmodReal(a,range),range); }
35c4762a1bSJed Brown 
36c4762a1bSJed Brown /* --------------------------------- Advection ----------------------------------- */
37c4762a1bSJed Brown 
38c4762a1bSJed Brown typedef struct {
39c4762a1bSJed Brown   PetscReal a[2];                  /* advective velocity */
40c4762a1bSJed Brown } AdvectCtx;
41c4762a1bSJed Brown 
42c4762a1bSJed Brown static PetscErrorCode PhysicsRiemann_Advect(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed,PetscReal x,PetscReal xmin, PetscReal xmax)
43c4762a1bSJed Brown {
44c4762a1bSJed Brown   AdvectCtx *ctx = (AdvectCtx*)vctx;
45c4762a1bSJed Brown   PetscReal *speed;
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   PetscFunctionBeginUser;
48c4762a1bSJed Brown   speed = ctx->a;
49c4762a1bSJed Brown   if (x==0 || x == xmin || x == xmax) flux[0] = PetscMax(0,(speed[0]+speed[1])/2.0)*uL[0] + PetscMin(0,(speed[0]+speed[1])/2.0)*uR[0]; /* if condition need to be changed base on different problem, '0' is the discontinuous point of a */
50c4762a1bSJed Brown   else if (x<0) flux[0] = PetscMax(0,speed[0])*uL[0] + PetscMin(0,speed[0])*uR[0];  /* else if condition need to be changed based on diferent problem, 'x = 0' is discontinuous point of a */
51c4762a1bSJed Brown   else flux[0] = PetscMax(0,speed[1])*uL[0] + PetscMin(0,speed[1])*uR[0];
52c4762a1bSJed Brown   *maxspeed = *speed;
53c4762a1bSJed Brown   PetscFunctionReturn(0);
54c4762a1bSJed Brown }
55c4762a1bSJed Brown 
56c4762a1bSJed Brown static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds,PetscReal x)
57c4762a1bSJed Brown {
58c4762a1bSJed Brown   AdvectCtx *ctx = (AdvectCtx*)vctx;
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   PetscFunctionBeginUser;
61c4762a1bSJed Brown   X[0]      = 1.;
62c4762a1bSJed Brown   Xi[0]     = 1.;
63c4762a1bSJed Brown   if (x<0) speeds[0] = ctx->a[0];  /* x is discontinuous point of a */
64c4762a1bSJed Brown   else    speeds[0] = ctx->a[1];
65c4762a1bSJed Brown   PetscFunctionReturn(0);
66c4762a1bSJed Brown }
67c4762a1bSJed Brown 
68c4762a1bSJed Brown static PetscErrorCode PhysicsSample_Advect(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
69c4762a1bSJed Brown {
70c4762a1bSJed Brown   AdvectCtx *ctx = (AdvectCtx*)vctx;
71c4762a1bSJed Brown   PetscReal *a    = ctx->a,x0;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   PetscFunctionBeginUser;
74c4762a1bSJed Brown   if (x<0){   /* x is cell center */
75c4762a1bSJed Brown     switch (bctype) {
76c4762a1bSJed Brown       case FVBC_OUTFLOW:  x0 = x-a[0]*t; break;
77c4762a1bSJed Brown       case FVBC_PERIODIC: x0 = RangeMod(x-a[0]*t,xmin,xmax); break;
78c4762a1bSJed Brown       default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown BCType");
79c4762a1bSJed Brown     }
80c4762a1bSJed Brown     switch (initial) {
81c4762a1bSJed Brown       case 0: u[0] = (x0 < 0) ? 1 : -1; break;
82c4762a1bSJed Brown       case 1: u[0] = (x0 < 0) ? -1 : 1; break;
83c4762a1bSJed Brown       case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break;
84c4762a1bSJed Brown       case 3: u[0] = PetscSinReal(2*PETSC_PI*x0); break;
85c4762a1bSJed Brown       case 4: u[0] = PetscAbs(x0); break;
86c4762a1bSJed Brown       case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2*PETSC_PI*x0)); break;
87c4762a1bSJed Brown       case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2-x0 : 0)); break;
88c4762a1bSJed Brown       case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI*x0),10.0);break;
89c4762a1bSJed Brown       default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition");
90c4762a1bSJed Brown     }
91c4762a1bSJed Brown   }
92c4762a1bSJed Brown   else{
93c4762a1bSJed Brown     switch (bctype) {
94c4762a1bSJed Brown       case FVBC_OUTFLOW:  x0 = x-a[1]*t; break;
95c4762a1bSJed Brown       case FVBC_PERIODIC: x0 = RangeMod(x-a[1]*t,xmin,xmax); break;
96c4762a1bSJed Brown       default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown BCType");
97c4762a1bSJed Brown     }
98c4762a1bSJed Brown     switch (initial) {
99c4762a1bSJed Brown       case 0: u[0] = (x0 < 0) ? 1 : -1; break;
100c4762a1bSJed Brown       case 1: u[0] = (x0 < 0) ? -1 : 1; break;
101c4762a1bSJed Brown       case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break;
102c4762a1bSJed Brown       case 3: u[0] = PetscSinReal(2*PETSC_PI*x0); break;
103c4762a1bSJed Brown       case 4: u[0] = PetscAbs(x0); break;
104c4762a1bSJed Brown       case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2*PETSC_PI*x0)); break;
105c4762a1bSJed Brown       case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2-x0 : 0)); break;
106c4762a1bSJed Brown       case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI*x0),10.0);break;
107c4762a1bSJed Brown       default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition");
108c4762a1bSJed Brown     }
109c4762a1bSJed Brown   }
110c4762a1bSJed Brown   PetscFunctionReturn(0);
111c4762a1bSJed Brown }
112c4762a1bSJed Brown 
113c4762a1bSJed Brown static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx)
114c4762a1bSJed Brown {
115c4762a1bSJed Brown   PetscErrorCode ierr;
116c4762a1bSJed Brown   AdvectCtx      *user;
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   PetscFunctionBeginUser;
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNew(&user));
120c4762a1bSJed Brown   ctx->physics.sample         = PhysicsSample_Advect;
121c4762a1bSJed Brown   ctx->physics.riemann        = PhysicsRiemann_Advect;
122c4762a1bSJed Brown   ctx->physics.characteristic = PhysicsCharacteristic_Advect;
123c4762a1bSJed Brown   ctx->physics.destroy        = PhysicsDestroy_SimpleFree;
124c4762a1bSJed Brown   ctx->physics.user           = user;
125c4762a1bSJed Brown   ctx->physics.dof            = 1;
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrallocpy("u",&ctx->physics.fieldname[0]));
127c4762a1bSJed Brown   user->a[0] = 1;
128c4762a1bSJed Brown   user->a[1] = 1;
129c4762a1bSJed Brown   ierr = PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");CHKERRQ(ierr);
130c4762a1bSJed Brown   {
1315f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-physics_advect_a1","Speed1","",user->a[0],&user->a[0],NULL));
1325f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-physics_advect_a2","Speed2","",user->a[1],&user->a[1],NULL));
133c4762a1bSJed Brown   }
134c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
135c4762a1bSJed Brown   PetscFunctionReturn(0);
136c4762a1bSJed Brown }
137c4762a1bSJed Brown 
138c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow parts ----------------------------------- */
139c4762a1bSJed Brown 
140c4762a1bSJed Brown PetscErrorCode FVRHSFunctionslow(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
141c4762a1bSJed Brown {
142c4762a1bSJed Brown   FVCtx          *ctx = (FVCtx*)vctx;
143c4762a1bSJed Brown   PetscInt       i,j,k,Mx,dof,xs,xm,len_slow;
144c4762a1bSJed Brown   PetscReal      hx,cfl_idt = 0;
145c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
146c4762a1bSJed Brown   Vec            Xloc;
147c4762a1bSJed Brown   DM             da;
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   PetscFunctionBeginUser;
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&Xloc));
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
153c4762a1bSJed Brown   hx   = (ctx->xmax-ctx->xmin)/Mx;
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(F));
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,Xloc,&x));
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope));
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetSize(ctx->iss,&len_slow));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
164c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
165c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
166c4762a1bSJed Brown     }
167c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
168c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
169c4762a1bSJed Brown     }
170c4762a1bSJed Brown   }
171c4762a1bSJed Brown   for (i=xs-1; i<xs+xm+1; i++) {
172c4762a1bSJed Brown     struct _LimitInfo info;
173c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
174c4762a1bSJed Brown     if (i < len_slow+1) {
175c4762a1bSJed Brown       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
1765f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds,ctx->xmin+hx*i));
177c4762a1bSJed Brown       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
1785f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof));
179c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
180c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
181c4762a1bSJed Brown       for (j=0; j<dof; j++) {
182c4762a1bSJed Brown         PetscScalar jmpL,jmpR;
183c4762a1bSJed Brown         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
184c4762a1bSJed Brown         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
185c4762a1bSJed Brown         for (k=0; k<dof; k++) {
186c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
187c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
188c4762a1bSJed Brown         }
189c4762a1bSJed Brown       }
190c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
191c4762a1bSJed Brown       info.m  = dof;
192c4762a1bSJed Brown       info.hx = hx;
193c4762a1bSJed Brown       (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope);
194c4762a1bSJed Brown       for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */
195c4762a1bSJed Brown       for (j=0; j<dof; j++) {
196c4762a1bSJed Brown         PetscScalar tmp = 0;
197c4762a1bSJed Brown         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
198c4762a1bSJed Brown         slope[i*dof+j] = tmp;
199c4762a1bSJed Brown       }
200c4762a1bSJed Brown     }
201c4762a1bSJed Brown   }
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
204c4762a1bSJed Brown     PetscReal   maxspeed;
205c4762a1bSJed Brown     PetscScalar *uL,*uR;
206c4762a1bSJed Brown     if (i < len_slow) { /* slow parts can be changed based on a */
207c4762a1bSJed Brown       uL = &ctx->uLR[0];
208c4762a1bSJed Brown       uR = &ctx->uLR[dof];
209c4762a1bSJed Brown       for (j=0; j<dof; j++) {
210c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
211c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
212c4762a1bSJed Brown       }
2135f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax));
214c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
215c4762a1bSJed Brown       if (i > xs) {
216c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hx;
217c4762a1bSJed Brown       }
218c4762a1bSJed Brown       if (i < xs+xm) {
219c4762a1bSJed Brown         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hx;
220c4762a1bSJed Brown       }
221c4762a1bSJed Brown     }
222c4762a1bSJed Brown     if (i == len_slow) { /* slow parts can be changed based on a */
223c4762a1bSJed Brown       uL = &ctx->uLR[0];
224c4762a1bSJed Brown       uR = &ctx->uLR[dof];
225c4762a1bSJed Brown       for (j=0; j<dof; j++) {
226c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
227c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
228c4762a1bSJed Brown       }
2295f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax));
230c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
231c4762a1bSJed Brown       if (i > xs) {
232c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hx;
233c4762a1bSJed Brown       }
234c4762a1bSJed Brown     }
235c4762a1bSJed Brown   }
2365f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x));
2375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
2385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope));
2395f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&Xloc));
240c4762a1bSJed Brown   PetscFunctionReturn(0);
241c4762a1bSJed Brown }
242c4762a1bSJed Brown 
243c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for fast  parts ----------------------------------- */
244c4762a1bSJed Brown 
245c4762a1bSJed Brown PetscErrorCode FVRHSFunctionfast(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
246c4762a1bSJed Brown {
247c4762a1bSJed Brown   FVCtx          *ctx = (FVCtx*)vctx;
248c4762a1bSJed Brown   PetscInt       i,j,k,Mx,dof,xs,xm,len_slow;
249c4762a1bSJed Brown   PetscReal      hx,cfl_idt = 0;
250c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
251c4762a1bSJed Brown   Vec            Xloc;
252c4762a1bSJed Brown   DM             da;
253c4762a1bSJed Brown 
254c4762a1bSJed Brown   PetscFunctionBeginUser;
2555f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
2565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&Xloc));
2575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
258c4762a1bSJed Brown   hx   = (ctx->xmax-ctx->xmin)/Mx;
2595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));
2605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc));
2615f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(F));
2625f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,Xloc,&x));
2635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
2645f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope));
2655f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
2665f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetSize(ctx->iss,&len_slow));
267c4762a1bSJed Brown 
268c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
269c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
270c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
271c4762a1bSJed Brown     }
272c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
273c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
274c4762a1bSJed Brown     }
275c4762a1bSJed Brown   }
276c4762a1bSJed Brown   for (i=xs-1; i<xs+xm+1; i++) {
277c4762a1bSJed Brown     struct _LimitInfo info;
278c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
279c4762a1bSJed Brown     if (i > len_slow-2) {
2805f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds,ctx->xmin+hx*i));
2815f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof));
282c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
283c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
284c4762a1bSJed Brown       for (j=0; j<dof; j++) {
285c4762a1bSJed Brown         PetscScalar jmpL,jmpR;
286c4762a1bSJed Brown         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
287c4762a1bSJed Brown         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
288c4762a1bSJed Brown         for (k=0; k<dof; k++) {
289c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
290c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
291c4762a1bSJed Brown         }
292c4762a1bSJed Brown       }
293c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
294c4762a1bSJed Brown       info.m  = dof;
295c4762a1bSJed Brown       info.hx = hx;
296c4762a1bSJed Brown       (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope);
297c4762a1bSJed Brown       for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */
298c4762a1bSJed Brown       for (j=0; j<dof; j++) {
299c4762a1bSJed Brown         PetscScalar tmp = 0;
300c4762a1bSJed Brown         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
301c4762a1bSJed Brown         slope[i*dof+j] = tmp;
302c4762a1bSJed Brown       }
303c4762a1bSJed Brown     }
304c4762a1bSJed Brown   }
305c4762a1bSJed Brown 
306c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
307c4762a1bSJed Brown     PetscReal   maxspeed;
308c4762a1bSJed Brown     PetscScalar *uL,*uR;
309c4762a1bSJed Brown     if (i > len_slow) {
310c4762a1bSJed Brown       uL = &ctx->uLR[0];
311c4762a1bSJed Brown       uR = &ctx->uLR[dof];
312c4762a1bSJed Brown       for (j=0; j<dof; j++) {
313c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
314c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
315c4762a1bSJed Brown       }
3165f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax));
317c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
318c4762a1bSJed Brown       if (i > xs) {
319c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow-1)*dof+j] -= ctx->flux[j]/hx;
320c4762a1bSJed Brown       }
321c4762a1bSJed Brown       if (i < xs+xm) {
322c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow)*dof+j] += ctx->flux[j]/hx;
323c4762a1bSJed Brown       }
324c4762a1bSJed Brown     }
325c4762a1bSJed Brown     if (i == len_slow) {
326c4762a1bSJed Brown       uL = &ctx->uLR[0];
327c4762a1bSJed Brown       uR = &ctx->uLR[dof];
328c4762a1bSJed Brown       for (j=0; j<dof; j++) {
329c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
330c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
331c4762a1bSJed Brown       }
3325f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax));
333c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
334c4762a1bSJed Brown       if (i < xs+xm) {
335c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow)*dof+j] += ctx->flux[j]/hx;
336c4762a1bSJed Brown       }
337c4762a1bSJed Brown     }
338c4762a1bSJed Brown   }
3395f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x));
3405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
3415f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope));
3425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&Xloc));
343c4762a1bSJed Brown   PetscFunctionReturn(0);
344c4762a1bSJed Brown }
345c4762a1bSJed Brown 
346c4762a1bSJed Brown PetscErrorCode FVRHSFunctionslow2(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
347c4762a1bSJed Brown {
348c4762a1bSJed Brown   FVCtx          *ctx = (FVCtx*)vctx;
349c4762a1bSJed Brown   PetscInt       i,j,k,Mx,dof,xs,xm,len_slow1,len_slow2;
350c4762a1bSJed Brown   PetscReal      hx,cfl_idt = 0;
351c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
352c4762a1bSJed Brown   Vec            Xloc;
353c4762a1bSJed Brown   DM             da;
354c4762a1bSJed Brown 
355c4762a1bSJed Brown   PetscFunctionBeginUser;
3565f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
3575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&Xloc));
3585f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
359c4762a1bSJed Brown   hx   = (ctx->xmax-ctx->xmin)/Mx;
3605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));
3615f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc));
3625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(F));
3635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,Xloc,&x));
3645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
3655f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope));
3665f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
3675f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetSize(ctx->iss,&len_slow1));
3685f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetSize(ctx->iss2,&len_slow2));
369c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
370c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
371c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
372c4762a1bSJed Brown     }
373c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
374c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
375c4762a1bSJed Brown     }
376c4762a1bSJed Brown   }
377c4762a1bSJed Brown   for (i=xs-1; i<xs+xm+1; i++) {
378c4762a1bSJed Brown     struct _LimitInfo info;
379c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
380c4762a1bSJed Brown     if (i < len_slow1+len_slow2+1 && i > len_slow1-2) {
3815f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds,ctx->xmin+hx*i));
3825f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof));
383c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
384c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
385c4762a1bSJed Brown       for (j=0; j<dof; j++) {
386c4762a1bSJed Brown         PetscScalar jmpL,jmpR;
387c4762a1bSJed Brown         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
388c4762a1bSJed Brown         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
389c4762a1bSJed Brown         for (k=0; k<dof; k++) {
390c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
391c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
392c4762a1bSJed Brown         }
393c4762a1bSJed Brown       }
394c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
395c4762a1bSJed Brown       info.m  = dof;
396c4762a1bSJed Brown       info.hx = hx;
397c4762a1bSJed Brown       (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope);
398c4762a1bSJed Brown       for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */
399c4762a1bSJed Brown       for (j=0; j<dof; j++) {
400c4762a1bSJed Brown         PetscScalar tmp = 0;
401c4762a1bSJed Brown         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
402c4762a1bSJed Brown         slope[i*dof+j] = tmp;
403c4762a1bSJed Brown       }
404c4762a1bSJed Brown     }
405c4762a1bSJed Brown   }
406c4762a1bSJed Brown 
407c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
408c4762a1bSJed Brown     PetscReal   maxspeed;
409c4762a1bSJed Brown     PetscScalar *uL,*uR;
410c4762a1bSJed Brown     if (i < len_slow1+len_slow2 && i > len_slow1) { /* slow parts can be changed based on a */
411c4762a1bSJed Brown       uL = &ctx->uLR[0];
412c4762a1bSJed Brown       uR = &ctx->uLR[dof];
413c4762a1bSJed Brown       for (j=0; j<dof; j++) {
414c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
415c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
416c4762a1bSJed Brown       }
4175f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax));
418c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
419c4762a1bSJed Brown       if (i > xs) {
420c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow1-1)*dof+j] -= ctx->flux[j]/hx;
421c4762a1bSJed Brown       }
422c4762a1bSJed Brown       if (i < xs+xm) {
423c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow1)*dof+j] += ctx->flux[j]/hx;
424c4762a1bSJed Brown       }
425c4762a1bSJed Brown     }
426c4762a1bSJed Brown     if (i == len_slow1+len_slow2) { /* slow parts can be changed based on a */
427c4762a1bSJed Brown       uL = &ctx->uLR[0];
428c4762a1bSJed Brown       uR = &ctx->uLR[dof];
429c4762a1bSJed Brown       for (j=0; j<dof; j++) {
430c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
431c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
432c4762a1bSJed Brown       }
4335f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax));
434c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
435c4762a1bSJed Brown       if (i > xs) {
436c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow1-1)*dof+j] -= ctx->flux[j]/hx;
437c4762a1bSJed Brown       }
438c4762a1bSJed Brown     }
439c4762a1bSJed Brown     if (i == len_slow1) { /* slow parts can be changed based on a */
440c4762a1bSJed Brown       uL = &ctx->uLR[0];
441c4762a1bSJed Brown       uR = &ctx->uLR[dof];
442c4762a1bSJed Brown       for (j=0; j<dof; j++) {
443c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
444c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
445c4762a1bSJed Brown       }
4465f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax));
447c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
448c4762a1bSJed Brown       if (i < xs+xm) {
449c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow1)*dof+j] += ctx->flux[j]/hx;
450c4762a1bSJed Brown       }
451c4762a1bSJed Brown     }
452c4762a1bSJed Brown   }
453c4762a1bSJed Brown 
4545f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x));
4555f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
4565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope));
4575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&Xloc));
458c4762a1bSJed Brown   PetscFunctionReturn(0);
459c4762a1bSJed Brown }
460c4762a1bSJed Brown 
461c4762a1bSJed Brown PetscErrorCode FVRHSFunctionfast2(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
462c4762a1bSJed Brown {
463c4762a1bSJed Brown   FVCtx          *ctx = (FVCtx*)vctx;
464c4762a1bSJed Brown   PetscInt       i,j,k,Mx,dof,xs,xm,len_slow1,len_slow2;
465c4762a1bSJed Brown   PetscReal      hx,cfl_idt = 0;
466c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
467c4762a1bSJed Brown   Vec            Xloc;
468c4762a1bSJed Brown   DM             da;
469c4762a1bSJed Brown 
470c4762a1bSJed Brown   PetscFunctionBeginUser;
4715f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
4725f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&Xloc));
4735f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
474c4762a1bSJed Brown   hx   = (ctx->xmax-ctx->xmin)/Mx;
4755f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));
4765f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc));
4775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(F));
4785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,Xloc,&x));
4795f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
4805f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope));
4815f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
4825f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetSize(ctx->iss,&len_slow1));
4835f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetSize(ctx->iss2,&len_slow2));
484c4762a1bSJed Brown 
485c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
486c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
487c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
488c4762a1bSJed Brown     }
489c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
490c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
491c4762a1bSJed Brown     }
492c4762a1bSJed Brown   }
493c4762a1bSJed Brown   for (i=xs-1; i<xs+xm+1; i++) {
494c4762a1bSJed Brown     struct _LimitInfo info;
495c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
496c4762a1bSJed Brown     if (i > len_slow1+len_slow2-2) {
4975f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds,ctx->xmin+hx*i));
4985f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof));
499c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
500c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
501c4762a1bSJed Brown       for (j=0; j<dof; j++) {
502c4762a1bSJed Brown         PetscScalar jmpL,jmpR;
503c4762a1bSJed Brown         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
504c4762a1bSJed Brown         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
505c4762a1bSJed Brown         for (k=0; k<dof; k++) {
506c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
507c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
508c4762a1bSJed Brown         }
509c4762a1bSJed Brown       }
510c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
511c4762a1bSJed Brown       info.m  = dof;
512c4762a1bSJed Brown       info.hx = hx;
513c4762a1bSJed Brown       (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope);
514c4762a1bSJed Brown       for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */
515c4762a1bSJed Brown       for (j=0; j<dof; j++) {
516c4762a1bSJed Brown         PetscScalar tmp = 0;
517c4762a1bSJed Brown         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
518c4762a1bSJed Brown         slope[i*dof+j] = tmp;
519c4762a1bSJed Brown       }
520c4762a1bSJed Brown     }
521c4762a1bSJed Brown   }
522c4762a1bSJed Brown 
523c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
524c4762a1bSJed Brown     PetscReal   maxspeed;
525c4762a1bSJed Brown     PetscScalar *uL,*uR;
526c4762a1bSJed Brown     if (i > len_slow1+len_slow2) {
527c4762a1bSJed Brown       uL = &ctx->uLR[0];
528c4762a1bSJed Brown       uR = &ctx->uLR[dof];
529c4762a1bSJed Brown       for (j=0; j<dof; j++) {
530c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
531c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
532c4762a1bSJed Brown       }
5335f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax));
534c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
535c4762a1bSJed Brown       if (i > xs) {
536c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow1-len_slow2-1)*dof+j] -= ctx->flux[j]/hx;
537c4762a1bSJed Brown       }
538c4762a1bSJed Brown       if (i < xs+xm) {
539c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow1-len_slow2)*dof+j] += ctx->flux[j]/hx;
540c4762a1bSJed Brown       }
541c4762a1bSJed Brown     }
542c4762a1bSJed Brown     if (i == len_slow1+len_slow2) {
543c4762a1bSJed Brown       uL = &ctx->uLR[0];
544c4762a1bSJed Brown       uR = &ctx->uLR[dof];
545c4762a1bSJed Brown       for (j=0; j<dof; j++) {
546c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
547c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
548c4762a1bSJed Brown       }
5495f80ce2aSJacob Faibussowitsch       CHKERRQ((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax));
550c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
551c4762a1bSJed Brown       if (i < xs+xm) {
552c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow1-len_slow2)*dof+j] += ctx->flux[j]/hx;
553c4762a1bSJed Brown       }
554c4762a1bSJed Brown     }
555c4762a1bSJed Brown   }
5565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x));
5575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
5585f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope));
5595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&Xloc));
560c4762a1bSJed Brown   PetscFunctionReturn(0);
561c4762a1bSJed Brown }
562c4762a1bSJed Brown 
563c4762a1bSJed Brown int main(int argc,char *argv[])
564c4762a1bSJed Brown {
565c4762a1bSJed Brown   char              lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m";
566c4762a1bSJed Brown   PetscFunctionList limiters   = 0,physics = 0;
567c4762a1bSJed Brown   MPI_Comm          comm;
568c4762a1bSJed Brown   TS                ts;
569c4762a1bSJed Brown   DM                da;
570c4762a1bSJed Brown   Vec               X,X0,R;
571c4762a1bSJed Brown   FVCtx             ctx;
572c4762a1bSJed Brown   PetscInt          i,k,dof,xs,xm,Mx,draw = 0,*index_slow,*index_fast,islow = 0,ifast = 0;
573c4762a1bSJed Brown   PetscBool         view_final = PETSC_FALSE;
574c4762a1bSJed Brown   PetscReal         ptime;
575c4762a1bSJed Brown   PetscErrorCode    ierr;
576c4762a1bSJed Brown 
577*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,0,help));
578c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
5795f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMemzero(&ctx,sizeof(ctx)));
580c4762a1bSJed Brown 
581c4762a1bSJed Brown   /* Register limiters to be available on the command line */
5825f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"upwind"              ,Limit_Upwind));
5835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"lax-wendroff"        ,Limit_LaxWendroff));
5845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"beam-warming"        ,Limit_BeamWarming));
5855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"fromm"               ,Limit_Fromm));
5865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"minmod"              ,Limit_Minmod));
5875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"superbee"            ,Limit_Superbee));
5885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"mc"                  ,Limit_MC));
5895f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"vanleer"             ,Limit_VanLeer));
5905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"vanalbada"           ,Limit_VanAlbada));
5915f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"vanalbadatvd"        ,Limit_VanAlbadaTVD));
5925f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"koren"               ,Limit_Koren));
5935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"korensym"            ,Limit_KorenSym));
5945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"koren3"              ,Limit_Koren3));
5955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"cada-torrilhon2"     ,Limit_CadaTorrilhon2));
5965f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"cada-torrilhon3-r0p1",Limit_CadaTorrilhon3R0p1));
5975f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"cada-torrilhon3-r1"  ,Limit_CadaTorrilhon3R1));
5985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"cada-torrilhon3-r10" ,Limit_CadaTorrilhon3R10));
5995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&limiters,"cada-torrilhon3-r100",Limit_CadaTorrilhon3R100));
600c4762a1bSJed Brown 
601c4762a1bSJed Brown   /* Register physical models to be available on the command line */
6025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&physics,"advect"          ,PhysicsCreate_Advect));
603c4762a1bSJed Brown 
604c4762a1bSJed Brown   ctx.comm = comm;
605c4762a1bSJed Brown   ctx.cfl  = 0.9;
606c4762a1bSJed Brown   ctx.bctype = FVBC_PERIODIC;
607c4762a1bSJed Brown   ctx.xmin = -1.0;
608c4762a1bSJed Brown   ctx.xmax = 1.0;
609c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm,NULL,"Finite Volume solver options","");CHKERRQ(ierr);
6105f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL));
6115f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL));
6125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsFList("-limit","Name of flux limiter to use","",limiters,lname,lname,sizeof(lname),NULL));
6135f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL));
6145f80ce2aSJacob 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));
6155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL));
6165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL));
6175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL));
6185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL));
6195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-recursive_split","Split the domain recursively","",ctx.recursive,&ctx.recursive,NULL));
620c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
621c4762a1bSJed Brown 
622c4762a1bSJed Brown   /* Choose the limiter from the list of registered limiters */
6235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListFind(limiters,lname,&ctx.limit));
6243c633725SBarry Smith   PetscCheck(ctx.limit,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Limiter '%s' not found",lname);
625c4762a1bSJed Brown 
626c4762a1bSJed Brown   /* Choose the physics from the list of registered models */
627c4762a1bSJed Brown   {
628c4762a1bSJed Brown     PetscErrorCode (*r)(FVCtx*);
6295f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFunctionListFind(physics,physname,&r));
6303c633725SBarry Smith     PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Physics '%s' not found",physname);
631c4762a1bSJed Brown     /* Create the physics, will set the number of fields and their names */
6325f80ce2aSJacob Faibussowitsch     CHKERRQ((*r)(&ctx));
633c4762a1bSJed Brown   }
634c4762a1bSJed Brown 
635c4762a1bSJed Brown   /* Create a DMDA to manage the parallel grid */
6365f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics.dof,2,NULL,&da));
6375f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
6385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
639c4762a1bSJed Brown   /* Inform the DMDA of the field names provided by the physics. */
640c4762a1bSJed Brown   /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
641c4762a1bSJed Brown   for (i=0; i<ctx.physics.dof; i++) {
6425f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDASetFieldName(da,i,ctx.physics.fieldname[i]));
643c4762a1bSJed Brown   }
6445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
6455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
646c4762a1bSJed Brown 
647c4762a1bSJed Brown   /* Set coordinates of cell centers */
6485f80ce2aSJacob 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));
649c4762a1bSJed Brown 
650c4762a1bSJed Brown   /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
6515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc4(dof*dof,&ctx.R,dof*dof,&ctx.Rinv,2*dof,&ctx.cjmpLR,1*dof,&ctx.cslope));
6525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3(2*dof,&ctx.uLR,dof,&ctx.flux,dof,&ctx.speeds));
653c4762a1bSJed Brown 
654c4762a1bSJed Brown   /* Create a vector to store the solution and to save the initial state */
6555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&X));
6565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(X,&X0));
6575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(X,&R));
658c4762a1bSJed Brown 
659c4762a1bSJed Brown   /*-------------------------------- create index for slow parts and fast parts ----------------------------------------*/
6605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(xm*dof,&index_slow));
6615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(xm*dof,&index_fast));
662c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
663c4762a1bSJed Brown     if (ctx.xmin+i*(ctx.xmax-ctx.xmin)/(PetscReal)Mx+0.5*(ctx.xmax-ctx.xmin)/(PetscReal)Mx < 0)
664c4762a1bSJed Brown       for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k;
665c4762a1bSJed Brown     else
666c4762a1bSJed Brown       for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k;
667c4762a1bSJed Brown   }  /* this step need to be changed based on discontinuous point of a */
6685f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss));
6695f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf));
670c4762a1bSJed Brown 
671c4762a1bSJed Brown   /* Create a time-stepping object */
6725f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(comm,&ts));
6735f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts,da));
6745f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,R,FVRHSFunction,&ctx));
6755f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetIS(ts,"slow",ctx.iss));
6765f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetIS(ts,"fast",ctx.isf));
6775f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetRHSFunction(ts,"slow",NULL,FVRHSFunctionslow,&ctx));
6785f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetRHSFunction(ts,"fast",NULL,FVRHSFunctionfast,&ctx));
679c4762a1bSJed Brown 
680c4762a1bSJed Brown   if (ctx.recursive) {
681c4762a1bSJed Brown     TS subts;
682c4762a1bSJed Brown     islow = 0;
683c4762a1bSJed Brown     ifast = 0;
684c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
685c4762a1bSJed Brown       PetscReal coord = ctx.xmin+i*(ctx.xmax-ctx.xmin)/(PetscReal)Mx+0.5*(ctx.xmax-ctx.xmin)/(PetscReal)Mx;
686c4762a1bSJed Brown       if (coord >= 0 && coord < ctx.xmin+(ctx.xmax-ctx.xmin)*3/4.)
687c4762a1bSJed Brown         for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k;
688c4762a1bSJed Brown       if (coord >= ctx.xmin+(ctx.xmax-ctx.xmin)*3/4.)
689c4762a1bSJed Brown         for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k;
690c4762a1bSJed Brown     }  /* this step need to be changed based on discontinuous point of a */
6915f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss2));
6925f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf2));
693c4762a1bSJed Brown 
6945f80ce2aSJacob Faibussowitsch     CHKERRQ(TSRHSSplitGetSubTS(ts,"fast",&subts));
6955f80ce2aSJacob Faibussowitsch     CHKERRQ(TSRHSSplitSetIS(subts,"slow",ctx.iss2));
6965f80ce2aSJacob Faibussowitsch     CHKERRQ(TSRHSSplitSetIS(subts,"fast",ctx.isf2));
6975f80ce2aSJacob Faibussowitsch     CHKERRQ(TSRHSSplitSetRHSFunction(subts,"slow",NULL,FVRHSFunctionslow2,&ctx));
6985f80ce2aSJacob Faibussowitsch     CHKERRQ(TSRHSSplitSetRHSFunction(subts,"fast",NULL,FVRHSFunctionfast2,&ctx));
699c4762a1bSJed Brown   }
700c4762a1bSJed Brown 
7015f80ce2aSJacob Faibussowitsch   /*CHKERRQ(TSSetType(ts,TSSSP));*/
7025f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSMPRK));
7035f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,10));
7045f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
705c4762a1bSJed Brown 
706c4762a1bSJed Brown   /* Compute initial conditions and starting time step */
7075f80ce2aSJacob Faibussowitsch   CHKERRQ(FVSample(&ctx,da,0,X0));
7085f80ce2aSJacob Faibussowitsch   CHKERRQ(FVRHSFunction(ts,0,X0,X,(void*)&ctx)); /* Initial function evaluation, only used to determine max speed */
7095f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(X0,X));                        /* The function value was not used so we set X=X0 again */
7105f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt));
7115f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts)); /* Take runtime options */
7125f80ce2aSJacob Faibussowitsch   CHKERRQ(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD));
713c4762a1bSJed Brown   {
714c4762a1bSJed Brown     PetscInt    steps;
715c4762a1bSJed Brown     PetscScalar mass_initial,mass_final,mass_difference;
716c4762a1bSJed Brown 
7175f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSolve(ts,X));
7185f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetSolveTime(ts,&ptime));
7195f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetStepNumber(ts,&steps));
7205f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm,"Final time %g, steps %D\n",(double)ptime,steps));
721c4762a1bSJed Brown     /* calculate the total mass at initial time and final time */
722c4762a1bSJed Brown     mass_initial = 0.0;
723c4762a1bSJed Brown     mass_final   = 0.0;
7245f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSum(X0,&mass_initial));
7255f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSum(X,&mass_final));
726c4762a1bSJed Brown     mass_difference = (ctx.xmax-ctx.xmin)/(PetscScalar)Mx*(mass_final - mass_initial);
7275f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm,"Mass difference %g\n",(double)mass_difference));
728c4762a1bSJed Brown     if (ctx.simulation) {
729c4762a1bSJed Brown       PetscViewer  fd;
730c4762a1bSJed Brown       char         filename[PETSC_MAX_PATH_LEN] = "binaryoutput";
731c4762a1bSJed Brown       Vec          XR;
732c4762a1bSJed Brown       PetscReal    nrm1,nrmsup;
733c4762a1bSJed Brown       PetscBool    flg;
734d8185827SBarry Smith 
7355f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg));
7363c633725SBarry Smith       PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
7375f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&fd));
7385f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDuplicate(X0,&XR));
7395f80ce2aSJacob Faibussowitsch       CHKERRQ(VecLoad(XR,fd));
7405f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerDestroy(&fd));
7415f80ce2aSJacob Faibussowitsch       CHKERRQ(VecAYPX(XR,-1,X));
7425f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(XR,NORM_1,&nrm1));
7435f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(XR,NORM_INFINITY,&nrmsup));
744c4762a1bSJed Brown       nrm1 /= Mx;
7455f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm,"Error ||x-x_e||_1 %g  ||x-x_e||_sup %g\n",(double)nrm1,(double)nrmsup));
7465f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&XR));
747c4762a1bSJed Brown     }
748c4762a1bSJed Brown   }
749c4762a1bSJed Brown 
7505f80ce2aSJacob Faibussowitsch   CHKERRQ(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD));
7515f80ce2aSJacob Faibussowitsch   if (draw & 0x1) CHKERRQ(VecView(X0,PETSC_VIEWER_DRAW_WORLD));
7525f80ce2aSJacob Faibussowitsch   if (draw & 0x2) CHKERRQ(VecView(X,PETSC_VIEWER_DRAW_WORLD));
753c4762a1bSJed Brown   if (draw & 0x4) {
754c4762a1bSJed Brown     Vec Y;
7555f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(X,&Y));
7565f80ce2aSJacob Faibussowitsch     CHKERRQ(FVSample(&ctx,da,ptime,Y));
7575f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAYPX(Y,-1,X));
7585f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(Y,PETSC_VIEWER_DRAW_WORLD));
7595f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&Y));
760c4762a1bSJed Brown   }
761c4762a1bSJed Brown 
762c4762a1bSJed Brown   if (view_final) {
763c4762a1bSJed Brown     PetscViewer viewer;
7645f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer));
7655f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
7665f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(X,viewer));
7675f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerPopFormat(viewer));
7685f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&viewer));
769c4762a1bSJed Brown   }
770c4762a1bSJed Brown 
771c4762a1bSJed Brown   /* Clean up */
7725f80ce2aSJacob Faibussowitsch   CHKERRQ((*ctx.physics.destroy)(ctx.physics.user));
7735f80ce2aSJacob Faibussowitsch   for (i=0; i<ctx.physics.dof; i++) CHKERRQ(PetscFree(ctx.physics.fieldname[i]));
7745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree4(ctx.R,ctx.Rinv,ctx.cjmpLR,ctx.cslope));
7755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(ctx.uLR,ctx.flux,ctx.speeds));
7765f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&ctx.iss));
7775f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&ctx.isf));
7785f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&ctx.iss2));
7795f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&ctx.isf2));
7805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X));
7815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X0));
7825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&R));
7835f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
7845f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
7855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(index_slow));
7865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(index_fast));
7875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListDestroy(&limiters));
7885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListDestroy(&physics));
789*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
790*b122ec5aSJacob Faibussowitsch   return 0;
791c4762a1bSJed Brown }
792c4762a1bSJed Brown 
793c4762a1bSJed Brown /*TEST
794c4762a1bSJed Brown     build:
795f56ea12dSJed Brown       requires: !complex
796c4762a1bSJed Brown       depends: finitevolume1d.c
797c4762a1bSJed Brown 
798c4762a1bSJed Brown     test:
799c4762a1bSJed Brown       args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 0
800c4762a1bSJed Brown 
801c4762a1bSJed Brown     test:
802c4762a1bSJed Brown       suffix: 2
803c4762a1bSJed Brown       args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 1
804c4762a1bSJed Brown       output_file: output/ex5_1.out
805c4762a1bSJed Brown 
806c4762a1bSJed Brown     test:
807c4762a1bSJed Brown       suffix: 3
808c4762a1bSJed Brown       args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 0
809c4762a1bSJed Brown 
810c4762a1bSJed Brown     test:
811c4762a1bSJed Brown       suffix: 4
812c4762a1bSJed Brown       args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 1
813c4762a1bSJed Brown       output_file: output/ex5_3.out
814c4762a1bSJed Brown 
815c4762a1bSJed Brown     test:
816c4762a1bSJed Brown       suffix: 5
817c4762a1bSJed Brown       args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 0 -recursive_split
818c4762a1bSJed Brown 
819c4762a1bSJed Brown     test:
820c4762a1bSJed Brown       suffix: 6
821c4762a1bSJed Brown       args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 1 -recursive_split
822c4762a1bSJed Brown TEST*/
823