xref: /petsc/src/ts/tutorials/multirate/ex5.c (revision f56ea12d3d7e20aed2e24c9530b477864b138a71) !
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 
34c4762a1bSJed Brown PETSC_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;
119c4762a1bSJed Brown   ierr = PetscNew(&user);CHKERRQ(ierr);
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;
126c4762a1bSJed Brown   ierr = PetscStrallocpy("u",&ctx->physics.fieldname[0]);CHKERRQ(ierr);
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   {
131c4762a1bSJed Brown     ierr = PetscOptionsReal("-physics_advect_a1","Speed1","",user->a[0],&user->a[0],NULL);CHKERRQ(ierr);
132c4762a1bSJed Brown     ierr = PetscOptionsReal("-physics_advect_a2","Speed2","",user->a[1],&user->a[1],NULL);CHKERRQ(ierr);
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   PetscErrorCode ierr;
144c4762a1bSJed Brown   PetscInt       i,j,k,Mx,dof,xs,xm,len_slow;
145c4762a1bSJed Brown   PetscReal      hx,cfl_idt = 0;
146c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
147c4762a1bSJed Brown   Vec            Xloc;
148c4762a1bSJed Brown   DM             da;
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   PetscFunctionBeginUser;
151c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
152c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr);
153c4762a1bSJed Brown   ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr);
154c4762a1bSJed Brown   hx   = (ctx->xmax-ctx->xmin)/Mx;
155c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
156c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
157c4762a1bSJed Brown   ierr = VecZeroEntries(F);CHKERRQ(ierr);
158c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr);
159c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
160c4762a1bSJed Brown   ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
161c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr);
162c4762a1bSJed Brown   ierr = ISGetSize(ctx->iss,&len_slow);CHKERRQ(ierr);
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
165c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
166c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
167c4762a1bSJed Brown     }
168c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
169c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
170c4762a1bSJed Brown     }
171c4762a1bSJed Brown   }
172c4762a1bSJed Brown   for (i=xs-1; i<xs+xm+1; i++) {
173c4762a1bSJed Brown     struct _LimitInfo info;
174c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
175c4762a1bSJed Brown     if (i < len_slow+1) {
176c4762a1bSJed Brown       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
177c4762a1bSJed Brown       ierr = (*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds,ctx->xmin+hx*i);CHKERRQ(ierr);
178c4762a1bSJed Brown       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
179c4762a1bSJed Brown       ierr  = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr);
180c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
181c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
182c4762a1bSJed Brown       for (j=0; j<dof; j++) {
183c4762a1bSJed Brown         PetscScalar jmpL,jmpR;
184c4762a1bSJed Brown         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
185c4762a1bSJed Brown         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
186c4762a1bSJed Brown         for (k=0; k<dof; k++) {
187c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
188c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
189c4762a1bSJed Brown         }
190c4762a1bSJed Brown       }
191c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
192c4762a1bSJed Brown       info.m  = dof;
193c4762a1bSJed Brown       info.hx = hx;
194c4762a1bSJed Brown       (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope);
195c4762a1bSJed Brown       for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */
196c4762a1bSJed Brown       for (j=0; j<dof; j++) {
197c4762a1bSJed Brown         PetscScalar tmp = 0;
198c4762a1bSJed Brown         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
199c4762a1bSJed Brown         slope[i*dof+j] = tmp;
200c4762a1bSJed Brown       }
201c4762a1bSJed Brown     }
202c4762a1bSJed Brown   }
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
205c4762a1bSJed Brown     PetscReal   maxspeed;
206c4762a1bSJed Brown     PetscScalar *uL,*uR;
207c4762a1bSJed Brown     if (i < len_slow) { /* slow parts can be changed based on a */
208c4762a1bSJed Brown       uL = &ctx->uLR[0];
209c4762a1bSJed Brown       uR = &ctx->uLR[dof];
210c4762a1bSJed Brown       for (j=0; j<dof; j++) {
211c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
212c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
213c4762a1bSJed Brown       }
214c4762a1bSJed Brown       ierr    = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr);
215c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
216c4762a1bSJed Brown       if (i > xs) {
217c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hx;
218c4762a1bSJed Brown       }
219c4762a1bSJed Brown       if (i < xs+xm) {
220c4762a1bSJed Brown         for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hx;
221c4762a1bSJed Brown       }
222c4762a1bSJed Brown     }
223c4762a1bSJed Brown     if (i == len_slow) { /* slow parts can be changed based on a */
224c4762a1bSJed Brown       uL = &ctx->uLR[0];
225c4762a1bSJed Brown       uR = &ctx->uLR[dof];
226c4762a1bSJed Brown       for (j=0; j<dof; j++) {
227c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
228c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
229c4762a1bSJed Brown       }
230c4762a1bSJed Brown       ierr    = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr);
231c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
232c4762a1bSJed Brown       if (i > xs) {
233c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hx;
234c4762a1bSJed Brown       }
235c4762a1bSJed Brown     }
236c4762a1bSJed Brown   }
237c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr);
238c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
239c4762a1bSJed Brown   ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
240c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr);
241c4762a1bSJed Brown   PetscFunctionReturn(0);
242c4762a1bSJed Brown }
243c4762a1bSJed Brown 
244c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for fast  parts ----------------------------------- */
245c4762a1bSJed Brown 
246c4762a1bSJed Brown PetscErrorCode FVRHSFunctionfast(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
247c4762a1bSJed Brown {
248c4762a1bSJed Brown   FVCtx          *ctx = (FVCtx*)vctx;
249c4762a1bSJed Brown   PetscErrorCode ierr;
250c4762a1bSJed Brown   PetscInt       i,j,k,Mx,dof,xs,xm,len_slow;
251c4762a1bSJed Brown   PetscReal      hx,cfl_idt = 0;
252c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
253c4762a1bSJed Brown   Vec            Xloc;
254c4762a1bSJed Brown   DM             da;
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   PetscFunctionBeginUser;
257c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
258c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr);
259c4762a1bSJed Brown   ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr);
260c4762a1bSJed Brown   hx   = (ctx->xmax-ctx->xmin)/Mx;
261c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
262c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
263c4762a1bSJed Brown   ierr = VecZeroEntries(F);CHKERRQ(ierr);
264c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr);
265c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
266c4762a1bSJed Brown   ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
267c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr);
268c4762a1bSJed Brown   ierr = ISGetSize(ctx->iss,&len_slow);CHKERRQ(ierr);
269c4762a1bSJed Brown 
270c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
271c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
272c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
273c4762a1bSJed Brown     }
274c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
275c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
276c4762a1bSJed Brown     }
277c4762a1bSJed Brown   }
278c4762a1bSJed Brown   for (i=xs-1; i<xs+xm+1; i++) {
279c4762a1bSJed Brown     struct _LimitInfo info;
280c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
281c4762a1bSJed Brown     if (i > len_slow-2) {
282c4762a1bSJed Brown       ierr = (*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds,ctx->xmin+hx*i);CHKERRQ(ierr);
283c4762a1bSJed Brown       ierr  = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr);
284c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
285c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
286c4762a1bSJed Brown       for (j=0; j<dof; j++) {
287c4762a1bSJed Brown         PetscScalar jmpL,jmpR;
288c4762a1bSJed Brown         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
289c4762a1bSJed Brown         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
290c4762a1bSJed Brown         for (k=0; k<dof; k++) {
291c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
292c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
293c4762a1bSJed Brown         }
294c4762a1bSJed Brown       }
295c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
296c4762a1bSJed Brown       info.m  = dof;
297c4762a1bSJed Brown       info.hx = hx;
298c4762a1bSJed Brown       (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope);
299c4762a1bSJed Brown       for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */
300c4762a1bSJed Brown       for (j=0; j<dof; j++) {
301c4762a1bSJed Brown         PetscScalar tmp = 0;
302c4762a1bSJed Brown         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
303c4762a1bSJed Brown         slope[i*dof+j] = tmp;
304c4762a1bSJed Brown       }
305c4762a1bSJed Brown     }
306c4762a1bSJed Brown   }
307c4762a1bSJed Brown 
308c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
309c4762a1bSJed Brown     PetscReal   maxspeed;
310c4762a1bSJed Brown     PetscScalar *uL,*uR;
311c4762a1bSJed Brown     if (i > len_slow) {
312c4762a1bSJed Brown       uL = &ctx->uLR[0];
313c4762a1bSJed Brown       uR = &ctx->uLR[dof];
314c4762a1bSJed Brown       for (j=0; j<dof; j++) {
315c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
316c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
317c4762a1bSJed Brown       }
318c4762a1bSJed Brown       ierr    = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr);
319c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
320c4762a1bSJed Brown       if (i > xs) {
321c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow-1)*dof+j] -= ctx->flux[j]/hx;
322c4762a1bSJed Brown       }
323c4762a1bSJed Brown       if (i < xs+xm) {
324c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow)*dof+j] += ctx->flux[j]/hx;
325c4762a1bSJed Brown       }
326c4762a1bSJed Brown     }
327c4762a1bSJed Brown     if (i == len_slow) {
328c4762a1bSJed Brown       uL = &ctx->uLR[0];
329c4762a1bSJed Brown       uR = &ctx->uLR[dof];
330c4762a1bSJed Brown       for (j=0; j<dof; j++) {
331c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
332c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
333c4762a1bSJed Brown       }
334c4762a1bSJed Brown       ierr    = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr);
335c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
336c4762a1bSJed Brown       if (i < xs+xm) {
337c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow)*dof+j] += ctx->flux[j]/hx;
338c4762a1bSJed Brown       }
339c4762a1bSJed Brown     }
340c4762a1bSJed Brown   }
341c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr);
342c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
343c4762a1bSJed Brown   ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
344c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr);
345c4762a1bSJed Brown   PetscFunctionReturn(0);
346c4762a1bSJed Brown }
347c4762a1bSJed Brown 
348c4762a1bSJed Brown PetscErrorCode FVRHSFunctionslow2(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
349c4762a1bSJed Brown {
350c4762a1bSJed Brown   FVCtx          *ctx = (FVCtx*)vctx;
351c4762a1bSJed Brown   PetscErrorCode ierr;
352c4762a1bSJed Brown   PetscInt       i,j,k,Mx,dof,xs,xm,len_slow1,len_slow2;
353c4762a1bSJed Brown   PetscReal      hx,cfl_idt = 0;
354c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
355c4762a1bSJed Brown   Vec            Xloc;
356c4762a1bSJed Brown   DM             da;
357c4762a1bSJed Brown 
358c4762a1bSJed Brown   PetscFunctionBeginUser;
359c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
360c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr);
361c4762a1bSJed Brown   ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr);
362c4762a1bSJed Brown   hx   = (ctx->xmax-ctx->xmin)/Mx;
363c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
364c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
365c4762a1bSJed Brown   ierr = VecZeroEntries(F);CHKERRQ(ierr);
366c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr);
367c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
368c4762a1bSJed Brown   ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
369c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr);
370c4762a1bSJed Brown   ierr = ISGetSize(ctx->iss,&len_slow1);CHKERRQ(ierr);
371c4762a1bSJed Brown   ierr = ISGetSize(ctx->iss2,&len_slow2);CHKERRQ(ierr);
372c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
373c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
374c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
375c4762a1bSJed Brown     }
376c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
377c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
378c4762a1bSJed Brown     }
379c4762a1bSJed Brown   }
380c4762a1bSJed Brown   for (i=xs-1; i<xs+xm+1; i++) {
381c4762a1bSJed Brown     struct _LimitInfo info;
382c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
383c4762a1bSJed Brown     if (i < len_slow1+len_slow2+1 && i > len_slow1-2) {
384c4762a1bSJed Brown       ierr = (*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds,ctx->xmin+hx*i);CHKERRQ(ierr);
385c4762a1bSJed Brown       ierr  = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr);
386c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
387c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
388c4762a1bSJed Brown       for (j=0; j<dof; j++) {
389c4762a1bSJed Brown         PetscScalar jmpL,jmpR;
390c4762a1bSJed Brown         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
391c4762a1bSJed Brown         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
392c4762a1bSJed Brown         for (k=0; k<dof; k++) {
393c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
394c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
395c4762a1bSJed Brown         }
396c4762a1bSJed Brown       }
397c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
398c4762a1bSJed Brown       info.m  = dof;
399c4762a1bSJed Brown       info.hx = hx;
400c4762a1bSJed Brown       (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope);
401c4762a1bSJed Brown       for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */
402c4762a1bSJed Brown       for (j=0; j<dof; j++) {
403c4762a1bSJed Brown         PetscScalar tmp = 0;
404c4762a1bSJed Brown         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
405c4762a1bSJed Brown         slope[i*dof+j] = tmp;
406c4762a1bSJed Brown       }
407c4762a1bSJed Brown     }
408c4762a1bSJed Brown   }
409c4762a1bSJed Brown 
410c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
411c4762a1bSJed Brown     PetscReal   maxspeed;
412c4762a1bSJed Brown     PetscScalar *uL,*uR;
413c4762a1bSJed Brown     if (i < len_slow1+len_slow2 && i > len_slow1) { /* slow parts can be changed based on a */
414c4762a1bSJed Brown       uL = &ctx->uLR[0];
415c4762a1bSJed Brown       uR = &ctx->uLR[dof];
416c4762a1bSJed Brown       for (j=0; j<dof; j++) {
417c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
418c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
419c4762a1bSJed Brown       }
420c4762a1bSJed Brown       ierr    = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr);
421c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
422c4762a1bSJed Brown       if (i > xs) {
423c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow1-1)*dof+j] -= ctx->flux[j]/hx;
424c4762a1bSJed Brown       }
425c4762a1bSJed Brown       if (i < xs+xm) {
426c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow1)*dof+j] += ctx->flux[j]/hx;
427c4762a1bSJed Brown       }
428c4762a1bSJed Brown     }
429c4762a1bSJed Brown     if (i == len_slow1+len_slow2) { /* slow parts can be changed based on a */
430c4762a1bSJed Brown       uL = &ctx->uLR[0];
431c4762a1bSJed Brown       uR = &ctx->uLR[dof];
432c4762a1bSJed Brown       for (j=0; j<dof; j++) {
433c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
434c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
435c4762a1bSJed Brown       }
436c4762a1bSJed Brown       ierr    = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr);
437c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
438c4762a1bSJed Brown       if (i > xs) {
439c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow1-1)*dof+j] -= ctx->flux[j]/hx;
440c4762a1bSJed Brown       }
441c4762a1bSJed Brown     }
442c4762a1bSJed Brown     if (i == len_slow1) { /* slow parts can be changed based on a */
443c4762a1bSJed Brown       uL = &ctx->uLR[0];
444c4762a1bSJed Brown       uR = &ctx->uLR[dof];
445c4762a1bSJed Brown       for (j=0; j<dof; j++) {
446c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
447c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
448c4762a1bSJed Brown       }
449c4762a1bSJed Brown       ierr    = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr);
450c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
451c4762a1bSJed Brown       if (i < xs+xm) {
452c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow1)*dof+j] += ctx->flux[j]/hx;
453c4762a1bSJed Brown       }
454c4762a1bSJed Brown     }
455c4762a1bSJed Brown   }
456c4762a1bSJed Brown 
457c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr);
458c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
459c4762a1bSJed Brown   ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
460c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr);
461c4762a1bSJed Brown   PetscFunctionReturn(0);
462c4762a1bSJed Brown }
463c4762a1bSJed Brown 
464c4762a1bSJed Brown PetscErrorCode FVRHSFunctionfast2(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
465c4762a1bSJed Brown {
466c4762a1bSJed Brown   FVCtx          *ctx = (FVCtx*)vctx;
467c4762a1bSJed Brown   PetscErrorCode ierr;
468c4762a1bSJed Brown   PetscInt       i,j,k,Mx,dof,xs,xm,len_slow1,len_slow2;
469c4762a1bSJed Brown   PetscReal      hx,cfl_idt = 0;
470c4762a1bSJed Brown   PetscScalar    *x,*f,*slope;
471c4762a1bSJed Brown   Vec            Xloc;
472c4762a1bSJed Brown   DM             da;
473c4762a1bSJed Brown 
474c4762a1bSJed Brown   PetscFunctionBeginUser;
475c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
476c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr);
477c4762a1bSJed Brown   ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr);
478c4762a1bSJed Brown   hx   = (ctx->xmax-ctx->xmin)/Mx;
479c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
480c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
481c4762a1bSJed Brown   ierr = VecZeroEntries(F);CHKERRQ(ierr);
482c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr);
483c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
484c4762a1bSJed Brown   ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
485c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr);
486c4762a1bSJed Brown   ierr = ISGetSize(ctx->iss,&len_slow1);CHKERRQ(ierr);
487c4762a1bSJed Brown   ierr = ISGetSize(ctx->iss2,&len_slow2);CHKERRQ(ierr);
488c4762a1bSJed Brown 
489c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
490c4762a1bSJed Brown     for (i=xs-2; i<0; i++) {
491c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
492c4762a1bSJed Brown     }
493c4762a1bSJed Brown     for (i=Mx; i<xs+xm+2; i++) {
494c4762a1bSJed Brown       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
495c4762a1bSJed Brown     }
496c4762a1bSJed Brown   }
497c4762a1bSJed Brown   for (i=xs-1; i<xs+xm+1; i++) {
498c4762a1bSJed Brown     struct _LimitInfo info;
499c4762a1bSJed Brown     PetscScalar       *cjmpL,*cjmpR;
500c4762a1bSJed Brown     if (i > len_slow1+len_slow2-2) {
501c4762a1bSJed Brown       ierr = (*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds,ctx->xmin+hx*i);CHKERRQ(ierr);
502c4762a1bSJed Brown       ierr  = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr);
503c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
504c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
505c4762a1bSJed Brown       for (j=0; j<dof; j++) {
506c4762a1bSJed Brown         PetscScalar jmpL,jmpR;
507c4762a1bSJed Brown         jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j];
508c4762a1bSJed Brown         jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j];
509c4762a1bSJed Brown         for (k=0; k<dof; k++) {
510c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL;
511c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR;
512c4762a1bSJed Brown         }
513c4762a1bSJed Brown       }
514c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
515c4762a1bSJed Brown       info.m  = dof;
516c4762a1bSJed Brown       info.hx = hx;
517c4762a1bSJed Brown       (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope);
518c4762a1bSJed Brown       for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */
519c4762a1bSJed Brown       for (j=0; j<dof; j++) {
520c4762a1bSJed Brown         PetscScalar tmp = 0;
521c4762a1bSJed Brown         for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k];
522c4762a1bSJed Brown         slope[i*dof+j] = tmp;
523c4762a1bSJed Brown       }
524c4762a1bSJed Brown     }
525c4762a1bSJed Brown   }
526c4762a1bSJed Brown 
527c4762a1bSJed Brown   for (i=xs; i<xs+xm+1; i++) {
528c4762a1bSJed Brown     PetscReal   maxspeed;
529c4762a1bSJed Brown     PetscScalar *uL,*uR;
530c4762a1bSJed Brown     if (i > len_slow1+len_slow2) {
531c4762a1bSJed Brown       uL = &ctx->uLR[0];
532c4762a1bSJed Brown       uR = &ctx->uLR[dof];
533c4762a1bSJed Brown       for (j=0; j<dof; j++) {
534c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
535c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
536c4762a1bSJed Brown       }
537c4762a1bSJed Brown       ierr    = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr);
538c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
539c4762a1bSJed Brown       if (i > xs) {
540c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow1-len_slow2-1)*dof+j] -= ctx->flux[j]/hx;
541c4762a1bSJed Brown       }
542c4762a1bSJed Brown       if (i < xs+xm) {
543c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow1-len_slow2)*dof+j] += ctx->flux[j]/hx;
544c4762a1bSJed Brown       }
545c4762a1bSJed Brown     }
546c4762a1bSJed Brown     if (i == len_slow1+len_slow2) {
547c4762a1bSJed Brown       uL = &ctx->uLR[0];
548c4762a1bSJed Brown       uR = &ctx->uLR[dof];
549c4762a1bSJed Brown       for (j=0; j<dof; j++) {
550c4762a1bSJed Brown         uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2;
551c4762a1bSJed Brown         uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2;
552c4762a1bSJed Brown       }
553c4762a1bSJed Brown       ierr    = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr);
554c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
555c4762a1bSJed Brown       if (i < xs+xm) {
556c4762a1bSJed Brown         for (j=0; j<dof; j++) f[(i-len_slow1-len_slow2)*dof+j] += ctx->flux[j]/hx;
557c4762a1bSJed Brown       }
558c4762a1bSJed Brown     }
559c4762a1bSJed Brown   }
560c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr);
561c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
562c4762a1bSJed Brown   ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr);
563c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr);
564c4762a1bSJed Brown   PetscFunctionReturn(0);
565c4762a1bSJed Brown }
566c4762a1bSJed Brown 
567c4762a1bSJed Brown int main(int argc,char *argv[])
568c4762a1bSJed Brown {
569c4762a1bSJed Brown   char              lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m";
570c4762a1bSJed Brown   PetscFunctionList limiters   = 0,physics = 0;
571c4762a1bSJed Brown   MPI_Comm          comm;
572c4762a1bSJed Brown   TS                ts;
573c4762a1bSJed Brown   DM                da;
574c4762a1bSJed Brown   Vec               X,X0,R;
575c4762a1bSJed Brown   FVCtx             ctx;
576c4762a1bSJed Brown   PetscInt          i,k,dof,xs,xm,Mx,draw = 0,*index_slow,*index_fast,islow = 0,ifast = 0;
577c4762a1bSJed Brown   PetscBool         view_final = PETSC_FALSE;
578c4762a1bSJed Brown   PetscReal         ptime;
579c4762a1bSJed Brown   PetscErrorCode    ierr;
580c4762a1bSJed Brown 
581c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
582c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
583c4762a1bSJed Brown   ierr = PetscMemzero(&ctx,sizeof(ctx));CHKERRQ(ierr);
584c4762a1bSJed Brown 
585c4762a1bSJed Brown   /* Register limiters to be available on the command line */
586c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"upwind"              ,Limit_Upwind);CHKERRQ(ierr);
587c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"lax-wendroff"        ,Limit_LaxWendroff);CHKERRQ(ierr);
588c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"beam-warming"        ,Limit_BeamWarming);CHKERRQ(ierr);
589c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"fromm"               ,Limit_Fromm);CHKERRQ(ierr);
590c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"minmod"              ,Limit_Minmod);CHKERRQ(ierr);
591c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"superbee"            ,Limit_Superbee);CHKERRQ(ierr);
592c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"mc"                  ,Limit_MC);CHKERRQ(ierr);
593c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"vanleer"             ,Limit_VanLeer);CHKERRQ(ierr);
594c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"vanalbada"           ,Limit_VanAlbada);CHKERRQ(ierr);
595c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"vanalbadatvd"        ,Limit_VanAlbadaTVD);CHKERRQ(ierr);
596c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"koren"               ,Limit_Koren);CHKERRQ(ierr);
597c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"korensym"            ,Limit_KorenSym);CHKERRQ(ierr);
598c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"koren3"              ,Limit_Koren3);CHKERRQ(ierr);
599c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"cada-torrilhon2"     ,Limit_CadaTorrilhon2);CHKERRQ(ierr);
600c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"cada-torrilhon3-r0p1",Limit_CadaTorrilhon3R0p1);CHKERRQ(ierr);
601c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"cada-torrilhon3-r1"  ,Limit_CadaTorrilhon3R1);CHKERRQ(ierr);
602c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"cada-torrilhon3-r10" ,Limit_CadaTorrilhon3R10);CHKERRQ(ierr);
603c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&limiters,"cada-torrilhon3-r100",Limit_CadaTorrilhon3R100);CHKERRQ(ierr);
604c4762a1bSJed Brown 
605c4762a1bSJed Brown   /* Register physical models to be available on the command line */
606c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&physics,"advect"          ,PhysicsCreate_Advect);CHKERRQ(ierr);
607c4762a1bSJed Brown 
608c4762a1bSJed Brown   ctx.comm = comm;
609c4762a1bSJed Brown   ctx.cfl  = 0.9;
610c4762a1bSJed Brown   ctx.bctype = FVBC_PERIODIC;
611c4762a1bSJed Brown   ctx.xmin = -1.0;
612c4762a1bSJed Brown   ctx.xmax = 1.0;
613c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm,NULL,"Finite Volume solver options","");CHKERRQ(ierr);
614c4762a1bSJed Brown   ierr = PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL);CHKERRQ(ierr);
615c4762a1bSJed Brown   ierr = PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL);CHKERRQ(ierr);
616c4762a1bSJed Brown   ierr = PetscOptionsFList("-limit","Name of flux limiter to use","",limiters,lname,lname,sizeof(lname),NULL);CHKERRQ(ierr);
617c4762a1bSJed Brown   ierr = PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL);CHKERRQ(ierr);
618c4762a1bSJed Brown   ierr = PetscOptionsString("-view_final","Write final solution in ASCII MATLAB format to given file name","",final_fname,final_fname,sizeof(final_fname),&view_final);CHKERRQ(ierr);
619c4762a1bSJed Brown   ierr = PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL);CHKERRQ(ierr);
620c4762a1bSJed Brown   ierr = PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL);CHKERRQ(ierr);
621c4762a1bSJed Brown   ierr = PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL);CHKERRQ(ierr);
622c4762a1bSJed Brown   ierr = PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL);CHKERRQ(ierr);
623c4762a1bSJed Brown   ierr = PetscOptionsBool("-recursive_split","Split the domain recursively","",ctx.recursive,&ctx.recursive,NULL);CHKERRQ(ierr);
624c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
625c4762a1bSJed Brown 
626c4762a1bSJed Brown   /* Choose the limiter from the list of registered limiters */
627c4762a1bSJed Brown   ierr = PetscFunctionListFind(limiters,lname,&ctx.limit);CHKERRQ(ierr);
628c4762a1bSJed Brown   if (!ctx.limit) SETERRQ1(PETSC_COMM_SELF,1,"Limiter '%s' not found",lname);
629c4762a1bSJed Brown 
630c4762a1bSJed Brown   /* Choose the physics from the list of registered models */
631c4762a1bSJed Brown   {
632c4762a1bSJed Brown     PetscErrorCode (*r)(FVCtx*);
633c4762a1bSJed Brown     ierr = PetscFunctionListFind(physics,physname,&r);CHKERRQ(ierr);
634c4762a1bSJed Brown     if (!r) SETERRQ1(PETSC_COMM_SELF,1,"Physics '%s' not found",physname);
635c4762a1bSJed Brown     /* Create the physics, will set the number of fields and their names */
636c4762a1bSJed Brown     ierr = (*r)(&ctx);CHKERRQ(ierr);
637c4762a1bSJed Brown   }
638c4762a1bSJed Brown 
639c4762a1bSJed Brown   /* Create a DMDA to manage the parallel grid */
640c4762a1bSJed Brown   ierr = DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics.dof,2,NULL,&da);CHKERRQ(ierr);
641c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
642c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
643c4762a1bSJed Brown   /* Inform the DMDA of the field names provided by the physics. */
644c4762a1bSJed Brown   /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
645c4762a1bSJed Brown   for (i=0; i<ctx.physics.dof; i++) {
646c4762a1bSJed Brown     ierr = DMDASetFieldName(da,i,ctx.physics.fieldname[i]);CHKERRQ(ierr);
647c4762a1bSJed Brown   }
648c4762a1bSJed Brown   ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr);
649c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr);
650c4762a1bSJed Brown 
651c4762a1bSJed Brown   /* Set coordinates of cell centers */
652c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(da,ctx.xmin+0.5*(ctx.xmax-ctx.xmin)/Mx,ctx.xmax+0.5*(ctx.xmax-ctx.xmin)/Mx,0,0,0,0);CHKERRQ(ierr);
653c4762a1bSJed Brown 
654c4762a1bSJed Brown   /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
655c4762a1bSJed Brown   ierr = PetscMalloc4(dof*dof,&ctx.R,dof*dof,&ctx.Rinv,2*dof,&ctx.cjmpLR,1*dof,&ctx.cslope);CHKERRQ(ierr);
656c4762a1bSJed Brown   ierr = PetscMalloc3(2*dof,&ctx.uLR,dof,&ctx.flux,dof,&ctx.speeds);CHKERRQ(ierr);
657c4762a1bSJed Brown 
658c4762a1bSJed Brown   /* Create a vector to store the solution and to save the initial state */
659c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&X);CHKERRQ(ierr);
660c4762a1bSJed Brown   ierr = VecDuplicate(X,&X0);CHKERRQ(ierr);
661c4762a1bSJed Brown   ierr = VecDuplicate(X,&R);CHKERRQ(ierr);
662c4762a1bSJed Brown 
663c4762a1bSJed Brown   /*-------------------------------- create index for slow parts and fast parts ----------------------------------------*/
664c4762a1bSJed Brown   ierr = PetscMalloc1(xm*dof,&index_slow);CHKERRQ(ierr);
665c4762a1bSJed Brown   ierr = PetscMalloc1(xm*dof,&index_fast);CHKERRQ(ierr);
666c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
667c4762a1bSJed Brown     if (ctx.xmin+i*(ctx.xmax-ctx.xmin)/(PetscReal)Mx+0.5*(ctx.xmax-ctx.xmin)/(PetscReal)Mx < 0)
668c4762a1bSJed Brown       for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k;
669c4762a1bSJed Brown     else
670c4762a1bSJed Brown       for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k;
671c4762a1bSJed Brown   }  /* this step need to be changed based on discontinuous point of a */
672c4762a1bSJed Brown   ierr = ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss);CHKERRQ(ierr);
673c4762a1bSJed Brown   ierr = ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf);CHKERRQ(ierr);
674c4762a1bSJed Brown 
675c4762a1bSJed Brown   /* Create a time-stepping object */
676c4762a1bSJed Brown   ierr = TSCreate(comm,&ts);CHKERRQ(ierr);
677c4762a1bSJed Brown   ierr = TSSetDM(ts,da);CHKERRQ(ierr);
678c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,R,FVRHSFunction,&ctx);CHKERRQ(ierr);
679c4762a1bSJed Brown   ierr = TSRHSSplitSetIS(ts,"slow",ctx.iss);CHKERRQ(ierr);
680c4762a1bSJed Brown   ierr = TSRHSSplitSetIS(ts,"fast",ctx.isf);CHKERRQ(ierr);
681c4762a1bSJed Brown   ierr = TSRHSSplitSetRHSFunction(ts,"slow",NULL,FVRHSFunctionslow,&ctx);CHKERRQ(ierr);
682c4762a1bSJed Brown   ierr = TSRHSSplitSetRHSFunction(ts,"fast",NULL,FVRHSFunctionfast,&ctx);CHKERRQ(ierr);
683c4762a1bSJed Brown 
684c4762a1bSJed Brown   if (ctx.recursive) {
685c4762a1bSJed Brown     TS subts;
686c4762a1bSJed Brown     islow = 0;
687c4762a1bSJed Brown     ifast = 0;
688c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
689c4762a1bSJed Brown       PetscReal coord = ctx.xmin+i*(ctx.xmax-ctx.xmin)/(PetscReal)Mx+0.5*(ctx.xmax-ctx.xmin)/(PetscReal)Mx;
690c4762a1bSJed Brown       if (coord >= 0 && coord < ctx.xmin+(ctx.xmax-ctx.xmin)*3/4.)
691c4762a1bSJed Brown         for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k;
692c4762a1bSJed Brown       if (coord >= ctx.xmin+(ctx.xmax-ctx.xmin)*3/4.)
693c4762a1bSJed Brown         for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k;
694c4762a1bSJed Brown     }  /* this step need to be changed based on discontinuous point of a */
695c4762a1bSJed Brown     ierr = ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss2);CHKERRQ(ierr);
696c4762a1bSJed Brown     ierr = ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf2);CHKERRQ(ierr);
697c4762a1bSJed Brown 
698c4762a1bSJed Brown     ierr = TSRHSSplitGetSubTS(ts,"fast",&subts);CHKERRQ(ierr);
699c4762a1bSJed Brown     ierr = TSRHSSplitSetIS(subts,"slow",ctx.iss2);CHKERRQ(ierr);
700c4762a1bSJed Brown     ierr = TSRHSSplitSetIS(subts,"fast",ctx.isf2);CHKERRQ(ierr);
701c4762a1bSJed Brown     ierr = TSRHSSplitSetRHSFunction(subts,"slow",NULL,FVRHSFunctionslow2,&ctx);CHKERRQ(ierr);
702c4762a1bSJed Brown     ierr = TSRHSSplitSetRHSFunction(subts,"fast",NULL,FVRHSFunctionfast2,&ctx);CHKERRQ(ierr);
703c4762a1bSJed Brown   }
704c4762a1bSJed Brown 
705c4762a1bSJed Brown   /*ierr = TSSetType(ts,TSSSP);CHKERRQ(ierr);*/
706c4762a1bSJed Brown   ierr = TSSetType(ts,TSMPRK);CHKERRQ(ierr);
707c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,10);CHKERRQ(ierr);
708c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
709c4762a1bSJed Brown 
710c4762a1bSJed Brown   /* Compute initial conditions and starting time step */
711c4762a1bSJed Brown   ierr = FVSample(&ctx,da,0,X0);CHKERRQ(ierr);
712c4762a1bSJed Brown   ierr = FVRHSFunction(ts,0,X0,X,(void*)&ctx);CHKERRQ(ierr); /* Initial function evaluation, only used to determine max speed */
713c4762a1bSJed Brown   ierr = VecCopy(X0,X);CHKERRQ(ierr);                        /* The function value was not used so we set X=X0 again */
714c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt);CHKERRQ(ierr);
715c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr); /* Take runtime options */
716c4762a1bSJed Brown   ierr = SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
717c4762a1bSJed Brown   {
718c4762a1bSJed Brown     PetscInt    steps;
719c4762a1bSJed Brown     PetscScalar mass_initial,mass_final,mass_difference;
720c4762a1bSJed Brown 
721c4762a1bSJed Brown     ierr = TSSolve(ts,X);CHKERRQ(ierr);
722c4762a1bSJed Brown     ierr = TSGetSolveTime(ts,&ptime);CHKERRQ(ierr);
723c4762a1bSJed Brown     ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
724c4762a1bSJed Brown     ierr = PetscPrintf(comm,"Final time %g, steps %D\n",(double)ptime,steps);CHKERRQ(ierr);
725c4762a1bSJed Brown     /* calculate the total mass at initial time and final time */
726c4762a1bSJed Brown     mass_initial = 0.0;
727c4762a1bSJed Brown     mass_final   = 0.0;
728c4762a1bSJed Brown     ierr = VecSum(X0,&mass_initial);CHKERRQ(ierr);
729c4762a1bSJed Brown     ierr = VecSum(X,&mass_final);CHKERRQ(ierr);
730c4762a1bSJed Brown     mass_difference = (ctx.xmax-ctx.xmin)/(PetscScalar)Mx*(mass_final - mass_initial);
731c4762a1bSJed Brown     ierr = PetscPrintf(comm,"Mass difference %g\n",(double)mass_difference);CHKERRQ(ierr);
732c4762a1bSJed Brown     if (ctx.simulation) {
733c4762a1bSJed Brown       PetscViewer  fd;
734c4762a1bSJed Brown       char         filename[PETSC_MAX_PATH_LEN] = "binaryoutput";
735c4762a1bSJed Brown       Vec          XR;
736c4762a1bSJed Brown       PetscReal    nrm1,nrmsup;
737c4762a1bSJed Brown       PetscBool    flg;
738c4762a1bSJed Brown       ierr = PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg);CHKERRQ(ierr);
739c4762a1bSJed Brown       if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
740c4762a1bSJed Brown       ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&fd);CHKERRQ(ierr);
741c4762a1bSJed Brown       ierr = VecDuplicate(X0,&XR);CHKERRQ(ierr);
742c4762a1bSJed Brown       ierr = VecLoad(XR,fd);CHKERRQ(ierr);
743c4762a1bSJed Brown       ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
744c4762a1bSJed Brown       ierr = VecAYPX(XR,-1,X);CHKERRQ(ierr);
745c4762a1bSJed Brown       ierr = VecNorm(XR,NORM_1,&nrm1);CHKERRQ(ierr);
746c4762a1bSJed Brown       ierr = VecNorm(XR,NORM_INFINITY,&nrmsup);CHKERRQ(ierr);
747c4762a1bSJed Brown       nrm1 /= Mx;
748c4762a1bSJed Brown       ierr = PetscPrintf(comm,"Error ||x-x_e||_1 %g  ||x-x_e||_sup %g\n",(double)nrm1,(double)nrmsup);CHKERRQ(ierr);
749c4762a1bSJed Brown       ierr = VecDestroy(&XR);CHKERRQ(ierr);
750c4762a1bSJed Brown     }
751c4762a1bSJed Brown   }
752c4762a1bSJed Brown 
753c4762a1bSJed Brown   ierr = SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
754c4762a1bSJed Brown   if (draw & 0x1) {ierr = VecView(X0,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);}
755c4762a1bSJed Brown   if (draw & 0x2) {ierr = VecView(X,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);}
756c4762a1bSJed Brown   if (draw & 0x4) {
757c4762a1bSJed Brown     Vec Y;
758c4762a1bSJed Brown     ierr = VecDuplicate(X,&Y);CHKERRQ(ierr);
759c4762a1bSJed Brown     ierr = FVSample(&ctx,da,ptime,Y);CHKERRQ(ierr);
760c4762a1bSJed Brown     ierr = VecAYPX(Y,-1,X);CHKERRQ(ierr);
761c4762a1bSJed Brown     ierr = VecView(Y,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);
762c4762a1bSJed Brown     ierr = VecDestroy(&Y);CHKERRQ(ierr);
763c4762a1bSJed Brown   }
764c4762a1bSJed Brown 
765c4762a1bSJed Brown   if (view_final) {
766c4762a1bSJed Brown     PetscViewer viewer;
767c4762a1bSJed Brown     ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer);CHKERRQ(ierr);
768c4762a1bSJed Brown     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
769c4762a1bSJed Brown     ierr = VecView(X,viewer);CHKERRQ(ierr);
770c4762a1bSJed Brown     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
771c4762a1bSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
772c4762a1bSJed Brown   }
773c4762a1bSJed Brown 
774c4762a1bSJed Brown   /* Clean up */
775c4762a1bSJed Brown   ierr = (*ctx.physics.destroy)(ctx.physics.user);CHKERRQ(ierr);
776c4762a1bSJed Brown   for (i=0; i<ctx.physics.dof; i++) {ierr = PetscFree(ctx.physics.fieldname[i]);CHKERRQ(ierr);}
777c4762a1bSJed Brown   ierr = PetscFree4(ctx.R,ctx.Rinv,ctx.cjmpLR,ctx.cslope);CHKERRQ(ierr);
778c4762a1bSJed Brown   ierr = PetscFree3(ctx.uLR,ctx.flux,ctx.speeds);CHKERRQ(ierr);
779c4762a1bSJed Brown   ierr = ISDestroy(&ctx.iss);CHKERRQ(ierr);
780c4762a1bSJed Brown   ierr = ISDestroy(&ctx.isf);CHKERRQ(ierr);
781c4762a1bSJed Brown   ierr = ISDestroy(&ctx.iss2);CHKERRQ(ierr);
782c4762a1bSJed Brown   ierr = ISDestroy(&ctx.isf2);CHKERRQ(ierr);
783c4762a1bSJed Brown   ierr = VecDestroy(&X);CHKERRQ(ierr);
784c4762a1bSJed Brown   ierr = VecDestroy(&X0);CHKERRQ(ierr);
785c4762a1bSJed Brown   ierr = VecDestroy(&R);CHKERRQ(ierr);
786c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
787c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
788c4762a1bSJed Brown   ierr = PetscFree(index_slow);CHKERRQ(ierr);
789c4762a1bSJed Brown   ierr = PetscFree(index_fast);CHKERRQ(ierr);
790c4762a1bSJed Brown   ierr = PetscFunctionListDestroy(&limiters);CHKERRQ(ierr);
791c4762a1bSJed Brown   ierr = PetscFunctionListDestroy(&physics);CHKERRQ(ierr);
792c4762a1bSJed Brown   ierr = PetscFinalize();
793c4762a1bSJed Brown   return ierr;
794c4762a1bSJed Brown }
795c4762a1bSJed Brown 
796c4762a1bSJed Brown /*TEST
797c4762a1bSJed Brown     build:
798*f56ea12dSJed Brown       requires: !complex
799c4762a1bSJed Brown       depends: finitevolume1d.c
800c4762a1bSJed Brown 
801c4762a1bSJed Brown     test:
802c4762a1bSJed 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
803c4762a1bSJed Brown 
804c4762a1bSJed Brown     test:
805c4762a1bSJed Brown       suffix: 2
806c4762a1bSJed 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
807c4762a1bSJed Brown       output_file: output/ex5_1.out
808c4762a1bSJed Brown 
809c4762a1bSJed Brown     test:
810c4762a1bSJed Brown       suffix: 3
811c4762a1bSJed 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
812c4762a1bSJed Brown 
813c4762a1bSJed Brown     test:
814c4762a1bSJed Brown       suffix: 4
815c4762a1bSJed 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
816c4762a1bSJed Brown       output_file: output/ex5_3.out
817c4762a1bSJed Brown 
818c4762a1bSJed Brown     test:
819c4762a1bSJed Brown       suffix: 5
820c4762a1bSJed 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
821c4762a1bSJed Brown 
822c4762a1bSJed Brown     test:
823c4762a1bSJed Brown       suffix: 6
824c4762a1bSJed 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
825c4762a1bSJed Brown TEST*/
826