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" 7c4762a1bSJed Brown " with x belongs to (xmin,xmax), the number of total mesh points is Mx and the ratio between the meshsize of corse\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 PetscErrorCode ierr; 120c4762a1bSJed Brown AdvectCtx *user; 121c4762a1bSJed Brown 122c4762a1bSJed Brown PetscFunctionBeginUser; 1239566063dSJacob Faibussowitsch PetscCall(PetscNew(&user)); 124c4762a1bSJed Brown ctx->physics.sample = PhysicsSample_Advect; 125c4762a1bSJed Brown ctx->physics.flux = PhysicsFlux_Advect; 126c4762a1bSJed Brown ctx->physics.destroy = PhysicsDestroy_SimpleFree; 127c4762a1bSJed Brown ctx->physics.user = user; 128c4762a1bSJed Brown ctx->physics.dof = 1; 1299566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("u",&ctx->physics.fieldname[0])); 130c4762a1bSJed Brown user->a = 1; 1319566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");PetscCall(ierr); 132c4762a1bSJed Brown { 1339566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-physics_advect_a","Speed","",user->a,&user->a,NULL)); 134c4762a1bSJed Brown } 1359566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 136c4762a1bSJed Brown PetscFunctionReturn(0); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver ----------------------------------- */ 140c4762a1bSJed Brown 141c4762a1bSJed Brown static PetscErrorCode FVRHSFunction(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 142c4762a1bSJed Brown { 143c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 144c4762a1bSJed Brown PetscInt i,j,Mx,dof,xs,xm,sf = ctx->sf,fs = ctx->fs; 145c4762a1bSJed Brown PetscReal hf,hs,cfl_idt = 0; 146c4762a1bSJed Brown PetscScalar *x,*f,*r,*min,*alpha,*gamma; 147c4762a1bSJed Brown Vec Xloc; 148c4762a1bSJed Brown DM da; 149c4762a1bSJed Brown 150c4762a1bSJed Brown PetscFunctionBeginUser; 1519566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 1529566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); /* Xloc contains ghost points */ 1539566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0)); /* Mx is the number of center points */ 154c4762a1bSJed Brown hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; 155c4762a1bSJed Brown hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 1569566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); /* X is solution vector which does not contain ghost points */ 1579566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 1589566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); /* F is the right hand side function corresponds to center points */ 1599566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 1609566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,F,&f)); 1619566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 1629566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dof,&r,dof,&min,dof,&alpha,dof,&gamma)); 163c4762a1bSJed Brown 164c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 165c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 166c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 167c4762a1bSJed Brown } 168c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 169c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 170c4762a1bSJed Brown } 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 174c4762a1bSJed Brown PetscReal maxspeed; 175c4762a1bSJed Brown PetscScalar *u; 176c4762a1bSJed Brown if (i < sf || i > fs+1) { 177c4762a1bSJed Brown u = &ctx->u[0]; 178c4762a1bSJed Brown alpha[0] = 1.0/6.0; 179c4762a1bSJed Brown gamma[0] = 1.0/3.0; 180c4762a1bSJed Brown for (j=0; j<dof; j++) { 181c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 182c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 183c4762a1bSJed 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]); 184c4762a1bSJed Brown } 1859566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 186c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hs)); 187c4762a1bSJed Brown if (i > xs) { 188c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hs; 189c4762a1bSJed Brown } 190c4762a1bSJed Brown if (i < xs+xm) { 191c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hs; 192c4762a1bSJed Brown } 193c4762a1bSJed Brown } else if (i == sf) { 194c4762a1bSJed Brown u = &ctx->u[0]; 195c4762a1bSJed Brown alpha[0] = hs*hf/(hs+hs)/(hs+hs+hf); 196c4762a1bSJed Brown gamma[0] = hs*(hs+hs)/(hs+hf)/(hs+hs+hf); 197c4762a1bSJed Brown for (j=0; j<dof; j++) { 198c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 199c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 200c4762a1bSJed 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]); 201c4762a1bSJed Brown } 2029566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 203c4762a1bSJed Brown if (i > xs) { 204c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hs; 205c4762a1bSJed Brown } 206c4762a1bSJed Brown if (i < xs+xm) { 207c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hf; 208c4762a1bSJed Brown } 209c4762a1bSJed Brown } else if (i == sf+1) { 210c4762a1bSJed Brown u = &ctx->u[0]; 211c4762a1bSJed Brown alpha[0] = hf*hf/(hs+hf)/(hs+hf+hf); 212c4762a1bSJed Brown gamma[0] = hf*(hs+hf)/(hf+hf)/(hs+hf+hf); 213c4762a1bSJed Brown for (j=0; j<dof; j++) { 214c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 215c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 216c4762a1bSJed 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]); 217c4762a1bSJed Brown } 2189566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 219c4762a1bSJed Brown if (i > xs) { 220c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hf; 221c4762a1bSJed Brown } 222c4762a1bSJed Brown if (i < xs+xm) { 223c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hf; 224c4762a1bSJed Brown } 225c4762a1bSJed Brown } else if (i > sf+1 && i < fs) { 226c4762a1bSJed Brown u = &ctx->u[0]; 227c4762a1bSJed Brown alpha[0] = 1.0/6.0; 228c4762a1bSJed Brown gamma[0] = 1.0/3.0; 229c4762a1bSJed Brown for (j=0; j<dof; j++) { 230c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 231c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 232c4762a1bSJed 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]); 233c4762a1bSJed Brown } 2349566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 235c4762a1bSJed Brown if (i > xs) { 236c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hf; 237c4762a1bSJed Brown } 238c4762a1bSJed Brown if (i < xs+xm) { 239c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hf; 240c4762a1bSJed Brown } 241c4762a1bSJed Brown } else if (i == fs) { 242c4762a1bSJed Brown u = &ctx->u[0]; 243c4762a1bSJed Brown alpha[0] = hf*hs/(hf+hf)/(hf+hf+hs); 244c4762a1bSJed Brown gamma[0] = hf*(hf+hf)/(hf+hs)/(hf+hf+hs); 245c4762a1bSJed Brown for (j=0; j<dof; j++) { 246c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 247c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 248c4762a1bSJed 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]); 249c4762a1bSJed Brown } 2509566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 251c4762a1bSJed Brown if (i > xs) { 252c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hf; 253c4762a1bSJed Brown } 254c4762a1bSJed Brown if (i < xs+xm) { 255c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hs; 256c4762a1bSJed Brown } 257c4762a1bSJed Brown } else if (i == fs+1) { 258c4762a1bSJed Brown u = &ctx->u[0]; 259c4762a1bSJed Brown alpha[0] = hs*hs/(hf+hs)/(hf+hs+hs); 260c4762a1bSJed Brown gamma[0] = hs*(hf+hs)/(hs+hs)/(hf+hs+hs); 261c4762a1bSJed Brown for (j=0; j<dof; j++) { 262c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 263c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 264c4762a1bSJed 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]); 265c4762a1bSJed Brown } 2669566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 267c4762a1bSJed Brown if (i > xs) { 268c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hs; 269c4762a1bSJed Brown } 270c4762a1bSJed Brown if (i < xs+xm) { 271c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hs; 272c4762a1bSJed Brown } 273c4762a1bSJed Brown } 274c4762a1bSJed Brown } 2759566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 2769566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,F,&f)); 2779566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 2789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da))); 279c4762a1bSJed Brown if (0) { 280c2a16db8SHong Zhang /* We need a way to inform the TS of a CFL constraint, this is a debugging fragment */ 281c4762a1bSJed Brown PetscReal dt,tnow; 2829566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&dt)); 2839566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts,&tnow)); 284c4762a1bSJed Brown if (dt > 0.5/ctx->cfl_idt) { 2859566063dSJacob 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))); 286c4762a1bSJed Brown } 287c4762a1bSJed Brown } 2889566063dSJacob Faibussowitsch PetscCall(PetscFree4(r,min,alpha,gamma)); 289c4762a1bSJed Brown PetscFunctionReturn(0); 290c4762a1bSJed Brown } 291c4762a1bSJed Brown 292c4762a1bSJed Brown static PetscErrorCode FVRHSFunctionslow(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 293c4762a1bSJed Brown { 294c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 295c4762a1bSJed Brown PetscInt i,j,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs; 296c4762a1bSJed Brown PetscReal hf,hs; 297c4762a1bSJed Brown PetscScalar *x,*f,*r,*min,*alpha,*gamma; 298c4762a1bSJed Brown Vec Xloc; 299c4762a1bSJed Brown DM da; 300c4762a1bSJed Brown 301c4762a1bSJed Brown PetscFunctionBeginUser; 3029566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 3039566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); /* Xloc contains ghost points */ 3049566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0)); /* Mx is the number of center points */ 305c4762a1bSJed Brown hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; 306c4762a1bSJed Brown hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 3079566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); /* X is solution vector which does not contain ghost points */ 3089566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd (da,X,INSERT_VALUES,Xloc)); 3099566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); /* F is the right hand side function corresponds to center points */ 3109566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 3119566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 3129566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 3139566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dof,&r,dof,&min,dof,&alpha,dof,&gamma)); 314c4762a1bSJed Brown 315c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 316c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 317c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 318c4762a1bSJed Brown } 319c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 320c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 321c4762a1bSJed Brown } 322c4762a1bSJed Brown } 323c4762a1bSJed Brown 324c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 325c4762a1bSJed Brown PetscReal maxspeed; 326c4762a1bSJed Brown PetscScalar *u; 327c4762a1bSJed Brown if (i < sf) { 328c4762a1bSJed Brown u = &ctx->u[0]; 329c4762a1bSJed Brown alpha[0] = 1.0/6.0; 330c4762a1bSJed Brown gamma[0] = 1.0/3.0; 331c4762a1bSJed Brown for (j=0; j<dof; j++) { 332c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 333c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 334c4762a1bSJed 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]); 335c4762a1bSJed Brown } 3369566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 337c4762a1bSJed Brown if (i > xs) { 338c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hs; 339c4762a1bSJed Brown } 340c4762a1bSJed Brown if (i < xs+xm) { 341c2a16db8SHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hs; 342c4762a1bSJed Brown islow++; 343c4762a1bSJed Brown } 344c4762a1bSJed Brown } else if (i == sf) { 345c4762a1bSJed Brown u = &ctx->u[0]; 346c4762a1bSJed Brown alpha[0] = hs*hf/(hs+hs)/(hs+hs+hf); 347c4762a1bSJed Brown gamma[0] = hs*(hs+hs)/(hs+hf)/(hs+hs+hf); 348c4762a1bSJed Brown for (j=0; j<dof; j++) { 349c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 350c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 351c4762a1bSJed 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]); 352c4762a1bSJed Brown } 3539566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 354c2a16db8SHong Zhang if (i > xs) { 355c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hs; 356c4762a1bSJed Brown } 357c4762a1bSJed Brown } else if (i == fs) { 358c4762a1bSJed Brown u = &ctx->u[0]; 359c4762a1bSJed Brown alpha[0] = hf*hs/(hf+hf)/(hf+hf+hs); 360c4762a1bSJed Brown gamma[0] = hf*(hf+hf)/(hf+hs)/(hf+hf+hs); 361c4762a1bSJed Brown for (j=0; j<dof; j++) { 362c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 363c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 364c4762a1bSJed 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]); 365c4762a1bSJed Brown } 3669566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 367c4762a1bSJed Brown if (i < xs+xm) { 368c2a16db8SHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hs; 369c4762a1bSJed Brown islow++; 370c4762a1bSJed Brown } 371c4762a1bSJed Brown } else if (i == fs+1) { 372c4762a1bSJed Brown u = &ctx->u[0]; 373c4762a1bSJed Brown alpha[0] = hs*hs/(hf+hs)/(hf+hs+hs); 374c4762a1bSJed Brown gamma[0] = hs*(hf+hs)/(hs+hs)/(hf+hs+hs); 375c4762a1bSJed Brown for (j=0; j<dof; j++) { 376c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 377c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 378c4762a1bSJed 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]); 379c4762a1bSJed Brown } 3809566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 381c4762a1bSJed Brown if (i > xs) { 382c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hs; 383c4762a1bSJed Brown } 384c4762a1bSJed Brown if (i < xs+xm) { 385c2a16db8SHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hs; 386c2a16db8SHong Zhang islow++; 387c4762a1bSJed Brown } 388c4762a1bSJed Brown } else if (i > fs+1) { 389c4762a1bSJed Brown u = &ctx->u[0]; 390c4762a1bSJed Brown alpha[0] = 1.0/6.0; 391c4762a1bSJed Brown gamma[0] = 1.0/3.0; 392c4762a1bSJed Brown for (j=0; j<dof; j++) { 393c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 394c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 395c4762a1bSJed 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]); 396c4762a1bSJed Brown } 3979566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 398c4762a1bSJed Brown if (i > xs) { 399c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hs; 400c4762a1bSJed Brown } 401c4762a1bSJed Brown if (i < xs+xm) { 402c2a16db8SHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hs; 403c2a16db8SHong Zhang islow++; 404c4762a1bSJed Brown } 405c4762a1bSJed Brown } 406c4762a1bSJed Brown } 4079566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 4089566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 4099566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 4109566063dSJacob Faibussowitsch PetscCall(PetscFree4(r,min,alpha,gamma)); 411c4762a1bSJed Brown PetscFunctionReturn(0); 412c4762a1bSJed Brown } 413c4762a1bSJed Brown 414c4762a1bSJed Brown static PetscErrorCode FVRHSFunctionfast(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 415c4762a1bSJed Brown { 416c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 417c4762a1bSJed Brown PetscInt i,j,Mx,dof,xs,xm,ifast = 0,sf = ctx->sf,fs = ctx->fs; 418c4762a1bSJed Brown PetscReal hf,hs; 419c4762a1bSJed Brown PetscScalar *x,*f,*r,*min,*alpha,*gamma; 420c4762a1bSJed Brown Vec Xloc; 421c4762a1bSJed Brown DM da; 422c4762a1bSJed Brown 423c4762a1bSJed Brown PetscFunctionBeginUser; 4249566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 4259566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); /* Xloc contains ghost points */ 4269566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0)); /* Mx is the number of center points */ 427c4762a1bSJed Brown hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; 428c4762a1bSJed Brown hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 4299566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); /* X is solution vector which does not contain ghost points */ 4309566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 4319566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); /* F is the right hand side function corresponds to center points */ 4329566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 4339566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 4349566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 4359566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dof,&r,dof,&min,dof,&alpha,dof,&gamma)); 436c4762a1bSJed Brown 437c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 438c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 439c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 440c4762a1bSJed Brown } 441c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 442c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 443c4762a1bSJed Brown } 444c4762a1bSJed Brown } 445c4762a1bSJed Brown 446c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 447c4762a1bSJed Brown PetscReal maxspeed; 448c4762a1bSJed Brown PetscScalar *u; 449c4762a1bSJed Brown if (i == sf) { 450c4762a1bSJed Brown u = &ctx->u[0]; 451c4762a1bSJed Brown alpha[0] = hs*hf/(hs+hs)/(hs+hs+hf); 452c4762a1bSJed Brown gamma[0] = hs*(hs+hs)/(hs+hf)/(hs+hs+hf); 453c4762a1bSJed Brown for (j=0; j<dof; j++) { 454c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 455c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 456c4762a1bSJed 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]); 457c4762a1bSJed Brown } 4589566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 459c4762a1bSJed Brown if (i < xs+xm) { 460c2a16db8SHong Zhang for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hf; 461c4762a1bSJed Brown ifast++; 462c4762a1bSJed Brown } 463c4762a1bSJed Brown } else if (i == sf+1) { 464c4762a1bSJed Brown u = &ctx->u[0]; 465c4762a1bSJed Brown alpha[0] = hf*hf/(hs+hf)/(hs+hf+hf); 466c4762a1bSJed Brown gamma[0] = hf*(hs+hf)/(hf+hf)/(hs+hf+hf); 467c4762a1bSJed Brown for (j=0; j<dof; j++) { 468c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 469c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 470c4762a1bSJed 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]); 471c4762a1bSJed Brown } 4729566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 473c4762a1bSJed Brown if (i > xs) { 474c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hf; 475c4762a1bSJed Brown } 476c4762a1bSJed Brown if (i < xs+xm) { 477c2a16db8SHong Zhang for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hf; 478c4762a1bSJed Brown ifast++; 479c4762a1bSJed Brown } 480c4762a1bSJed Brown } else if (i > sf+1 && i < fs) { 481c4762a1bSJed Brown u = &ctx->u[0]; 482c4762a1bSJed Brown alpha[0] = 1.0/6.0; 483c4762a1bSJed Brown gamma[0] = 1.0/3.0; 484c4762a1bSJed Brown for (j=0; j<dof; j++) { 485c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 486c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 487c4762a1bSJed 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]); 488c4762a1bSJed Brown } 4899566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 490c4762a1bSJed Brown if (i > xs) { 491c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hf; 492c4762a1bSJed Brown } 493c4762a1bSJed Brown if (i < xs+xm) { 494c2a16db8SHong Zhang for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hf; 495c4762a1bSJed Brown ifast++; 496c4762a1bSJed Brown } 497c4762a1bSJed Brown } else if (i == fs) { 498c4762a1bSJed Brown u = &ctx->u[0]; 499c4762a1bSJed Brown alpha[0] = hf*hs/(hf+hf)/(hf+hf+hs); 500c4762a1bSJed Brown gamma[0] = hf*(hf+hf)/(hf+hs)/(hf+hf+hs); 501c4762a1bSJed Brown for (j=0; j<dof; j++) { 502c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 503c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 504c4762a1bSJed 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]); 505c4762a1bSJed Brown } 5069566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed)); 507c4762a1bSJed Brown if (i > xs) { 508c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hf; 509c4762a1bSJed Brown } 510c4762a1bSJed Brown } 511c4762a1bSJed Brown } 5129566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 5139566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 5149566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 5159566063dSJacob Faibussowitsch PetscCall(PetscFree4(r,min,alpha,gamma)); 516c4762a1bSJed Brown PetscFunctionReturn(0); 517c4762a1bSJed Brown } 518c4762a1bSJed Brown 519c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */ 520c4762a1bSJed Brown 521c4762a1bSJed Brown PetscErrorCode FVSample(FVCtx *ctx,DM da,PetscReal time,Vec U) 522c4762a1bSJed Brown { 523c4762a1bSJed Brown PetscScalar *u,*uj,xj,xi; 524c4762a1bSJed Brown PetscInt i,j,k,dof,xs,xm,Mx,count_slow,count_fast; 525c4762a1bSJed Brown const PetscInt N=200; 526c4762a1bSJed Brown 527c4762a1bSJed Brown PetscFunctionBeginUser; 5283c633725SBarry Smith PetscCheck(ctx->physics.sample,PETSC_COMM_SELF,PETSC_ERR_SUP,"Physics has not provided a sampling function"); 5299566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0)); 5309566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 5319566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,U,&u)); 5329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof,&uj)); 533c4762a1bSJed Brown const PetscReal hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; 534c4762a1bSJed Brown const PetscReal hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 535c4762a1bSJed Brown count_slow = Mx/(1+ctx->hratio); 536c4762a1bSJed Brown count_fast = Mx-count_slow; 537c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 538c4762a1bSJed Brown if (i*hs+0.5*hs<(ctx->xmax-ctx->xmin)*0.25) { 539c4762a1bSJed Brown xi = ctx->xmin+0.5*hs+i*hs; 540c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 541c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 542c4762a1bSJed Brown for (j=0; j<N+1; j++) { 543c4762a1bSJed Brown xj = xi+hs*(j-N/2)/(PetscReal)N; 5449566063dSJacob Faibussowitsch PetscCall((*ctx->physics.sample)(ctx->physics.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj)); 545c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 546c4762a1bSJed Brown } 547c4762a1bSJed Brown } else if ((ctx->xmax-ctx->xmin)*0.25+(i-count_slow/2)*hf+0.5*hf<(ctx->xmax-ctx->xmin)*0.75) { 548c4762a1bSJed Brown xi = ctx->xmin+(ctx->xmax-ctx->xmin)*0.25+0.5*hf+(i-count_slow/2)*hf; 549c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 550c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 551c4762a1bSJed Brown for (j=0; j<N+1; j++) { 552c4762a1bSJed Brown xj = xi+hf*(j-N/2)/(PetscReal)N; 5539566063dSJacob Faibussowitsch PetscCall((*ctx->physics.sample)(ctx->physics.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj)); 554c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 555c4762a1bSJed Brown } 556c4762a1bSJed Brown } else { 557c4762a1bSJed Brown xi = ctx->xmin+(ctx->xmax-ctx->xmin)*0.75+0.5*hs+(i-count_slow/2-count_fast)*hs; 558c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 559c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 560c4762a1bSJed Brown for (j=0; j<N+1; j++) { 561c4762a1bSJed Brown xj = xi+hs*(j-N/2)/(PetscReal)N; 5629566063dSJacob Faibussowitsch PetscCall((*ctx->physics.sample)(ctx->physics.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj)); 563c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 564c4762a1bSJed Brown } 565c4762a1bSJed Brown } 566c4762a1bSJed Brown } 5679566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,U,&u)); 5689566063dSJacob Faibussowitsch PetscCall(PetscFree(uj)); 569c4762a1bSJed Brown PetscFunctionReturn(0); 570c4762a1bSJed Brown } 571c4762a1bSJed Brown 572c4762a1bSJed Brown static PetscErrorCode SolutionStatsView(DM da,Vec X,PetscViewer viewer) 573c4762a1bSJed Brown { 574c4762a1bSJed Brown PetscReal xmin,xmax; 575c4762a1bSJed Brown PetscScalar sum,tvsum,tvgsum; 576c4762a1bSJed Brown const PetscScalar *x; 577c4762a1bSJed Brown PetscInt imin,imax,Mx,i,j,xs,xm,dof; 578c4762a1bSJed Brown Vec Xloc; 579c4762a1bSJed Brown PetscBool iascii; 580c4762a1bSJed Brown 581c4762a1bSJed Brown PetscFunctionBeginUser; 5829566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 583c4762a1bSJed Brown if (iascii) { 584c4762a1bSJed Brown /* PETSc lacks a function to compute total variation norm (difficult in multiple dimensions), we do it here */ 5859566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); 5869566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 5879566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 5889566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,Xloc,(void*)&x)); 5899566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 5909566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0)); 591c4762a1bSJed Brown tvsum = 0; 592c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 593c4762a1bSJed Brown for (j=0; j<dof; j++) tvsum += PetscAbsScalar(x[i*dof+j]-x[(i-1)*dof+j]); 594c4762a1bSJed Brown } 5959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&tvsum,&tvgsum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)da))); 5969566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,Xloc,(void*)&x)); 5979566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 598c4762a1bSJed Brown 5999566063dSJacob Faibussowitsch PetscCall(VecMin(X,&imin,&xmin)); 6009566063dSJacob Faibussowitsch PetscCall(VecMax(X,&imax,&xmax)); 6019566063dSJacob Faibussowitsch PetscCall(VecSum(X,&sum)); 6029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Solution range [%g,%g] with minimum at %D, mean %g, ||x||_TV %g\n",(double)xmin,(double)xmax,imin,(double)(sum/Mx),(double)(tvgsum/Mx))); 603c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type not supported"); 604c4762a1bSJed Brown PetscFunctionReturn(0); 605c4762a1bSJed Brown } 606c4762a1bSJed Brown 607c4762a1bSJed Brown static PetscErrorCode SolutionErrorNorms(FVCtx *ctx,DM da,PetscReal t,Vec X,PetscReal *nrm1) 608c4762a1bSJed Brown { 609c4762a1bSJed Brown Vec Y; 610c4762a1bSJed Brown PetscInt i,Mx,count_slow=0,count_fast=0; 611c4762a1bSJed Brown const PetscScalar *ptr_X,*ptr_Y; 612c4762a1bSJed Brown 613c4762a1bSJed Brown PetscFunctionBeginUser; 6149566063dSJacob Faibussowitsch PetscCall(VecGetSize(X,&Mx)); 6159566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&Y)); 6169566063dSJacob Faibussowitsch PetscCall(FVSample(ctx,da,t,Y)); 617c4762a1bSJed Brown const PetscReal hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; 618c4762a1bSJed Brown const PetscReal hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 619c4762a1bSJed Brown count_slow = (PetscReal)Mx/(1.0+ctx->hratio); 620c4762a1bSJed Brown count_fast = Mx-count_slow; 6219566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&ptr_X)); 6229566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&ptr_Y)); 623c4762a1bSJed Brown for (i=0; i<Mx; i++) { 624c4762a1bSJed Brown if (i < count_slow/2 || i > count_slow/2+count_fast-1) *nrm1 += hs*PetscAbs(ptr_X[i]-ptr_Y[i]); 625c4762a1bSJed Brown else *nrm1 += hf*PetscAbs(ptr_X[i]-ptr_Y[i]); 626c4762a1bSJed Brown } 6279566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&ptr_X)); 6289566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&ptr_Y)); 6299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 630c4762a1bSJed Brown PetscFunctionReturn(0); 631c4762a1bSJed Brown } 632c4762a1bSJed Brown 633c4762a1bSJed Brown int main(int argc,char *argv[]) 634c4762a1bSJed Brown { 635c4762a1bSJed Brown char physname[256] = "advect",final_fname[256] = "solution.m"; 636c4762a1bSJed Brown PetscFunctionList physics = 0; 637c4762a1bSJed Brown MPI_Comm comm; 638c4762a1bSJed Brown TS ts; 639c4762a1bSJed Brown DM da; 640c4762a1bSJed Brown Vec X,X0,R; 641c4762a1bSJed Brown FVCtx ctx; 642c4762a1bSJed Brown PetscInt i,k,dof,xs,xm,Mx,draw = 0,count_slow,count_fast,islow = 0,ifast = 0,*index_slow,*index_fast; 643c4762a1bSJed Brown PetscBool view_final = PETSC_FALSE; 644c4762a1bSJed Brown PetscReal ptime; 645c4762a1bSJed Brown PetscErrorCode ierr; 646c4762a1bSJed Brown 6479566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,0,help)); 648c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 6499566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&ctx,sizeof(ctx))); 650c4762a1bSJed Brown 651c4762a1bSJed Brown /* Register physical models to be available on the command line */ 6529566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&physics,"advect",PhysicsCreate_Advect)); 653c4762a1bSJed Brown 654c4762a1bSJed Brown ctx.comm = comm; 655c4762a1bSJed Brown ctx.cfl = 0.9; 656c4762a1bSJed Brown ctx.bctype = FVBC_PERIODIC; 657c4762a1bSJed Brown ctx.xmin = -1.0; 658c4762a1bSJed Brown ctx.xmax = 1.0; 6599566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(comm,NULL,"Finite Volume solver options","");PetscCall(ierr); 6609566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL)); 6619566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL)); 6629566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL)); 6639566063dSJacob 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)); 6649566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL)); 6659566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-exact","Compare errors with exact solution","",ctx.exact,&ctx.exact,NULL)); 6669566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL)); 6679566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL)); 6689566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL)); 6699566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-hratio","Spacing ratio","",ctx.hratio,&ctx.hratio,NULL)); 6709566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 671c4762a1bSJed Brown 672c4762a1bSJed Brown /* Choose the physics from the list of registered models */ 673c4762a1bSJed Brown { 674c4762a1bSJed Brown PetscErrorCode (*r)(FVCtx*); 6759566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(physics,physname,&r)); 6763c633725SBarry Smith PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Physics '%s' not found",physname); 677c4762a1bSJed Brown /* Create the physics, will set the number of fields and their names */ 6789566063dSJacob Faibussowitsch PetscCall((*r)(&ctx)); 679c4762a1bSJed Brown } 680c4762a1bSJed Brown 681c4762a1bSJed Brown /* Create a DMDA to manage the parallel grid */ 6829566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics.dof,2,NULL,&da)); 6839566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 6849566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 685c4762a1bSJed Brown /* Inform the DMDA of the field names provided by the physics. */ 686c4762a1bSJed Brown /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ 687c4762a1bSJed Brown for (i=0; i<ctx.physics.dof; i++) { 6889566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da,i,ctx.physics.fieldname[i])); 689c4762a1bSJed Brown } 6909566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0)); 6919566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 692c4762a1bSJed Brown 693c4762a1bSJed Brown /* Set coordinates of cell centers */ 6949566063dSJacob 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)); 695c4762a1bSJed Brown 696c4762a1bSJed Brown /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */ 6979566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(dof,&ctx.u,dof,&ctx.flux,dof,&ctx.speeds)); 698c4762a1bSJed Brown 699c4762a1bSJed Brown /* Create a vector to store the solution and to save the initial state */ 7009566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da,&X)); 7019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&X0)); 7029566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&R)); 703c4762a1bSJed Brown 704c4762a1bSJed Brown /* create index for slow parts and fast parts*/ 705c4762a1bSJed Brown count_slow = Mx/(1+ctx.hratio); 706*08401ef6SPierre 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"); 707c4762a1bSJed Brown count_fast = Mx-count_slow; 708c4762a1bSJed Brown ctx.sf = count_slow/2; 709c4762a1bSJed Brown ctx.fs = ctx.sf + count_fast; 7109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm*dof,&index_slow)); 7119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm*dof,&index_fast)); 712c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 713c4762a1bSJed Brown if (i < count_slow/2 || i > count_slow/2+count_fast-1) 714c4762a1bSJed Brown for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k; 715c4762a1bSJed Brown else 716c4762a1bSJed Brown for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k; 717c4762a1bSJed Brown } 7189566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss)); 7199566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf)); 720c4762a1bSJed Brown 721c4762a1bSJed Brown /* Create a time-stepping object */ 7229566063dSJacob Faibussowitsch PetscCall(TSCreate(comm,&ts)); 7239566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts,da)); 7249566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,R,FVRHSFunction,&ctx)); 7259566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"slow",ctx.iss)); 7269566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"fast",ctx.isf)); 7279566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"slow",NULL,FVRHSFunctionslow,&ctx)); 7289566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"fast",NULL,FVRHSFunctionfast,&ctx)); 729c4762a1bSJed Brown 7309566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSMPRK)); 7319566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,10)); 7329566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 733c4762a1bSJed Brown 734c4762a1bSJed Brown /* Compute initial conditions and starting time step */ 7359566063dSJacob Faibussowitsch PetscCall(FVSample(&ctx,da,0,X0)); 7369566063dSJacob Faibussowitsch PetscCall(FVRHSFunction(ts,0,X0,X,(void*)&ctx)); /* Initial function evaluation, only used to determine max speed */ 7379566063dSJacob Faibussowitsch PetscCall(VecCopy(X0,X)); /* The function value was not used so we set X=X0 again */ 7389566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt)); 7399566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); /* Take runtime options */ 7409566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD)); 741c4762a1bSJed Brown { 742c4762a1bSJed Brown PetscInt steps; 743c4762a1bSJed Brown PetscScalar mass_initial,mass_final,mass_difference,mass_differenceg; 744c4762a1bSJed Brown const PetscScalar *ptr_X,*ptr_X0; 745c2a16db8SHong Zhang const PetscReal hs = (ctx.xmax-ctx.xmin)/2.0/count_slow; 746c2a16db8SHong Zhang const PetscReal hf = (ctx.xmax-ctx.xmin)/2.0/count_fast; 7479566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,X)); 7489566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&ptime)); 7499566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&steps)); 750c4762a1bSJed Brown /* calculate the total mass at initial time and final time */ 751c4762a1bSJed Brown mass_initial = 0.0; 752c4762a1bSJed Brown mass_final = 0.0; 7539566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,X0,(void*)&ptr_X0)); 7549566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,X,(void*)&ptr_X)); 755c2a16db8SHong Zhang for (i=xs; i<xs+xm; i++) { 756c2a16db8SHong Zhang if (i < ctx.sf || i > ctx.fs-1) { 757c2a16db8SHong Zhang for (k=0; k<dof; k++) { 758c2a16db8SHong Zhang mass_initial = mass_initial+hs*ptr_X0[i*dof+k]; 759c2a16db8SHong Zhang mass_final = mass_final+hs*ptr_X[i*dof+k]; 760c2a16db8SHong Zhang } 761c4762a1bSJed Brown } else { 762c2a16db8SHong Zhang for (k=0; k<dof; k++) { 763c2a16db8SHong Zhang mass_initial = mass_initial+hf*ptr_X0[i*dof+k]; 764c2a16db8SHong Zhang mass_final = mass_final+hf*ptr_X[i*dof+k]; 765c2a16db8SHong Zhang } 766c4762a1bSJed Brown } 767c4762a1bSJed Brown } 7689566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,X0,(void*)&ptr_X0)); 7699566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,X,(void*)&ptr_X)); 770c4762a1bSJed Brown mass_difference = mass_final-mass_initial; 7719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&mass_difference,&mass_differenceg,1,MPIU_SCALAR,MPIU_SUM,comm)); 7729566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Mass difference %g\n",(double)mass_differenceg)); 7739566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Final time %g, steps %D\n",(double)ptime,steps)); 774c4762a1bSJed Brown if (ctx.exact) { 775c4762a1bSJed Brown PetscReal nrm1 = 0; 7769566063dSJacob Faibussowitsch PetscCall(SolutionErrorNorms(&ctx,da,ptime,X,&nrm1)); 7779566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1)); 778c4762a1bSJed Brown } 779c4762a1bSJed Brown if (ctx.simulation) { 780c4762a1bSJed Brown PetscReal nrm1 = 0; 781c4762a1bSJed Brown PetscViewer fd; 782c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = "binaryoutput"; 783c4762a1bSJed Brown Vec XR; 784c4762a1bSJed Brown PetscBool flg; 785c4762a1bSJed Brown const PetscScalar *ptr_XR; 7869566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg)); 7873c633725SBarry Smith PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 7889566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&fd)); 7899566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X0,&XR)); 7909566063dSJacob Faibussowitsch PetscCall(VecLoad(XR,fd)); 7919566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 7929566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&ptr_X)); 7939566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XR,&ptr_XR)); 794c4762a1bSJed Brown for (i=0; i<Mx; i++) { 795c4762a1bSJed Brown if (i < count_slow/2 || i > count_slow/2+count_fast-1) nrm1 = nrm1 + hs*PetscAbs(ptr_X[i]-ptr_XR[i]); 796c4762a1bSJed Brown else nrm1 = nrm1 + hf*PetscAbs(ptr_X[i]-ptr_XR[i]); 797c4762a1bSJed Brown } 7989566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&ptr_X)); 7999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XR,&ptr_XR)); 8009566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1)); 8019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XR)); 802c4762a1bSJed Brown } 803c4762a1bSJed Brown } 804c4762a1bSJed Brown 8059566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD)); 8069566063dSJacob Faibussowitsch if (draw & 0x1) PetscCall(VecView(X0,PETSC_VIEWER_DRAW_WORLD)); 8079566063dSJacob Faibussowitsch if (draw & 0x2) PetscCall(VecView(X,PETSC_VIEWER_DRAW_WORLD)); 808c4762a1bSJed Brown if (draw & 0x4) { 809c4762a1bSJed Brown Vec Y; 8109566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&Y)); 8119566063dSJacob Faibussowitsch PetscCall(FVSample(&ctx,da,ptime,Y)); 8129566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y,-1,X)); 8139566063dSJacob Faibussowitsch PetscCall(VecView(Y,PETSC_VIEWER_DRAW_WORLD)); 8149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 815c4762a1bSJed Brown } 816c4762a1bSJed Brown 817c4762a1bSJed Brown if (view_final) { 818c4762a1bSJed Brown PetscViewer viewer; 8199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer)); 8209566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB)); 8219566063dSJacob Faibussowitsch PetscCall(VecView(X,viewer)); 8229566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 8239566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 824c4762a1bSJed Brown } 825c4762a1bSJed Brown 826c4762a1bSJed Brown /* Clean up */ 8279566063dSJacob Faibussowitsch PetscCall((*ctx.physics.destroy)(ctx.physics.user)); 8289566063dSJacob Faibussowitsch for (i=0; i<ctx.physics.dof; i++) PetscCall(PetscFree(ctx.physics.fieldname[i])); 8299566063dSJacob Faibussowitsch PetscCall(PetscFree3(ctx.u,ctx.flux,ctx.speeds)); 8309566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.iss)); 8319566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.isf)); 8329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 8339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X0)); 8349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&R)); 8359566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 8369566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 8379566063dSJacob Faibussowitsch PetscCall(PetscFree(index_slow)); 8389566063dSJacob Faibussowitsch PetscCall(PetscFree(index_fast)); 8399566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&physics)); 8409566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 841b122ec5aSJacob Faibussowitsch return 0; 842c4762a1bSJed Brown } 843c4762a1bSJed Brown 844c4762a1bSJed Brown /*TEST 845c4762a1bSJed Brown 846c4762a1bSJed Brown build: 847f56ea12dSJed Brown requires: !complex 848c4762a1bSJed Brown 849c4762a1bSJed Brown test: 850c4762a1bSJed 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 851c4762a1bSJed Brown 852c4762a1bSJed Brown test: 853c4762a1bSJed Brown suffix: 2 854c4762a1bSJed 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 855c4762a1bSJed Brown output_file: output/ex7_1.out 856c4762a1bSJed Brown 857c4762a1bSJed Brown test: 858c4762a1bSJed Brown suffix: 3 859c4762a1bSJed 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 860c4762a1bSJed Brown 861c4762a1bSJed Brown test: 862c4762a1bSJed Brown suffix: 4 863c4762a1bSJed 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 864c4762a1bSJed Brown output_file: output/ex7_3.out 865c2a16db8SHong Zhang 866c2a16db8SHong Zhang test: 867c2a16db8SHong Zhang suffix: 5 868c2a16db8SHong Zhang nsize: 2 869c2a16db8SHong 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 870c2a16db8SHong Zhang output_file: output/ex7_3.out 871c4762a1bSJed Brown TEST*/ 872