1c4762a1bSJed Brown static const char help[] = "1D periodic Finite Volume solver by a particular slope limiter with semidiscrete time stepping.\n" 2c4762a1bSJed Brown " advection - Constant coefficient scalar advection\n" 3c4762a1bSJed Brown " u_t + (a*u)_x = 0\n" 4c4762a1bSJed Brown " for this toy problem, we choose different meshsizes for different sub-domains, say\n" 5c4762a1bSJed Brown " hxs = (xmax - xmin)/2.0*(hratio+1.0)/Mx, \n" 6c4762a1bSJed Brown " hxf = (xmax - xmin)/2.0*(1.0+1.0/hratio)/Mx, \n" 76aad120cSJose E. Roman " with x belongs to (xmin,xmax), the number of total mesh points is Mx and the ratio between the meshsize of coarse\n\n" 8c4762a1bSJed Brown " grids and fine grids is hratio.\n" 9c4762a1bSJed Brown " exact - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n" 10c4762a1bSJed Brown " the states across shocks and rarefactions\n" 11c4762a1bSJed Brown " simulation - use reference solution which is generated by smaller time step size to be true solution,\n" 12c4762a1bSJed Brown " also the reference solution should be generated by user and stored in a binary file.\n" 13c4762a1bSJed Brown " characteristic - Limit the characteristic variables, this is usually preferred (default)\n" 14c4762a1bSJed Brown "Several initial conditions can be chosen with -initial N\n\n" 15c4762a1bSJed Brown "The problem size should be set with -da_grid_x M\n\n" 16c4762a1bSJed Brown "This script choose the slope limiter by biased second-order upwind procedure which is proposed by Van Leer in 1994\n" 17c4762a1bSJed Brown " 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" 18c4762a1bSJed Brown " 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" 19c4762a1bSJed Brown " r(k+1/2) = (u(x_(k+1))-u(x_k))/(u(x_k)-u(x_(k-1))) \n" 20c4762a1bSJed Brown " alpha(k+1/2) = (h_k*h_(k+1))/(h_(k-1)+h_k)/(h_(k-1)+h_k+h_(k+1)) \n" 21c4762a1bSJed Brown " 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"; 22c4762a1bSJed Brown 23c4762a1bSJed Brown #include <petscts.h> 24c4762a1bSJed Brown #include <petscdm.h> 25c4762a1bSJed Brown #include <petscdmda.h> 26c4762a1bSJed Brown #include <petscdraw.h> 27c4762a1bSJed Brown #include <petscmath.h> 28c4762a1bSJed Brown 299fbee547SJacob Faibussowitsch static inline PetscReal RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax) { PetscReal range = xmax-xmin; return xmin +PetscFmodReal(range+PetscFmodReal(a,range),range); } 30c4762a1bSJed Brown 31c4762a1bSJed Brown /* --------------------------------- Finite Volume data structures ----------------------------------- */ 32c4762a1bSJed Brown 33c4762a1bSJed Brown typedef enum {FVBC_PERIODIC, FVBC_OUTFLOW} FVBCType; 34c4762a1bSJed Brown static const char *FVBCTypes[] = {"PERIODIC","OUTFLOW","FVBCType","FVBC_",0}; 35c4762a1bSJed Brown 36c4762a1bSJed Brown typedef struct { 37c4762a1bSJed Brown PetscErrorCode (*sample)(void*,PetscInt,FVBCType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*); 38c4762a1bSJed Brown PetscErrorCode (*flux)(void*,const PetscScalar*,PetscScalar*,PetscReal*); 39c4762a1bSJed Brown PetscErrorCode (*destroy)(void*); 40c4762a1bSJed Brown void *user; 41c4762a1bSJed Brown PetscInt dof; 42c4762a1bSJed Brown char *fieldname[16]; 43c4762a1bSJed Brown } PhysicsCtx; 44c4762a1bSJed Brown 45c4762a1bSJed Brown typedef struct { 46c4762a1bSJed Brown PhysicsCtx physics; 47c4762a1bSJed Brown MPI_Comm comm; 48c4762a1bSJed Brown char prefix[256]; 49c4762a1bSJed Brown 50c4762a1bSJed Brown /* Local work arrays */ 51c4762a1bSJed Brown PetscScalar *flux; /* Flux across interface */ 52c4762a1bSJed Brown PetscReal *speeds; /* Speeds of each wave */ 53c4762a1bSJed Brown PetscScalar *u; /* value at face */ 54c4762a1bSJed Brown 55c4762a1bSJed Brown PetscReal cfl_idt; /* Max allowable value of 1/Delta t */ 56c4762a1bSJed Brown PetscReal cfl; 57c4762a1bSJed Brown PetscReal xmin,xmax; 58c4762a1bSJed Brown PetscInt initial; 59c4762a1bSJed Brown PetscBool exact; 60c4762a1bSJed Brown PetscBool simulation; 61c4762a1bSJed Brown FVBCType bctype; 62c4762a1bSJed Brown PetscInt hratio; /* hratio = hslow/hfast */ 63c4762a1bSJed Brown IS isf,iss; 64c4762a1bSJed Brown PetscInt sf,fs; /* slow-fast and fast-slow interfaces */ 65c4762a1bSJed Brown } FVCtx; 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* --------------------------------- Physics ----------------------------------- */ 68c4762a1bSJed Brown static PetscErrorCode PhysicsDestroy_SimpleFree(void *vctx) 69c4762a1bSJed Brown { 70c4762a1bSJed Brown PetscFunctionBeginUser; 719566063dSJacob Faibussowitsch PetscCall(PetscFree(vctx)); 72c4762a1bSJed Brown PetscFunctionReturn(0); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 75c4762a1bSJed Brown /* --------------------------------- Advection ----------------------------------- */ 76c4762a1bSJed Brown typedef struct { 77c4762a1bSJed Brown PetscReal a; /* advective velocity */ 78c4762a1bSJed Brown } AdvectCtx; 79c4762a1bSJed Brown 80c4762a1bSJed Brown static PetscErrorCode PhysicsFlux_Advect(void *vctx,const PetscScalar *u,PetscScalar *flux,PetscReal *maxspeed) 81c4762a1bSJed Brown { 82c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx*)vctx; 83c4762a1bSJed Brown PetscReal speed; 84c4762a1bSJed Brown 85c4762a1bSJed Brown PetscFunctionBeginUser; 86c4762a1bSJed Brown speed = ctx->a; 87c4762a1bSJed Brown flux[0] = speed*u[0]; 88c4762a1bSJed Brown *maxspeed = speed; 89c4762a1bSJed Brown PetscFunctionReturn(0); 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92c4762a1bSJed Brown static PetscErrorCode PhysicsSample_Advect(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u) 93c4762a1bSJed Brown { 94c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx*)vctx; 95c4762a1bSJed Brown PetscReal a = ctx->a,x0; 96c4762a1bSJed Brown 97c4762a1bSJed Brown PetscFunctionBeginUser; 98c4762a1bSJed Brown switch (bctype) { 99c4762a1bSJed Brown case FVBC_OUTFLOW: x0 = x-a*t; break; 100c4762a1bSJed Brown case FVBC_PERIODIC: x0 = RangeMod(x-a*t,xmin,xmax); break; 101c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown BCType"); 102c4762a1bSJed Brown } 103c4762a1bSJed Brown switch (initial) { 104c4762a1bSJed Brown case 0: u[0] = (x0 < 0) ? 1 : -1; break; 105c4762a1bSJed Brown case 1: u[0] = (x0 < 0) ? -1 : 1; break; 106c4762a1bSJed Brown case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break; 107c4762a1bSJed Brown case 3: u[0] = PetscSinReal(2*PETSC_PI*x0); break; 108c4762a1bSJed Brown case 4: u[0] = PetscAbs(x0); break; 109c4762a1bSJed Brown case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2*PETSC_PI*x0)); break; 110c4762a1bSJed Brown case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2-x0 : 0)); break; 111c4762a1bSJed Brown case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI*x0),10.0);break; 112c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition"); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown PetscFunctionReturn(0); 115c4762a1bSJed Brown } 116c4762a1bSJed Brown 117c4762a1bSJed Brown static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx) 118c4762a1bSJed Brown { 119c4762a1bSJed Brown AdvectCtx *user; 120c4762a1bSJed Brown 121c4762a1bSJed Brown PetscFunctionBeginUser; 1229566063dSJacob Faibussowitsch PetscCall(PetscNew(&user)); 123c4762a1bSJed Brown ctx->physics.sample = PhysicsSample_Advect; 124c4762a1bSJed Brown ctx->physics.flux = PhysicsFlux_Advect; 125c4762a1bSJed Brown ctx->physics.destroy = PhysicsDestroy_SimpleFree; 126c4762a1bSJed Brown ctx->physics.user = user; 127c4762a1bSJed Brown ctx->physics.dof = 1; 1289566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("u",&ctx->physics.fieldname[0])); 129c4762a1bSJed Brown user->a = 1; 130d0609cedSBarry Smith PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection",""); 131c4762a1bSJed Brown { 1329566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-physics_advect_a","Speed","",user->a,&user->a,NULL)); 133c4762a1bSJed Brown } 134d0609cedSBarry Smith PetscOptionsEnd(); 135c4762a1bSJed Brown PetscFunctionReturn(0); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver ----------------------------------- */ 139c4762a1bSJed Brown 140c4762a1bSJed Brown static PetscErrorCode FVRHSFunction(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 141c4762a1bSJed Brown { 142c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 143c4762a1bSJed Brown PetscInt i,j,Mx,dof,xs,xm,sf = ctx->sf,fs = ctx->fs; 144c4762a1bSJed Brown PetscReal hf,hs,cfl_idt = 0; 145c4762a1bSJed Brown PetscScalar *x,*f,*r,*min,*alpha,*gamma; 146c4762a1bSJed Brown Vec Xloc; 147c4762a1bSJed Brown DM da; 148c4762a1bSJed Brown 149c4762a1bSJed Brown PetscFunctionBeginUser; 1509566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 1519566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); /* Xloc contains ghost points */ 1529566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0)); /* Mx is the number of center points */ 153c4762a1bSJed Brown hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; 154c4762a1bSJed Brown hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 1559566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); /* X is solution vector which does not contain ghost points */ 1569566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 1579566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); /* F is the right hand side function corresponds to center points */ 1589566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 1599566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,F,&f)); 1609566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 1619566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dof,&r,dof,&min,dof,&alpha,dof,&gamma)); 162c4762a1bSJed Brown 163c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 164c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 165c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 166c4762a1bSJed Brown } 167c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 168c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 169c4762a1bSJed Brown } 170c4762a1bSJed Brown } 171c4762a1bSJed Brown 172c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 173c4762a1bSJed Brown PetscReal maxspeed; 174c4762a1bSJed Brown PetscScalar *u; 175c4762a1bSJed Brown if (i < sf || i > fs+1) { 176c4762a1bSJed Brown u = &ctx->u[0]; 177c4762a1bSJed Brown alpha[0] = 1.0/6.0; 178c4762a1bSJed Brown gamma[0] = 1.0/3.0; 179c4762a1bSJed Brown for (j=0; j<dof; j++) { 180c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 181c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 182c4762a1bSJed Brown u[j] = x[(i-1)*dof+j]+PetscMax(0,PetscMin(min[j],alpha[0]+gamma[0]*r[j]))*(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 183c4762a1bSJed Brown } 1849566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 185c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hs)); 186c4762a1bSJed Brown if (i > xs) { 187c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hs; 188c4762a1bSJed Brown } 189c4762a1bSJed Brown if (i < xs+xm) { 190c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hs; 191c4762a1bSJed Brown } 192c4762a1bSJed Brown } else if (i == sf) { 193c4762a1bSJed Brown u = &ctx->u[0]; 194c4762a1bSJed Brown alpha[0] = hs*hf/(hs+hs)/(hs+hs+hf); 195c4762a1bSJed Brown gamma[0] = hs*(hs+hs)/(hs+hf)/(hs+hs+hf); 196c4762a1bSJed Brown for (j=0; j<dof; j++) { 197c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 198c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 199c4762a1bSJed Brown u[j] = x[(i-1)*dof+j]+PetscMax(0,PetscMin(min[j],alpha[0]+gamma[0]*r[j]))*(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 200c4762a1bSJed Brown } 2019566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 202c4762a1bSJed Brown if (i > xs) { 203c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hs; 204c4762a1bSJed Brown } 205c4762a1bSJed Brown if (i < xs+xm) { 206c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hf; 207c4762a1bSJed Brown } 208c4762a1bSJed Brown } else if (i == sf+1) { 209c4762a1bSJed Brown u = &ctx->u[0]; 210c4762a1bSJed Brown alpha[0] = hf*hf/(hs+hf)/(hs+hf+hf); 211c4762a1bSJed Brown gamma[0] = hf*(hs+hf)/(hf+hf)/(hs+hf+hf); 212c4762a1bSJed Brown for (j=0; j<dof; j++) { 213c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 214c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 215c4762a1bSJed Brown u[j] = x[(i-1)*dof+j]+PetscMax(0,PetscMin(min[j],alpha[0]+gamma[0]*r[j]))*(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 216c4762a1bSJed Brown } 2179566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 218c4762a1bSJed Brown if (i > xs) { 219c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hf; 220c4762a1bSJed Brown } 221c4762a1bSJed Brown if (i < xs+xm) { 222c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hf; 223c4762a1bSJed Brown } 224c4762a1bSJed Brown } else if (i > sf+1 && i < fs) { 225c4762a1bSJed Brown u = &ctx->u[0]; 226c4762a1bSJed Brown alpha[0] = 1.0/6.0; 227c4762a1bSJed Brown gamma[0] = 1.0/3.0; 228c4762a1bSJed Brown for (j=0; j<dof; j++) { 229c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 230c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 231c4762a1bSJed Brown u[j] = x[(i-1)*dof+j]+PetscMax(0,PetscMin(min[j],alpha[0]+gamma[0]*r[j]))*(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 232c4762a1bSJed Brown } 2339566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 234c4762a1bSJed Brown if (i > xs) { 235c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hf; 236c4762a1bSJed Brown } 237c4762a1bSJed Brown if (i < xs+xm) { 238c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hf; 239c4762a1bSJed Brown } 240c4762a1bSJed Brown } else if (i == fs) { 241c4762a1bSJed Brown u = &ctx->u[0]; 242c4762a1bSJed Brown alpha[0] = hf*hs/(hf+hf)/(hf+hf+hs); 243c4762a1bSJed Brown gamma[0] = hf*(hf+hf)/(hf+hs)/(hf+hf+hs); 244c4762a1bSJed Brown for (j=0; j<dof; j++) { 245c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 246c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 247c4762a1bSJed Brown u[j] = x[(i-1)*dof+j]+PetscMax(0,PetscMin(min[j],alpha[0]+gamma[0]*r[j]))*(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 248c4762a1bSJed Brown } 2499566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 250c4762a1bSJed Brown if (i > xs) { 251c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hf; 252c4762a1bSJed Brown } 253c4762a1bSJed Brown if (i < xs+xm) { 254c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hs; 255c4762a1bSJed Brown } 256c4762a1bSJed Brown } else if (i == fs+1) { 257c4762a1bSJed Brown u = &ctx->u[0]; 258c4762a1bSJed Brown alpha[0] = hs*hs/(hf+hs)/(hf+hs+hs); 259c4762a1bSJed Brown gamma[0] = hs*(hf+hs)/(hs+hs)/(hf+hs+hs); 260c4762a1bSJed Brown for (j=0; j<dof; j++) { 261c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 262c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 263c4762a1bSJed Brown u[j] = x[(i-1)*dof+j]+PetscMax(0,PetscMin(min[j],alpha[0]+gamma[0]*r[j]))*(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 264c4762a1bSJed Brown } 2659566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 266c4762a1bSJed Brown if (i > xs) { 267c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hs; 268c4762a1bSJed Brown } 269c4762a1bSJed Brown if (i < xs+xm) { 270c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hs; 271c4762a1bSJed Brown } 272c4762a1bSJed Brown } 273c4762a1bSJed Brown } 2749566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 2759566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,F,&f)); 2769566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 2779566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da))); 278c4762a1bSJed Brown if (0) { 279c2a16db8SHong Zhang /* We need a way to inform the TS of a CFL constraint, this is a debugging fragment */ 280c4762a1bSJed Brown PetscReal dt,tnow; 2819566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&dt)); 2829566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts,&tnow)); 283c4762a1bSJed Brown if (dt > 0.5/ctx->cfl_idt) { 2849566063dSJacob Faibussowitsch PetscCall(PetscPrintf(ctx->comm,"Stability constraint exceeded at t=%g, dt %g > %g\n",(double)tnow,(double)dt,(double)(0.5/ctx->cfl_idt))); 285c4762a1bSJed Brown } 286c4762a1bSJed Brown } 2879566063dSJacob Faibussowitsch PetscCall(PetscFree4(r,min,alpha,gamma)); 288c4762a1bSJed Brown PetscFunctionReturn(0); 289c4762a1bSJed Brown } 290c4762a1bSJed Brown 291c4762a1bSJed Brown static PetscErrorCode FVRHSFunctionslow(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 292c4762a1bSJed Brown { 293c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 294c4762a1bSJed Brown PetscInt i,j,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs; 295c4762a1bSJed Brown PetscReal hf,hs; 296c4762a1bSJed Brown PetscScalar *x,*f,*r,*min,*alpha,*gamma; 297c4762a1bSJed Brown Vec Xloc; 298c4762a1bSJed Brown DM da; 299c4762a1bSJed Brown 300c4762a1bSJed Brown PetscFunctionBeginUser; 3019566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 3029566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); /* Xloc contains ghost points */ 3039566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0)); /* Mx is the number of center points */ 304c4762a1bSJed Brown hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; 305c4762a1bSJed Brown hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 3069566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); /* X is solution vector which does not contain ghost points */ 3079566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd (da,X,INSERT_VALUES,Xloc)); 3089566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); /* F is the right hand side function corresponds to center points */ 3099566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 3109566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 3119566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 3129566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dof,&r,dof,&min,dof,&alpha,dof,&gamma)); 313c4762a1bSJed Brown 314c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 315c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 316c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 317c4762a1bSJed Brown } 318c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 319c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 320c4762a1bSJed Brown } 321c4762a1bSJed Brown } 322c4762a1bSJed Brown 323c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 324c4762a1bSJed Brown PetscReal maxspeed; 325c4762a1bSJed Brown PetscScalar *u; 326c4762a1bSJed Brown if (i < sf) { 327c4762a1bSJed Brown u = &ctx->u[0]; 328c4762a1bSJed Brown alpha[0] = 1.0/6.0; 329c4762a1bSJed Brown gamma[0] = 1.0/3.0; 330c4762a1bSJed Brown for (j=0; j<dof; j++) { 331c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 332c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 333c4762a1bSJed Brown u[j] = x[(i-1)*dof+j]+PetscMax(0,PetscMin(min[j],alpha[0]+gamma[0]*r[j]))*(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 334c4762a1bSJed Brown } 3359566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 336c4762a1bSJed Brown if (i > xs) { 337c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hs; 338c4762a1bSJed Brown } 339c4762a1bSJed Brown if (i < xs+xm) { 340c2a16db8SHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hs; 341c4762a1bSJed Brown islow++; 342c4762a1bSJed Brown } 343c4762a1bSJed Brown } else if (i == sf) { 344c4762a1bSJed Brown u = &ctx->u[0]; 345c4762a1bSJed Brown alpha[0] = hs*hf/(hs+hs)/(hs+hs+hf); 346c4762a1bSJed Brown gamma[0] = hs*(hs+hs)/(hs+hf)/(hs+hs+hf); 347c4762a1bSJed Brown for (j=0; j<dof; j++) { 348c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 349c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 350c4762a1bSJed Brown u[j] = x[(i-1)*dof+j]+PetscMax(0,PetscMin(min[j],alpha[0]+gamma[0]*r[j]))*(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 351c4762a1bSJed Brown } 3529566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 353c2a16db8SHong Zhang if (i > xs) { 354c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hs; 355c4762a1bSJed Brown } 356c4762a1bSJed Brown } else if (i == fs) { 357c4762a1bSJed Brown u = &ctx->u[0]; 358c4762a1bSJed Brown alpha[0] = hf*hs/(hf+hf)/(hf+hf+hs); 359c4762a1bSJed Brown gamma[0] = hf*(hf+hf)/(hf+hs)/(hf+hf+hs); 360c4762a1bSJed Brown for (j=0; j<dof; j++) { 361c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 362c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 363c4762a1bSJed Brown u[j] = x[(i-1)*dof+j]+PetscMax(0,PetscMin(min[j],alpha[0]+gamma[0]*r[j]))*(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 364c4762a1bSJed Brown } 3659566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 366c4762a1bSJed Brown if (i < xs+xm) { 367c2a16db8SHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hs; 368c4762a1bSJed Brown islow++; 369c4762a1bSJed Brown } 370c4762a1bSJed Brown } else if (i == fs+1) { 371c4762a1bSJed Brown u = &ctx->u[0]; 372c4762a1bSJed Brown alpha[0] = hs*hs/(hf+hs)/(hf+hs+hs); 373c4762a1bSJed Brown gamma[0] = hs*(hf+hs)/(hs+hs)/(hf+hs+hs); 374c4762a1bSJed Brown for (j=0; j<dof; j++) { 375c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 376c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 377c4762a1bSJed Brown u[j] = x[(i-1)*dof+j]+PetscMax(0,PetscMin(min[j],alpha[0]+gamma[0]*r[j]))*(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 378c4762a1bSJed Brown } 3799566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 380c4762a1bSJed Brown if (i > xs) { 381c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hs; 382c4762a1bSJed Brown } 383c4762a1bSJed Brown if (i < xs+xm) { 384c2a16db8SHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hs; 385c2a16db8SHong Zhang islow++; 386c4762a1bSJed Brown } 387c4762a1bSJed Brown } else if (i > fs+1) { 388c4762a1bSJed Brown u = &ctx->u[0]; 389c4762a1bSJed Brown alpha[0] = 1.0/6.0; 390c4762a1bSJed Brown gamma[0] = 1.0/3.0; 391c4762a1bSJed Brown for (j=0; j<dof; j++) { 392c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 393c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 394c4762a1bSJed Brown u[j] = x[(i-1)*dof+j]+PetscMax(0,PetscMin(min[j],alpha[0]+gamma[0]*r[j]))*(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 395c4762a1bSJed Brown } 3969566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 397c4762a1bSJed Brown if (i > xs) { 398c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hs; 399c4762a1bSJed Brown } 400c4762a1bSJed Brown if (i < xs+xm) { 401c2a16db8SHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hs; 402c2a16db8SHong Zhang islow++; 403c4762a1bSJed Brown } 404c4762a1bSJed Brown } 405c4762a1bSJed Brown } 4069566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 4079566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 4089566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 4099566063dSJacob Faibussowitsch PetscCall(PetscFree4(r,min,alpha,gamma)); 410c4762a1bSJed Brown PetscFunctionReturn(0); 411c4762a1bSJed Brown } 412c4762a1bSJed Brown 413c4762a1bSJed Brown static PetscErrorCode FVRHSFunctionfast(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 414c4762a1bSJed Brown { 415c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 416c4762a1bSJed Brown PetscInt i,j,Mx,dof,xs,xm,ifast = 0,sf = ctx->sf,fs = ctx->fs; 417c4762a1bSJed Brown PetscReal hf,hs; 418c4762a1bSJed Brown PetscScalar *x,*f,*r,*min,*alpha,*gamma; 419c4762a1bSJed Brown Vec Xloc; 420c4762a1bSJed Brown DM da; 421c4762a1bSJed Brown 422c4762a1bSJed Brown PetscFunctionBeginUser; 4239566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 4249566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); /* Xloc contains ghost points */ 4259566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0)); /* Mx is the number of center points */ 426c4762a1bSJed Brown hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; 427c4762a1bSJed Brown hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 4289566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); /* X is solution vector which does not contain ghost points */ 4299566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 4309566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); /* F is the right hand side function corresponds to center points */ 4319566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 4329566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 4339566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 4349566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dof,&r,dof,&min,dof,&alpha,dof,&gamma)); 435c4762a1bSJed Brown 436c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 437c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 438c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 439c4762a1bSJed Brown } 440c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 441c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 442c4762a1bSJed Brown } 443c4762a1bSJed Brown } 444c4762a1bSJed Brown 445c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 446c4762a1bSJed Brown PetscReal maxspeed; 447c4762a1bSJed Brown PetscScalar *u; 448c4762a1bSJed Brown if (i == sf) { 449c4762a1bSJed Brown u = &ctx->u[0]; 450c4762a1bSJed Brown alpha[0] = hs*hf/(hs+hs)/(hs+hs+hf); 451c4762a1bSJed Brown gamma[0] = hs*(hs+hs)/(hs+hf)/(hs+hs+hf); 452c4762a1bSJed Brown for (j=0; j<dof; j++) { 453c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 454c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 455c4762a1bSJed Brown u[j] = x[(i-1)*dof+j]+PetscMax(0,PetscMin(min[j],alpha[0]+gamma[0]*r[j]))*(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 456c4762a1bSJed Brown } 4579566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 458c4762a1bSJed Brown if (i < xs+xm) { 459c2a16db8SHong Zhang for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hf; 460c4762a1bSJed Brown ifast++; 461c4762a1bSJed Brown } 462c4762a1bSJed Brown } else if (i == sf+1) { 463c4762a1bSJed Brown u = &ctx->u[0]; 464c4762a1bSJed Brown alpha[0] = hf*hf/(hs+hf)/(hs+hf+hf); 465c4762a1bSJed Brown gamma[0] = hf*(hs+hf)/(hf+hf)/(hs+hf+hf); 466c4762a1bSJed Brown for (j=0; j<dof; j++) { 467c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 468c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 469c4762a1bSJed Brown u[j] = x[(i-1)*dof+j]+PetscMax(0,PetscMin(min[j],alpha[0]+gamma[0]*r[j]))*(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 470c4762a1bSJed Brown } 4719566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 472c4762a1bSJed Brown if (i > xs) { 473c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hf; 474c4762a1bSJed Brown } 475c4762a1bSJed Brown if (i < xs+xm) { 476c2a16db8SHong Zhang for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hf; 477c4762a1bSJed Brown ifast++; 478c4762a1bSJed Brown } 479c4762a1bSJed Brown } else if (i > sf+1 && i < fs) { 480c4762a1bSJed Brown u = &ctx->u[0]; 481c4762a1bSJed Brown alpha[0] = 1.0/6.0; 482c4762a1bSJed Brown gamma[0] = 1.0/3.0; 483c4762a1bSJed Brown for (j=0; j<dof; j++) { 484c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 485c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 486c4762a1bSJed Brown u[j] = x[(i-1)*dof+j]+PetscMax(0,PetscMin(min[j],alpha[0]+gamma[0]*r[j]))*(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 487c4762a1bSJed Brown } 4889566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 489c4762a1bSJed Brown if (i > xs) { 490c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hf; 491c4762a1bSJed Brown } 492c4762a1bSJed Brown if (i < xs+xm) { 493c2a16db8SHong Zhang for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hf; 494c4762a1bSJed Brown ifast++; 495c4762a1bSJed Brown } 496c4762a1bSJed Brown } else if (i == fs) { 497c4762a1bSJed Brown u = &ctx->u[0]; 498c4762a1bSJed Brown alpha[0] = hf*hs/(hf+hf)/(hf+hf+hs); 499c4762a1bSJed Brown gamma[0] = hf*(hf+hf)/(hf+hs)/(hf+hf+hs); 500c4762a1bSJed Brown for (j=0; j<dof; j++) { 501c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 502c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 503c4762a1bSJed Brown u[j] = x[(i-1)*dof+j]+PetscMax(0,PetscMin(min[j],alpha[0]+gamma[0]*r[j]))*(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 504c4762a1bSJed Brown } 5059566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 506c4762a1bSJed Brown if (i > xs) { 507c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hf; 508c4762a1bSJed Brown } 509c4762a1bSJed Brown } 510c4762a1bSJed Brown } 5119566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 5129566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 5139566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 5149566063dSJacob Faibussowitsch PetscCall(PetscFree4(r,min,alpha,gamma)); 515c4762a1bSJed Brown PetscFunctionReturn(0); 516c4762a1bSJed Brown } 517c4762a1bSJed Brown 518c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */ 519c4762a1bSJed Brown 520c4762a1bSJed Brown PetscErrorCode FVSample(FVCtx *ctx,DM da,PetscReal time,Vec U) 521c4762a1bSJed Brown { 522c4762a1bSJed Brown PetscScalar *u,*uj,xj,xi; 523c4762a1bSJed Brown PetscInt i,j,k,dof,xs,xm,Mx,count_slow,count_fast; 524c4762a1bSJed Brown const PetscInt N=200; 525c4762a1bSJed Brown 526c4762a1bSJed Brown PetscFunctionBeginUser; 5273c633725SBarry Smith PetscCheck(ctx->physics.sample,PETSC_COMM_SELF,PETSC_ERR_SUP,"Physics has not provided a sampling function"); 5289566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0)); 5299566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 5309566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,U,&u)); 5319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof,&uj)); 532c4762a1bSJed Brown const PetscReal hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; 533c4762a1bSJed Brown const PetscReal hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 534c4762a1bSJed Brown count_slow = Mx/(1+ctx->hratio); 535c4762a1bSJed Brown count_fast = Mx-count_slow; 536c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 537c4762a1bSJed Brown if (i*hs+0.5*hs<(ctx->xmax-ctx->xmin)*0.25) { 538c4762a1bSJed Brown xi = ctx->xmin+0.5*hs+i*hs; 539c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 540c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 541c4762a1bSJed Brown for (j=0; j<N+1; j++) { 542c4762a1bSJed Brown xj = xi+hs*(j-N/2)/(PetscReal)N; 5439566063dSJacob Faibussowitsch PetscCall((*ctx->physics.sample)(ctx->physics.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj)); 544c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 545c4762a1bSJed Brown } 546c4762a1bSJed Brown } else if ((ctx->xmax-ctx->xmin)*0.25+(i-count_slow/2)*hf+0.5*hf<(ctx->xmax-ctx->xmin)*0.75) { 547c4762a1bSJed Brown xi = ctx->xmin+(ctx->xmax-ctx->xmin)*0.25+0.5*hf+(i-count_slow/2)*hf; 548c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 549c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 550c4762a1bSJed Brown for (j=0; j<N+1; j++) { 551c4762a1bSJed Brown xj = xi+hf*(j-N/2)/(PetscReal)N; 5529566063dSJacob Faibussowitsch PetscCall((*ctx->physics.sample)(ctx->physics.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj)); 553c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 554c4762a1bSJed Brown } 555c4762a1bSJed Brown } else { 556c4762a1bSJed Brown xi = ctx->xmin+(ctx->xmax-ctx->xmin)*0.75+0.5*hs+(i-count_slow/2-count_fast)*hs; 557c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 558c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 559c4762a1bSJed Brown for (j=0; j<N+1; j++) { 560c4762a1bSJed Brown xj = xi+hs*(j-N/2)/(PetscReal)N; 5619566063dSJacob Faibussowitsch PetscCall((*ctx->physics.sample)(ctx->physics.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj)); 562c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 563c4762a1bSJed Brown } 564c4762a1bSJed Brown } 565c4762a1bSJed Brown } 5669566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,U,&u)); 5679566063dSJacob Faibussowitsch PetscCall(PetscFree(uj)); 568c4762a1bSJed Brown PetscFunctionReturn(0); 569c4762a1bSJed Brown } 570c4762a1bSJed Brown 571c4762a1bSJed Brown static PetscErrorCode SolutionStatsView(DM da,Vec X,PetscViewer viewer) 572c4762a1bSJed Brown { 573c4762a1bSJed Brown PetscReal xmin,xmax; 574c4762a1bSJed Brown PetscScalar sum,tvsum,tvgsum; 575c4762a1bSJed Brown const PetscScalar *x; 576c4762a1bSJed Brown PetscInt imin,imax,Mx,i,j,xs,xm,dof; 577c4762a1bSJed Brown Vec Xloc; 578c4762a1bSJed Brown PetscBool iascii; 579c4762a1bSJed Brown 580c4762a1bSJed Brown PetscFunctionBeginUser; 5819566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 582c4762a1bSJed Brown if (iascii) { 583c4762a1bSJed Brown /* PETSc lacks a function to compute total variation norm (difficult in multiple dimensions), we do it here */ 5849566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); 5859566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 5869566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 5879566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,Xloc,(void*)&x)); 5889566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 5899566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0)); 590c4762a1bSJed Brown tvsum = 0; 591c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 592c4762a1bSJed Brown for (j=0; j<dof; j++) tvsum += PetscAbsScalar(x[i*dof+j]-x[(i-1)*dof+j]); 593c4762a1bSJed Brown } 5949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&tvsum,&tvgsum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)da))); 5959566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,Xloc,(void*)&x)); 5969566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 597c4762a1bSJed Brown 5989566063dSJacob Faibussowitsch PetscCall(VecMin(X,&imin,&xmin)); 5999566063dSJacob Faibussowitsch PetscCall(VecMax(X,&imax,&xmax)); 6009566063dSJacob Faibussowitsch PetscCall(VecSum(X,&sum)); 60163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Solution range [%g,%g] with minimum at %" PetscInt_FMT ", mean %g, ||x||_TV %g\n",(double)xmin,(double)xmax,imin,(double)(sum/Mx),(double)(tvgsum/Mx))); 602c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type not supported"); 603c4762a1bSJed Brown PetscFunctionReturn(0); 604c4762a1bSJed Brown } 605c4762a1bSJed Brown 606c4762a1bSJed Brown static PetscErrorCode SolutionErrorNorms(FVCtx *ctx,DM da,PetscReal t,Vec X,PetscReal *nrm1) 607c4762a1bSJed Brown { 608c4762a1bSJed Brown Vec Y; 609c4762a1bSJed Brown PetscInt i,Mx,count_slow=0,count_fast=0; 610c4762a1bSJed Brown const PetscScalar *ptr_X,*ptr_Y; 611c4762a1bSJed Brown 612c4762a1bSJed Brown PetscFunctionBeginUser; 6139566063dSJacob Faibussowitsch PetscCall(VecGetSize(X,&Mx)); 6149566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&Y)); 6159566063dSJacob Faibussowitsch PetscCall(FVSample(ctx,da,t,Y)); 616c4762a1bSJed Brown const PetscReal hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; 617c4762a1bSJed Brown const PetscReal hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 618c4762a1bSJed Brown count_slow = (PetscReal)Mx/(1.0+ctx->hratio); 619c4762a1bSJed Brown count_fast = Mx-count_slow; 6209566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&ptr_X)); 6219566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&ptr_Y)); 622c4762a1bSJed Brown for (i=0; i<Mx; i++) { 623c4762a1bSJed Brown if (i < count_slow/2 || i > count_slow/2+count_fast-1) *nrm1 += hs*PetscAbs(ptr_X[i]-ptr_Y[i]); 624c4762a1bSJed Brown else *nrm1 += hf*PetscAbs(ptr_X[i]-ptr_Y[i]); 625c4762a1bSJed Brown } 6269566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&ptr_X)); 6279566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&ptr_Y)); 6289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 629c4762a1bSJed Brown PetscFunctionReturn(0); 630c4762a1bSJed Brown } 631c4762a1bSJed Brown 632c4762a1bSJed Brown int main(int argc,char *argv[]) 633c4762a1bSJed Brown { 634c4762a1bSJed Brown char physname[256] = "advect",final_fname[256] = "solution.m"; 635c4762a1bSJed Brown PetscFunctionList physics = 0; 636c4762a1bSJed Brown MPI_Comm comm; 637c4762a1bSJed Brown TS ts; 638c4762a1bSJed Brown DM da; 639c4762a1bSJed Brown Vec X,X0,R; 640c4762a1bSJed Brown FVCtx ctx; 641c4762a1bSJed Brown PetscInt i,k,dof,xs,xm,Mx,draw = 0,count_slow,count_fast,islow = 0,ifast = 0,*index_slow,*index_fast; 642c4762a1bSJed Brown PetscBool view_final = PETSC_FALSE; 643c4762a1bSJed Brown PetscReal ptime; 644c4762a1bSJed Brown 645*327415f7SBarry Smith PetscFunctionBeginUser; 6469566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,0,help)); 647c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 6489566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&ctx,sizeof(ctx))); 649c4762a1bSJed Brown 650c4762a1bSJed Brown /* Register physical models to be available on the command line */ 6519566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&physics,"advect",PhysicsCreate_Advect)); 652c4762a1bSJed Brown 653c4762a1bSJed Brown ctx.comm = comm; 654c4762a1bSJed Brown ctx.cfl = 0.9; 655c4762a1bSJed Brown ctx.bctype = FVBC_PERIODIC; 656c4762a1bSJed Brown ctx.xmin = -1.0; 657c4762a1bSJed Brown ctx.xmax = 1.0; 658d0609cedSBarry Smith PetscOptionsBegin(comm,NULL,"Finite Volume solver options",""); 6599566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL)); 6609566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL)); 6619566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL)); 6629566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-view_final","Write final solution in ASCII MATLAB format to given file name","",final_fname,final_fname,sizeof(final_fname),&view_final)); 6639566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL)); 6649566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-exact","Compare errors with exact solution","",ctx.exact,&ctx.exact,NULL)); 6659566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL)); 6669566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL)); 6679566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL)); 6689566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-hratio","Spacing ratio","",ctx.hratio,&ctx.hratio,NULL)); 669d0609cedSBarry Smith PetscOptionsEnd(); 670c4762a1bSJed Brown 671c4762a1bSJed Brown /* Choose the physics from the list of registered models */ 672c4762a1bSJed Brown { 673c4762a1bSJed Brown PetscErrorCode (*r)(FVCtx*); 6749566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(physics,physname,&r)); 6753c633725SBarry Smith PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Physics '%s' not found",physname); 676c4762a1bSJed Brown /* Create the physics, will set the number of fields and their names */ 6779566063dSJacob Faibussowitsch PetscCall((*r)(&ctx)); 678c4762a1bSJed Brown } 679c4762a1bSJed Brown 680c4762a1bSJed Brown /* Create a DMDA to manage the parallel grid */ 6819566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics.dof,2,NULL,&da)); 6829566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 6839566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 684c4762a1bSJed Brown /* Inform the DMDA of the field names provided by the physics. */ 685c4762a1bSJed Brown /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ 686c4762a1bSJed Brown for (i=0; i<ctx.physics.dof; i++) { 6879566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da,i,ctx.physics.fieldname[i])); 688c4762a1bSJed Brown } 6899566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0)); 6909566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 691c4762a1bSJed Brown 692c4762a1bSJed Brown /* Set coordinates of cell centers */ 6939566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da,ctx.xmin+0.5*(ctx.xmax-ctx.xmin)/Mx,ctx.xmax+0.5*(ctx.xmax-ctx.xmin)/Mx,0,0,0,0)); 694c4762a1bSJed Brown 695c4762a1bSJed Brown /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */ 6969566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(dof,&ctx.u,dof,&ctx.flux,dof,&ctx.speeds)); 697c4762a1bSJed Brown 698c4762a1bSJed Brown /* Create a vector to store the solution and to save the initial state */ 6999566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da,&X)); 7009566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&X0)); 7019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&R)); 702c4762a1bSJed Brown 703c4762a1bSJed Brown /* create index for slow parts and fast parts*/ 704c4762a1bSJed Brown count_slow = Mx/(1+ctx.hratio); 70508401ef6SPierre Jolivet PetscCheck(count_slow%2 == 0,PETSC_COMM_WORLD,PETSC_ERR_USER,"Please adjust grid size Mx (-da_grid_x) and hratio (-hratio) so that Mx/(1+hartio) is even"); 706c4762a1bSJed Brown count_fast = Mx-count_slow; 707c4762a1bSJed Brown ctx.sf = count_slow/2; 708c4762a1bSJed Brown ctx.fs = ctx.sf + count_fast; 7099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm*dof,&index_slow)); 7109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm*dof,&index_fast)); 711c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 712c4762a1bSJed Brown if (i < count_slow/2 || i > count_slow/2+count_fast-1) 713c4762a1bSJed Brown for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k; 714c4762a1bSJed Brown else 715c4762a1bSJed Brown for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k; 716c4762a1bSJed Brown } 7179566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss)); 7189566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf)); 719c4762a1bSJed Brown 720c4762a1bSJed Brown /* Create a time-stepping object */ 7219566063dSJacob Faibussowitsch PetscCall(TSCreate(comm,&ts)); 7229566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts,da)); 7239566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,R,FVRHSFunction,&ctx)); 7249566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"slow",ctx.iss)); 7259566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"fast",ctx.isf)); 7269566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"slow",NULL,FVRHSFunctionslow,&ctx)); 7279566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"fast",NULL,FVRHSFunctionfast,&ctx)); 728c4762a1bSJed Brown 7299566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSMPRK)); 7309566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,10)); 7319566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 732c4762a1bSJed Brown 733c4762a1bSJed Brown /* Compute initial conditions and starting time step */ 7349566063dSJacob Faibussowitsch PetscCall(FVSample(&ctx,da,0,X0)); 7359566063dSJacob Faibussowitsch PetscCall(FVRHSFunction(ts,0,X0,X,(void*)&ctx)); /* Initial function evaluation, only used to determine max speed */ 7369566063dSJacob Faibussowitsch PetscCall(VecCopy(X0,X)); /* The function value was not used so we set X=X0 again */ 7379566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt)); 7389566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); /* Take runtime options */ 7399566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD)); 740c4762a1bSJed Brown { 741c4762a1bSJed Brown PetscInt steps; 742c4762a1bSJed Brown PetscScalar mass_initial,mass_final,mass_difference,mass_differenceg; 743c4762a1bSJed Brown const PetscScalar *ptr_X,*ptr_X0; 744c2a16db8SHong Zhang const PetscReal hs = (ctx.xmax-ctx.xmin)/2.0/count_slow; 745c2a16db8SHong Zhang const PetscReal hf = (ctx.xmax-ctx.xmin)/2.0/count_fast; 7469566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,X)); 7479566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&ptime)); 7489566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&steps)); 749c4762a1bSJed Brown /* calculate the total mass at initial time and final time */ 750c4762a1bSJed Brown mass_initial = 0.0; 751c4762a1bSJed Brown mass_final = 0.0; 7529566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,X0,(void*)&ptr_X0)); 7539566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,X,(void*)&ptr_X)); 754c2a16db8SHong Zhang for (i=xs; i<xs+xm; i++) { 755c2a16db8SHong Zhang if (i < ctx.sf || i > ctx.fs-1) { 756c2a16db8SHong Zhang for (k=0; k<dof; k++) { 757c2a16db8SHong Zhang mass_initial = mass_initial+hs*ptr_X0[i*dof+k]; 758c2a16db8SHong Zhang mass_final = mass_final+hs*ptr_X[i*dof+k]; 759c2a16db8SHong Zhang } 760c4762a1bSJed Brown } else { 761c2a16db8SHong Zhang for (k=0; k<dof; k++) { 762c2a16db8SHong Zhang mass_initial = mass_initial+hf*ptr_X0[i*dof+k]; 763c2a16db8SHong Zhang mass_final = mass_final+hf*ptr_X[i*dof+k]; 764c2a16db8SHong Zhang } 765c4762a1bSJed Brown } 766c4762a1bSJed Brown } 7679566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,X0,(void*)&ptr_X0)); 7689566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,X,(void*)&ptr_X)); 769c4762a1bSJed Brown mass_difference = mass_final-mass_initial; 7709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&mass_difference,&mass_differenceg,1,MPIU_SCALAR,MPIU_SUM,comm)); 7719566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Mass difference %g\n",(double)mass_differenceg)); 77263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Final time %g, steps %" PetscInt_FMT "\n",(double)ptime,steps)); 773c4762a1bSJed Brown if (ctx.exact) { 774c4762a1bSJed Brown PetscReal nrm1 = 0; 7759566063dSJacob Faibussowitsch PetscCall(SolutionErrorNorms(&ctx,da,ptime,X,&nrm1)); 7769566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1)); 777c4762a1bSJed Brown } 778c4762a1bSJed Brown if (ctx.simulation) { 779c4762a1bSJed Brown PetscReal nrm1 = 0; 780c4762a1bSJed Brown PetscViewer fd; 781c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = "binaryoutput"; 782c4762a1bSJed Brown Vec XR; 783c4762a1bSJed Brown PetscBool flg; 784c4762a1bSJed Brown const PetscScalar *ptr_XR; 7859566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg)); 7863c633725SBarry Smith PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 7879566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&fd)); 7889566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X0,&XR)); 7899566063dSJacob Faibussowitsch PetscCall(VecLoad(XR,fd)); 7909566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 7919566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&ptr_X)); 7929566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XR,&ptr_XR)); 793c4762a1bSJed Brown for (i=0; i<Mx; i++) { 794c4762a1bSJed Brown if (i < count_slow/2 || i > count_slow/2+count_fast-1) nrm1 = nrm1 + hs*PetscAbs(ptr_X[i]-ptr_XR[i]); 795c4762a1bSJed Brown else nrm1 = nrm1 + hf*PetscAbs(ptr_X[i]-ptr_XR[i]); 796c4762a1bSJed Brown } 7979566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&ptr_X)); 7989566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XR,&ptr_XR)); 7999566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1)); 8009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XR)); 801c4762a1bSJed Brown } 802c4762a1bSJed Brown } 803c4762a1bSJed Brown 8049566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD)); 8059566063dSJacob Faibussowitsch if (draw & 0x1) PetscCall(VecView(X0,PETSC_VIEWER_DRAW_WORLD)); 8069566063dSJacob Faibussowitsch if (draw & 0x2) PetscCall(VecView(X,PETSC_VIEWER_DRAW_WORLD)); 807c4762a1bSJed Brown if (draw & 0x4) { 808c4762a1bSJed Brown Vec Y; 8099566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&Y)); 8109566063dSJacob Faibussowitsch PetscCall(FVSample(&ctx,da,ptime,Y)); 8119566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y,-1,X)); 8129566063dSJacob Faibussowitsch PetscCall(VecView(Y,PETSC_VIEWER_DRAW_WORLD)); 8139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 814c4762a1bSJed Brown } 815c4762a1bSJed Brown 816c4762a1bSJed Brown if (view_final) { 817c4762a1bSJed Brown PetscViewer viewer; 8189566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer)); 8199566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB)); 8209566063dSJacob Faibussowitsch PetscCall(VecView(X,viewer)); 8219566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 8229566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 823c4762a1bSJed Brown } 824c4762a1bSJed Brown 825c4762a1bSJed Brown /* Clean up */ 8269566063dSJacob Faibussowitsch PetscCall((*ctx.physics.destroy)(ctx.physics.user)); 8279566063dSJacob Faibussowitsch for (i=0; i<ctx.physics.dof; i++) PetscCall(PetscFree(ctx.physics.fieldname[i])); 8289566063dSJacob Faibussowitsch PetscCall(PetscFree3(ctx.u,ctx.flux,ctx.speeds)); 8299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.iss)); 8309566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.isf)); 8319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 8329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X0)); 8339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&R)); 8349566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 8359566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 8369566063dSJacob Faibussowitsch PetscCall(PetscFree(index_slow)); 8379566063dSJacob Faibussowitsch PetscCall(PetscFree(index_fast)); 8389566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&physics)); 8399566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 840b122ec5aSJacob Faibussowitsch return 0; 841c4762a1bSJed Brown } 842c4762a1bSJed Brown 843c4762a1bSJed Brown /*TEST 844c4762a1bSJed Brown 845c4762a1bSJed Brown build: 846f56ea12dSJed Brown requires: !complex 847c4762a1bSJed Brown 848c4762a1bSJed Brown test: 849c4762a1bSJed Brown args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 0 850c4762a1bSJed Brown 851c4762a1bSJed Brown test: 852c4762a1bSJed Brown suffix: 2 853c4762a1bSJed Brown args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 1 854c4762a1bSJed Brown output_file: output/ex7_1.out 855c4762a1bSJed Brown 856c4762a1bSJed Brown test: 857c4762a1bSJed Brown suffix: 3 858c4762a1bSJed Brown args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 0 859c4762a1bSJed Brown 860c4762a1bSJed Brown test: 861c4762a1bSJed Brown suffix: 4 862c4762a1bSJed Brown args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 1 863c4762a1bSJed Brown output_file: output/ex7_3.out 864c2a16db8SHong Zhang 865c2a16db8SHong Zhang test: 866c2a16db8SHong Zhang suffix: 5 867c2a16db8SHong Zhang nsize: 2 868c2a16db8SHong Zhang args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 1 869c2a16db8SHong Zhang output_file: output/ex7_3.out 870c4762a1bSJed Brown TEST*/ 871