static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n" " advection - Constant coefficient scalar advection\n" " u_t + (a*u)_x = 0\n" " for this toy problem, we choose different meshsizes for different sub-domains (slow-medium-fast-medium-slow), \n" " the meshsize ratio between two adjacient sub-domains is controlled with -hratio,\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" "Several initial conditions can be chosen with -initial N\n\n" "The problem size should be set with -da_grid_x M\n\n"; #include #include #include #include #include "finitevolume1d.h" PETSC_STATIC_INLINE PetscReal RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax) { PetscReal range = xmax-xmin; return xmin +PetscFmodReal(range+PetscFmodReal(a,range),range); } /* --------------------------------- 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); } PetscErrorCode FVSample_3WaySplit(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,hm,hf; PetscFunctionBeginUser; if (!ctx->physics2.sample2) SETERRQ(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)/8.0/ctx->sm; hm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); hf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); for (i=xs; ism) { 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; kmf) { xi = ctx->xmin+ctx->sm*hs+0.5*hm+(i-ctx->sm)*hm; /* 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; kfm) { xi = ctx->xmin+ctx->sm*hs+(ctx->mf-ctx->sm)*hm+0.5*hf+(i-ctx->mf)*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; kms) { xi = ctx->xmin+ctx->sm*hs+(ctx->mf-ctx->sm)*hm+(ctx->fm-ctx->mf)*hf+0.5*hm+(i-ctx->fm)*hm; /* 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->sm*hs+(ctx->mf-ctx->sm)*hm+(ctx->fm-ctx->mf)*hf+(ctx->ms-ctx->fm)*hm+0.5*hs+(i-ctx->ms)*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)/8.0/ctx->sm; hm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); hf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); ierr = VecDuplicate(X,&Y);CHKERRQ(ierr); ierr = FVSample_3WaySplit(ctx,da,t,Y);CHKERRQ(ierr); ierr = VecGetArrayRead(X,&ptr_X);CHKERRQ(ierr); ierr = VecGetArrayRead(Y,&ptr_Y);CHKERRQ(ierr); for (i=0;ism || i > ctx->ms - 1) *nrm1 += hs*PetscAbs(ptr_X[i]-ptr_Y[i]); else if (i < ctx->mf || i > ctx->fm - 1) *nrm1 += hm*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_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) { FVCtx *ctx = (FVCtx*)vctx; PetscErrorCode ierr; PetscInt i,j,k,Mx,dof,xs,xm,sm = ctx->sm,mf = ctx->mf,fm = ctx->fm,ms = ctx->ms; PetscReal hxf,hxm,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)/8.0/ctx->sm; hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 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; jphysics2.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.hxm = hxm; info.hxf = hxf; (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,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 < sm || i > ms) { /* 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 == sm) { /* interface between slow and medium component */ 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]/hxm; } } else if (i == ms) { /* interface between medium and slow regions */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxm; } if (i < xs+xm) { for (j=0; jflux[j]/hxs; } } else if (i < mf || i > fm) { /* medium 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]/hxm; } if (i < xs+xm) { for (j=0; jflux[j]/hxm; } } else if (i == mf) { /* interface between medium and fast regions */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxm; } if (i < xs+xm) { for (j=0; jflux[j]/hxf; } } else if (i == fm) { /* interface between fast and medium regions */ 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]/hxm; } } else { /* fast 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]/hxf; } if (i < xs+xm) { for (j=0; jflux[j]/hxf; } } } 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) { 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); } } PetscFunctionReturn(0); } /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */ PetscErrorCode FVRHSFunctionslow_3WaySplit(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,sm = ctx->sm,ms = ctx->ms,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth; PetscReal hxs,hxm,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)/8.0/ctx->sm; hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 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; j ms+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.hxm = hxm; info.hxf = hxf; (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,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 < sm-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 == sm-lsbwidth) { /* interface between slow and medium regions */ 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 == ms+rsbwidth) { /* interface between medium and slow regions */ 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 > ms+rsbwidth) { /* 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++; } } } 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_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) { FVCtx *ctx = (FVCtx*)vctx; PetscErrorCode ierr; PetscInt i,j,k,Mx,dof,xs,xm,islowbuffer = 0,sm = ctx->sm,ms = ctx->ms,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth; PetscReal hxs,hxm,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)/8.0/ctx->sm; hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 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; j sm-lsbwidth-2 && i < sm+1) || (i > ms-2 && i < ms+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.hxm = hxm; info.hxf = hxf; (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,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 == sm-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; islowbuffer++; } } if (i > sm-lsbwidth && i < sm) { 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; islowbuffer++; } } if (i == sm) { /* interface between the slow region and the medium 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 == ms) { /* interface between the medium 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; islowbuffer++; } } if (i > ms && i < ms+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; islowbuffer++; } } if (i == ms+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 medium components ----------------------------------- */ PetscErrorCode FVRHSFunctionmedium_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) { FVCtx *ctx = (FVCtx*)vctx; PetscErrorCode ierr; PetscInt i,j,k,Mx,dof,xs,xm,imedium = 0,sm = ctx->sm,mf = ctx->mf,fm = ctx->fm,ms = ctx->ms,lmbwidth = ctx->lmbwidth,rmbwidth = ctx->rmbwidth; PetscReal hxs,hxm,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)/8.0/ctx->sm; hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 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; j sm-2 && i < mf-lmbwidth+1) || (i > fm+rmbwidth-2 && i < ms+1)) { /* 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.hxm = hxm; info.hxf = hxf; (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,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 == sm) { /* interface between slow and medium regions */ 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]/hxm; imedium++; } } if (i > sm && i < mf-lmbwidth) { /* medium 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]/hxm; } if (i < xs+xm) { for (j=0; jflux[j]/hxm; imedium++; } } if (i == mf-lmbwidth) { /* interface between medium and fast regions */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxm; } } if (i == fm+rmbwidth) { /* interface between fast and medium regions */ 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]/hxm; imedium++; } } if (i > fm+rmbwidth && i < ms) { /* medium 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]/hxm; } if (i < xs+xm) { for (j=0; jflux[j]/hxm; imedium++; } } if (i == ms) { /* interface between medium and slow regions */ for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxm; } } } 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); } PetscErrorCode FVRHSFunctionmediumbuffer_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) { FVCtx *ctx = (FVCtx*)vctx; PetscErrorCode ierr; PetscInt i,j,k,Mx,dof,xs,xm,imediumbuffer = 0,mf = ctx->mf,fm = ctx->fm,lmbwidth = ctx->lmbwidth,rmbwidth = ctx->rmbwidth; PetscReal hxs,hxm,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)/8.0/ctx->sm; hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 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; j mf-lmbwidth-2 && i < mf+1) || (i > fm-2 && i < fm+rmbwidth+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.hxm = hxm; info.hxf = hxf; (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,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 == mf-lmbwidth) { 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]/hxm; imediumbuffer++; } } if (i > mf-lmbwidth && i < mf) { for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxm; } if (i < xs+xm) { for (j=0; jflux[j]/hxm; imediumbuffer++; } } if (i == mf) { /* interface between the medium 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]/hxm; } } if (i == fm) { /* interface between the fast region and the medium 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]/hxm; imediumbuffer++; } } if (i > fm && i < fm+rmbwidth) { for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxm; } if (i < xs+xm) { for (j=0; jflux[j]/hxm; imediumbuffer++; } } if (i == fm+rmbwidth) { for (j=0; jphysics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hxm; } } } 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_3WaySplit(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,mf = ctx->mf,fm = ctx->fm; PetscReal hxs,hxm,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)/8.0/ctx->sm; hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 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; j mf-2 && i < fm+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.hxm = hxm; info.hxf = hxf; (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,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 == mf) { /* interface between medium and fast regions */ 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 > mf && i < fm) { /* 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 == fm) { /* interface between fast and medium regions */ 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_medium,count_fast,islow = 0,imedium = 0,ifast = 0,*index_slow,*index_medium,*index_fast,islowbuffer = 0,*index_slowbuffer,imediumbuffer = 0,*index_mediumbuffer; PetscBool view_final = PETSC_FALSE; PetscReal ptime; 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" ,Limit3_Upwind);CHKERRQ(ierr); ierr = PetscFunctionListAdd(&limiters,"lax-wendroff" ,Limit3_LaxWendroff);CHKERRQ(ierr); ierr = PetscFunctionListAdd(&limiters,"beam-warming" ,Limit3_BeamWarming);CHKERRQ(ierr); ierr = PetscFunctionListAdd(&limiters,"fromm" ,Limit3_Fromm);CHKERRQ(ierr); ierr = PetscFunctionListAdd(&limiters,"minmod" ,Limit3_Minmod);CHKERRQ(ierr); ierr = PetscFunctionListAdd(&limiters,"superbee" ,Limit3_Superbee);CHKERRQ(ierr); ierr = PetscFunctionListAdd(&limiters,"mc" ,Limit3_MC);CHKERRQ(ierr); ierr = PetscFunctionListAdd(&limiters,"koren3" ,Limit3_Koren3);CHKERRQ(ierr); /* Register physical models to be available on the command line */ 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; 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 = 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.limit3);CHKERRQ(ierr); if (!ctx.limit3) SETERRQ1(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); if (!r) SETERRQ1(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; ia > 0) { ctx.lsbwidth = 2; ctx.rsbwidth = 4; ctx.lmbwidth = 2; ctx.rmbwidth = 4; } else { ctx.lsbwidth = 4; ctx.rsbwidth = 2; ctx.lmbwidth = 4; ctx.rmbwidth = 2; } for (i=xs; i ctx.ms+ctx.rsbwidth-1) for (k=0; k= ctx.sm-ctx.lsbwidth && i < ctx.sm) || (i > ctx.ms-1 && i <= ctx.ms+ctx.rsbwidth-1)) for (k=0; k ctx.fm+ctx.rmbwidth-1) for (k=0; k= ctx.mf-ctx.lmbwidth && i < ctx.mf) || (i > ctx.fm-1 && i <= ctx.fm+ctx.rmbwidth-1)) for (k=0; k ctx.ms-1) for (k=0; k ctx.fm-1) for (k=0; k ctx.ms-1) for (k=0; k