static const char help[] = "1D periodic Finite Volume solver by a particular slope limiter 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, say\n" " hxs = (xmax - xmin)/2.0*(hratio+1.0)/Mx, \n" " hxf = (xmax - xmin)/2.0*(1.0+1.0/hratio)/Mx, \n" " with x belongs to (xmin,xmax), the number of total mesh points is Mx and the ratio between the meshsize of corse\n\n" " grids and fine grids is 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" "This script choose the slope limiter by biased second-order upwind procedure which is proposed by Van Leer in 1994\n" " u(x_(k+1/2),t) = u(x_k,t) + phi(x_(k+1/2),t)*(u(x_k,t)-u(x_(k-1),t)) \n" " limiter phi(x_(k+1/2),t) = max(0,min(r(k+1/2),min(2,gamma(k+1/2)*r(k+1/2)+alpha(k+1/2)))) \n" " r(k+1/2) = (u(x_(k+1))-u(x_k))/(u(x_k)-u(x_(k-1))) \n" " alpha(k+1/2) = (h_k*h_(k+1))/(h_(k-1)+h_k)/(h_(k-1)+h_k+h_(k+1)) \n" " gamma(k+1/2) = h_k*(h_(k-1)+h_k)/(h_k+h_(k+1))/(h_(k-1)+h_k+h_(k+1)) \n"; #include #include #include #include #include PETSC_STATIC_INLINE PetscReal RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax) { PetscReal range = xmax-xmin; return xmin +PetscFmodReal(range+PetscFmodReal(a,range),range); } /* --------------------------------- Finite Volume data structures ----------------------------------- */ typedef enum {FVBC_PERIODIC, FVBC_OUTFLOW} FVBCType; static const char *FVBCTypes[] = {"PERIODIC","OUTFLOW","FVBCType","FVBC_",0}; typedef struct { PetscErrorCode (*sample)(void*,PetscInt,FVBCType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*); PetscErrorCode (*flux)(void*,const PetscScalar*,PetscScalar*,PetscReal*); PetscErrorCode (*destroy)(void*); void *user; PetscInt dof; char *fieldname[16]; } PhysicsCtx; typedef struct { PhysicsCtx physics; MPI_Comm comm; char prefix[256]; /* Local work arrays */ PetscScalar *flux; /* Flux across interface */ PetscReal *speeds; /* Speeds of each wave */ PetscScalar *u; /* value at face */ PetscReal cfl_idt; /* Max allowable value of 1/Delta t */ PetscReal cfl; PetscReal xmin,xmax; PetscInt initial; PetscBool exact; PetscBool simulation; FVBCType bctype; PetscInt hratio; /* hratio = hslow/hfast */ IS isf,iss; PetscInt sf,fs; /* slow-fast and fast-slow interfaces */ } FVCtx; /* --------------------------------- Physics ----------------------------------- */ static PetscErrorCode PhysicsDestroy_SimpleFree(void *vctx) { PetscErrorCode ierr; PetscFunctionBeginUser; ierr = PetscFree(vctx);CHKERRQ(ierr); PetscFunctionReturn(0); } /* --------------------------------- Advection ----------------------------------- */ typedef struct { PetscReal a; /* advective velocity */ } AdvectCtx; static PetscErrorCode PhysicsFlux_Advect(void *vctx,const PetscScalar *u,PetscScalar *flux,PetscReal *maxspeed) { AdvectCtx *ctx = (AdvectCtx*)vctx; PetscReal speed; PetscFunctionBeginUser; speed = ctx->a; flux[0] = speed*u[0]; *maxspeed = speed; 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->physics.sample = PhysicsSample_Advect; ctx->physics.flux = PhysicsFlux_Advect; ctx->physics.destroy = PhysicsDestroy_SimpleFree; ctx->physics.user = user; ctx->physics.dof = 1; ierr = PetscStrallocpy("u",&ctx->physics.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); } /* --------------------------------- Finite Volume Solver ----------------------------------- */ static PetscErrorCode FVRHSFunction(TS ts,PetscReal time,Vec X,Vec F,void *vctx) { FVCtx *ctx = (FVCtx*)vctx; PetscErrorCode ierr; PetscInt i,j,Mx,dof,xs,xm,sf = ctx->sf,fs = ctx->fs; PetscReal hf,hs,cfl_idt = 0; PetscScalar *x,*f,*r,*min,*alpha,*gamma; 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 */ hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 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 = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); ierr = PetscMalloc4(dof,&r,dof,&min,dof,&alpha,dof,&gamma);CHKERRQ(ierr); if (ctx->bctype == FVBC_OUTFLOW) { for (i=xs-2; i<0; i++) { for (j=0; j fs+1) { u = &ctx->u[0]; alpha[0] = 1.0/6.0; gamma[0] = 1.0/3.0; for (j=0; jphysics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hs)); if (i > xs) { for (j=0; jflux[j]/hs; } if (i < xs+xm) { for (j=0; jflux[j]/hs; } } else if (i == sf) { u = &ctx->u[0]; alpha[0] = hs*hf/(hs+hs)/(hs+hs+hf); gamma[0] = hs*(hs+hs)/(hs+hf)/(hs+hs+hf); for (j=0; jphysics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hs; } if (i < xs+xm) { for (j=0; jflux[j]/hf; } } else if (i == sf+1) { u = &ctx->u[0]; alpha[0] = hf*hf/(hs+hf)/(hs+hf+hf); gamma[0] = hf*(hs+hf)/(hf+hf)/(hs+hf+hf); for (j=0; jphysics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hf; } if (i < xs+xm) { for (j=0; jflux[j]/hf; } } else if (i > sf+1 && i < fs) { u = &ctx->u[0]; alpha[0] = 1.0/6.0; gamma[0] = 1.0/3.0; for (j=0; jphysics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hf; } if (i < xs+xm) { for (j=0; jflux[j]/hf; } } else if (i == fs) { u = &ctx->u[0]; alpha[0] = hf*hs/(hf+hf)/(hf+hf+hs); gamma[0] = hf*(hf+hf)/(hf+hs)/(hf+hf+hs); for (j=0; jphysics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hf; } if (i < xs+xm) { for (j=0; jflux[j]/hs; } } else if (i == fs+1) { u = &ctx->u[0]; alpha[0] = hs*hs/(hf+hs)/(hf+hs+hs); gamma[0] = hs*(hf+hs)/(hs+hs)/(hf+hs+hs); for (j=0; jphysics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hs; } if (i < xs+xm) { for (j=0; jflux[j]/hs; } } } ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); ierr = DMDAVecRestoreArray(da,F,&f);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 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); } } ierr = PetscFree4(r,min,alpha,gamma);CHKERRQ(ierr); PetscFunctionReturn(0); } static PetscErrorCode FVRHSFunctionslow(TS ts,PetscReal time,Vec X,Vec F,void *vctx) { FVCtx *ctx = (FVCtx*)vctx; PetscErrorCode ierr; PetscInt i,j,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs; PetscReal hf,hs; PetscScalar *x,*f,*r,*min,*alpha,*gamma; 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 */ hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 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 = VecGetArray(F,&f);CHKERRQ(ierr); ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); ierr = PetscMalloc4(dof,&r,dof,&min,dof,&alpha,dof,&gamma);CHKERRQ(ierr); if (ctx->bctype == FVBC_OUTFLOW) { for (i=xs-2; i<0; i++) { for (j=0; ju[0]; alpha[0] = 1.0/6.0; gamma[0] = 1.0/3.0; for (j=0; jphysics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hs; } if (i < xs+xm) { for (j=0; jflux[j]/hs; islow++; } } else if (i == sf) { u = &ctx->u[0]; alpha[0] = hs*hf/(hs+hs)/(hs+hs+hf); gamma[0] = hs*(hs+hs)/(hs+hf)/(hs+hs+hf); for (j=0; jphysics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hs; } } else if (i == fs) { u = &ctx->u[0]; alpha[0] = hf*hs/(hf+hf)/(hf+hf+hs); gamma[0] = hf*(hf+hf)/(hf+hs)/(hf+hf+hs); for (j=0; jphysics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i < xs+xm) { for (j=0; jflux[j]/hs; islow++; } } else if (i == fs+1) { u = &ctx->u[0]; alpha[0] = hs*hs/(hf+hs)/(hf+hs+hs); gamma[0] = hs*(hf+hs)/(hs+hs)/(hf+hs+hs); for (j=0; jphysics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hs; } if (i < xs+xm) { for (j=0; jflux[j]/hs; islow++; } } else if (i > fs+1) { u = &ctx->u[0]; alpha[0] = 1.0/6.0; gamma[0] = 1.0/3.0; for (j=0; jphysics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hs; } if (i < xs+xm) { for (j=0; jflux[j]/hs; islow++; } } } ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); ierr = PetscFree4(r,min,alpha,gamma);CHKERRQ(ierr); PetscFunctionReturn(0); } static PetscErrorCode FVRHSFunctionfast(TS ts,PetscReal time,Vec X,Vec F,void *vctx) { FVCtx *ctx = (FVCtx*)vctx; PetscErrorCode ierr; PetscInt i,j,Mx,dof,xs,xm,ifast = 0,sf = ctx->sf,fs = ctx->fs; PetscReal hf,hs; PetscScalar *x,*f,*r,*min,*alpha,*gamma; 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 */ hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 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 = VecGetArray(F,&f);CHKERRQ(ierr); ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); ierr = PetscMalloc4(dof,&r,dof,&min,dof,&alpha,dof,&gamma);CHKERRQ(ierr); if (ctx->bctype == FVBC_OUTFLOW) { for (i=xs-2; i<0; i++) { for (j=0; ju[0]; alpha[0] = hs*hf/(hs+hs)/(hs+hs+hf); gamma[0] = hs*(hs+hs)/(hs+hf)/(hs+hs+hf); for (j=0; jphysics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i < xs+xm) { for (j=0; jflux[j]/hf; ifast++; } } else if (i == sf+1) { u = &ctx->u[0]; alpha[0] = hf*hf/(hs+hf)/(hs+hf+hf); gamma[0] = hf*(hs+hf)/(hf+hf)/(hs+hf+hf); for (j=0; jphysics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hf; } if (i < xs+xm) { for (j=0; jflux[j]/hf; ifast++; } } else if (i > sf+1 && i < fs) { u = &ctx->u[0]; alpha[0] = 1.0/6.0; gamma[0] = 1.0/3.0; for (j=0; jphysics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hf; } if (i < xs+xm) { for (j=0; jflux[j]/hf; ifast++; } } else if (i == fs) { u = &ctx->u[0]; alpha[0] = hf*hs/(hf+hf)/(hf+hf+hs); gamma[0] = hf*(hf+hf)/(hf+hs)/(hf+hf+hs); for (j=0; jphysics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); if (i > xs) { for (j=0; jflux[j]/hf; } } } ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); ierr = PetscFree4(r,min,alpha,gamma);CHKERRQ(ierr); PetscFunctionReturn(0); } /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */ PetscErrorCode FVSample(FVCtx *ctx,DM da,PetscReal time,Vec U) { PetscErrorCode ierr; PetscScalar *u,*uj,xj,xi; PetscInt i,j,k,dof,xs,xm,Mx,count_slow,count_fast; const PetscInt N=200; PetscFunctionBeginUser; if (!ctx->physics.sample) 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); const PetscReal hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; const PetscReal hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; count_slow = Mx/(1+ctx->hratio); count_fast = Mx-count_slow; for (i=xs; ixmax-ctx->xmin)*0.25) { xi = ctx->xmin+0.5*hs+i*hs; /* Integrate over cell i using trapezoid rule with N points. */ for (k=0; kphysics.sample)(ctx->physics.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); for (k=0; kxmax-ctx->xmin)*0.25+(i-count_slow/2)*hf+0.5*hf<(ctx->xmax-ctx->xmin)*0.75) { xi = ctx->xmin+(ctx->xmax-ctx->xmin)*0.25+0.5*hf+(i-count_slow/2)*hf; /* Integrate over cell i using trapezoid rule with N points. */ for (k=0; kphysics.sample)(ctx->physics.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); for (k=0; kxmin+(ctx->xmax-ctx->xmin)*0.75+0.5*hs+(i-count_slow/2-count_fast)*hs; /* Integrate over cell i using trapezoid rule with N points. */ for (k=0; kphysics.sample)(ctx->physics.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); for (k=0; kxmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; const PetscReal hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; count_slow = (PetscReal)Mx/(1.0+ctx->hratio); count_fast = Mx-count_slow; ierr = VecGetArrayRead(X,&ptr_X);CHKERRQ(ierr); ierr = VecGetArrayRead(Y,&ptr_Y);CHKERRQ(ierr); for (i=0; i count_slow/2+count_fast-1) *nrm1 += hs*PetscAbs(ptr_X[i]-ptr_Y[i]); else *nrm1 += hf*PetscAbs(ptr_X[i]-ptr_Y[i]); } ierr = VecRestoreArrayRead(X,&ptr_X);CHKERRQ(ierr); ierr = VecRestoreArrayRead(Y,&ptr_Y);CHKERRQ(ierr); ierr = VecDestroy(&Y);CHKERRQ(ierr); PetscFunctionReturn(0); } int main(int argc,char *argv[]) { char physname[256] = "advect",final_fname[256] = "solution.m"; PetscFunctionList physics = 0; MPI_Comm comm; TS ts; DM da; Vec X,X0,R; FVCtx ctx; PetscInt i,k,dof,xs,xm,Mx,draw = 0,count_slow,count_fast,islow = 0,ifast = 0,*index_slow,*index_fast; 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 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 = 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 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.physics.dof,2,NULL,&da);CHKERRQ(ierr); ierr = DMSetFromOptions(da);CHKERRQ(ierr); ierr = DMSetUp(da);CHKERRQ(ierr); /* Inform the DMDA of the field names provided by the physics. */ /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ for (i=0; i count_slow/2+count_fast-1) for (k=0; k ctx.fs-1) { for (k=0; k count_slow/2+count_fast-1) nrm1 = nrm1 + hs*PetscAbs(ptr_X[i]-ptr_XR[i]); else nrm1 = nrm1 + hf*PetscAbs(ptr_X[i]-ptr_XR[i]); } ierr = VecRestoreArrayRead(X,&ptr_X);CHKERRQ(ierr); ierr = VecRestoreArrayRead(XR,&ptr_XR);CHKERRQ(ierr); ierr = PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1);CHKERRQ(ierr); ierr = VecDestroy(&XR);CHKERRQ(ierr); } } ierr = SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); if (draw & 0x1) { ierr = VecView(X0,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr); } if (draw & 0x2) { ierr = VecView(X,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr); } if (draw & 0x4) { Vec Y; ierr = VecDuplicate(X,&Y);CHKERRQ(ierr); ierr = FVSample(&ctx,da,ptime,Y);CHKERRQ(ierr); ierr = VecAYPX(Y,-1,X);CHKERRQ(ierr); ierr = VecView(Y,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr); ierr = VecDestroy(&Y);CHKERRQ(ierr); } if (view_final) { PetscViewer viewer; ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer);CHKERRQ(ierr); ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); ierr = VecView(X,viewer);CHKERRQ(ierr); ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); } /* Clean up */ ierr = (*ctx.physics.destroy)(ctx.physics.user);CHKERRQ(ierr); for (i=0; i