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 AdvectCtx *user; 116c4762a1bSJed Brown 117c4762a1bSJed Brown PetscFunctionBeginUser; 1189566063dSJacob Faibussowitsch PetscCall(PetscNew(&user)); 119c4762a1bSJed Brown ctx->physics.sample = PhysicsSample_Advect; 120c4762a1bSJed Brown ctx->physics.riemann = PhysicsRiemann_Advect; 121c4762a1bSJed Brown ctx->physics.characteristic = PhysicsCharacteristic_Advect; 122c4762a1bSJed Brown ctx->physics.destroy = PhysicsDestroy_SimpleFree; 123c4762a1bSJed Brown ctx->physics.user = user; 124c4762a1bSJed Brown ctx->physics.dof = 1; 1259566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("u",&ctx->physics.fieldname[0])); 126c4762a1bSJed Brown user->a[0] = 1; 127c4762a1bSJed Brown user->a[1] = 1; 128d0609cedSBarry Smith PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection",""); 129c4762a1bSJed Brown { 1309566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-physics_advect_a1","Speed1","",user->a[0],&user->a[0],NULL)); 1319566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-physics_advect_a2","Speed2","",user->a[1],&user->a[1],NULL)); 132c4762a1bSJed Brown } 133d0609cedSBarry Smith PetscOptionsEnd(); 134c4762a1bSJed Brown PetscFunctionReturn(0); 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow parts ----------------------------------- */ 138c4762a1bSJed Brown 139c4762a1bSJed Brown PetscErrorCode FVRHSFunctionslow(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 140c4762a1bSJed Brown { 141c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 142c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,len_slow; 143c4762a1bSJed Brown PetscReal hx,cfl_idt = 0; 144c4762a1bSJed Brown PetscScalar *x,*f,*slope; 145c4762a1bSJed Brown Vec Xloc; 146c4762a1bSJed Brown DM da; 147c4762a1bSJed Brown 148c4762a1bSJed Brown PetscFunctionBeginUser; 1499566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 1509566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); 1519566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 152c4762a1bSJed Brown hx = (ctx->xmax-ctx->xmin)/Mx; 1539566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 1549566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 1559566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 1569566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 1579566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 1589566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da,PETSC_TRUE,&slope)); 1599566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 1609566063dSJacob Faibussowitsch PetscCall(ISGetSize(ctx->iss,&len_slow)); 161c4762a1bSJed Brown 162c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 163c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 164c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 165c4762a1bSJed Brown } 166c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 167c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 168c4762a1bSJed Brown } 169c4762a1bSJed Brown } 170c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 171c4762a1bSJed Brown struct _LimitInfo info; 172c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 173c4762a1bSJed Brown if (i < len_slow+1) { 174c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 1759566063dSJacob Faibussowitsch PetscCall((*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds,ctx->xmin+hx*i)); 176c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 1779566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR,2*dof)); 178c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 179c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 180c4762a1bSJed Brown for (j=0; j<dof; j++) { 181c4762a1bSJed Brown PetscScalar jmpL,jmpR; 182c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 183c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 184c4762a1bSJed Brown for (k=0; k<dof; k++) { 185c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 186c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 187c4762a1bSJed Brown } 188c4762a1bSJed Brown } 189c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 190c4762a1bSJed Brown info.m = dof; 191c4762a1bSJed Brown info.hx = hx; 192c4762a1bSJed Brown (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope); 193c4762a1bSJed Brown for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 194c4762a1bSJed Brown for (j=0; j<dof; j++) { 195c4762a1bSJed Brown PetscScalar tmp = 0; 196c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 197c4762a1bSJed Brown slope[i*dof+j] = tmp; 198c4762a1bSJed Brown } 199c4762a1bSJed Brown } 200c4762a1bSJed Brown } 201c4762a1bSJed Brown 202c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 203c4762a1bSJed Brown PetscReal maxspeed; 204c4762a1bSJed Brown PetscScalar *uL,*uR; 205c4762a1bSJed Brown if (i < len_slow) { /* slow parts can be changed based on a */ 206c4762a1bSJed Brown uL = &ctx->uLR[0]; 207c4762a1bSJed Brown uR = &ctx->uLR[dof]; 208c4762a1bSJed Brown for (j=0; j<dof; j++) { 209c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 210c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 211c4762a1bSJed Brown } 2129566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax)); 213c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 214c4762a1bSJed Brown if (i > xs) { 215c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hx; 216c4762a1bSJed Brown } 217c4762a1bSJed Brown if (i < xs+xm) { 218c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hx; 219c4762a1bSJed Brown } 220c4762a1bSJed Brown } 221c4762a1bSJed Brown if (i == len_slow) { /* slow parts can be changed based on a */ 222c4762a1bSJed Brown uL = &ctx->uLR[0]; 223c4762a1bSJed Brown uR = &ctx->uLR[dof]; 224c4762a1bSJed Brown for (j=0; j<dof; j++) { 225c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 226c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 227c4762a1bSJed Brown } 2289566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax)); 229c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 230c4762a1bSJed Brown if (i > xs) { 231c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hx; 232c4762a1bSJed Brown } 233c4762a1bSJed Brown } 234c4762a1bSJed Brown } 2359566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 2369566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 2379566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da,PETSC_TRUE,&slope)); 2389566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 239c4762a1bSJed Brown PetscFunctionReturn(0); 240c4762a1bSJed Brown } 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for fast parts ----------------------------------- */ 243c4762a1bSJed Brown 244c4762a1bSJed Brown PetscErrorCode FVRHSFunctionfast(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 245c4762a1bSJed Brown { 246c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 247c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,len_slow; 248c4762a1bSJed Brown PetscReal hx,cfl_idt = 0; 249c4762a1bSJed Brown PetscScalar *x,*f,*slope; 250c4762a1bSJed Brown Vec Xloc; 251c4762a1bSJed Brown DM da; 252c4762a1bSJed Brown 253c4762a1bSJed Brown PetscFunctionBeginUser; 2549566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 2559566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); 2569566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 257c4762a1bSJed Brown hx = (ctx->xmax-ctx->xmin)/Mx; 2589566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 2599566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 2609566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 2619566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 2629566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 2639566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da,PETSC_TRUE,&slope)); 2649566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 2659566063dSJacob Faibussowitsch PetscCall(ISGetSize(ctx->iss,&len_slow)); 266c4762a1bSJed Brown 267c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 268c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 269c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 270c4762a1bSJed Brown } 271c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 272c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 273c4762a1bSJed Brown } 274c4762a1bSJed Brown } 275c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 276c4762a1bSJed Brown struct _LimitInfo info; 277c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 278c4762a1bSJed Brown if (i > len_slow-2) { 2799566063dSJacob Faibussowitsch PetscCall((*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds,ctx->xmin+hx*i)); 2809566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR,2*dof)); 281c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 282c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 283c4762a1bSJed Brown for (j=0; j<dof; j++) { 284c4762a1bSJed Brown PetscScalar jmpL,jmpR; 285c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 286c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 287c4762a1bSJed Brown for (k=0; k<dof; k++) { 288c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 289c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 290c4762a1bSJed Brown } 291c4762a1bSJed Brown } 292c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 293c4762a1bSJed Brown info.m = dof; 294c4762a1bSJed Brown info.hx = hx; 295c4762a1bSJed Brown (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope); 296c4762a1bSJed Brown for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 297c4762a1bSJed Brown for (j=0; j<dof; j++) { 298c4762a1bSJed Brown PetscScalar tmp = 0; 299c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 300c4762a1bSJed Brown slope[i*dof+j] = tmp; 301c4762a1bSJed Brown } 302c4762a1bSJed Brown } 303c4762a1bSJed Brown } 304c4762a1bSJed Brown 305c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 306c4762a1bSJed Brown PetscReal maxspeed; 307c4762a1bSJed Brown PetscScalar *uL,*uR; 308c4762a1bSJed Brown if (i > len_slow) { 309c4762a1bSJed Brown uL = &ctx->uLR[0]; 310c4762a1bSJed Brown uR = &ctx->uLR[dof]; 311c4762a1bSJed Brown for (j=0; j<dof; j++) { 312c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 313c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 314c4762a1bSJed Brown } 3159566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax)); 316c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 317c4762a1bSJed Brown if (i > xs) { 318c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow-1)*dof+j] -= ctx->flux[j]/hx; 319c4762a1bSJed Brown } 320c4762a1bSJed Brown if (i < xs+xm) { 321c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow)*dof+j] += ctx->flux[j]/hx; 322c4762a1bSJed Brown } 323c4762a1bSJed Brown } 324c4762a1bSJed Brown if (i == len_slow) { 325c4762a1bSJed Brown uL = &ctx->uLR[0]; 326c4762a1bSJed Brown uR = &ctx->uLR[dof]; 327c4762a1bSJed Brown for (j=0; j<dof; j++) { 328c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 329c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 330c4762a1bSJed Brown } 3319566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax)); 332c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 333c4762a1bSJed Brown if (i < xs+xm) { 334c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow)*dof+j] += ctx->flux[j]/hx; 335c4762a1bSJed Brown } 336c4762a1bSJed Brown } 337c4762a1bSJed Brown } 3389566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 3399566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 3409566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da,PETSC_TRUE,&slope)); 3419566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 342c4762a1bSJed Brown PetscFunctionReturn(0); 343c4762a1bSJed Brown } 344c4762a1bSJed Brown 345c4762a1bSJed Brown PetscErrorCode FVRHSFunctionslow2(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 346c4762a1bSJed Brown { 347c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 348c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,len_slow1,len_slow2; 349c4762a1bSJed Brown PetscReal hx,cfl_idt = 0; 350c4762a1bSJed Brown PetscScalar *x,*f,*slope; 351c4762a1bSJed Brown Vec Xloc; 352c4762a1bSJed Brown DM da; 353c4762a1bSJed Brown 354c4762a1bSJed Brown PetscFunctionBeginUser; 3559566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 3569566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); 3579566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 358c4762a1bSJed Brown hx = (ctx->xmax-ctx->xmin)/Mx; 3599566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 3609566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 3619566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 3629566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 3639566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 3649566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da,PETSC_TRUE,&slope)); 3659566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 3669566063dSJacob Faibussowitsch PetscCall(ISGetSize(ctx->iss,&len_slow1)); 3679566063dSJacob Faibussowitsch PetscCall(ISGetSize(ctx->iss2,&len_slow2)); 368c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 369c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 370c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 371c4762a1bSJed Brown } 372c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 373c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 374c4762a1bSJed Brown } 375c4762a1bSJed Brown } 376c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 377c4762a1bSJed Brown struct _LimitInfo info; 378c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 379c4762a1bSJed Brown if (i < len_slow1+len_slow2+1 && i > len_slow1-2) { 3809566063dSJacob Faibussowitsch PetscCall((*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds,ctx->xmin+hx*i)); 3819566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR,2*dof)); 382c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 383c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 384c4762a1bSJed Brown for (j=0; j<dof; j++) { 385c4762a1bSJed Brown PetscScalar jmpL,jmpR; 386c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 387c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 388c4762a1bSJed Brown for (k=0; k<dof; k++) { 389c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 390c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 391c4762a1bSJed Brown } 392c4762a1bSJed Brown } 393c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 394c4762a1bSJed Brown info.m = dof; 395c4762a1bSJed Brown info.hx = hx; 396c4762a1bSJed Brown (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope); 397c4762a1bSJed Brown for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 398c4762a1bSJed Brown for (j=0; j<dof; j++) { 399c4762a1bSJed Brown PetscScalar tmp = 0; 400c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 401c4762a1bSJed Brown slope[i*dof+j] = tmp; 402c4762a1bSJed Brown } 403c4762a1bSJed Brown } 404c4762a1bSJed Brown } 405c4762a1bSJed Brown 406c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 407c4762a1bSJed Brown PetscReal maxspeed; 408c4762a1bSJed Brown PetscScalar *uL,*uR; 409c4762a1bSJed Brown if (i < len_slow1+len_slow2 && i > len_slow1) { /* slow parts can be changed based on a */ 410c4762a1bSJed Brown uL = &ctx->uLR[0]; 411c4762a1bSJed Brown uR = &ctx->uLR[dof]; 412c4762a1bSJed Brown for (j=0; j<dof; j++) { 413c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 414c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 415c4762a1bSJed Brown } 4169566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax)); 417c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 418c4762a1bSJed Brown if (i > xs) { 419c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow1-1)*dof+j] -= ctx->flux[j]/hx; 420c4762a1bSJed Brown } 421c4762a1bSJed Brown if (i < xs+xm) { 422c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow1)*dof+j] += ctx->flux[j]/hx; 423c4762a1bSJed Brown } 424c4762a1bSJed Brown } 425c4762a1bSJed Brown if (i == len_slow1+len_slow2) { /* slow parts can be changed based on a */ 426c4762a1bSJed Brown uL = &ctx->uLR[0]; 427c4762a1bSJed Brown uR = &ctx->uLR[dof]; 428c4762a1bSJed Brown for (j=0; j<dof; j++) { 429c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 430c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 431c4762a1bSJed Brown } 4329566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax)); 433c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 434c4762a1bSJed Brown if (i > xs) { 435c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow1-1)*dof+j] -= ctx->flux[j]/hx; 436c4762a1bSJed Brown } 437c4762a1bSJed Brown } 438c4762a1bSJed Brown if (i == len_slow1) { /* slow parts can be changed based on a */ 439c4762a1bSJed Brown uL = &ctx->uLR[0]; 440c4762a1bSJed Brown uR = &ctx->uLR[dof]; 441c4762a1bSJed Brown for (j=0; j<dof; j++) { 442c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 443c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 444c4762a1bSJed Brown } 4459566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax)); 446c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 447c4762a1bSJed Brown if (i < xs+xm) { 448c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow1)*dof+j] += ctx->flux[j]/hx; 449c4762a1bSJed Brown } 450c4762a1bSJed Brown } 451c4762a1bSJed Brown } 452c4762a1bSJed Brown 4539566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 4549566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 4559566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da,PETSC_TRUE,&slope)); 4569566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 457c4762a1bSJed Brown PetscFunctionReturn(0); 458c4762a1bSJed Brown } 459c4762a1bSJed Brown 460c4762a1bSJed Brown PetscErrorCode FVRHSFunctionfast2(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 461c4762a1bSJed Brown { 462c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 463c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,len_slow1,len_slow2; 464c4762a1bSJed Brown PetscReal hx,cfl_idt = 0; 465c4762a1bSJed Brown PetscScalar *x,*f,*slope; 466c4762a1bSJed Brown Vec Xloc; 467c4762a1bSJed Brown DM da; 468c4762a1bSJed Brown 469c4762a1bSJed Brown PetscFunctionBeginUser; 4709566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 4719566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); 4729566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 473c4762a1bSJed Brown hx = (ctx->xmax-ctx->xmin)/Mx; 4749566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 4759566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 4769566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 4779566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 4789566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 4799566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da,PETSC_TRUE,&slope)); 4809566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 4819566063dSJacob Faibussowitsch PetscCall(ISGetSize(ctx->iss,&len_slow1)); 4829566063dSJacob Faibussowitsch PetscCall(ISGetSize(ctx->iss2,&len_slow2)); 483c4762a1bSJed Brown 484c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 485c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 486c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 487c4762a1bSJed Brown } 488c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 489c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 490c4762a1bSJed Brown } 491c4762a1bSJed Brown } 492c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 493c4762a1bSJed Brown struct _LimitInfo info; 494c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 495c4762a1bSJed Brown if (i > len_slow1+len_slow2-2) { 4969566063dSJacob Faibussowitsch PetscCall((*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds,ctx->xmin+hx*i)); 4979566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR,2*dof)); 498c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 499c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 500c4762a1bSJed Brown for (j=0; j<dof; j++) { 501c4762a1bSJed Brown PetscScalar jmpL,jmpR; 502c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 503c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 504c4762a1bSJed Brown for (k=0; k<dof; k++) { 505c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 506c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 507c4762a1bSJed Brown } 508c4762a1bSJed Brown } 509c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 510c4762a1bSJed Brown info.m = dof; 511c4762a1bSJed Brown info.hx = hx; 512c4762a1bSJed Brown (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope); 513c4762a1bSJed Brown for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 514c4762a1bSJed Brown for (j=0; j<dof; j++) { 515c4762a1bSJed Brown PetscScalar tmp = 0; 516c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 517c4762a1bSJed Brown slope[i*dof+j] = tmp; 518c4762a1bSJed Brown } 519c4762a1bSJed Brown } 520c4762a1bSJed Brown } 521c4762a1bSJed Brown 522c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 523c4762a1bSJed Brown PetscReal maxspeed; 524c4762a1bSJed Brown PetscScalar *uL,*uR; 525c4762a1bSJed Brown if (i > len_slow1+len_slow2) { 526c4762a1bSJed Brown uL = &ctx->uLR[0]; 527c4762a1bSJed Brown uR = &ctx->uLR[dof]; 528c4762a1bSJed Brown for (j=0; j<dof; j++) { 529c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 530c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 531c4762a1bSJed Brown } 5329566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax)); 533c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 534c4762a1bSJed Brown if (i > xs) { 535c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow1-len_slow2-1)*dof+j] -= ctx->flux[j]/hx; 536c4762a1bSJed Brown } 537c4762a1bSJed Brown if (i < xs+xm) { 538c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow1-len_slow2)*dof+j] += ctx->flux[j]/hx; 539c4762a1bSJed Brown } 540c4762a1bSJed Brown } 541c4762a1bSJed Brown if (i == len_slow1+len_slow2) { 542c4762a1bSJed Brown uL = &ctx->uLR[0]; 543c4762a1bSJed Brown uR = &ctx->uLR[dof]; 544c4762a1bSJed Brown for (j=0; j<dof; j++) { 545c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 546c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 547c4762a1bSJed Brown } 5489566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax)); 549c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 550c4762a1bSJed Brown if (i < xs+xm) { 551c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow1-len_slow2)*dof+j] += ctx->flux[j]/hx; 552c4762a1bSJed Brown } 553c4762a1bSJed Brown } 554c4762a1bSJed Brown } 5559566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 5569566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 5579566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da,PETSC_TRUE,&slope)); 5589566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 559c4762a1bSJed Brown PetscFunctionReturn(0); 560c4762a1bSJed Brown } 561c4762a1bSJed Brown 562c4762a1bSJed Brown int main(int argc,char *argv[]) 563c4762a1bSJed Brown { 564c4762a1bSJed Brown char lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m"; 565c4762a1bSJed Brown PetscFunctionList limiters = 0,physics = 0; 566c4762a1bSJed Brown MPI_Comm comm; 567c4762a1bSJed Brown TS ts; 568c4762a1bSJed Brown DM da; 569c4762a1bSJed Brown Vec X,X0,R; 570c4762a1bSJed Brown FVCtx ctx; 571c4762a1bSJed Brown PetscInt i,k,dof,xs,xm,Mx,draw = 0,*index_slow,*index_fast,islow = 0,ifast = 0; 572c4762a1bSJed Brown PetscBool view_final = PETSC_FALSE; 573c4762a1bSJed Brown PetscReal ptime; 574c4762a1bSJed Brown 575*327415f7SBarry Smith PetscFunctionBeginUser; 5769566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,0,help)); 577c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 5789566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&ctx,sizeof(ctx))); 579c4762a1bSJed Brown 580c4762a1bSJed Brown /* Register limiters to be available on the command line */ 5819566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"upwind" ,Limit_Upwind)); 5829566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"lax-wendroff" ,Limit_LaxWendroff)); 5839566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"beam-warming" ,Limit_BeamWarming)); 5849566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"fromm" ,Limit_Fromm)); 5859566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"minmod" ,Limit_Minmod)); 5869566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"superbee" ,Limit_Superbee)); 5879566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"mc" ,Limit_MC)); 5889566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"vanleer" ,Limit_VanLeer)); 5899566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"vanalbada" ,Limit_VanAlbada)); 5909566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"vanalbadatvd" ,Limit_VanAlbadaTVD)); 5919566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"koren" ,Limit_Koren)); 5929566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"korensym" ,Limit_KorenSym)); 5939566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"koren3" ,Limit_Koren3)); 5949566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"cada-torrilhon2" ,Limit_CadaTorrilhon2)); 5959566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"cada-torrilhon3-r0p1",Limit_CadaTorrilhon3R0p1)); 5969566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"cada-torrilhon3-r1" ,Limit_CadaTorrilhon3R1)); 5979566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"cada-torrilhon3-r10" ,Limit_CadaTorrilhon3R10)); 5989566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"cada-torrilhon3-r100",Limit_CadaTorrilhon3R100)); 599c4762a1bSJed Brown 600c4762a1bSJed Brown /* Register physical models to be available on the command line */ 6019566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&physics,"advect" ,PhysicsCreate_Advect)); 602c4762a1bSJed Brown 603c4762a1bSJed Brown ctx.comm = comm; 604c4762a1bSJed Brown ctx.cfl = 0.9; 605c4762a1bSJed Brown ctx.bctype = FVBC_PERIODIC; 606c4762a1bSJed Brown ctx.xmin = -1.0; 607c4762a1bSJed Brown ctx.xmax = 1.0; 608d0609cedSBarry Smith PetscOptionsBegin(comm,NULL,"Finite Volume solver options",""); 6099566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL)); 6109566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL)); 6119566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-limit","Name of flux limiter to use","",limiters,lname,lname,sizeof(lname),NULL)); 6129566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL)); 6139566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-view_final","Write final solution in ASCII MATLAB format to given file name","",final_fname,final_fname,sizeof(final_fname),&view_final)); 6149566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL)); 6159566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL)); 6169566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL)); 6179566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL)); 6189566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-recursive_split","Split the domain recursively","",ctx.recursive,&ctx.recursive,NULL)); 619d0609cedSBarry Smith PetscOptionsEnd(); 620c4762a1bSJed Brown 621c4762a1bSJed Brown /* Choose the limiter from the list of registered limiters */ 6229566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(limiters,lname,&ctx.limit)); 6233c633725SBarry Smith PetscCheck(ctx.limit,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Limiter '%s' not found",lname); 624c4762a1bSJed Brown 625c4762a1bSJed Brown /* Choose the physics from the list of registered models */ 626c4762a1bSJed Brown { 627c4762a1bSJed Brown PetscErrorCode (*r)(FVCtx*); 6289566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(physics,physname,&r)); 6293c633725SBarry Smith PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Physics '%s' not found",physname); 630c4762a1bSJed Brown /* Create the physics, will set the number of fields and their names */ 6319566063dSJacob Faibussowitsch PetscCall((*r)(&ctx)); 632c4762a1bSJed Brown } 633c4762a1bSJed Brown 634c4762a1bSJed Brown /* Create a DMDA to manage the parallel grid */ 6359566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics.dof,2,NULL,&da)); 6369566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 6379566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 638c4762a1bSJed Brown /* Inform the DMDA of the field names provided by the physics. */ 639c4762a1bSJed Brown /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ 640c4762a1bSJed Brown for (i=0; i<ctx.physics.dof; i++) { 6419566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da,i,ctx.physics.fieldname[i])); 642c4762a1bSJed Brown } 6439566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 6449566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 645c4762a1bSJed Brown 646c4762a1bSJed Brown /* Set coordinates of cell centers */ 6479566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da,ctx.xmin+0.5*(ctx.xmax-ctx.xmin)/Mx,ctx.xmax+0.5*(ctx.xmax-ctx.xmin)/Mx,0,0,0,0)); 648c4762a1bSJed Brown 649c4762a1bSJed Brown /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */ 6509566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dof*dof,&ctx.R,dof*dof,&ctx.Rinv,2*dof,&ctx.cjmpLR,1*dof,&ctx.cslope)); 6519566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(2*dof,&ctx.uLR,dof,&ctx.flux,dof,&ctx.speeds)); 652c4762a1bSJed Brown 653c4762a1bSJed Brown /* Create a vector to store the solution and to save the initial state */ 6549566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da,&X)); 6559566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&X0)); 6569566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&R)); 657c4762a1bSJed Brown 658c4762a1bSJed Brown /*-------------------------------- create index for slow parts and fast parts ----------------------------------------*/ 6599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm*dof,&index_slow)); 6609566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm*dof,&index_fast)); 661c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 662c4762a1bSJed Brown if (ctx.xmin+i*(ctx.xmax-ctx.xmin)/(PetscReal)Mx+0.5*(ctx.xmax-ctx.xmin)/(PetscReal)Mx < 0) 663c4762a1bSJed Brown for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k; 664c4762a1bSJed Brown else 665c4762a1bSJed Brown for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k; 666c4762a1bSJed Brown } /* this step need to be changed based on discontinuous point of a */ 6679566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss)); 6689566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf)); 669c4762a1bSJed Brown 670c4762a1bSJed Brown /* Create a time-stepping object */ 6719566063dSJacob Faibussowitsch PetscCall(TSCreate(comm,&ts)); 6729566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts,da)); 6739566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,R,FVRHSFunction,&ctx)); 6749566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"slow",ctx.iss)); 6759566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"fast",ctx.isf)); 6769566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"slow",NULL,FVRHSFunctionslow,&ctx)); 6779566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"fast",NULL,FVRHSFunctionfast,&ctx)); 678c4762a1bSJed Brown 679c4762a1bSJed Brown if (ctx.recursive) { 680c4762a1bSJed Brown TS subts; 681c4762a1bSJed Brown islow = 0; 682c4762a1bSJed Brown ifast = 0; 683c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 684c4762a1bSJed Brown PetscReal coord = ctx.xmin+i*(ctx.xmax-ctx.xmin)/(PetscReal)Mx+0.5*(ctx.xmax-ctx.xmin)/(PetscReal)Mx; 685c4762a1bSJed Brown if (coord >= 0 && coord < ctx.xmin+(ctx.xmax-ctx.xmin)*3/4.) 686c4762a1bSJed Brown for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k; 687c4762a1bSJed Brown if (coord >= ctx.xmin+(ctx.xmax-ctx.xmin)*3/4.) 688c4762a1bSJed Brown for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k; 689c4762a1bSJed Brown } /* this step need to be changed based on discontinuous point of a */ 6909566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss2)); 6919566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf2)); 692c4762a1bSJed Brown 6939566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts,"fast",&subts)); 6949566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(subts,"slow",ctx.iss2)); 6959566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(subts,"fast",ctx.isf2)); 6969566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(subts,"slow",NULL,FVRHSFunctionslow2,&ctx)); 6979566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(subts,"fast",NULL,FVRHSFunctionfast2,&ctx)); 698c4762a1bSJed Brown } 699c4762a1bSJed Brown 7009566063dSJacob Faibussowitsch /*PetscCall(TSSetType(ts,TSSSP));*/ 7019566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSMPRK)); 7029566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,10)); 7039566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 704c4762a1bSJed Brown 705c4762a1bSJed Brown /* Compute initial conditions and starting time step */ 7069566063dSJacob Faibussowitsch PetscCall(FVSample(&ctx,da,0,X0)); 7079566063dSJacob Faibussowitsch PetscCall(FVRHSFunction(ts,0,X0,X,(void*)&ctx)); /* Initial function evaluation, only used to determine max speed */ 7089566063dSJacob Faibussowitsch PetscCall(VecCopy(X0,X)); /* The function value was not used so we set X=X0 again */ 7099566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt)); 7109566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); /* Take runtime options */ 7119566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD)); 712c4762a1bSJed Brown { 713c4762a1bSJed Brown PetscInt steps; 714c4762a1bSJed Brown PetscScalar mass_initial,mass_final,mass_difference; 715c4762a1bSJed Brown 7169566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,X)); 7179566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&ptime)); 7189566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&steps)); 71963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Final time %g, steps %" PetscInt_FMT "\n",(double)ptime,steps)); 720c4762a1bSJed Brown /* calculate the total mass at initial time and final time */ 721c4762a1bSJed Brown mass_initial = 0.0; 722c4762a1bSJed Brown mass_final = 0.0; 7239566063dSJacob Faibussowitsch PetscCall(VecSum(X0,&mass_initial)); 7249566063dSJacob Faibussowitsch PetscCall(VecSum(X,&mass_final)); 725c4762a1bSJed Brown mass_difference = (ctx.xmax-ctx.xmin)/(PetscScalar)Mx*(mass_final - mass_initial); 7269566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Mass difference %g\n",(double)mass_difference)); 727c4762a1bSJed Brown if (ctx.simulation) { 728c4762a1bSJed Brown PetscViewer fd; 729c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = "binaryoutput"; 730c4762a1bSJed Brown Vec XR; 731c4762a1bSJed Brown PetscReal nrm1,nrmsup; 732c4762a1bSJed Brown PetscBool flg; 733d8185827SBarry Smith 7349566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg)); 7353c633725SBarry Smith PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 7369566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&fd)); 7379566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X0,&XR)); 7389566063dSJacob Faibussowitsch PetscCall(VecLoad(XR,fd)); 7399566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 7409566063dSJacob Faibussowitsch PetscCall(VecAYPX(XR,-1,X)); 7419566063dSJacob Faibussowitsch PetscCall(VecNorm(XR,NORM_1,&nrm1)); 7429566063dSJacob Faibussowitsch PetscCall(VecNorm(XR,NORM_INFINITY,&nrmsup)); 743c4762a1bSJed Brown nrm1 /= Mx; 7449566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Error ||x-x_e||_1 %g ||x-x_e||_sup %g\n",(double)nrm1,(double)nrmsup)); 7459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XR)); 746c4762a1bSJed Brown } 747c4762a1bSJed Brown } 748c4762a1bSJed Brown 7499566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD)); 7509566063dSJacob Faibussowitsch if (draw & 0x1) PetscCall(VecView(X0,PETSC_VIEWER_DRAW_WORLD)); 7519566063dSJacob Faibussowitsch if (draw & 0x2) PetscCall(VecView(X,PETSC_VIEWER_DRAW_WORLD)); 752c4762a1bSJed Brown if (draw & 0x4) { 753c4762a1bSJed Brown Vec Y; 7549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&Y)); 7559566063dSJacob Faibussowitsch PetscCall(FVSample(&ctx,da,ptime,Y)); 7569566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y,-1,X)); 7579566063dSJacob Faibussowitsch PetscCall(VecView(Y,PETSC_VIEWER_DRAW_WORLD)); 7589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 759c4762a1bSJed Brown } 760c4762a1bSJed Brown 761c4762a1bSJed Brown if (view_final) { 762c4762a1bSJed Brown PetscViewer viewer; 7639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer)); 7649566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB)); 7659566063dSJacob Faibussowitsch PetscCall(VecView(X,viewer)); 7669566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 7679566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 768c4762a1bSJed Brown } 769c4762a1bSJed Brown 770c4762a1bSJed Brown /* Clean up */ 7719566063dSJacob Faibussowitsch PetscCall((*ctx.physics.destroy)(ctx.physics.user)); 7729566063dSJacob Faibussowitsch for (i=0; i<ctx.physics.dof; i++) PetscCall(PetscFree(ctx.physics.fieldname[i])); 7739566063dSJacob Faibussowitsch PetscCall(PetscFree4(ctx.R,ctx.Rinv,ctx.cjmpLR,ctx.cslope)); 7749566063dSJacob Faibussowitsch PetscCall(PetscFree3(ctx.uLR,ctx.flux,ctx.speeds)); 7759566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.iss)); 7769566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.isf)); 7779566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.iss2)); 7789566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.isf2)); 7799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 7809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X0)); 7819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&R)); 7829566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 7839566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 7849566063dSJacob Faibussowitsch PetscCall(PetscFree(index_slow)); 7859566063dSJacob Faibussowitsch PetscCall(PetscFree(index_fast)); 7869566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&limiters)); 7879566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&physics)); 7889566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 789b122ec5aSJacob Faibussowitsch return 0; 790c4762a1bSJed Brown } 791c4762a1bSJed Brown 792c4762a1bSJed Brown /*TEST 793c4762a1bSJed Brown build: 794f56ea12dSJed Brown requires: !complex 795c4762a1bSJed Brown depends: finitevolume1d.c 796c4762a1bSJed Brown 797c4762a1bSJed Brown test: 798c4762a1bSJed 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 799c4762a1bSJed Brown 800c4762a1bSJed Brown test: 801c4762a1bSJed Brown suffix: 2 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 1 803c4762a1bSJed Brown output_file: output/ex5_1.out 804c4762a1bSJed Brown 805c4762a1bSJed Brown test: 806c4762a1bSJed Brown suffix: 3 807c4762a1bSJed 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 808c4762a1bSJed Brown 809c4762a1bSJed Brown test: 810c4762a1bSJed Brown suffix: 4 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 1 812c4762a1bSJed Brown output_file: output/ex5_3.out 813c4762a1bSJed Brown 814c4762a1bSJed Brown test: 815c4762a1bSJed Brown suffix: 5 816c4762a1bSJed 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 817c4762a1bSJed Brown 818c4762a1bSJed Brown test: 819c4762a1bSJed Brown suffix: 6 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 1 -recursive_split 821c4762a1bSJed Brown TEST*/ 822