/* Note: -hratio is the ratio between mesh size of corse grids and fine grids */ static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n" " advect - Constant coefficient scalar advection\n" " u_t + (a*u)_x = 0\n" " shallow - 1D Shallow water equations (Saint Venant System)\n" " h_t + (q)_x = 0 \n" " q_t + (\frac{q^2}{h} + g/2*h^2)_x = 0 \n" " where, h(x,t) denotes the height of the water and q(x,t) the momentum.\n" " for this toy problem, we choose different meshsizes for different sub-domains (slow-fast-slow), say\n" " hxs = hratio*hxf \n" " where hxs and hxf are the grid spacings for coarse and fine grids respectively.\n" " exact - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n" " the states across shocks and rarefactions\n" " simulation - use reference solution which is generated by smaller time step size to be true solution,\n" " also the reference solution should be generated by user and stored in a binary file.\n" " characteristic - Limit the characteristic variables, this is usually preferred (default)\n" " bc_type - Boundary condition for the problem, options are: periodic, outflow, inflow " "Several problem descriptions (initial data, physics specific features, boundary data) can be chosen with -initial N\n\n" "The problem size should be set with -da_grid_x M\n\n"; /* Example: ./ex4 -da_grid_x 40 -initial 1 -hratio 1 -limit mc -ts_dt 0.01 -ts_max_time 7.0 -ts_type mprk -ts_mprk_type 2a22 -ts_monitor_draw_solution -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0 ./ex4 -da_grid_x 40 -initial 2 -hratio 1 -limit mc -ts_dt 0.01 -ts_max_time 2.5 -ts_type mprk -ts_mprk_type 2a22 -ts_monitor_draw_solution -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0 ./ex4 -da_grid_x 40 -initial 3 -hratio 1 -limit mc -ts_dt 0.01 -ts_max_time 4.0 -ts_type mprk -ts_mprk_type 2a22 -ts_monitor_draw_solution -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0 ./ex4 -da_grid_x 40 -initial 4 -hratio 1 -limit koren3 -ts_dt 0.01 -ts_max_time 4.0 -ts_type mprk -ts_mprk_type 2a22 -ts_monitor_draw_solution -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0 ./ex4 -da_grid_x 40 -initial 5 -hratio 1 -limit mc -ts_dt 0.01 -ts_max_time 5.0 -ts_type mprk -ts_mprk_type 2a22 -ts_monitor_draw_solution -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0 Contributed by: Aidan Hamilton */ #include #include #include #include #include "finitevolume1d.h" #include static inline PetscReal RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax) { PetscReal range = xmax-xmin; return xmin +PetscFmodReal(range+PetscFmodReal(a,range),range); } static inline PetscReal MaxAbs(PetscReal a,PetscReal b) { return (PetscAbs(a) > PetscAbs(b)) ? a : b; } /* --------------------------------- Advection ----------------------------------- */ typedef struct { PetscReal a; /* advective velocity */ } AdvectCtx; static PetscErrorCode PhysicsRiemann_Advect(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed) { AdvectCtx *ctx = (AdvectCtx*)vctx; PetscReal speed; PetscFunctionBeginUser; speed = ctx->a; flux[0] = PetscMax(0,speed)*uL[0] + PetscMin(0,speed)*uR[0]; *maxspeed = speed; PetscFunctionReturn(0); } static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds) { AdvectCtx *ctx = (AdvectCtx*)vctx; PetscFunctionBeginUser; X[0] = 1.; Xi[0] = 1.; speeds[0] = ctx->a; PetscFunctionReturn(0); } static PetscErrorCode PhysicsSample_Advect(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u) { AdvectCtx *ctx = (AdvectCtx*)vctx; PetscReal a = ctx->a,x0; PetscFunctionBeginUser; switch (bctype) { case FVBC_OUTFLOW: x0 = x-a*t; break; case FVBC_PERIODIC: x0 = RangeMod(x-a*t,xmin,xmax); break; default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown BCType"); } switch (initial) { case 0: u[0] = (x0 < 0) ? 1 : -1; break; case 1: u[0] = (x0 < 0) ? -1 : 1; break; case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break; case 3: u[0] = PetscSinReal(2*PETSC_PI*x0); break; case 4: u[0] = PetscAbs(x0); break; case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2*PETSC_PI*x0)); break; case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2-x0 : 0)); break; case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI*x0),10.0);break; default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition"); } PetscFunctionReturn(0); } static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx) { PetscErrorCode ierr; AdvectCtx *user; PetscFunctionBeginUser; ierr = PetscNew(&user);CHKERRQ(ierr); ctx->physics2.sample2 = PhysicsSample_Advect; ctx->physics2.riemann2 = PhysicsRiemann_Advect; ctx->physics2.characteristic2 = PhysicsCharacteristic_Advect; ctx->physics2.destroy = PhysicsDestroy_SimpleFree; ctx->physics2.user = user; ctx->physics2.dof = 1; ierr = PetscStrallocpy("u",&ctx->physics2.fieldname[0]);CHKERRQ(ierr); user->a = 1; ierr = PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");CHKERRQ(ierr); ierr = PetscOptionsReal("-physics_advect_a","Speed","",user->a,&user->a,NULL);CHKERRQ(ierr); ierr = PetscOptionsEnd();CHKERRQ(ierr); PetscFunctionReturn(0); } /* --------------------------------- Shallow Water ----------------------------------- */ typedef struct { PetscReal gravity; } ShallowCtx; static inline void ShallowFlux(ShallowCtx *phys,const PetscScalar *u,PetscScalar *f) { f[0] = u[1]; f[1] = PetscSqr(u[1])/u[0] + 0.5*phys->gravity*PetscSqr(u[0]); } static PetscErrorCode PhysicsRiemann_Shallow_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed) { ShallowCtx *phys = (ShallowCtx*)vctx; PetscScalar g = phys->gravity,ustar[2],cL,cR,c,cstar; struct {PetscScalar h,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]},star; PetscInt i; PetscFunctionBeginUser; PetscAssertFalse(!(L.h > 0 && R.h > 0),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative"); cL = PetscSqrtScalar(g*L.h); cR = PetscSqrtScalar(g*R.h); c = PetscMax(cL,cR); { /* Solve for star state */ const PetscInt maxits = 50; PetscScalar tmp,res,res0=0,h0,h = 0.5*(L.h + R.h); /* initial guess */ h0 = h; for (i=0; i 0 && PetscAbsScalar(h-h0) < PETSC_SQRT_MACHINE_EPSILON)) { star.h = h; star.u = L.u - fl; goto converged; } else if (i > 0 && PetscAbsScalar(res) >= PetscAbsScalar(res0)) { /* Line search */ h = 0.8*h0 + 0.2*h; continue; } /* Accept the last step and take another */ res0 = res; h0 = h; dfl = (L.h < h) ? 0.5/fl*0.5*g*(-L.h*L.h/(h*h) - 1 + 2*h/L.h) : PetscSqrtScalar(g/h); dfr = (R.h < h) ? 0.5/fr*0.5*g*(-R.h*R.h/(h*h) - 1 + 2*h/R.h) : PetscSqrtScalar(g/h); tmp = h - res/(dfr+dfl); if (tmp <= 0) h /= 2; /* Guard against Newton shooting off to a negative thickness */ else h = tmp; PetscAssertFalse(!((h > 0) && PetscIsNormalScalar(h)),PETSC_COMM_SELF,PETSC_ERR_FP,"non-normal iterate h=%g",(double)h); } SETERRQ(PETSC_COMM_SELF,1,"Newton iteration for star.h diverged after %D iterations",i); } converged: cstar = PetscSqrtScalar(g*star.h); if (L.u-cL < 0 && 0 < star.u-cstar) { /* 1-wave is sonic rarefaction */ PetscScalar ufan[2]; ufan[0] = 1/g*PetscSqr(L.u/3 + 2./3*cL); ufan[1] = PetscSqrtScalar(g*ufan[0])*ufan[0]; ShallowFlux(phys,ufan,flux); } else if (star.u+cstar < 0 && 0 < R.u+cR) { /* 2-wave is sonic rarefaction */ PetscScalar ufan[2]; ufan[0] = 1/g*PetscSqr(R.u/3 - 2./3*cR); ufan[1] = -PetscSqrtScalar(g*ufan[0])*ufan[0]; ShallowFlux(phys,ufan,flux); } else if ((L.h >= star.h && L.u-c >= 0) || (L.h 0)) { /* 1-wave is right-travelling shock (supersonic) */ ShallowFlux(phys,uL,flux); } else if ((star.h <= R.h && R.u+c <= 0) || (star.h>R.h && (R.h*R.u-star.h*star.h)/(R.h-star.h) < 0)) { /* 2-wave is left-travelling shock (supersonic) */ ShallowFlux(phys,uR,flux); } else { ustar[0] = star.h; ustar[1] = star.h*star.u; ShallowFlux(phys,ustar,flux); } *maxspeed = MaxAbs(MaxAbs(star.u-cstar,star.u+cstar),MaxAbs(L.u-cL,R.u+cR)); PetscFunctionReturn(0); } static PetscErrorCode PhysicsRiemann_Shallow_Rusanov(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed) { ShallowCtx *phys = (ShallowCtx*)vctx; PetscScalar g = phys->gravity,fL[2],fR[2],s; struct {PetscScalar h,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]}; PetscReal tol = 1e-6; PetscFunctionBeginUser; /* Positivity preserving modification*/ if (L.h < tol) L.u = 0.0; if (R.h < tol) R.u = 0.0; /*simple pos preserve limiter*/ if (L.h < 0) L.h = 0; if (R.h < 0) R.h = 0; ShallowFlux(phys,uL,fL); ShallowFlux(phys,uR,fR); s = PetscMax(PetscAbs(L.u)+PetscSqrtScalar(g*L.h),PetscAbs(R.u)+PetscSqrtScalar(g*R.h)); flux[0] = 0.5*(fL[0] + fR[0]) + 0.5*s*(L.h - R.h); flux[1] = 0.5*(fL[1] + fR[1]) + 0.5*s*(uL[1] - uR[1]); *maxspeed = s; PetscFunctionReturn(0); } static PetscErrorCode PhysicsCharacteristic_Conservative(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds) { PetscInt i,j; PetscFunctionBeginUser; for (i=0; igravity); if (u[0] < tol) { /*Use conservative variables*/ X[0*2+0] = 1; X[0*2+1] = 0; X[1*2+0] = 0; X[1*2+1] = 1; } else { speeds[0] = u[1]/u[0] - c; speeds[1] = u[1]/u[0] + c; X[0*2+0] = 1; X[0*2+1] = speeds[0]; X[1*2+0] = 1; X[1*2+1] = speeds[1]; } ierr = PetscArraycpy(Xi,X,4);CHKERRQ(ierr); ierr = PetscKernel_A_gets_inverse_A_2(Xi,0,PETSC_FALSE,NULL);CHKERRQ(ierr); PetscFunctionReturn(0); } static PetscErrorCode PhysicsSample_Shallow(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u) { PetscFunctionBeginUser; PetscAssertFalse(t > 0,PETSC_COMM_SELF,PETSC_ERR_SUP,"Exact solutions not implemented for t > 0"); switch (initial) { case 0: u[0] = (x < 0) ? 2 : 1; /* Standard Dam Break Problem */ u[1] = (x < 0) ? 0 : 0; break; case 1: u[0] = (x < 10) ? 1 : 0.1; /*The Next 5 problems are standard Riemann problem tests */ u[1] = (x < 10) ? 2.5 : 0; break; case 2: u[0] = (x < 25) ? 1 : 1; u[1] = (x < 25) ? -5 : 5; break; case 3: u[0] = (x < 20) ? 1 : 1e-12; u[1] = (x < 20) ? 0 : 0; break; case 4: u[0] = (x < 30) ? 1e-12 : 1; u[1] = (x < 30) ? 0 : 0; break; case 5: u[0] = (x < 25) ? 0.1 : 0.1; u[1] = (x < 25) ? -0.3 : 0.3; break; case 6: u[0] = 1+0.5*PetscSinReal(2*PETSC_PI*x); u[1] = 1*u[0]; break; case 7: if (x < -0.1) { u[0] = 1e-9; u[1] = 0.0; } else if (x < 0.1) { u[0] = 1.0; u[1] = 0.0; } else { u[0] = 1e-9; u[1] = 0.0; } break; case 8: if (x < -0.1) { u[0] = 2; u[1] = 0.0; } else if (x < 0.1) { u[0] = 3.0; u[1] = 0.0; } else { u[0] = 2; u[1] = 0.0; } break; default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition"); } PetscFunctionReturn(0); } /* Implements inflow conditions for the given initial conditions. Which conditions are actually enforced depends on on the results of PhysicsSetInflowType_Shallow. */ static PetscErrorCode PhysicsInflow_Shallow(void *vctx,PetscReal t,PetscReal x,PetscReal *u) { FVCtx *ctx = (FVCtx*)vctx; PetscFunctionBeginUser; if (ctx->bctype == FVBC_INFLOW) { switch (ctx->initial) { case 0: case 1: case 2: case 3: case 4: case 5: u[0] = 0; u[1] = 0.0; /* Left boundary conditions */ u[2] = 0; u[3] = 0.0; /* Right boundary conditions */ break; case 6: u[0] = 0; u[1] = 0.0; /* Left boundary conditions */ u[2] = 0; u[3] = 0.0; /* Right boundary conditions */ break; case 7: u[0] = 0; u[1] = 0.0; /* Left boundary conditions */ u[2] = 0; u[3] = 0.0; /* Right boundary conditions */ break; case 8: u[0] = 0; u[1] = 1.0; /* Left boundary conditions */ u[2] = 0; u[3] = -1.0;/* Right boundary conditions */ break; default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition"); } } PetscFunctionReturn(0); } /* Selects which boundary conditions are marked as inflow and which as outflow when FVBC_INFLOW is selected. */ static PetscErrorCode PhysicsSetInflowType_Shallow(FVCtx *ctx) { PetscFunctionBeginUser; switch (ctx->initial) { case 0: case 1: case 2: case 3: case 4: case 5: case 6: case 7: case 8: /* Fix left and right momentum, height is outflow */ ctx->physics2.bcinflowindex[0] = PETSC_FALSE; ctx->physics2.bcinflowindex[1] = PETSC_TRUE; ctx->physics2.bcinflowindex[2] = PETSC_FALSE; ctx->physics2.bcinflowindex[3] = PETSC_TRUE; break; default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition"); } PetscFunctionReturn(0); } static PetscErrorCode PhysicsCreate_Shallow(FVCtx *ctx) { PetscErrorCode ierr; ShallowCtx *user; PetscFunctionList rlist = 0,rclist = 0; char rname[256] = "rusanov",rcname[256] = "characteristic"; PetscFunctionBeginUser; ierr = PetscNew(&user);CHKERRQ(ierr); ctx->physics2.sample2 = PhysicsSample_Shallow; ctx->physics2.inflow = PhysicsInflow_Shallow; ctx->physics2.destroy = PhysicsDestroy_SimpleFree; ctx->physics2.riemann2 = PhysicsRiemann_Shallow_Rusanov; ctx->physics2.characteristic2 = PhysicsCharacteristic_Shallow; ctx->physics2.user = user; ctx->physics2.dof = 2; PetscMalloc1(2*(ctx->physics2.dof),&ctx->physics2.bcinflowindex); PhysicsSetInflowType_Shallow(ctx); ierr = PetscStrallocpy("height",&ctx->physics2.fieldname[0]);CHKERRQ(ierr); ierr = PetscStrallocpy("momentum",&ctx->physics2.fieldname[1]);CHKERRQ(ierr); user->gravity = 9.81; ierr = RiemannListAdd_2WaySplit(&rlist,"exact", PhysicsRiemann_Shallow_Exact);CHKERRQ(ierr); ierr = RiemannListAdd_2WaySplit(&rlist,"rusanov",PhysicsRiemann_Shallow_Rusanov);CHKERRQ(ierr); ierr = ReconstructListAdd_2WaySplit(&rclist,"characteristic",PhysicsCharacteristic_Shallow);CHKERRQ(ierr); ierr = ReconstructListAdd_2WaySplit(&rclist,"conservative",PhysicsCharacteristic_Conservative);CHKERRQ(ierr); ierr = PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for Shallow","");CHKERRQ(ierr); ierr = PetscOptionsReal("-physics_shallow_gravity","Gravity","",user->gravity,&user->gravity,NULL);CHKERRQ(ierr); ierr = PetscOptionsFList("-physics_shallow_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL);CHKERRQ(ierr); ierr = PetscOptionsFList("-physics_shallow_reconstruct","Reconstruction","",rclist,rcname,rcname,sizeof(rcname),NULL);CHKERRQ(ierr); ierr = PetscOptionsEnd();CHKERRQ(ierr); ierr = RiemannListFind_2WaySplit(rlist,rname,&ctx->physics2.riemann2);CHKERRQ(ierr); ierr = ReconstructListFind_2WaySplit(rclist,rcname,&ctx->physics2.characteristic2);CHKERRQ(ierr); ierr = PetscFunctionListDestroy(&rlist);CHKERRQ(ierr); ierr = PetscFunctionListDestroy(&rclist);CHKERRQ(ierr); PetscFunctionReturn(0); } PetscErrorCode FVSample_2WaySplit(FVCtx *ctx,DM da,PetscReal time,Vec U) { PetscErrorCode ierr; PetscScalar *u,*uj,xj,xi; PetscInt i,j,k,dof,xs,xm,Mx; const PetscInt N = 200; PetscReal hs,hf; PetscFunctionBeginUser; PetscAssertFalse(!ctx->physics2.sample2,PETSC_COMM_SELF,PETSC_ERR_SUP,"Physics has not provided a sampling function"); ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr); ierr = PetscMalloc1(dof,&uj);CHKERRQ(ierr); hs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; hf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); for (i=xs; isf) { xi = ctx->xmin+0.5*hs+i*hs; /* Integrate over cell i using trapezoid rule with N points. */ for (k=0; kphysics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); for (k=0; kfs) { xi = ctx->xmin+ctx->sf*hs+0.5*hf+(i-ctx->sf)*hf; /* Integrate over cell i using trapezoid rule with N points. */ for (k=0; kphysics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); for (k=0; kxmin+ctx->sf*hs+(ctx->fs-ctx->sf)*hf+0.5*hs+(i-ctx->fs)*hs; /* Integrate over cell i using trapezoid rule with N points. */ for (k=0; kphysics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); for (k=0; kxmax-ctx->xmin)*3.0/8.0/ctx->sf; hf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); ierr = VecGetArrayRead(X,&ptr_X);CHKERRQ(ierr); ierr = VecGetArrayRead(Y,&ptr_Y);CHKERRQ(ierr); for (i=0; isf || i > ctx->fs -1) *nrm1 += hs*PetscAbs(ptr_X[i]-ptr_Y[i]); else *nrm1 += hf*PetscAbs(ptr_X[i]-ptr_Y[i]); } ierr = VecRestoreArrayRead(X,&ptr_X);CHKERRQ(ierr); ierr = VecRestoreArrayRead(Y,&ptr_Y);CHKERRQ(ierr); ierr = VecDestroy(&Y);CHKERRQ(ierr); PetscFunctionReturn(0); } PetscErrorCode FVRHSFunction_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) { FVCtx *ctx = (FVCtx*)vctx; PetscErrorCode ierr; PetscInt i,j,k,Mx,dof,xs,xm,sf = ctx->sf,fs = ctx->fs; PetscReal hxf,hxs,cfl_idt = 0; PetscScalar *x,*f,*slope; Vec Xloc; DM da; PetscFunctionBeginUser; ierr = TSGetDM(ts,&da);CHKERRQ(ierr); ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); /* Xloc contains ghost points */ ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); /* Mx is the number of center points */ hxs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); /* X is solution vector which does not contain ghost points */ ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); ierr = VecZeroEntries(F);CHKERRQ(ierr); /* F is the right hand side function corresponds to center points */ ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr); ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); /* contains ghost points */ ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); if (ctx->bctype == FVBC_OUTFLOW) { for (i=xs-2; i<0; i++) { for (j=0; jbctype == FVBC_INFLOW) { /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253 pages 137-138 for the scheme. */ if (xs==0) { /* Left Boundary */ ctx->physics2.inflow(ctx,time,ctx->xmin,ctx->ub); for (j=0; jphysics2.bcinflowindex[j]) { for (i=-2; i<0; i++) x[i*dof+j] = 2.0*ctx->ub[j]-x[-(i+1)*dof+j]; } else { for (i=-2; i<0; i++) x[i*dof+j] = x[j]; /* Outflow */ } } } if (xs+xm==Mx) { /* Right Boundary */ ctx->physics2.inflow(ctx,time,ctx->xmax,ctx->ub); for (j=0; jphysics2.bcinflowindex[dof+j]) { for (i=Mx; iub[dof+j]-x[(2*Mx-(i+1))*dof+j]; } else { for (i=Mx; iphysics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); cjmpL = &ctx->cjmpLR[0]; cjmpR = &ctx->cjmpLR[dof]; for (j=0; jRinv[k+j*dof]*jmpL; cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; } } /* Apply limiter to the left and right characteristic jumps */ info.m = dof; info.hxs = hxs; info.hxf = hxf; (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope); for (j=0; jR[j+k*dof]*ctx->cslope[k]; slope[i*dof+j] = tmp; } } for (i=xs; iuLR[0]; uR = &ctx->uLR[dof]; if (i < sf) { /* slow region */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxs; } if (i < xs+xm) { for (j=0; jflux[j]/hxs; } } else if (i == sf) { /* interface between the slow region and the fast region */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxs; } if (i < xs+xm) { for (j=0; jflux[j]/hxf; } } else if (i > sf && i < fs) { /* fast region */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxf; } if (i < xs+xm) { for (j=0; jflux[j]/hxf; } } else if (i == fs) { /* interface between the fast region and the slow region */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxf; } if (i < xs+xm) { for (j=0; jflux[j]/hxs; } } else { /* slow region */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */ if (i > xs) { for (j=0; jflux[j]/hxs; } if (i < xs+xm) { for (j=0; jflux[j]/hxs; } } } ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr); ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); ierr = MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr); if (0) { /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */ PetscReal dt,tnow; ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); ierr = TSGetTime(ts,&tnow);CHKERRQ(ierr); if (dt > 0.5/ctx->cfl_idt) { if (1) { ierr = PetscPrintf(ctx->comm,"Stability constraint exceeded at t=%g, dt %g > %g\n",(double)tnow,(double)dt,(double)(0.5/ctx->cfl_idt));CHKERRQ(ierr); } else SETERRQ(PETSC_COMM_SELF,1,"Stability constraint exceeded, %g > %g",(double)dt,(double)(ctx->cfl/ctx->cfl_idt)); } } PetscFunctionReturn(0); } /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */ PetscErrorCode FVRHSFunctionslow_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) { FVCtx *ctx = (FVCtx*)vctx; PetscErrorCode ierr; PetscInt i,j,k,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth; PetscReal hxs,hxf,cfl_idt = 0; PetscScalar *x,*f,*slope; Vec Xloc; DM da; PetscFunctionBeginUser; ierr = TSGetDM(ts,&da);CHKERRQ(ierr); ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); hxs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); ierr = VecZeroEntries(F);CHKERRQ(ierr); ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); ierr = VecGetArray(F,&f);CHKERRQ(ierr); ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); if (ctx->bctype == FVBC_OUTFLOW) { for (i=xs-2; i<0; i++) { for (j=0; jbctype == FVBC_INFLOW) { /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253 pages 137-138 for the scheme. */ if (xs==0) { /* Left Boundary */ ctx->physics2.inflow(ctx,time,ctx->xmin,ctx->ub); for (j=0; jphysics2.bcinflowindex[j]==PETSC_TRUE) { for (i=-2; i<0; i++) x[i*dof+j] = 2.0*ctx->ub[j]-x[-(i+1)*dof+j]; } else { for (i=-2; i<0; i++) x[i*dof+j] = x[j]; /* Outflow */ } } } if (xs+xm==Mx) { /* Right Boundary */ ctx->physics2.inflow(ctx,time,ctx->xmax,ctx->ub); for (j=0; jphysics2.bcinflowindex[dof+j]==PETSC_TRUE) { for (i=Mx; iub[dof+j]-x[(2*Mx-(i+1))*dof+j]; } else { for (i=Mx; i fs+rsbwidth-2) { /* slow components and the first and last fast components */ /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); cjmpL = &ctx->cjmpLR[0]; cjmpR = &ctx->cjmpLR[dof]; for (j=0; jRinv[k+j*dof]*jmpL; cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; } } /* Apply limiter to the left and right characteristic jumps */ info.m = dof; info.hxs = hxs; info.hxf = hxf; (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope); for (j=0; jR[j+k*dof]*ctx->cslope[k]; slope[i*dof+j] = tmp; } } } for (i=xs; iuLR[0]; uR = &ctx->uLR[dof]; if (i < sf-lsbwidth) { /* slow region */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */ if (i > xs) { for (j=0; jflux[j]/hxs; } if (i < xs+xm) { for (j=0; jflux[j]/hxs; islow++; } } if (i == sf-lsbwidth) { /* interface between the slow region and the fast region */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxs; } } if (i == fs+rsbwidth) { /* slow region */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i < xs+xm) { for (j=0; jflux[j]/hxs; islow++; } } if (i > fs+rsbwidth) { /* slow region */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxs; } if (i < xs+xm) { for (j=0; jflux[j]/hxs; islow++; } } } ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); ierr = MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr); PetscFunctionReturn(0); } PetscErrorCode FVRHSFunctionslowbuffer_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) { FVCtx *ctx = (FVCtx*)vctx; PetscErrorCode ierr; PetscInt i,j,k,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth; PetscReal hxs,hxf; PetscScalar *x,*f,*slope; Vec Xloc; DM da; PetscFunctionBeginUser; ierr = TSGetDM(ts,&da);CHKERRQ(ierr); ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); hxs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); ierr = VecZeroEntries(F);CHKERRQ(ierr); ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); ierr = VecGetArray(F,&f);CHKERRQ(ierr); ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); if (ctx->bctype == FVBC_OUTFLOW) { for (i=xs-2; i<0; i++) { for (j=0; jbctype == FVBC_INFLOW) { /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253 pages 137-138 for the scheme. */ if (xs==0) { /* Left Boundary */ ctx->physics2.inflow(ctx,time,ctx->xmin,ctx->ub); for (j=0; jphysics2.bcinflowindex[j]==PETSC_TRUE) { for (i=-2; i<0; i++) x[i*dof+j] = 2.0*ctx->ub[j]-x[-(i+1)*dof+j]; } else { for (i=-2; i<0; i++) x[i*dof+j] = x[j]; /* Outflow */ } } } if (xs+xm==Mx) { /* Right Boundary */ ctx->physics2.inflow(ctx,time,ctx->xmax,ctx->ub); for (j=0; jphysics2.bcinflowindex[dof+j]==PETSC_TRUE) { for (i=Mx; iub[dof+j]-x[(2*Mx-(i+1))*dof+j]; } else { for (i=Mx; i sf-lsbwidth-2 && i < sf+1) || (i > fs-2 && i < fs+rsbwidth+1)) { /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); cjmpL = &ctx->cjmpLR[0]; cjmpR = &ctx->cjmpLR[dof]; for (j=0; jRinv[k+j*dof]*jmpL; cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; } } /* Apply limiter to the left and right characteristic jumps */ info.m = dof; info.hxs = hxs; info.hxf = hxf; (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope); for (j=0; jR[j+k*dof]*ctx->cslope[k]; slope[i*dof+j] = tmp; } } } for (i=xs; iuLR[0]; uR = &ctx->uLR[dof]; if (i == sf-lsbwidth) { for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i < xs+xm) { for (j=0; jflux[j]/hxs; islow++; } } if (i > sf-lsbwidth && i < sf) { for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxs; } if (i < xs+xm) { for (j=0; jflux[j]/hxs; islow++; } } if (i == sf) { /* interface between the slow region and the fast region */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxs; } } if (i == fs) { /* interface between the fast region and the slow region */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i < xs+xm) { for (j=0; jflux[j]/hxs; islow++; } } if (i > fs && i < fs+rsbwidth) { for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxs; } if (i < xs+xm) { for (j=0; jflux[j]/hxs; islow++; } } if (i == fs+rsbwidth) { for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxs; } } } ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); PetscFunctionReturn(0); } /* --------------------------------- Finite Volume Solver for fast parts ----------------------------------- */ PetscErrorCode FVRHSFunctionfast_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) { FVCtx *ctx = (FVCtx*)vctx; PetscErrorCode ierr; PetscInt i,j,k,Mx,dof,xs,xm,ifast = 0,sf = ctx->sf,fs = ctx->fs; PetscReal hxs,hxf; PetscScalar *x,*f,*slope; Vec Xloc; DM da; PetscFunctionBeginUser; ierr = TSGetDM(ts,&da);CHKERRQ(ierr); ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); hxs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); ierr = VecZeroEntries(F);CHKERRQ(ierr); ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); ierr = VecGetArray(F,&f);CHKERRQ(ierr); ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); if (ctx->bctype == FVBC_OUTFLOW) { for (i=xs-2; i<0; i++) { for (j=0; jbctype == FVBC_INFLOW) { /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253 pages 137-138 for the scheme. */ if (xs==0) { /* Left Boundary */ ctx->physics2.inflow(ctx,time,ctx->xmin,ctx->ub); for (j=0; jphysics2.bcinflowindex[j]==PETSC_TRUE) { for (i=-2; i<0; i++) x[i*dof+j] = 2.0*ctx->ub[j]-x[-(i+1)*dof+j]; } else { for (i=-2; i<0; i++) x[i*dof+j] = x[j]; /* Outflow */ } } } if (xs+xm==Mx) { /* Right Boundary */ ctx->physics2.inflow(ctx,time,ctx->xmax,ctx->ub); for (j=0; jphysics2.bcinflowindex[dof+j]==PETSC_TRUE) { for (i=Mx; iub[dof+j]-x[(2*Mx-(i+1))*dof+j]; } else { for (i=Mx; i sf-2 && i < fs+1) { ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); cjmpL = &ctx->cjmpLR[0]; cjmpR = &ctx->cjmpLR[dof]; for (j=0; jRinv[k+j*dof]*jmpL; cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; } } /* Apply limiter to the left and right characteristic jumps */ info.m = dof; info.hxs = hxs; info.hxf = hxf; (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope); for (j=0; jR[j+k*dof]*ctx->cslope[k]; slope[i*dof+j] = tmp; } } } for (i=xs; iuLR[0]; uR = &ctx->uLR[dof]; if (i == sf) { /* interface between the slow region and the fast region */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i < xs+xm) { for (j=0; jflux[j]/hxf; ifast++; } } if (i > sf && i < fs) { /* fast region */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxf; } if (i < xs+xm) { for (j=0; jflux[j]/hxf; ifast++; } } if (i == fs) { /* interface between the fast region and the slow region */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxf; } } } ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); PetscFunctionReturn(0); } int main(int argc,char *argv[]) { char lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m"; PetscFunctionList limiters = 0,physics = 0; MPI_Comm comm; TS ts; DM da; Vec X,X0,R; FVCtx ctx; PetscInt i,k,dof,xs,xm,Mx,draw = 0,count_slow,count_fast,islow = 0,ifast =0,islowbuffer = 0,*index_slow,*index_fast,*index_slowbuffer; PetscBool view_final = PETSC_FALSE; PetscReal ptime,maxtime; PetscErrorCode ierr; ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr; comm = PETSC_COMM_WORLD; ierr = PetscMemzero(&ctx,sizeof(ctx));CHKERRQ(ierr); /* Register limiters to be available on the command line */ ierr = PetscFunctionListAdd(&limiters,"upwind" ,Limit2_Upwind);CHKERRQ(ierr); ierr = PetscFunctionListAdd(&limiters,"lax-wendroff" ,Limit2_LaxWendroff);CHKERRQ(ierr); ierr = PetscFunctionListAdd(&limiters,"beam-warming" ,Limit2_BeamWarming);CHKERRQ(ierr); ierr = PetscFunctionListAdd(&limiters,"fromm" ,Limit2_Fromm);CHKERRQ(ierr); ierr = PetscFunctionListAdd(&limiters,"minmod" ,Limit2_Minmod);CHKERRQ(ierr); ierr = PetscFunctionListAdd(&limiters,"superbee" ,Limit2_Superbee);CHKERRQ(ierr); ierr = PetscFunctionListAdd(&limiters,"mc" ,Limit2_MC);CHKERRQ(ierr); ierr = PetscFunctionListAdd(&limiters,"koren3" ,Limit2_Koren3);CHKERRQ(ierr); /* Register physical models to be available on the command line */ ierr = PetscFunctionListAdd(&physics,"shallow" ,PhysicsCreate_Shallow);CHKERRQ(ierr); ierr = PetscFunctionListAdd(&physics,"advect" ,PhysicsCreate_Advect);CHKERRQ(ierr); ctx.comm = comm; ctx.cfl = 0.9; ctx.bctype = FVBC_PERIODIC; ctx.xmin = -1.0; ctx.xmax = 1.0; ctx.initial = 1; ctx.hratio = 2; maxtime = 10.0; ctx.simulation = PETSC_FALSE; ierr = PetscOptionsBegin(comm,NULL,"Finite Volume solver options","");CHKERRQ(ierr); ierr = PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL);CHKERRQ(ierr); ierr = PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL);CHKERRQ(ierr); ierr = PetscOptionsFList("-limit","Name of flux imiter to use","",limiters,lname,lname,sizeof(lname),NULL);CHKERRQ(ierr); ierr = PetscOptionsFList("-physics","Name of physics model to use","",physics,physname,physname,sizeof(physname),NULL);CHKERRQ(ierr); ierr = PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL);CHKERRQ(ierr); 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); ierr = PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL);CHKERRQ(ierr); ierr = PetscOptionsBool("-exact","Compare errors with exact solution","",ctx.exact,&ctx.exact,NULL);CHKERRQ(ierr); ierr = PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL);CHKERRQ(ierr); ierr = PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL);CHKERRQ(ierr); ierr = PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL);CHKERRQ(ierr); ierr = PetscOptionsInt("-hratio","Spacing ratio","",ctx.hratio,&ctx.hratio,NULL);CHKERRQ(ierr); ierr = PetscOptionsEnd();CHKERRQ(ierr); /* Choose the limiter from the list of registered limiters */ ierr = PetscFunctionListFind(limiters,lname,&ctx.limit2);CHKERRQ(ierr); PetscAssertFalse(!ctx.limit2,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Limiter '%s' not found",lname); /* Choose the physics from the list of registered models */ { PetscErrorCode (*r)(FVCtx*); ierr = PetscFunctionListFind(physics,physname,&r);CHKERRQ(ierr); PetscAssertFalse(!r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Physics '%s' not found",physname); /* Create the physics, will set the number of fields and their names */ ierr = (*r)(&ctx);CHKERRQ(ierr); } /* Create a DMDA to manage the parallel grid */ ierr = DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics2.dof,2,NULL,&da);CHKERRQ(ierr); ierr = DMSetFromOptions(da);CHKERRQ(ierr); ierr = DMSetUp(da);CHKERRQ(ierr); /* Inform the DMDA of the field names provided by the physics. */ /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ for (i=0; i ctx.fs+ctx.rsbwidth-1) for (k=0; k= ctx.sf-ctx.lsbwidth && i < ctx.sf) || (i > ctx.fs-1 && i <= ctx.fs+ctx.rsbwidth-1)) for (k=0; k ctx.fs-1) { for (k=0; k ctx.fs-1) for (k=0; k