184ff8c8bSHong Zhang /* 284ff8c8bSHong Zhang Note: 384ff8c8bSHong Zhang -hratio is the ratio between mesh size of corse grids and fine grids 484ff8c8bSHong Zhang */ 584ff8c8bSHong Zhang 684ff8c8bSHong Zhang static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n" 784ff8c8bSHong Zhang " advect - Constant coefficient scalar advection\n" 884ff8c8bSHong Zhang " u_t + (a*u)_x = 0\n" 984ff8c8bSHong Zhang " shallow - 1D Shallow water equations (Saint Venant System)\n" 1084ff8c8bSHong Zhang " h_t + (q)_x = 0 \n" 1184ff8c8bSHong Zhang " q_t + (\frac{q^2}{h} + g/2*h^2)_x = 0 \n" 1284ff8c8bSHong Zhang " where, h(x,t) denotes the height of the water and q(x,t) the momentum.\n" 1384ff8c8bSHong Zhang " for this toy problem, we choose different meshsizes for different sub-domains (slow-fast-slow), say\n" 1484ff8c8bSHong Zhang " hxs = hratio*hxf \n" 1584ff8c8bSHong Zhang " where hxs and hxf are the grid spacings for coarse and fine grids respectively.\n" 1684ff8c8bSHong Zhang " exact - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n" 1784ff8c8bSHong Zhang " the states across shocks and rarefactions\n" 1884ff8c8bSHong Zhang " simulation - use reference solution which is generated by smaller time step size to be true solution,\n" 1984ff8c8bSHong Zhang " also the reference solution should be generated by user and stored in a binary file.\n" 2084ff8c8bSHong Zhang " characteristic - Limit the characteristic variables, this is usually preferred (default)\n" 2184ff8c8bSHong Zhang " bc_type - Boundary condition for the problem, options are: periodic, outflow, inflow " 2284ff8c8bSHong Zhang "Several problem descriptions (initial data, physics specific features, boundary data) can be chosen with -initial N\n\n" 2384ff8c8bSHong Zhang "The problem size should be set with -da_grid_x M\n\n"; 2484ff8c8bSHong Zhang 2584ff8c8bSHong Zhang /* 2684ff8c8bSHong Zhang Example: 2784ff8c8bSHong Zhang ./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 2884ff8c8bSHong Zhang ./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 2984ff8c8bSHong Zhang ./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 3084ff8c8bSHong Zhang ./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 3184ff8c8bSHong Zhang ./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 3284ff8c8bSHong Zhang 3384ff8c8bSHong Zhang Contributed by: Aidan Hamilton <aidan@udel.edu> 3484ff8c8bSHong Zhang */ 3584ff8c8bSHong Zhang 3684ff8c8bSHong Zhang #include <petscts.h> 3784ff8c8bSHong Zhang #include <petscdm.h> 3884ff8c8bSHong Zhang #include <petscdmda.h> 3984ff8c8bSHong Zhang #include <petscdraw.h> 4084ff8c8bSHong Zhang #include "finitevolume1d.h" 4184ff8c8bSHong Zhang #include <petsc/private/kernels/blockinvert.h> 4284ff8c8bSHong Zhang 439fbee547SJacob Faibussowitsch static inline PetscReal RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax) { PetscReal range = xmax-xmin; return xmin +PetscFmodReal(range+PetscFmodReal(a,range),range); } 449fbee547SJacob Faibussowitsch static inline PetscReal MaxAbs(PetscReal a,PetscReal b) { return (PetscAbs(a) > PetscAbs(b)) ? a : b; } 4584ff8c8bSHong Zhang /* --------------------------------- Advection ----------------------------------- */ 4684ff8c8bSHong Zhang typedef struct { 4784ff8c8bSHong Zhang PetscReal a; /* advective velocity */ 4884ff8c8bSHong Zhang } AdvectCtx; 4984ff8c8bSHong Zhang 5084ff8c8bSHong Zhang static PetscErrorCode PhysicsRiemann_Advect(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed) 5184ff8c8bSHong Zhang { 5284ff8c8bSHong Zhang AdvectCtx *ctx = (AdvectCtx*)vctx; 5384ff8c8bSHong Zhang PetscReal speed; 5484ff8c8bSHong Zhang 5584ff8c8bSHong Zhang PetscFunctionBeginUser; 5684ff8c8bSHong Zhang speed = ctx->a; 5784ff8c8bSHong Zhang flux[0] = PetscMax(0,speed)*uL[0] + PetscMin(0,speed)*uR[0]; 5884ff8c8bSHong Zhang *maxspeed = speed; 5984ff8c8bSHong Zhang PetscFunctionReturn(0); 6084ff8c8bSHong Zhang } 6184ff8c8bSHong Zhang 6284ff8c8bSHong Zhang static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds) 6384ff8c8bSHong Zhang { 6484ff8c8bSHong Zhang AdvectCtx *ctx = (AdvectCtx*)vctx; 6584ff8c8bSHong Zhang 6684ff8c8bSHong Zhang PetscFunctionBeginUser; 6784ff8c8bSHong Zhang X[0] = 1.; 6884ff8c8bSHong Zhang Xi[0] = 1.; 6984ff8c8bSHong Zhang speeds[0] = ctx->a; 7084ff8c8bSHong Zhang PetscFunctionReturn(0); 7184ff8c8bSHong Zhang } 7284ff8c8bSHong Zhang 7384ff8c8bSHong Zhang static PetscErrorCode PhysicsSample_Advect(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u) 7484ff8c8bSHong Zhang { 7584ff8c8bSHong Zhang AdvectCtx *ctx = (AdvectCtx*)vctx; 7684ff8c8bSHong Zhang PetscReal a = ctx->a,x0; 7784ff8c8bSHong Zhang 7884ff8c8bSHong Zhang PetscFunctionBeginUser; 7984ff8c8bSHong Zhang switch (bctype) { 8084ff8c8bSHong Zhang case FVBC_OUTFLOW: x0 = x-a*t; break; 8184ff8c8bSHong Zhang case FVBC_PERIODIC: x0 = RangeMod(x-a*t,xmin,xmax); break; 8284ff8c8bSHong Zhang default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown BCType"); 8384ff8c8bSHong Zhang } 8484ff8c8bSHong Zhang switch (initial) { 8584ff8c8bSHong Zhang case 0: u[0] = (x0 < 0) ? 1 : -1; break; 8684ff8c8bSHong Zhang case 1: u[0] = (x0 < 0) ? -1 : 1; break; 8784ff8c8bSHong Zhang case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break; 8884ff8c8bSHong Zhang case 3: u[0] = PetscSinReal(2*PETSC_PI*x0); break; 8984ff8c8bSHong Zhang case 4: u[0] = PetscAbs(x0); break; 9084ff8c8bSHong Zhang case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2*PETSC_PI*x0)); break; 9184ff8c8bSHong Zhang case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2-x0 : 0)); break; 9284ff8c8bSHong Zhang case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI*x0),10.0);break; 9384ff8c8bSHong Zhang default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition"); 9484ff8c8bSHong Zhang } 9584ff8c8bSHong Zhang PetscFunctionReturn(0); 9684ff8c8bSHong Zhang } 9784ff8c8bSHong Zhang 9884ff8c8bSHong Zhang static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx) 9984ff8c8bSHong Zhang { 10084ff8c8bSHong Zhang PetscErrorCode ierr; 10184ff8c8bSHong Zhang AdvectCtx *user; 10284ff8c8bSHong Zhang 10384ff8c8bSHong Zhang PetscFunctionBeginUser; 1045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&user)); 10584ff8c8bSHong Zhang ctx->physics2.sample2 = PhysicsSample_Advect; 10684ff8c8bSHong Zhang ctx->physics2.riemann2 = PhysicsRiemann_Advect; 10784ff8c8bSHong Zhang ctx->physics2.characteristic2 = PhysicsCharacteristic_Advect; 10884ff8c8bSHong Zhang ctx->physics2.destroy = PhysicsDestroy_SimpleFree; 10984ff8c8bSHong Zhang ctx->physics2.user = user; 11084ff8c8bSHong Zhang ctx->physics2.dof = 1; 11184ff8c8bSHong Zhang 1125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrallocpy("u",&ctx->physics2.fieldname[0])); 11384ff8c8bSHong Zhang user->a = 1; 11484ff8c8bSHong Zhang ierr = PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");CHKERRQ(ierr); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-physics_advect_a","Speed","",user->a,&user->a,NULL)); 11684ff8c8bSHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 11784ff8c8bSHong Zhang PetscFunctionReturn(0); 11884ff8c8bSHong Zhang } 11984ff8c8bSHong Zhang /* --------------------------------- Shallow Water ----------------------------------- */ 12084ff8c8bSHong Zhang 12184ff8c8bSHong Zhang typedef struct { 12284ff8c8bSHong Zhang PetscReal gravity; 12384ff8c8bSHong Zhang } ShallowCtx; 12484ff8c8bSHong Zhang 1259fbee547SJacob Faibussowitsch static inline void ShallowFlux(ShallowCtx *phys,const PetscScalar *u,PetscScalar *f) 12684ff8c8bSHong Zhang { 12784ff8c8bSHong Zhang f[0] = u[1]; 12884ff8c8bSHong Zhang f[1] = PetscSqr(u[1])/u[0] + 0.5*phys->gravity*PetscSqr(u[0]); 12984ff8c8bSHong Zhang } 13084ff8c8bSHong Zhang 13184ff8c8bSHong Zhang static PetscErrorCode PhysicsRiemann_Shallow_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed) 13284ff8c8bSHong Zhang { 13384ff8c8bSHong Zhang ShallowCtx *phys = (ShallowCtx*)vctx; 13484ff8c8bSHong Zhang PetscScalar g = phys->gravity,ustar[2],cL,cR,c,cstar; 13584ff8c8bSHong Zhang struct {PetscScalar h,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]},star; 13684ff8c8bSHong Zhang PetscInt i; 13784ff8c8bSHong Zhang 13884ff8c8bSHong Zhang PetscFunctionBeginUser; 1393c633725SBarry Smith PetscCheck((L.h > 0 && R.h > 0),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative"); 14084ff8c8bSHong Zhang cL = PetscSqrtScalar(g*L.h); 14184ff8c8bSHong Zhang cR = PetscSqrtScalar(g*R.h); 14284ff8c8bSHong Zhang c = PetscMax(cL,cR); 14384ff8c8bSHong Zhang { 14484ff8c8bSHong Zhang /* Solve for star state */ 14584ff8c8bSHong Zhang const PetscInt maxits = 50; 14684ff8c8bSHong Zhang PetscScalar tmp,res,res0=0,h0,h = 0.5*(L.h + R.h); /* initial guess */ 14784ff8c8bSHong Zhang h0 = h; 14884ff8c8bSHong Zhang for (i=0; i<maxits; i++) { 14984ff8c8bSHong Zhang PetscScalar fr,fl,dfr,dfl; 15084ff8c8bSHong Zhang fl = (L.h < h) 15184ff8c8bSHong Zhang ? PetscSqrtScalar(0.5*g*(h*h - L.h*L.h)*(1/L.h - 1/h)) /* shock */ 15284ff8c8bSHong Zhang : 2*PetscSqrtScalar(g*h) - 2*PetscSqrtScalar(g*L.h); /* rarefaction */ 15384ff8c8bSHong Zhang fr = (R.h < h) 15484ff8c8bSHong Zhang ? PetscSqrtScalar(0.5*g*(h*h - R.h*R.h)*(1/R.h - 1/h)) /* shock */ 15584ff8c8bSHong Zhang : 2*PetscSqrtScalar(g*h) - 2*PetscSqrtScalar(g*R.h); /* rarefaction */ 15684ff8c8bSHong Zhang res = R.u - L.u + fr + fl; 1572c71b3e2SJacob Faibussowitsch PetscCheckFalse(PetscIsInfOrNanScalar(res),PETSC_COMM_SELF,PETSC_ERR_FP,"Infinity or Not-a-Number generated in computation"); 15884ff8c8bSHong Zhang if (PetscAbsScalar(res) < PETSC_SQRT_MACHINE_EPSILON || (i > 0 && PetscAbsScalar(h-h0) < PETSC_SQRT_MACHINE_EPSILON)) { 15984ff8c8bSHong Zhang star.h = h; 16084ff8c8bSHong Zhang star.u = L.u - fl; 16184ff8c8bSHong Zhang goto converged; 16284ff8c8bSHong Zhang } else if (i > 0 && PetscAbsScalar(res) >= PetscAbsScalar(res0)) { /* Line search */ 16384ff8c8bSHong Zhang h = 0.8*h0 + 0.2*h; 16484ff8c8bSHong Zhang continue; 16584ff8c8bSHong Zhang } 16684ff8c8bSHong Zhang /* Accept the last step and take another */ 16784ff8c8bSHong Zhang res0 = res; 16884ff8c8bSHong Zhang h0 = h; 16984ff8c8bSHong Zhang dfl = (L.h < h) ? 0.5/fl*0.5*g*(-L.h*L.h/(h*h) - 1 + 2*h/L.h) : PetscSqrtScalar(g/h); 17084ff8c8bSHong Zhang dfr = (R.h < h) ? 0.5/fr*0.5*g*(-R.h*R.h/(h*h) - 1 + 2*h/R.h) : PetscSqrtScalar(g/h); 17184ff8c8bSHong Zhang tmp = h - res/(dfr+dfl); 17284ff8c8bSHong Zhang if (tmp <= 0) h /= 2; /* Guard against Newton shooting off to a negative thickness */ 17384ff8c8bSHong Zhang else h = tmp; 1743c633725SBarry Smith PetscCheck(((h > 0) && PetscIsNormalScalar(h)),PETSC_COMM_SELF,PETSC_ERR_FP,"non-normal iterate h=%g",(double)h); 17584ff8c8bSHong Zhang } 17698921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,1,"Newton iteration for star.h diverged after %D iterations",i); 17784ff8c8bSHong Zhang } 17884ff8c8bSHong Zhang converged: 17984ff8c8bSHong Zhang cstar = PetscSqrtScalar(g*star.h); 18084ff8c8bSHong Zhang if (L.u-cL < 0 && 0 < star.u-cstar) { /* 1-wave is sonic rarefaction */ 18184ff8c8bSHong Zhang PetscScalar ufan[2]; 18284ff8c8bSHong Zhang ufan[0] = 1/g*PetscSqr(L.u/3 + 2./3*cL); 18384ff8c8bSHong Zhang ufan[1] = PetscSqrtScalar(g*ufan[0])*ufan[0]; 18484ff8c8bSHong Zhang ShallowFlux(phys,ufan,flux); 18584ff8c8bSHong Zhang } else if (star.u+cstar < 0 && 0 < R.u+cR) { /* 2-wave is sonic rarefaction */ 18684ff8c8bSHong Zhang PetscScalar ufan[2]; 18784ff8c8bSHong Zhang ufan[0] = 1/g*PetscSqr(R.u/3 - 2./3*cR); 18884ff8c8bSHong Zhang ufan[1] = -PetscSqrtScalar(g*ufan[0])*ufan[0]; 18984ff8c8bSHong Zhang ShallowFlux(phys,ufan,flux); 19084ff8c8bSHong Zhang } else if ((L.h >= star.h && L.u-c >= 0) || (L.h<star.h && (star.h*star.u-L.h*L.u)/(star.h-L.h) > 0)) { 19184ff8c8bSHong Zhang /* 1-wave is right-travelling shock (supersonic) */ 19284ff8c8bSHong Zhang ShallowFlux(phys,uL,flux); 19384ff8c8bSHong Zhang } 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)) { 19484ff8c8bSHong Zhang /* 2-wave is left-travelling shock (supersonic) */ 19584ff8c8bSHong Zhang ShallowFlux(phys,uR,flux); 19684ff8c8bSHong Zhang } else { 19784ff8c8bSHong Zhang ustar[0] = star.h; 19884ff8c8bSHong Zhang ustar[1] = star.h*star.u; 19984ff8c8bSHong Zhang ShallowFlux(phys,ustar,flux); 20084ff8c8bSHong Zhang } 20184ff8c8bSHong Zhang *maxspeed = MaxAbs(MaxAbs(star.u-cstar,star.u+cstar),MaxAbs(L.u-cL,R.u+cR)); 20284ff8c8bSHong Zhang PetscFunctionReturn(0); 20384ff8c8bSHong Zhang } 20484ff8c8bSHong Zhang 20584ff8c8bSHong Zhang static PetscErrorCode PhysicsRiemann_Shallow_Rusanov(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed) 20684ff8c8bSHong Zhang { 20784ff8c8bSHong Zhang ShallowCtx *phys = (ShallowCtx*)vctx; 20884ff8c8bSHong Zhang PetscScalar g = phys->gravity,fL[2],fR[2],s; 20984ff8c8bSHong Zhang struct {PetscScalar h,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]}; 21084ff8c8bSHong Zhang PetscReal tol = 1e-6; 21184ff8c8bSHong Zhang 21284ff8c8bSHong Zhang PetscFunctionBeginUser; 21384ff8c8bSHong Zhang /* Positivity preserving modification*/ 21484ff8c8bSHong Zhang if (L.h < tol) L.u = 0.0; 21584ff8c8bSHong Zhang if (R.h < tol) R.u = 0.0; 21684ff8c8bSHong Zhang 21784ff8c8bSHong Zhang /*simple pos preserve limiter*/ 21884ff8c8bSHong Zhang if (L.h < 0) L.h = 0; 21984ff8c8bSHong Zhang if (R.h < 0) R.h = 0; 22084ff8c8bSHong Zhang 22184ff8c8bSHong Zhang ShallowFlux(phys,uL,fL); 22284ff8c8bSHong Zhang ShallowFlux(phys,uR,fR); 22384ff8c8bSHong Zhang 22484ff8c8bSHong Zhang s = PetscMax(PetscAbs(L.u)+PetscSqrtScalar(g*L.h),PetscAbs(R.u)+PetscSqrtScalar(g*R.h)); 22584ff8c8bSHong Zhang flux[0] = 0.5*(fL[0] + fR[0]) + 0.5*s*(L.h - R.h); 22684ff8c8bSHong Zhang flux[1] = 0.5*(fL[1] + fR[1]) + 0.5*s*(uL[1] - uR[1]); 22784ff8c8bSHong Zhang *maxspeed = s; 22884ff8c8bSHong Zhang PetscFunctionReturn(0); 22984ff8c8bSHong Zhang } 23084ff8c8bSHong Zhang 23184ff8c8bSHong Zhang static PetscErrorCode PhysicsCharacteristic_Conservative(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds) 23284ff8c8bSHong Zhang { 23384ff8c8bSHong Zhang PetscInt i,j; 23484ff8c8bSHong Zhang 23584ff8c8bSHong Zhang PetscFunctionBeginUser; 23684ff8c8bSHong Zhang for (i=0; i<m; i++) { 23784ff8c8bSHong Zhang for (j=0; j<m; j++) Xi[i*m+j] = X[i*m+j] = (PetscScalar)(i==j); 23884ff8c8bSHong Zhang speeds[i] = PETSC_MAX_REAL; /* Indicates invalid */ 23984ff8c8bSHong Zhang } 24084ff8c8bSHong Zhang PetscFunctionReturn(0); 24184ff8c8bSHong Zhang } 24284ff8c8bSHong Zhang 24384ff8c8bSHong Zhang static PetscErrorCode PhysicsCharacteristic_Shallow(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds) 24484ff8c8bSHong Zhang { 24584ff8c8bSHong Zhang ShallowCtx *phys = (ShallowCtx*)vctx; 24684ff8c8bSHong Zhang PetscReal c; 24784ff8c8bSHong Zhang PetscReal tol = 1e-6; 24884ff8c8bSHong Zhang 24984ff8c8bSHong Zhang PetscFunctionBeginUser; 25084ff8c8bSHong Zhang c = PetscSqrtScalar(u[0]*phys->gravity); 25184ff8c8bSHong Zhang 25284ff8c8bSHong Zhang if (u[0] < tol) { /*Use conservative variables*/ 25384ff8c8bSHong Zhang X[0*2+0] = 1; 25484ff8c8bSHong Zhang X[0*2+1] = 0; 25584ff8c8bSHong Zhang X[1*2+0] = 0; 25684ff8c8bSHong Zhang X[1*2+1] = 1; 25784ff8c8bSHong Zhang } else { 25884ff8c8bSHong Zhang speeds[0] = u[1]/u[0] - c; 25984ff8c8bSHong Zhang speeds[1] = u[1]/u[0] + c; 26084ff8c8bSHong Zhang X[0*2+0] = 1; 26184ff8c8bSHong Zhang X[0*2+1] = speeds[0]; 26284ff8c8bSHong Zhang X[1*2+0] = 1; 26384ff8c8bSHong Zhang X[1*2+1] = speeds[1]; 26484ff8c8bSHong Zhang } 26584ff8c8bSHong Zhang 2665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(Xi,X,4)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscKernel_A_gets_inverse_A_2(Xi,0,PETSC_FALSE,NULL)); 26884ff8c8bSHong Zhang PetscFunctionReturn(0); 26984ff8c8bSHong Zhang } 27084ff8c8bSHong Zhang 27184ff8c8bSHong Zhang static PetscErrorCode PhysicsSample_Shallow(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u) 27284ff8c8bSHong Zhang { 27384ff8c8bSHong Zhang PetscFunctionBeginUser; 2742c71b3e2SJacob Faibussowitsch PetscCheckFalse(t > 0,PETSC_COMM_SELF,PETSC_ERR_SUP,"Exact solutions not implemented for t > 0"); 27584ff8c8bSHong Zhang switch (initial) { 27684ff8c8bSHong Zhang case 0: 27784ff8c8bSHong Zhang u[0] = (x < 0) ? 2 : 1; /* Standard Dam Break Problem */ 27884ff8c8bSHong Zhang u[1] = (x < 0) ? 0 : 0; 27984ff8c8bSHong Zhang break; 28084ff8c8bSHong Zhang case 1: 28184ff8c8bSHong Zhang u[0] = (x < 10) ? 1 : 0.1; /*The Next 5 problems are standard Riemann problem tests */ 28284ff8c8bSHong Zhang u[1] = (x < 10) ? 2.5 : 0; 28384ff8c8bSHong Zhang break; 28484ff8c8bSHong Zhang case 2: 28584ff8c8bSHong Zhang u[0] = (x < 25) ? 1 : 1; 28684ff8c8bSHong Zhang u[1] = (x < 25) ? -5 : 5; 28784ff8c8bSHong Zhang break; 28884ff8c8bSHong Zhang case 3: 28984ff8c8bSHong Zhang u[0] = (x < 20) ? 1 : 1e-12; 29084ff8c8bSHong Zhang u[1] = (x < 20) ? 0 : 0; 29184ff8c8bSHong Zhang break; 29284ff8c8bSHong Zhang case 4: 29384ff8c8bSHong Zhang u[0] = (x < 30) ? 1e-12 : 1; 29484ff8c8bSHong Zhang u[1] = (x < 30) ? 0 : 0; 29584ff8c8bSHong Zhang break; 29684ff8c8bSHong Zhang case 5: 29784ff8c8bSHong Zhang u[0] = (x < 25) ? 0.1 : 0.1; 29884ff8c8bSHong Zhang u[1] = (x < 25) ? -0.3 : 0.3; 29984ff8c8bSHong Zhang break; 30084ff8c8bSHong Zhang case 6: 30184ff8c8bSHong Zhang u[0] = 1+0.5*PetscSinReal(2*PETSC_PI*x); 30284ff8c8bSHong Zhang u[1] = 1*u[0]; 30384ff8c8bSHong Zhang break; 30484ff8c8bSHong Zhang case 7: 30584ff8c8bSHong Zhang if (x < -0.1) { 30684ff8c8bSHong Zhang u[0] = 1e-9; 30784ff8c8bSHong Zhang u[1] = 0.0; 30884ff8c8bSHong Zhang } else if (x < 0.1) { 30984ff8c8bSHong Zhang u[0] = 1.0; 31084ff8c8bSHong Zhang u[1] = 0.0; 31184ff8c8bSHong Zhang } else { 31284ff8c8bSHong Zhang u[0] = 1e-9; 31384ff8c8bSHong Zhang u[1] = 0.0; 31484ff8c8bSHong Zhang } 31584ff8c8bSHong Zhang break; 31684ff8c8bSHong Zhang case 8: 31784ff8c8bSHong Zhang if (x < -0.1) { 31884ff8c8bSHong Zhang u[0] = 2; 31984ff8c8bSHong Zhang u[1] = 0.0; 32084ff8c8bSHong Zhang } else if (x < 0.1) { 32184ff8c8bSHong Zhang u[0] = 3.0; 32284ff8c8bSHong Zhang u[1] = 0.0; 32384ff8c8bSHong Zhang } else { 32484ff8c8bSHong Zhang u[0] = 2; 32584ff8c8bSHong Zhang u[1] = 0.0; 32684ff8c8bSHong Zhang } 32784ff8c8bSHong Zhang break; 32884ff8c8bSHong Zhang default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition"); 32984ff8c8bSHong Zhang } 33084ff8c8bSHong Zhang PetscFunctionReturn(0); 33184ff8c8bSHong Zhang } 33284ff8c8bSHong Zhang 33384ff8c8bSHong Zhang /* Implements inflow conditions for the given initial conditions. Which conditions are actually enforced depends on 33484ff8c8bSHong Zhang on the results of PhysicsSetInflowType_Shallow. */ 33584ff8c8bSHong Zhang static PetscErrorCode PhysicsInflow_Shallow(void *vctx,PetscReal t,PetscReal x,PetscReal *u) 33684ff8c8bSHong Zhang { 33784ff8c8bSHong Zhang FVCtx *ctx = (FVCtx*)vctx; 33884ff8c8bSHong Zhang 33984ff8c8bSHong Zhang PetscFunctionBeginUser; 34084ff8c8bSHong Zhang if (ctx->bctype == FVBC_INFLOW) { 34184ff8c8bSHong Zhang switch (ctx->initial) { 34284ff8c8bSHong Zhang case 0: 34384ff8c8bSHong Zhang case 1: 34484ff8c8bSHong Zhang case 2: 34584ff8c8bSHong Zhang case 3: 34684ff8c8bSHong Zhang case 4: 34784ff8c8bSHong Zhang case 5: 34884ff8c8bSHong Zhang u[0] = 0; u[1] = 0.0; /* Left boundary conditions */ 34984ff8c8bSHong Zhang u[2] = 0; u[3] = 0.0; /* Right boundary conditions */ 35084ff8c8bSHong Zhang break; 35184ff8c8bSHong Zhang case 6: 35284ff8c8bSHong Zhang u[0] = 0; u[1] = 0.0; /* Left boundary conditions */ 35384ff8c8bSHong Zhang u[2] = 0; u[3] = 0.0; /* Right boundary conditions */ 35484ff8c8bSHong Zhang break; 35584ff8c8bSHong Zhang case 7: 35684ff8c8bSHong Zhang u[0] = 0; u[1] = 0.0; /* Left boundary conditions */ 35784ff8c8bSHong Zhang u[2] = 0; u[3] = 0.0; /* Right boundary conditions */ 35884ff8c8bSHong Zhang break; 35984ff8c8bSHong Zhang case 8: 36084ff8c8bSHong Zhang u[0] = 0; u[1] = 1.0; /* Left boundary conditions */ 36184ff8c8bSHong Zhang u[2] = 0; u[3] = -1.0;/* Right boundary conditions */ 36284ff8c8bSHong Zhang break; 36384ff8c8bSHong Zhang default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition"); 36484ff8c8bSHong Zhang } 36584ff8c8bSHong Zhang } 36684ff8c8bSHong Zhang PetscFunctionReturn(0); 36784ff8c8bSHong Zhang } 36884ff8c8bSHong Zhang 36984ff8c8bSHong Zhang /* Selects which boundary conditions are marked as inflow and which as outflow when FVBC_INFLOW is selected. */ 37084ff8c8bSHong Zhang static PetscErrorCode PhysicsSetInflowType_Shallow(FVCtx *ctx) 37184ff8c8bSHong Zhang { 37284ff8c8bSHong Zhang PetscFunctionBeginUser; 37384ff8c8bSHong Zhang switch (ctx->initial) { 37484ff8c8bSHong Zhang case 0: 37584ff8c8bSHong Zhang case 1: 37684ff8c8bSHong Zhang case 2: 37784ff8c8bSHong Zhang case 3: 37884ff8c8bSHong Zhang case 4: 37984ff8c8bSHong Zhang case 5: 38084ff8c8bSHong Zhang case 6: 38184ff8c8bSHong Zhang case 7: 38284ff8c8bSHong Zhang case 8: /* Fix left and right momentum, height is outflow */ 38384ff8c8bSHong Zhang ctx->physics2.bcinflowindex[0] = PETSC_FALSE; 38484ff8c8bSHong Zhang ctx->physics2.bcinflowindex[1] = PETSC_TRUE; 38584ff8c8bSHong Zhang ctx->physics2.bcinflowindex[2] = PETSC_FALSE; 38684ff8c8bSHong Zhang ctx->physics2.bcinflowindex[3] = PETSC_TRUE; 38784ff8c8bSHong Zhang break; 38884ff8c8bSHong Zhang default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition"); 38984ff8c8bSHong Zhang } 39084ff8c8bSHong Zhang PetscFunctionReturn(0); 39184ff8c8bSHong Zhang } 39284ff8c8bSHong Zhang 39384ff8c8bSHong Zhang static PetscErrorCode PhysicsCreate_Shallow(FVCtx *ctx) 39484ff8c8bSHong Zhang { 39584ff8c8bSHong Zhang PetscErrorCode ierr; 39684ff8c8bSHong Zhang ShallowCtx *user; 39784ff8c8bSHong Zhang PetscFunctionList rlist = 0,rclist = 0; 39884ff8c8bSHong Zhang char rname[256] = "rusanov",rcname[256] = "characteristic"; 39984ff8c8bSHong Zhang 40084ff8c8bSHong Zhang PetscFunctionBeginUser; 4015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&user)); 40284ff8c8bSHong Zhang ctx->physics2.sample2 = PhysicsSample_Shallow; 40384ff8c8bSHong Zhang ctx->physics2.inflow = PhysicsInflow_Shallow; 40484ff8c8bSHong Zhang ctx->physics2.destroy = PhysicsDestroy_SimpleFree; 40584ff8c8bSHong Zhang ctx->physics2.riemann2 = PhysicsRiemann_Shallow_Rusanov; 40684ff8c8bSHong Zhang ctx->physics2.characteristic2 = PhysicsCharacteristic_Shallow; 40784ff8c8bSHong Zhang ctx->physics2.user = user; 40884ff8c8bSHong Zhang ctx->physics2.dof = 2; 40984ff8c8bSHong Zhang 41084ff8c8bSHong Zhang PetscMalloc1(2*(ctx->physics2.dof),&ctx->physics2.bcinflowindex); 41184ff8c8bSHong Zhang PhysicsSetInflowType_Shallow(ctx); 41284ff8c8bSHong Zhang 4135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrallocpy("height",&ctx->physics2.fieldname[0])); 4145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrallocpy("momentum",&ctx->physics2.fieldname[1])); 41584ff8c8bSHong Zhang 41684ff8c8bSHong Zhang user->gravity = 9.81; 41784ff8c8bSHong Zhang 4185f80ce2aSJacob Faibussowitsch CHKERRQ(RiemannListAdd_2WaySplit(&rlist,"exact", PhysicsRiemann_Shallow_Exact)); 4195f80ce2aSJacob Faibussowitsch CHKERRQ(RiemannListAdd_2WaySplit(&rlist,"rusanov",PhysicsRiemann_Shallow_Rusanov)); 4205f80ce2aSJacob Faibussowitsch CHKERRQ(ReconstructListAdd_2WaySplit(&rclist,"characteristic",PhysicsCharacteristic_Shallow)); 4215f80ce2aSJacob Faibussowitsch CHKERRQ(ReconstructListAdd_2WaySplit(&rclist,"conservative",PhysicsCharacteristic_Conservative)); 42284ff8c8bSHong Zhang ierr = PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for Shallow","");CHKERRQ(ierr); 4235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-physics_shallow_gravity","Gravity","",user->gravity,&user->gravity,NULL)); 4245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsFList("-physics_shallow_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL)); 4255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsFList("-physics_shallow_reconstruct","Reconstruction","",rclist,rcname,rcname,sizeof(rcname),NULL)); 42684ff8c8bSHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 4275f80ce2aSJacob Faibussowitsch CHKERRQ(RiemannListFind_2WaySplit(rlist,rname,&ctx->physics2.riemann2)); 4285f80ce2aSJacob Faibussowitsch CHKERRQ(ReconstructListFind_2WaySplit(rclist,rcname,&ctx->physics2.characteristic2)); 4295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListDestroy(&rlist)); 4305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListDestroy(&rclist)); 43184ff8c8bSHong Zhang PetscFunctionReturn(0); 43284ff8c8bSHong Zhang } 43384ff8c8bSHong Zhang 43484ff8c8bSHong Zhang PetscErrorCode FVSample_2WaySplit(FVCtx *ctx,DM da,PetscReal time,Vec U) 43584ff8c8bSHong Zhang { 43684ff8c8bSHong Zhang PetscScalar *u,*uj,xj,xi; 43784ff8c8bSHong Zhang PetscInt i,j,k,dof,xs,xm,Mx; 43884ff8c8bSHong Zhang const PetscInt N = 200; 43984ff8c8bSHong Zhang PetscReal hs,hf; 44084ff8c8bSHong Zhang 44184ff8c8bSHong Zhang PetscFunctionBeginUser; 4423c633725SBarry Smith PetscCheck(ctx->physics2.sample2,PETSC_COMM_SELF,PETSC_ERR_SUP,"Physics has not provided a sampling function"); 4435f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 4445f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 4455f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,U,&u)); 4465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(dof,&uj)); 44784ff8c8bSHong Zhang hs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 44884ff8c8bSHong Zhang hf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 44984ff8c8bSHong Zhang for (i=xs; i<xs+xm; i++) { 45084ff8c8bSHong Zhang if (i < ctx->sf) { 45184ff8c8bSHong Zhang xi = ctx->xmin+0.5*hs+i*hs; 45284ff8c8bSHong Zhang /* Integrate over cell i using trapezoid rule with N points. */ 45384ff8c8bSHong Zhang for (k=0; k<dof; k++) u[i*dof+k] = 0; 45484ff8c8bSHong Zhang for (j=0; j<N+1; j++) { 45584ff8c8bSHong Zhang xj = xi+hs*(j-N/2)/(PetscReal)N; 4565f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj)); 45784ff8c8bSHong Zhang for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 45884ff8c8bSHong Zhang } 45984ff8c8bSHong Zhang } else if (i < ctx->fs) { 46084ff8c8bSHong Zhang xi = ctx->xmin+ctx->sf*hs+0.5*hf+(i-ctx->sf)*hf; 46184ff8c8bSHong Zhang /* Integrate over cell i using trapezoid rule with N points. */ 46284ff8c8bSHong Zhang for (k=0; k<dof; k++) u[i*dof+k] = 0; 46384ff8c8bSHong Zhang for (j=0; j<N+1; j++) { 46484ff8c8bSHong Zhang xj = xi+hf*(j-N/2)/(PetscReal)N; 4655f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj)); 46684ff8c8bSHong Zhang for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 46784ff8c8bSHong Zhang } 46884ff8c8bSHong Zhang } else { 46984ff8c8bSHong Zhang xi = ctx->xmin+ctx->sf*hs+(ctx->fs-ctx->sf)*hf+0.5*hs+(i-ctx->fs)*hs; 47084ff8c8bSHong Zhang /* Integrate over cell i using trapezoid rule with N points. */ 47184ff8c8bSHong Zhang for (k=0; k<dof; k++) u[i*dof+k] = 0; 47284ff8c8bSHong Zhang for (j=0; j<N+1; j++) { 47384ff8c8bSHong Zhang xj = xi+hs*(j-N/2)/(PetscReal)N; 4745f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj)); 47584ff8c8bSHong Zhang for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 47684ff8c8bSHong Zhang } 47784ff8c8bSHong Zhang } 47884ff8c8bSHong Zhang } 4795f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,U,&u)); 4805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(uj)); 48184ff8c8bSHong Zhang PetscFunctionReturn(0); 48284ff8c8bSHong Zhang } 48384ff8c8bSHong Zhang 48484ff8c8bSHong Zhang static PetscErrorCode SolutionErrorNorms_2WaySplit(FVCtx *ctx,DM da,PetscReal t,Vec X,PetscReal *nrm1) 48584ff8c8bSHong Zhang { 48684ff8c8bSHong Zhang Vec Y; 48784ff8c8bSHong Zhang PetscInt i,Mx; 48884ff8c8bSHong Zhang const PetscScalar *ptr_X,*ptr_Y; 48984ff8c8bSHong Zhang PetscReal hs,hf; 49084ff8c8bSHong Zhang 49184ff8c8bSHong Zhang PetscFunctionBeginUser; 4925f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(X,&Mx)); 4935f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(X,&Y)); 4945f80ce2aSJacob Faibussowitsch CHKERRQ(FVSample_2WaySplit(ctx,da,t,Y)); 49584ff8c8bSHong Zhang hs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 49684ff8c8bSHong Zhang hf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 4975f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X,&ptr_X)); 4985f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(Y,&ptr_Y)); 49984ff8c8bSHong Zhang for (i=0; i<Mx; i++) { 50084ff8c8bSHong Zhang if (i < ctx->sf || i > ctx->fs -1) *nrm1 += hs*PetscAbs(ptr_X[i]-ptr_Y[i]); 50184ff8c8bSHong Zhang else *nrm1 += hf*PetscAbs(ptr_X[i]-ptr_Y[i]); 50284ff8c8bSHong Zhang } 5035f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X,&ptr_X)); 5045f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(Y,&ptr_Y)); 5055f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Y)); 50684ff8c8bSHong Zhang PetscFunctionReturn(0); 50784ff8c8bSHong Zhang } 50884ff8c8bSHong Zhang 50984ff8c8bSHong Zhang PetscErrorCode FVRHSFunction_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 51084ff8c8bSHong Zhang { 51184ff8c8bSHong Zhang FVCtx *ctx = (FVCtx*)vctx; 51284ff8c8bSHong Zhang PetscInt i,j,k,Mx,dof,xs,xm,sf = ctx->sf,fs = ctx->fs; 51384ff8c8bSHong Zhang PetscReal hxf,hxs,cfl_idt = 0; 51484ff8c8bSHong Zhang PetscScalar *x,*f,*slope; 51584ff8c8bSHong Zhang Vec Xloc; 51684ff8c8bSHong Zhang DM da; 51784ff8c8bSHong Zhang 51884ff8c8bSHong Zhang PetscFunctionBeginUser; 5195f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 5205f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&Xloc)); /* Xloc contains ghost points */ 5215f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); /* Mx is the number of center points */ 52284ff8c8bSHong Zhang hxs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 52384ff8c8bSHong Zhang hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 5245f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); /* X is solution vector which does not contain ghost points */ 5255f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 52684ff8c8bSHong Zhang 5275f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(F)); /* F is the right hand side function corresponds to center points */ 52884ff8c8bSHong Zhang 5295f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,Xloc,&x)); 5305f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,F,&f)); 5315f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope)); /* contains ghost points */ 5325f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 53384ff8c8bSHong Zhang 53484ff8c8bSHong Zhang if (ctx->bctype == FVBC_OUTFLOW) { 53584ff8c8bSHong Zhang for (i=xs-2; i<0; i++) { 53684ff8c8bSHong Zhang for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 53784ff8c8bSHong Zhang } 53884ff8c8bSHong Zhang for (i=Mx; i<xs+xm+2; i++) { 53984ff8c8bSHong Zhang for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 54084ff8c8bSHong Zhang } 54184ff8c8bSHong Zhang } 54284ff8c8bSHong Zhang 54384ff8c8bSHong Zhang if (ctx->bctype == FVBC_INFLOW) { 54484ff8c8bSHong Zhang /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253 54584ff8c8bSHong Zhang pages 137-138 for the scheme. */ 54684ff8c8bSHong Zhang if (xs==0) { /* Left Boundary */ 54784ff8c8bSHong Zhang ctx->physics2.inflow(ctx,time,ctx->xmin,ctx->ub); 54884ff8c8bSHong Zhang for (j=0; j<dof; j++) { 54984ff8c8bSHong Zhang if (ctx->physics2.bcinflowindex[j]) { 55084ff8c8bSHong Zhang for (i=-2; i<0; i++) x[i*dof+j] = 2.0*ctx->ub[j]-x[-(i+1)*dof+j]; 55184ff8c8bSHong Zhang } else { 55284ff8c8bSHong Zhang for (i=-2; i<0; i++) x[i*dof+j] = x[j]; /* Outflow */ 55384ff8c8bSHong Zhang } 55484ff8c8bSHong Zhang } 55584ff8c8bSHong Zhang } 55684ff8c8bSHong Zhang if (xs+xm==Mx) { /* Right Boundary */ 55784ff8c8bSHong Zhang ctx->physics2.inflow(ctx,time,ctx->xmax,ctx->ub); 55884ff8c8bSHong Zhang for (j=0; j<dof; j++) { 55984ff8c8bSHong Zhang if (ctx->physics2.bcinflowindex[dof+j]) { 56084ff8c8bSHong Zhang for (i=Mx; i<Mx+2; i++) x[i*dof+j] = 2.0*ctx->ub[dof+j]-x[(2*Mx-(i+1))*dof+j]; 56184ff8c8bSHong Zhang } else { 56284ff8c8bSHong Zhang for (i=Mx; i<Mx+2; i++) x[i*dof+j] = x[(Mx-1)*dof+j]; /* Outflow */ 56384ff8c8bSHong Zhang } 56484ff8c8bSHong Zhang } 56584ff8c8bSHong Zhang } 56684ff8c8bSHong Zhang } 56784ff8c8bSHong Zhang 56884ff8c8bSHong Zhang for (i=xs-1; i<xs+xm+1; i++) { 56984ff8c8bSHong Zhang struct _LimitInfo info; 57084ff8c8bSHong Zhang PetscScalar *cjmpL,*cjmpR; 57184ff8c8bSHong Zhang /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 5725f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds)); 57384ff8c8bSHong Zhang /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 5745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof)); 57584ff8c8bSHong Zhang cjmpL = &ctx->cjmpLR[0]; 57684ff8c8bSHong Zhang cjmpR = &ctx->cjmpLR[dof]; 57784ff8c8bSHong Zhang for (j=0; j<dof; j++) { 57884ff8c8bSHong Zhang PetscScalar jmpL,jmpR; 57984ff8c8bSHong Zhang jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 58084ff8c8bSHong Zhang jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 58184ff8c8bSHong Zhang for (k=0; k<dof; k++) { 58284ff8c8bSHong Zhang cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 58384ff8c8bSHong Zhang cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 58484ff8c8bSHong Zhang } 58584ff8c8bSHong Zhang } 58684ff8c8bSHong Zhang /* Apply limiter to the left and right characteristic jumps */ 58784ff8c8bSHong Zhang info.m = dof; 58884ff8c8bSHong Zhang info.hxs = hxs; 58984ff8c8bSHong Zhang info.hxf = hxf; 59084ff8c8bSHong Zhang (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope); 59184ff8c8bSHong Zhang for (j=0; j<dof; j++) { 59284ff8c8bSHong Zhang PetscScalar tmp = 0; 59384ff8c8bSHong Zhang for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 59484ff8c8bSHong Zhang slope[i*dof+j] = tmp; 59584ff8c8bSHong Zhang } 59684ff8c8bSHong Zhang } 59784ff8c8bSHong Zhang 59884ff8c8bSHong Zhang for (i=xs; i<xs+xm+1; i++) { 59984ff8c8bSHong Zhang PetscReal maxspeed; 60084ff8c8bSHong Zhang PetscScalar *uL,*uR; 60184ff8c8bSHong Zhang uL = &ctx->uLR[0]; 60284ff8c8bSHong Zhang uR = &ctx->uLR[dof]; 60384ff8c8bSHong Zhang if (i < sf) { /* slow region */ 60484ff8c8bSHong Zhang for (j=0; j<dof; j++) { 60584ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 60684ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 60784ff8c8bSHong Zhang } 6085f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 60984ff8c8bSHong Zhang if (i > xs) { 61084ff8c8bSHong Zhang for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs; 61184ff8c8bSHong Zhang } 61284ff8c8bSHong Zhang if (i < xs+xm) { 61384ff8c8bSHong Zhang for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs; 61484ff8c8bSHong Zhang } 61584ff8c8bSHong Zhang } else if (i == sf) { /* interface between the slow region and the fast region */ 61684ff8c8bSHong Zhang for (j=0; j<dof; j++) { 61784ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 61884ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 61984ff8c8bSHong Zhang } 6205f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 62184ff8c8bSHong Zhang if (i > xs) { 62284ff8c8bSHong Zhang for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs; 62384ff8c8bSHong Zhang } 62484ff8c8bSHong Zhang if (i < xs+xm) { 62584ff8c8bSHong Zhang for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxf; 62684ff8c8bSHong Zhang } 62784ff8c8bSHong Zhang } else if (i > sf && i < fs) { /* fast region */ 62884ff8c8bSHong Zhang for (j=0; j<dof; j++) { 62984ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 63084ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 63184ff8c8bSHong Zhang } 6325f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 63384ff8c8bSHong Zhang if (i > xs) { 63484ff8c8bSHong Zhang for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf; 63584ff8c8bSHong Zhang } 63684ff8c8bSHong Zhang if (i < xs+xm) { 63784ff8c8bSHong Zhang for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxf; 63884ff8c8bSHong Zhang } 63984ff8c8bSHong Zhang } else if (i == fs) { /* interface between the fast region and the slow region */ 64084ff8c8bSHong Zhang for (j=0; j<dof; j++) { 64184ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 64284ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 64384ff8c8bSHong Zhang } 6445f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 64584ff8c8bSHong Zhang if (i > xs) { 64684ff8c8bSHong Zhang for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf; 64784ff8c8bSHong Zhang } 64884ff8c8bSHong Zhang if (i < xs+xm) { 64984ff8c8bSHong Zhang for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs; 65084ff8c8bSHong Zhang } 65184ff8c8bSHong Zhang } else { /* slow region */ 65284ff8c8bSHong Zhang for (j=0; j<dof; j++) { 65384ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 65484ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 65584ff8c8bSHong Zhang } 6565f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 65784ff8c8bSHong Zhang cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */ 65884ff8c8bSHong Zhang if (i > xs) { 65984ff8c8bSHong Zhang for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs; 66084ff8c8bSHong Zhang } 66184ff8c8bSHong Zhang if (i < xs+xm) { 66284ff8c8bSHong Zhang for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs; 66384ff8c8bSHong Zhang } 66484ff8c8bSHong Zhang } 66584ff8c8bSHong Zhang } 6665f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x)); 6675f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,F,&f)); 6685f80ce2aSJacob Faibussowitsch CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope)); 6695f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&Xloc)); 6705f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da))); 67184ff8c8bSHong Zhang if (0) { 67284ff8c8bSHong Zhang /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */ 67384ff8c8bSHong Zhang PetscReal dt,tnow; 6745f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTimeStep(ts,&dt)); 6755f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTime(ts,&tnow)); 67684ff8c8bSHong Zhang if (dt > 0.5/ctx->cfl_idt) { 67784ff8c8bSHong Zhang if (1) { 6785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(ctx->comm,"Stability constraint exceeded at t=%g, dt %g > %g\n",(double)tnow,(double)dt,(double)(0.5/ctx->cfl_idt))); 67998921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,1,"Stability constraint exceeded, %g > %g",(double)dt,(double)(ctx->cfl/ctx->cfl_idt)); 68084ff8c8bSHong Zhang } 68184ff8c8bSHong Zhang } 68284ff8c8bSHong Zhang PetscFunctionReturn(0); 68384ff8c8bSHong Zhang } 68484ff8c8bSHong Zhang 68584ff8c8bSHong Zhang /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */ 68684ff8c8bSHong Zhang PetscErrorCode FVRHSFunctionslow_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 68784ff8c8bSHong Zhang { 68884ff8c8bSHong Zhang FVCtx *ctx = (FVCtx*)vctx; 68984ff8c8bSHong Zhang PetscInt i,j,k,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth; 69084ff8c8bSHong Zhang PetscReal hxs,hxf,cfl_idt = 0; 69184ff8c8bSHong Zhang PetscScalar *x,*f,*slope; 69284ff8c8bSHong Zhang Vec Xloc; 69384ff8c8bSHong Zhang DM da; 69484ff8c8bSHong Zhang 69584ff8c8bSHong Zhang PetscFunctionBeginUser; 6965f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 6975f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&Xloc)); 6985f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 69984ff8c8bSHong Zhang hxs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 70084ff8c8bSHong Zhang hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 7015f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 7025f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 7035f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(F)); 7045f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,Xloc,&x)); 7055f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(F,&f)); 7065f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope)); 7075f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 70884ff8c8bSHong Zhang 70984ff8c8bSHong Zhang if (ctx->bctype == FVBC_OUTFLOW) { 71084ff8c8bSHong Zhang for (i=xs-2; i<0; i++) { 71184ff8c8bSHong Zhang for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 71284ff8c8bSHong Zhang } 71384ff8c8bSHong Zhang for (i=Mx; i<xs+xm+2; i++) { 71484ff8c8bSHong Zhang for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 71584ff8c8bSHong Zhang } 71684ff8c8bSHong Zhang } 71784ff8c8bSHong Zhang if (ctx->bctype == FVBC_INFLOW) { 71884ff8c8bSHong Zhang /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253 71984ff8c8bSHong Zhang pages 137-138 for the scheme. */ 72084ff8c8bSHong Zhang if (xs==0) { /* Left Boundary */ 72184ff8c8bSHong Zhang ctx->physics2.inflow(ctx,time,ctx->xmin,ctx->ub); 72284ff8c8bSHong Zhang for (j=0; j<dof; j++) { 72384ff8c8bSHong Zhang if (ctx->physics2.bcinflowindex[j]==PETSC_TRUE) { 72484ff8c8bSHong Zhang for (i=-2; i<0; i++) x[i*dof+j] = 2.0*ctx->ub[j]-x[-(i+1)*dof+j]; 72584ff8c8bSHong Zhang } else { 72684ff8c8bSHong Zhang for (i=-2; i<0; i++) x[i*dof+j] = x[j]; /* Outflow */ 72784ff8c8bSHong Zhang } 72884ff8c8bSHong Zhang } 72984ff8c8bSHong Zhang } 73084ff8c8bSHong Zhang if (xs+xm==Mx) { /* Right Boundary */ 73184ff8c8bSHong Zhang ctx->physics2.inflow(ctx,time,ctx->xmax,ctx->ub); 73284ff8c8bSHong Zhang for (j=0; j<dof; j++) { 73384ff8c8bSHong Zhang if (ctx->physics2.bcinflowindex[dof+j]==PETSC_TRUE) { 73484ff8c8bSHong Zhang for (i=Mx; i<Mx+2; i++) x[i*dof+j] = 2.0*ctx->ub[dof+j]-x[(2*Mx-(i+1))*dof+j]; 73584ff8c8bSHong Zhang } else { 73684ff8c8bSHong Zhang for (i=Mx; i<Mx+2; i++) x[i*dof+j] = x[(Mx-1)*dof+j]; /* Outflow */ 73784ff8c8bSHong Zhang } 73884ff8c8bSHong Zhang } 73984ff8c8bSHong Zhang } 74084ff8c8bSHong Zhang } 74184ff8c8bSHong Zhang for (i=xs-1; i<xs+xm+1; i++) { 74284ff8c8bSHong Zhang struct _LimitInfo info; 74384ff8c8bSHong Zhang PetscScalar *cjmpL,*cjmpR; 74484ff8c8bSHong Zhang if (i < sf-lsbwidth+1 || i > fs+rsbwidth-2) { /* slow components and the first and last fast components */ 74584ff8c8bSHong Zhang /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 7465f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds)); 74784ff8c8bSHong Zhang /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 7485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof)); 74984ff8c8bSHong Zhang cjmpL = &ctx->cjmpLR[0]; 75084ff8c8bSHong Zhang cjmpR = &ctx->cjmpLR[dof]; 75184ff8c8bSHong Zhang for (j=0; j<dof; j++) { 75284ff8c8bSHong Zhang PetscScalar jmpL,jmpR; 75384ff8c8bSHong Zhang jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 75484ff8c8bSHong Zhang jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 75584ff8c8bSHong Zhang for (k=0; k<dof; k++) { 75684ff8c8bSHong Zhang cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 75784ff8c8bSHong Zhang cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 75884ff8c8bSHong Zhang } 75984ff8c8bSHong Zhang } 76084ff8c8bSHong Zhang /* Apply limiter to the left and right characteristic jumps */ 76184ff8c8bSHong Zhang info.m = dof; 76284ff8c8bSHong Zhang info.hxs = hxs; 76384ff8c8bSHong Zhang info.hxf = hxf; 76484ff8c8bSHong Zhang (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope); 76584ff8c8bSHong Zhang for (j=0; j<dof; j++) { 76684ff8c8bSHong Zhang PetscScalar tmp = 0; 76784ff8c8bSHong Zhang for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 76884ff8c8bSHong Zhang slope[i*dof+j] = tmp; 76984ff8c8bSHong Zhang } 77084ff8c8bSHong Zhang } 77184ff8c8bSHong Zhang } 77284ff8c8bSHong Zhang 77384ff8c8bSHong Zhang for (i=xs; i<xs+xm+1; i++) { 77484ff8c8bSHong Zhang PetscReal maxspeed; 77584ff8c8bSHong Zhang PetscScalar *uL,*uR; 77684ff8c8bSHong Zhang uL = &ctx->uLR[0]; 77784ff8c8bSHong Zhang uR = &ctx->uLR[dof]; 77884ff8c8bSHong Zhang if (i < sf-lsbwidth) { /* slow region */ 77984ff8c8bSHong Zhang for (j=0; j<dof; j++) { 78084ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 78184ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 78284ff8c8bSHong Zhang } 7835f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 78484ff8c8bSHong Zhang cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */ 78584ff8c8bSHong Zhang if (i > xs) { 78684ff8c8bSHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 78784ff8c8bSHong Zhang } 78884ff8c8bSHong Zhang if (i < xs+xm) { 78984ff8c8bSHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 79084ff8c8bSHong Zhang islow++; 79184ff8c8bSHong Zhang } 79284ff8c8bSHong Zhang } 79384ff8c8bSHong Zhang if (i == sf-lsbwidth) { /* interface between the slow region and the fast region */ 79484ff8c8bSHong Zhang for (j=0; j<dof; j++) { 79584ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 79684ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 79784ff8c8bSHong Zhang } 7985f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 79984ff8c8bSHong Zhang if (i > xs) { 80084ff8c8bSHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 80184ff8c8bSHong Zhang } 80284ff8c8bSHong Zhang } 80384ff8c8bSHong Zhang if (i == fs+rsbwidth) { /* slow region */ 80484ff8c8bSHong Zhang for (j=0; j<dof; j++) { 80584ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 80684ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 80784ff8c8bSHong Zhang } 8085f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 80984ff8c8bSHong Zhang if (i < xs+xm) { 81084ff8c8bSHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 81184ff8c8bSHong Zhang islow++; 81284ff8c8bSHong Zhang } 81384ff8c8bSHong Zhang } 81484ff8c8bSHong Zhang if (i > fs+rsbwidth) { /* slow region */ 81584ff8c8bSHong Zhang for (j=0; j<dof; j++) { 81684ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 81784ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 81884ff8c8bSHong Zhang } 8195f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 82084ff8c8bSHong Zhang if (i > xs) { 82184ff8c8bSHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 82284ff8c8bSHong Zhang } 82384ff8c8bSHong Zhang if (i < xs+xm) { 82484ff8c8bSHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 82584ff8c8bSHong Zhang islow++; 82684ff8c8bSHong Zhang } 82784ff8c8bSHong Zhang } 82884ff8c8bSHong Zhang } 8295f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x)); 8305f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(F,&f)); 8315f80ce2aSJacob Faibussowitsch CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope)); 8325f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&Xloc)); 8335f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da))); 83484ff8c8bSHong Zhang PetscFunctionReturn(0); 83584ff8c8bSHong Zhang } 83684ff8c8bSHong Zhang 83784ff8c8bSHong Zhang PetscErrorCode FVRHSFunctionslowbuffer_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 83884ff8c8bSHong Zhang { 83984ff8c8bSHong Zhang FVCtx *ctx = (FVCtx*)vctx; 84084ff8c8bSHong Zhang PetscInt i,j,k,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth; 84184ff8c8bSHong Zhang PetscReal hxs,hxf; 84284ff8c8bSHong Zhang PetscScalar *x,*f,*slope; 84384ff8c8bSHong Zhang Vec Xloc; 84484ff8c8bSHong Zhang DM da; 84584ff8c8bSHong Zhang 84684ff8c8bSHong Zhang PetscFunctionBeginUser; 8475f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 8485f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&Xloc)); 8495f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 85084ff8c8bSHong Zhang hxs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 85184ff8c8bSHong Zhang hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 8525f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 8535f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 8545f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(F)); 8555f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,Xloc,&x)); 8565f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(F,&f)); 8575f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope)); 8585f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 85984ff8c8bSHong Zhang 86084ff8c8bSHong Zhang if (ctx->bctype == FVBC_OUTFLOW) { 86184ff8c8bSHong Zhang for (i=xs-2; i<0; i++) { 86284ff8c8bSHong Zhang for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 86384ff8c8bSHong Zhang } 86484ff8c8bSHong Zhang for (i=Mx; i<xs+xm+2; i++) { 86584ff8c8bSHong Zhang for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 86684ff8c8bSHong Zhang } 86784ff8c8bSHong Zhang } 86884ff8c8bSHong Zhang if (ctx->bctype == FVBC_INFLOW) { 86984ff8c8bSHong Zhang /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253 87084ff8c8bSHong Zhang pages 137-138 for the scheme. */ 87184ff8c8bSHong Zhang if (xs==0) { /* Left Boundary */ 87284ff8c8bSHong Zhang ctx->physics2.inflow(ctx,time,ctx->xmin,ctx->ub); 87384ff8c8bSHong Zhang for (j=0; j<dof; j++) { 87484ff8c8bSHong Zhang if (ctx->physics2.bcinflowindex[j]==PETSC_TRUE) { 87584ff8c8bSHong Zhang for (i=-2; i<0; i++) x[i*dof+j] = 2.0*ctx->ub[j]-x[-(i+1)*dof+j]; 87684ff8c8bSHong Zhang } else { 87784ff8c8bSHong Zhang for (i=-2; i<0; i++) x[i*dof+j] = x[j]; /* Outflow */ 87884ff8c8bSHong Zhang } 87984ff8c8bSHong Zhang } 88084ff8c8bSHong Zhang } 88184ff8c8bSHong Zhang if (xs+xm==Mx) { /* Right Boundary */ 88284ff8c8bSHong Zhang ctx->physics2.inflow(ctx,time,ctx->xmax,ctx->ub); 88384ff8c8bSHong Zhang for (j=0; j<dof; j++) { 88484ff8c8bSHong Zhang if (ctx->physics2.bcinflowindex[dof+j]==PETSC_TRUE) { 88584ff8c8bSHong Zhang for (i=Mx; i<Mx+2; i++) x[i*dof+j] = 2.0*ctx->ub[dof+j]-x[(2*Mx-(i+1))*dof+j]; 88684ff8c8bSHong Zhang } else { 88784ff8c8bSHong Zhang for (i=Mx; i<Mx+2; i++) x[i*dof+j] = x[(Mx-1)*dof+j]; /* Outflow */ 88884ff8c8bSHong Zhang } 88984ff8c8bSHong Zhang } 89084ff8c8bSHong Zhang } 89184ff8c8bSHong Zhang } 89284ff8c8bSHong Zhang for (i=xs-1; i<xs+xm+1; i++) { 89384ff8c8bSHong Zhang struct _LimitInfo info; 89484ff8c8bSHong Zhang PetscScalar *cjmpL,*cjmpR; 89584ff8c8bSHong Zhang if ((i > sf-lsbwidth-2 && i < sf+1) || (i > fs-2 && i < fs+rsbwidth+1)) { 89684ff8c8bSHong Zhang /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 8975f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds)); 89884ff8c8bSHong Zhang /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 8995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof)); 90084ff8c8bSHong Zhang cjmpL = &ctx->cjmpLR[0]; 90184ff8c8bSHong Zhang cjmpR = &ctx->cjmpLR[dof]; 90284ff8c8bSHong Zhang for (j=0; j<dof; j++) { 90384ff8c8bSHong Zhang PetscScalar jmpL,jmpR; 90484ff8c8bSHong Zhang jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 90584ff8c8bSHong Zhang jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 90684ff8c8bSHong Zhang for (k=0; k<dof; k++) { 90784ff8c8bSHong Zhang cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 90884ff8c8bSHong Zhang cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 90984ff8c8bSHong Zhang } 91084ff8c8bSHong Zhang } 91184ff8c8bSHong Zhang /* Apply limiter to the left and right characteristic jumps */ 91284ff8c8bSHong Zhang info.m = dof; 91384ff8c8bSHong Zhang info.hxs = hxs; 91484ff8c8bSHong Zhang info.hxf = hxf; 91584ff8c8bSHong Zhang (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope); 91684ff8c8bSHong Zhang for (j=0; j<dof; j++) { 91784ff8c8bSHong Zhang PetscScalar tmp = 0; 91884ff8c8bSHong Zhang for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 91984ff8c8bSHong Zhang slope[i*dof+j] = tmp; 92084ff8c8bSHong Zhang } 92184ff8c8bSHong Zhang } 92284ff8c8bSHong Zhang } 92384ff8c8bSHong Zhang 92484ff8c8bSHong Zhang for (i=xs; i<xs+xm+1; i++) { 92584ff8c8bSHong Zhang PetscReal maxspeed; 92684ff8c8bSHong Zhang PetscScalar *uL,*uR; 92784ff8c8bSHong Zhang uL = &ctx->uLR[0]; 92884ff8c8bSHong Zhang uR = &ctx->uLR[dof]; 92984ff8c8bSHong Zhang if (i == sf-lsbwidth) { 93084ff8c8bSHong Zhang for (j=0; j<dof; j++) { 93184ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 93284ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 93384ff8c8bSHong Zhang } 9345f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 93584ff8c8bSHong Zhang if (i < xs+xm) { 93684ff8c8bSHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 93784ff8c8bSHong Zhang islow++; 93884ff8c8bSHong Zhang } 93984ff8c8bSHong Zhang } 94084ff8c8bSHong Zhang if (i > sf-lsbwidth && i < sf) { 94184ff8c8bSHong Zhang for (j=0; j<dof; j++) { 94284ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 94384ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 94484ff8c8bSHong Zhang } 9455f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 94684ff8c8bSHong Zhang if (i > xs) { 94784ff8c8bSHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 94884ff8c8bSHong Zhang } 94984ff8c8bSHong Zhang if (i < xs+xm) { 95084ff8c8bSHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 95184ff8c8bSHong Zhang islow++; 95284ff8c8bSHong Zhang } 95384ff8c8bSHong Zhang } 95484ff8c8bSHong Zhang if (i == sf) { /* interface between the slow region and the fast region */ 95584ff8c8bSHong Zhang for (j=0; j<dof; j++) { 95684ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 95784ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 95884ff8c8bSHong Zhang } 9595f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 96084ff8c8bSHong Zhang if (i > xs) { 96184ff8c8bSHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 96284ff8c8bSHong Zhang } 96384ff8c8bSHong Zhang } 96484ff8c8bSHong Zhang if (i == fs) { /* interface between the fast region and the slow region */ 96584ff8c8bSHong Zhang for (j=0; j<dof; j++) { 96684ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 96784ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 96884ff8c8bSHong Zhang } 9695f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 97084ff8c8bSHong Zhang if (i < xs+xm) { 97184ff8c8bSHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 97284ff8c8bSHong Zhang islow++; 97384ff8c8bSHong Zhang } 97484ff8c8bSHong Zhang } 97584ff8c8bSHong Zhang if (i > fs && i < fs+rsbwidth) { 97684ff8c8bSHong Zhang for (j=0; j<dof; j++) { 97784ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 97884ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 97984ff8c8bSHong Zhang } 9805f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 98184ff8c8bSHong Zhang if (i > xs) { 98284ff8c8bSHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 98384ff8c8bSHong Zhang } 98484ff8c8bSHong Zhang if (i < xs+xm) { 98584ff8c8bSHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 98684ff8c8bSHong Zhang islow++; 98784ff8c8bSHong Zhang } 98884ff8c8bSHong Zhang } 98984ff8c8bSHong Zhang if (i == fs+rsbwidth) { 99084ff8c8bSHong Zhang for (j=0; j<dof; j++) { 99184ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 99284ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 99384ff8c8bSHong Zhang } 9945f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 99584ff8c8bSHong Zhang if (i > xs) { 99684ff8c8bSHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 99784ff8c8bSHong Zhang } 99884ff8c8bSHong Zhang } 99984ff8c8bSHong Zhang } 10005f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x)); 10015f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(F,&f)); 10025f80ce2aSJacob Faibussowitsch CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope)); 10035f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&Xloc)); 100484ff8c8bSHong Zhang PetscFunctionReturn(0); 100584ff8c8bSHong Zhang } 100684ff8c8bSHong Zhang 100784ff8c8bSHong Zhang /* --------------------------------- Finite Volume Solver for fast parts ----------------------------------- */ 100884ff8c8bSHong Zhang PetscErrorCode FVRHSFunctionfast_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 100984ff8c8bSHong Zhang { 101084ff8c8bSHong Zhang FVCtx *ctx = (FVCtx*)vctx; 101184ff8c8bSHong Zhang PetscInt i,j,k,Mx,dof,xs,xm,ifast = 0,sf = ctx->sf,fs = ctx->fs; 101284ff8c8bSHong Zhang PetscReal hxs,hxf; 101384ff8c8bSHong Zhang PetscScalar *x,*f,*slope; 101484ff8c8bSHong Zhang Vec Xloc; 101584ff8c8bSHong Zhang DM da; 101684ff8c8bSHong Zhang 101784ff8c8bSHong Zhang PetscFunctionBeginUser; 10185f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 10195f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&Xloc)); 10205f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 102184ff8c8bSHong Zhang hxs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 102284ff8c8bSHong Zhang hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 10235f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 10245f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 10255f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(F)); 10265f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,Xloc,&x)); 10275f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(F,&f)); 10285f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetArray(da,PETSC_TRUE,&slope)); 10295f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 103084ff8c8bSHong Zhang 103184ff8c8bSHong Zhang if (ctx->bctype == FVBC_OUTFLOW) { 103284ff8c8bSHong Zhang for (i=xs-2; i<0; i++) { 103384ff8c8bSHong Zhang for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 103484ff8c8bSHong Zhang } 103584ff8c8bSHong Zhang for (i=Mx; i<xs+xm+2; i++) { 103684ff8c8bSHong Zhang for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 103784ff8c8bSHong Zhang } 103884ff8c8bSHong Zhang } 103984ff8c8bSHong Zhang if (ctx->bctype == FVBC_INFLOW) { 104084ff8c8bSHong Zhang /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253 104184ff8c8bSHong Zhang pages 137-138 for the scheme. */ 104284ff8c8bSHong Zhang if (xs==0) { /* Left Boundary */ 104384ff8c8bSHong Zhang ctx->physics2.inflow(ctx,time,ctx->xmin,ctx->ub); 104484ff8c8bSHong Zhang for (j=0; j<dof; j++) { 104584ff8c8bSHong Zhang if (ctx->physics2.bcinflowindex[j]==PETSC_TRUE) { 104684ff8c8bSHong Zhang for (i=-2; i<0; i++) x[i*dof+j] = 2.0*ctx->ub[j]-x[-(i+1)*dof+j]; 104784ff8c8bSHong Zhang } else { 104884ff8c8bSHong Zhang for (i=-2; i<0; i++) x[i*dof+j] = x[j]; /* Outflow */ 104984ff8c8bSHong Zhang } 105084ff8c8bSHong Zhang } 105184ff8c8bSHong Zhang } 105284ff8c8bSHong Zhang if (xs+xm==Mx) { /* Right Boundary */ 105384ff8c8bSHong Zhang ctx->physics2.inflow(ctx,time,ctx->xmax,ctx->ub); 105484ff8c8bSHong Zhang for (j=0; j<dof; j++) { 105584ff8c8bSHong Zhang if (ctx->physics2.bcinflowindex[dof+j]==PETSC_TRUE) { 105684ff8c8bSHong Zhang for (i=Mx; i<Mx+2; i++) x[i*dof+j] = 2.0*ctx->ub[dof+j]-x[(2*Mx-(i+1))*dof+j]; 105784ff8c8bSHong Zhang } else { 105884ff8c8bSHong Zhang for (i=Mx; i<Mx+2; i++) x[i*dof+j] = x[(Mx-1)*dof+j]; /* Outflow */ 105984ff8c8bSHong Zhang } 106084ff8c8bSHong Zhang } 106184ff8c8bSHong Zhang } 106284ff8c8bSHong Zhang } 106384ff8c8bSHong Zhang for (i=xs-1; i<xs+xm+1; i++) { 106484ff8c8bSHong Zhang struct _LimitInfo info; 106584ff8c8bSHong Zhang PetscScalar *cjmpL,*cjmpR; 106684ff8c8bSHong Zhang if (i > sf-2 && i < fs+1) { 10675f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds)); 10685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(ctx->cjmpLR,2*dof)); 106984ff8c8bSHong Zhang cjmpL = &ctx->cjmpLR[0]; 107084ff8c8bSHong Zhang cjmpR = &ctx->cjmpLR[dof]; 107184ff8c8bSHong Zhang for (j=0; j<dof; j++) { 107284ff8c8bSHong Zhang PetscScalar jmpL,jmpR; 107384ff8c8bSHong Zhang jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 107484ff8c8bSHong Zhang jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 107584ff8c8bSHong Zhang for (k=0; k<dof; k++) { 107684ff8c8bSHong Zhang cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 107784ff8c8bSHong Zhang cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 107884ff8c8bSHong Zhang } 107984ff8c8bSHong Zhang } 108084ff8c8bSHong Zhang /* Apply limiter to the left and right characteristic jumps */ 108184ff8c8bSHong Zhang info.m = dof; 108284ff8c8bSHong Zhang info.hxs = hxs; 108384ff8c8bSHong Zhang info.hxf = hxf; 108484ff8c8bSHong Zhang (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope); 108584ff8c8bSHong Zhang for (j=0; j<dof; j++) { 108684ff8c8bSHong Zhang PetscScalar tmp = 0; 108784ff8c8bSHong Zhang for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 108884ff8c8bSHong Zhang slope[i*dof+j] = tmp; 108984ff8c8bSHong Zhang } 109084ff8c8bSHong Zhang } 109184ff8c8bSHong Zhang } 109284ff8c8bSHong Zhang 109384ff8c8bSHong Zhang for (i=xs; i<xs+xm+1; i++) { 109484ff8c8bSHong Zhang PetscReal maxspeed; 109584ff8c8bSHong Zhang PetscScalar *uL,*uR; 109684ff8c8bSHong Zhang uL = &ctx->uLR[0]; 109784ff8c8bSHong Zhang uR = &ctx->uLR[dof]; 109884ff8c8bSHong Zhang if (i == sf) { /* interface between the slow region and the fast region */ 109984ff8c8bSHong Zhang for (j=0; j<dof; j++) { 110084ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 110184ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 110284ff8c8bSHong Zhang } 11035f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 110484ff8c8bSHong Zhang if (i < xs+xm) { 110584ff8c8bSHong Zhang for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf; 110684ff8c8bSHong Zhang ifast++; 110784ff8c8bSHong Zhang } 110884ff8c8bSHong Zhang } 110984ff8c8bSHong Zhang if (i > sf && i < fs) { /* fast region */ 111084ff8c8bSHong Zhang for (j=0; j<dof; j++) { 111184ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 111284ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 111384ff8c8bSHong Zhang } 11145f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 111584ff8c8bSHong Zhang if (i > xs) { 111684ff8c8bSHong Zhang for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf; 111784ff8c8bSHong Zhang } 111884ff8c8bSHong Zhang if (i < xs+xm) { 111984ff8c8bSHong Zhang for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf; 112084ff8c8bSHong Zhang ifast++; 112184ff8c8bSHong Zhang } 112284ff8c8bSHong Zhang } 112384ff8c8bSHong Zhang if (i == fs) { /* interface between the fast region and the slow region */ 112484ff8c8bSHong Zhang for (j=0; j<dof; j++) { 112584ff8c8bSHong Zhang uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 112684ff8c8bSHong Zhang uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 112784ff8c8bSHong Zhang } 11285f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 112984ff8c8bSHong Zhang if (i > xs) { 113084ff8c8bSHong Zhang for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf; 113184ff8c8bSHong Zhang } 113284ff8c8bSHong Zhang } 113384ff8c8bSHong Zhang } 11345f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,Xloc,&x)); 11355f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(F,&f)); 11365f80ce2aSJacob Faibussowitsch CHKERRQ(DMDARestoreArray(da,PETSC_TRUE,&slope)); 11375f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&Xloc)); 113884ff8c8bSHong Zhang PetscFunctionReturn(0); 113984ff8c8bSHong Zhang } 114084ff8c8bSHong Zhang 114184ff8c8bSHong Zhang int main(int argc,char *argv[]) 114284ff8c8bSHong Zhang { 114384ff8c8bSHong Zhang char lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m"; 114484ff8c8bSHong Zhang PetscFunctionList limiters = 0,physics = 0; 114584ff8c8bSHong Zhang MPI_Comm comm; 114684ff8c8bSHong Zhang TS ts; 114784ff8c8bSHong Zhang DM da; 114884ff8c8bSHong Zhang Vec X,X0,R; 114984ff8c8bSHong Zhang FVCtx ctx; 115084ff8c8bSHong Zhang 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; 115184ff8c8bSHong Zhang PetscBool view_final = PETSC_FALSE; 115284ff8c8bSHong Zhang PetscReal ptime,maxtime; 115384ff8c8bSHong Zhang PetscErrorCode ierr; 115484ff8c8bSHong Zhang 1155*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,0,help)); 115684ff8c8bSHong Zhang comm = PETSC_COMM_WORLD; 11575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemzero(&ctx,sizeof(ctx))); 115884ff8c8bSHong Zhang 115984ff8c8bSHong Zhang /* Register limiters to be available on the command line */ 11605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListAdd(&limiters,"upwind" ,Limit2_Upwind)); 11615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListAdd(&limiters,"lax-wendroff" ,Limit2_LaxWendroff)); 11625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListAdd(&limiters,"beam-warming" ,Limit2_BeamWarming)); 11635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListAdd(&limiters,"fromm" ,Limit2_Fromm)); 11645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListAdd(&limiters,"minmod" ,Limit2_Minmod)); 11655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListAdd(&limiters,"superbee" ,Limit2_Superbee)); 11665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListAdd(&limiters,"mc" ,Limit2_MC)); 11675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListAdd(&limiters,"koren3" ,Limit2_Koren3)); 116884ff8c8bSHong Zhang 116984ff8c8bSHong Zhang /* Register physical models to be available on the command line */ 11705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListAdd(&physics,"shallow" ,PhysicsCreate_Shallow)); 11715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListAdd(&physics,"advect" ,PhysicsCreate_Advect)); 117284ff8c8bSHong Zhang 117384ff8c8bSHong Zhang ctx.comm = comm; 117484ff8c8bSHong Zhang ctx.cfl = 0.9; 117584ff8c8bSHong Zhang ctx.bctype = FVBC_PERIODIC; 117684ff8c8bSHong Zhang ctx.xmin = -1.0; 117784ff8c8bSHong Zhang ctx.xmax = 1.0; 117884ff8c8bSHong Zhang ctx.initial = 1; 117984ff8c8bSHong Zhang ctx.hratio = 2; 118084ff8c8bSHong Zhang maxtime = 10.0; 118184ff8c8bSHong Zhang ctx.simulation = PETSC_FALSE; 118284ff8c8bSHong Zhang ierr = PetscOptionsBegin(comm,NULL,"Finite Volume solver options","");CHKERRQ(ierr); 11835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL)); 11845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL)); 11855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsFList("-limit","Name of flux imiter to use","",limiters,lname,lname,sizeof(lname),NULL)); 11865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsFList("-physics","Name of physics model to use","",physics,physname,physname,sizeof(physname),NULL)); 11875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL)); 11885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsString("-view_final","Write final solution in ASCII MATLAB format to given file name","",final_fname,final_fname,sizeof(final_fname),&view_final)); 11895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL)); 11905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-exact","Compare errors with exact solution","",ctx.exact,&ctx.exact,NULL)); 11915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL)); 11925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL)); 11935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL)); 11945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-hratio","Spacing ratio","",ctx.hratio,&ctx.hratio,NULL)); 119584ff8c8bSHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 119684ff8c8bSHong Zhang 119784ff8c8bSHong Zhang /* Choose the limiter from the list of registered limiters */ 11985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListFind(limiters,lname,&ctx.limit2)); 11993c633725SBarry Smith PetscCheck(ctx.limit2,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Limiter '%s' not found",lname); 120084ff8c8bSHong Zhang 120184ff8c8bSHong Zhang /* Choose the physics from the list of registered models */ 120284ff8c8bSHong Zhang { 120384ff8c8bSHong Zhang PetscErrorCode (*r)(FVCtx*); 12045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListFind(physics,physname,&r)); 12053c633725SBarry Smith PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Physics '%s' not found",physname); 120684ff8c8bSHong Zhang /* Create the physics, will set the number of fields and their names */ 12075f80ce2aSJacob Faibussowitsch CHKERRQ((*r)(&ctx)); 120884ff8c8bSHong Zhang } 120984ff8c8bSHong Zhang 121084ff8c8bSHong Zhang /* Create a DMDA to manage the parallel grid */ 12115f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics2.dof,2,NULL,&da)); 12125f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 12135f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 121484ff8c8bSHong Zhang /* Inform the DMDA of the field names provided by the physics. */ 121584ff8c8bSHong Zhang /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ 121684ff8c8bSHong Zhang for (i=0; i<ctx.physics2.dof; i++) { 12175f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetFieldName(da,i,ctx.physics2.fieldname[i])); 121884ff8c8bSHong Zhang } 12195f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 12205f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 122184ff8c8bSHong Zhang 122284ff8c8bSHong Zhang /* Set coordinates of cell centers */ 12235f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(da,ctx.xmin+0.5*(ctx.xmax-ctx.xmin)/Mx,ctx.xmax+0.5*(ctx.xmax-ctx.xmin)/Mx,0,0,0,0)); 122484ff8c8bSHong Zhang 122584ff8c8bSHong Zhang /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */ 12265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc4(dof*dof,&ctx.R,dof*dof,&ctx.Rinv,2*dof,&ctx.cjmpLR,1*dof,&ctx.cslope)); 12275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc3(2*dof,&ctx.uLR,dof,&ctx.flux,dof,&ctx.speeds)); 12285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(2*dof,&ctx.ub)); 122984ff8c8bSHong Zhang 123084ff8c8bSHong Zhang /* Create a vector to store the solution and to save the initial state */ 12315f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&X)); 12325f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(X,&X0)); 12335f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(X,&R)); 123484ff8c8bSHong Zhang 123584ff8c8bSHong Zhang /* create index for slow parts and fast parts, 123684ff8c8bSHong Zhang count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */ 123784ff8c8bSHong Zhang count_slow = Mx/(1.0+ctx.hratio/3.0); 12382c71b3e2SJacob Faibussowitsch PetscCheckFalse(count_slow%2,PETSC_COMM_WORLD,PETSC_ERR_USER,"Please adjust grid size Mx (-da_grid_x) and hratio (-hratio) so that Mx/(1+hratio/3) is even"); 123984ff8c8bSHong Zhang count_fast = Mx-count_slow; 124084ff8c8bSHong Zhang ctx.sf = count_slow/2; 124184ff8c8bSHong Zhang ctx.fs = ctx.sf+count_fast; 12425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(xm*dof,&index_slow)); 12435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(xm*dof,&index_fast)); 12445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(8*dof,&index_slowbuffer)); 124584ff8c8bSHong Zhang ctx.lsbwidth = 4; 124684ff8c8bSHong Zhang ctx.rsbwidth = 4; 124784ff8c8bSHong Zhang 124884ff8c8bSHong Zhang for (i=xs; i<xs+xm; i++) { 124984ff8c8bSHong Zhang if (i < ctx.sf-ctx.lsbwidth || i > ctx.fs+ctx.rsbwidth-1) 125084ff8c8bSHong Zhang for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k; 125184ff8c8bSHong Zhang else if ((i >= ctx.sf-ctx.lsbwidth && i < ctx.sf) || (i > ctx.fs-1 && i <= ctx.fs+ctx.rsbwidth-1)) 125284ff8c8bSHong Zhang for (k=0; k<dof; k++) index_slowbuffer[islowbuffer++] = i*dof+k; 125384ff8c8bSHong Zhang else 125484ff8c8bSHong Zhang for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k; 125584ff8c8bSHong Zhang } 12565f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss)); 12575f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf)); 12585f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,islowbuffer,index_slowbuffer,PETSC_COPY_VALUES,&ctx.issb)); 125984ff8c8bSHong Zhang 126084ff8c8bSHong Zhang /* Create a time-stepping object */ 12615f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(comm,&ts)); 12625f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts,da)); 12635f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,R,FVRHSFunction_2WaySplit,&ctx)); 12645f80ce2aSJacob Faibussowitsch CHKERRQ(TSRHSSplitSetIS(ts,"slow",ctx.iss)); 12655f80ce2aSJacob Faibussowitsch CHKERRQ(TSRHSSplitSetIS(ts,"slowbuffer",ctx.issb)); 12665f80ce2aSJacob Faibussowitsch CHKERRQ(TSRHSSplitSetIS(ts,"fast",ctx.isf)); 12675f80ce2aSJacob Faibussowitsch CHKERRQ(TSRHSSplitSetRHSFunction(ts,"slow",NULL,FVRHSFunctionslow_2WaySplit,&ctx)); 12685f80ce2aSJacob Faibussowitsch CHKERRQ(TSRHSSplitSetRHSFunction(ts,"fast",NULL,FVRHSFunctionfast_2WaySplit,&ctx)); 12695f80ce2aSJacob Faibussowitsch CHKERRQ(TSRHSSplitSetRHSFunction(ts,"slowbuffer",NULL,FVRHSFunctionslowbuffer_2WaySplit,&ctx)); 127084ff8c8bSHong Zhang 12715f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSMPRK)); 12725f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,maxtime)); 12735f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 127484ff8c8bSHong Zhang 127584ff8c8bSHong Zhang /* Compute initial conditions and starting time step */ 12765f80ce2aSJacob Faibussowitsch CHKERRQ(FVSample_2WaySplit(&ctx,da,0,X0)); 12775f80ce2aSJacob Faibussowitsch CHKERRQ(FVRHSFunction_2WaySplit(ts,0,X0,X,(void*)&ctx)); /* Initial function evaluation, only used to determine max speed */ 12785f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(X0,X)); /* The function value was not used so we set X=X0 again */ 12795f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt)); 12805f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); /* Take runtime options */ 12815f80ce2aSJacob Faibussowitsch CHKERRQ(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD)); 128284ff8c8bSHong Zhang { 128384ff8c8bSHong Zhang PetscInt steps; 128484ff8c8bSHong Zhang PetscScalar mass_initial,mass_final,mass_difference,mass_differenceg; 128584ff8c8bSHong Zhang const PetscScalar *ptr_X,*ptr_X0; 128684ff8c8bSHong Zhang const PetscReal hs = (ctx.xmax-ctx.xmin)*3.0/4.0/count_slow; 128784ff8c8bSHong Zhang const PetscReal hf = (ctx.xmax-ctx.xmin)/4.0/count_fast; 128884ff8c8bSHong Zhang 12895f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,X)); 12905f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolveTime(ts,&ptime)); 12915f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetStepNumber(ts,&steps)); 129284ff8c8bSHong Zhang /* calculate the total mass at initial time and final time */ 129384ff8c8bSHong Zhang mass_initial = 0.0; 129484ff8c8bSHong Zhang mass_final = 0.0; 12955f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,X0,(void*)&ptr_X0)); 12965f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,X,(void*)&ptr_X)); 129784ff8c8bSHong Zhang for (i=xs;i<xs+xm;i++) { 129884ff8c8bSHong Zhang if (i < ctx.sf || i > ctx.fs-1) { 129984ff8c8bSHong Zhang for (k=0; k<dof; k++) { 130084ff8c8bSHong Zhang mass_initial = mass_initial + hs*ptr_X0[i*dof+k]; 130184ff8c8bSHong Zhang mass_final = mass_final + hs*ptr_X[i*dof+k]; 130284ff8c8bSHong Zhang } 130384ff8c8bSHong Zhang } else { 130484ff8c8bSHong Zhang for (k=0; k<dof; k++) { 130584ff8c8bSHong Zhang mass_initial = mass_initial + hf*ptr_X0[i*dof+k]; 130684ff8c8bSHong Zhang mass_final = mass_final + hf*ptr_X[i*dof+k]; 130784ff8c8bSHong Zhang } 130884ff8c8bSHong Zhang } 130984ff8c8bSHong Zhang } 13105f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,X0,(void*)&ptr_X0)); 13115f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,X,(void*)&ptr_X)); 131284ff8c8bSHong Zhang mass_difference = mass_final - mass_initial; 13135f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(&mass_difference,&mass_differenceg,1,MPIU_SCALAR,MPIU_SUM,comm)); 13145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(comm,"Mass difference %g\n",(double)mass_differenceg)); 13155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(comm,"Final time %g, steps %D\n",(double)ptime,steps)); 13165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(comm,"Maximum allowable stepsize according to CFL %g\n",(double)(1.0/ctx.cfl_idt))); 131784ff8c8bSHong Zhang if (ctx.exact) { 131884ff8c8bSHong Zhang PetscReal nrm1=0; 13195f80ce2aSJacob Faibussowitsch CHKERRQ(SolutionErrorNorms_2WaySplit(&ctx,da,ptime,X,&nrm1)); 13205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1)); 132184ff8c8bSHong Zhang } 132284ff8c8bSHong Zhang if (ctx.simulation) { 132384ff8c8bSHong Zhang PetscReal nrm1=0; 132484ff8c8bSHong Zhang PetscViewer fd; 132584ff8c8bSHong Zhang char filename[PETSC_MAX_PATH_LEN] = "binaryoutput"; 132684ff8c8bSHong Zhang Vec XR; 132784ff8c8bSHong Zhang PetscBool flg; 132884ff8c8bSHong Zhang const PetscScalar *ptr_XR; 13295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg)); 13303c633725SBarry Smith PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 13315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&fd)); 13325f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(X0,&XR)); 13335f80ce2aSJacob Faibussowitsch CHKERRQ(VecLoad(XR,fd)); 13345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&fd)); 13355f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X,&ptr_X)); 13365f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(XR,&ptr_XR)); 133784ff8c8bSHong Zhang for (i=xs;i<xs+xm;i++) { 133884ff8c8bSHong Zhang if (i < ctx.sf || i > ctx.fs-1) 133984ff8c8bSHong Zhang for (k=0; k<dof; k++) nrm1 = nrm1 + hs*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]); 134084ff8c8bSHong Zhang else 134184ff8c8bSHong Zhang for (k=0; k<dof; k++) nrm1 = nrm1 + hf*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]); 134284ff8c8bSHong Zhang } 13435f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X,&ptr_X)); 13445f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(XR,&ptr_XR)); 13455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1)); 13465f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&XR)); 134784ff8c8bSHong Zhang } 134884ff8c8bSHong Zhang } 134984ff8c8bSHong Zhang 13505f80ce2aSJacob Faibussowitsch CHKERRQ(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD)); 13515f80ce2aSJacob Faibussowitsch if (draw & 0x1) CHKERRQ(VecView(X0,PETSC_VIEWER_DRAW_WORLD)); 13525f80ce2aSJacob Faibussowitsch if (draw & 0x2) CHKERRQ(VecView(X,PETSC_VIEWER_DRAW_WORLD)); 135384ff8c8bSHong Zhang if (draw & 0x4) { 135484ff8c8bSHong Zhang Vec Y; 13555f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(X,&Y)); 13565f80ce2aSJacob Faibussowitsch CHKERRQ(FVSample_2WaySplit(&ctx,da,ptime,Y)); 13575f80ce2aSJacob Faibussowitsch CHKERRQ(VecAYPX(Y,-1,X)); 13585f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(Y,PETSC_VIEWER_DRAW_WORLD)); 13595f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Y)); 136084ff8c8bSHong Zhang } 136184ff8c8bSHong Zhang 136284ff8c8bSHong Zhang if (view_final) { 136384ff8c8bSHong Zhang PetscViewer viewer; 13645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer)); 13655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB)); 13665f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(X,viewer)); 13675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(viewer)); 13685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 136984ff8c8bSHong Zhang } 137084ff8c8bSHong Zhang 137184ff8c8bSHong Zhang /* Clean up */ 13725f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx.physics2.destroy)(ctx.physics2.user)); 13735f80ce2aSJacob Faibussowitsch for (i=0; i<ctx.physics2.dof; i++) CHKERRQ(PetscFree(ctx.physics2.fieldname[i])); 13745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ctx.physics2.bcinflowindex)); 137584ff8c8bSHong Zhang ierr = PetscFree(ctx.ub); 13765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree4(ctx.R,ctx.Rinv,ctx.cjmpLR,ctx.cslope)); 13775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree3(ctx.uLR,ctx.flux,ctx.speeds)); 13785f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X)); 13795f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X0)); 13805f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&R)); 13815f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 13825f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 13835f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&ctx.iss)); 13845f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&ctx.isf)); 13855f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&ctx.issb)); 13865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(index_slow)); 13875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(index_fast)); 13885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(index_slowbuffer)); 13895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListDestroy(&limiters)); 13905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFunctionListDestroy(&physics)); 1391*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 1392*b122ec5aSJacob Faibussowitsch return 0; 139384ff8c8bSHong Zhang } 139484ff8c8bSHong Zhang 139584ff8c8bSHong Zhang /*TEST 139684ff8c8bSHong Zhang 139784ff8c8bSHong Zhang build: 139884ff8c8bSHong Zhang requires: !complex !single 139984ff8c8bSHong Zhang depends: finitevolume1d.c 140084ff8c8bSHong Zhang 140184ff8c8bSHong Zhang test: 140284ff8c8bSHong Zhang suffix: 1 140384ff8c8bSHong Zhang args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 140484ff8c8bSHong Zhang output_file: output/ex4_1.out 140584ff8c8bSHong Zhang 140684ff8c8bSHong Zhang test: 140784ff8c8bSHong Zhang suffix: 2 140884ff8c8bSHong Zhang args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 1 140984ff8c8bSHong Zhang output_file: output/ex4_1.out 141084ff8c8bSHong Zhang 141184ff8c8bSHong Zhang test: 141284ff8c8bSHong Zhang suffix: 3 141384ff8c8bSHong Zhang args: -da_grid_x 40 -initial 1 -hratio 2 -limit mc -ts_dt 0.1 -ts_max_steps 24 -ts_max_time 7.0 -ts_type mprk -ts_mprk_type 2a22 -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0 141484ff8c8bSHong Zhang output_file: output/ex4_3.out 141584ff8c8bSHong Zhang 141684ff8c8bSHong Zhang test: 141784ff8c8bSHong Zhang suffix: 4 141884ff8c8bSHong Zhang nsize: 2 141984ff8c8bSHong Zhang args: -da_grid_x 40 -initial 1 -hratio 2 -limit mc -ts_dt 0.1 -ts_max_steps 24 -ts_max_time 7.0 -ts_type mprk -ts_mprk_type 2a22 -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 1 142084ff8c8bSHong Zhang output_file: output/ex4_3.out 142184ff8c8bSHong Zhang 142284ff8c8bSHong Zhang test: 142384ff8c8bSHong Zhang suffix: 5 142484ff8c8bSHong Zhang nsize: 4 142584ff8c8bSHong Zhang args: args: -da_grid_x 40 -initial 1 -hratio 2 -limit mc -ts_dt 0.1 -ts_max_steps 24 -ts_max_time 7.0 -ts_type mprk -ts_mprk_type 2a22 -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 1 142684ff8c8bSHong Zhang output_file: output/ex4_3.out 142784ff8c8bSHong Zhang TEST*/ 1428