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 PetscErrorCode ierr; 71c4762a1bSJed Brown 72c4762a1bSJed Brown PetscFunctionBeginUser; 73c4762a1bSJed Brown ierr = PetscFree(vctx);CHKERRQ(ierr); 74c4762a1bSJed Brown PetscFunctionReturn(0); 75c4762a1bSJed Brown } 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* --------------------------------- Advection ----------------------------------- */ 78c4762a1bSJed Brown typedef struct { 79c4762a1bSJed Brown PetscReal a; /* advective velocity */ 80c4762a1bSJed Brown } AdvectCtx; 81c4762a1bSJed Brown 82c4762a1bSJed Brown static PetscErrorCode PhysicsFlux_Advect(void *vctx,const PetscScalar *u,PetscScalar *flux,PetscReal *maxspeed) 83c4762a1bSJed Brown { 84c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx*)vctx; 85c4762a1bSJed Brown PetscReal speed; 86c4762a1bSJed Brown 87c4762a1bSJed Brown PetscFunctionBeginUser; 88c4762a1bSJed Brown speed = ctx->a; 89c4762a1bSJed Brown flux[0] = speed*u[0]; 90c4762a1bSJed Brown *maxspeed = speed; 91c4762a1bSJed Brown PetscFunctionReturn(0); 92c4762a1bSJed Brown } 93c4762a1bSJed Brown 94c4762a1bSJed Brown static PetscErrorCode PhysicsSample_Advect(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u) 95c4762a1bSJed Brown { 96c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx*)vctx; 97c4762a1bSJed Brown PetscReal a = ctx->a,x0; 98c4762a1bSJed Brown 99c4762a1bSJed Brown PetscFunctionBeginUser; 100c4762a1bSJed Brown switch (bctype) { 101c4762a1bSJed Brown case FVBC_OUTFLOW: x0 = x-a*t; break; 102c4762a1bSJed Brown case FVBC_PERIODIC: x0 = RangeMod(x-a*t,xmin,xmax); break; 103c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown BCType"); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown switch (initial) { 106c4762a1bSJed Brown case 0: u[0] = (x0 < 0) ? 1 : -1; break; 107c4762a1bSJed Brown case 1: u[0] = (x0 < 0) ? -1 : 1; break; 108c4762a1bSJed Brown case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break; 109c4762a1bSJed Brown case 3: u[0] = PetscSinReal(2*PETSC_PI*x0); break; 110c4762a1bSJed Brown case 4: u[0] = PetscAbs(x0); break; 111c4762a1bSJed Brown case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2*PETSC_PI*x0)); break; 112c4762a1bSJed Brown case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2-x0 : 0)); break; 113c4762a1bSJed Brown case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI*x0),10.0);break; 114c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition"); 115c4762a1bSJed Brown } 116c4762a1bSJed Brown PetscFunctionReturn(0); 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 119c4762a1bSJed Brown static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx) 120c4762a1bSJed Brown { 121c4762a1bSJed Brown PetscErrorCode ierr; 122c4762a1bSJed Brown AdvectCtx *user; 123c4762a1bSJed Brown 124c4762a1bSJed Brown PetscFunctionBeginUser; 125c4762a1bSJed Brown ierr = PetscNew(&user);CHKERRQ(ierr); 126c4762a1bSJed Brown ctx->physics.sample = PhysicsSample_Advect; 127c4762a1bSJed Brown ctx->physics.flux = PhysicsFlux_Advect; 128c4762a1bSJed Brown ctx->physics.destroy = PhysicsDestroy_SimpleFree; 129c4762a1bSJed Brown ctx->physics.user = user; 130c4762a1bSJed Brown ctx->physics.dof = 1; 131c4762a1bSJed Brown ierr = PetscStrallocpy("u",&ctx->physics.fieldname[0]);CHKERRQ(ierr); 132c4762a1bSJed Brown user->a = 1; 133c4762a1bSJed Brown ierr = PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");CHKERRQ(ierr); 134c4762a1bSJed Brown { 135c4762a1bSJed Brown ierr = PetscOptionsReal("-physics_advect_a","Speed","",user->a,&user->a,NULL);CHKERRQ(ierr); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 138c4762a1bSJed Brown PetscFunctionReturn(0); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver ----------------------------------- */ 142c4762a1bSJed Brown 143c4762a1bSJed Brown static PetscErrorCode FVRHSFunction(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 144c4762a1bSJed Brown { 145c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 146c4762a1bSJed Brown PetscErrorCode ierr; 147c4762a1bSJed Brown PetscInt i,j,Mx,dof,xs,xm,sf = ctx->sf,fs = ctx->fs; 148c4762a1bSJed Brown PetscReal hf,hs,cfl_idt = 0; 149c4762a1bSJed Brown PetscScalar *x,*f,*r,*min,*alpha,*gamma; 150c4762a1bSJed Brown Vec Xloc; 151c4762a1bSJed Brown DM da; 152c4762a1bSJed Brown 153c4762a1bSJed Brown PetscFunctionBeginUser; 154c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); /* Xloc contains ghost points */ 156c4762a1bSJed Brown ierr = DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); /* Mx is the number of center points */ 157c4762a1bSJed Brown hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; 158c4762a1bSJed Brown hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 159c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); /* X is solution vector which does not contain ghost points */ 160c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 161c4762a1bSJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); /* F is the right hand side function corresponds to center points */ 162c4762a1bSJed Brown ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr); 164c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = PetscMalloc4(dof,&r,dof,&min,dof,&alpha,dof,&gamma);CHKERRQ(ierr); 166c4762a1bSJed Brown 167c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 168c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 169c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 170c4762a1bSJed Brown } 171c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 172c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 173c4762a1bSJed Brown } 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 177c4762a1bSJed Brown PetscReal maxspeed; 178c4762a1bSJed Brown PetscScalar *u; 179c4762a1bSJed Brown if (i < sf || i > fs+1) { 180c4762a1bSJed Brown u = &ctx->u[0]; 181c4762a1bSJed Brown alpha[0] = 1.0/6.0; 182c4762a1bSJed Brown gamma[0] = 1.0/3.0; 183c4762a1bSJed Brown for (j=0; j<dof; j++) { 184c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 185c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 186c4762a1bSJed 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]); 187c4762a1bSJed Brown } 188c4762a1bSJed Brown ierr = (*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); 189c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hs)); 190c4762a1bSJed Brown if (i > xs) { 191c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hs; 192c4762a1bSJed Brown } 193c4762a1bSJed Brown if (i < xs+xm) { 194c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hs; 195c4762a1bSJed Brown } 196c4762a1bSJed Brown } else if (i == sf) { 197c4762a1bSJed Brown u = &ctx->u[0]; 198c4762a1bSJed Brown alpha[0] = hs*hf/(hs+hs)/(hs+hs+hf); 199c4762a1bSJed Brown gamma[0] = hs*(hs+hs)/(hs+hf)/(hs+hs+hf); 200c4762a1bSJed Brown for (j=0; j<dof; j++) { 201c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 202c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 203c4762a1bSJed 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]); 204c4762a1bSJed Brown } 205c4762a1bSJed Brown ierr = (*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); 206c4762a1bSJed Brown if (i > xs) { 207c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hs; 208c4762a1bSJed Brown } 209c4762a1bSJed Brown if (i < xs+xm) { 210c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hf; 211c4762a1bSJed Brown } 212c4762a1bSJed Brown } else if (i == sf+1) { 213c4762a1bSJed Brown u = &ctx->u[0]; 214c4762a1bSJed Brown alpha[0] = hf*hf/(hs+hf)/(hs+hf+hf); 215c4762a1bSJed Brown gamma[0] = hf*(hs+hf)/(hf+hf)/(hs+hf+hf); 216c4762a1bSJed Brown for (j=0; j<dof; j++) { 217c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 218c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 219c4762a1bSJed 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]); 220c4762a1bSJed Brown } 221c4762a1bSJed Brown ierr = (*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); 222c4762a1bSJed Brown if (i > xs) { 223c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hf; 224c4762a1bSJed Brown } 225c4762a1bSJed Brown if (i < xs+xm) { 226c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hf; 227c4762a1bSJed Brown } 228c4762a1bSJed Brown } else if (i > sf+1 && i < fs) { 229c4762a1bSJed Brown u = &ctx->u[0]; 230c4762a1bSJed Brown alpha[0] = 1.0/6.0; 231c4762a1bSJed Brown gamma[0] = 1.0/3.0; 232c4762a1bSJed Brown for (j=0; j<dof; j++) { 233c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 234c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 235c4762a1bSJed 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]); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown ierr = (*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); 238c4762a1bSJed Brown if (i > xs) { 239c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hf; 240c4762a1bSJed Brown } 241c4762a1bSJed Brown if (i < xs+xm) { 242c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hf; 243c4762a1bSJed Brown } 244c4762a1bSJed Brown } else if (i == fs) { 245c4762a1bSJed Brown u = &ctx->u[0]; 246c4762a1bSJed Brown alpha[0] = hf*hs/(hf+hf)/(hf+hf+hs); 247c4762a1bSJed Brown gamma[0] = hf*(hf+hf)/(hf+hs)/(hf+hf+hs); 248c4762a1bSJed Brown for (j=0; j<dof; j++) { 249c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 250c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 251c4762a1bSJed 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]); 252c4762a1bSJed Brown } 253c4762a1bSJed Brown ierr = (*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); 254c4762a1bSJed Brown if (i > xs) { 255c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hf; 256c4762a1bSJed Brown } 257c4762a1bSJed Brown if (i < xs+xm) { 258c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hs; 259c4762a1bSJed Brown } 260c4762a1bSJed Brown } else if (i == fs+1) { 261c4762a1bSJed Brown u = &ctx->u[0]; 262c4762a1bSJed Brown alpha[0] = hs*hs/(hf+hs)/(hf+hs+hs); 263c4762a1bSJed Brown gamma[0] = hs*(hf+hs)/(hs+hs)/(hf+hs+hs); 264c4762a1bSJed Brown for (j=0; j<dof; j++) { 265c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 266c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 267c4762a1bSJed 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]); 268c4762a1bSJed Brown } 269c4762a1bSJed Brown ierr = (*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); 270c4762a1bSJed Brown if (i > xs) { 271c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hs; 272c4762a1bSJed Brown } 273c4762a1bSJed Brown if (i < xs+xm) { 274c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hs; 275c4762a1bSJed Brown } 276c4762a1bSJed Brown } 277c4762a1bSJed Brown } 278c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 279c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr); 280c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 281ffc4695bSBarry Smith ierr = MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr); 282c4762a1bSJed Brown if (0) { 283c2a16db8SHong Zhang /* We need a way to inform the TS of a CFL constraint, this is a debugging fragment */ 284c4762a1bSJed Brown PetscReal dt,tnow; 285c4762a1bSJed Brown ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 286c4762a1bSJed Brown ierr = TSGetTime(ts,&tnow);CHKERRQ(ierr); 287c4762a1bSJed Brown if (dt > 0.5/ctx->cfl_idt) { 288c4762a1bSJed Brown ierr = PetscPrintf(ctx->comm,"Stability constraint exceeded at t=%g, dt %g > %g\n",(double)tnow,(double)dt,(double)(0.5/ctx->cfl_idt));CHKERRQ(ierr); 289c4762a1bSJed Brown } 290c4762a1bSJed Brown } 291c4762a1bSJed Brown ierr = PetscFree4(r,min,alpha,gamma);CHKERRQ(ierr); 292c4762a1bSJed Brown PetscFunctionReturn(0); 293c4762a1bSJed Brown } 294c4762a1bSJed Brown 295c4762a1bSJed Brown static PetscErrorCode FVRHSFunctionslow(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 296c4762a1bSJed Brown { 297c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 298c4762a1bSJed Brown PetscErrorCode ierr; 299c4762a1bSJed Brown PetscInt i,j,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs; 300c4762a1bSJed Brown PetscReal hf,hs; 301c4762a1bSJed Brown PetscScalar *x,*f,*r,*min,*alpha,*gamma; 302c4762a1bSJed Brown Vec Xloc; 303c4762a1bSJed Brown DM da; 304c4762a1bSJed Brown 305c4762a1bSJed Brown PetscFunctionBeginUser; 306c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 307c4762a1bSJed Brown ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); /* Xloc contains ghost points */ 308c4762a1bSJed Brown ierr = DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); /* Mx is the number of center points */ 309c4762a1bSJed Brown hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; 310c4762a1bSJed Brown hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 311c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); /* X is solution vector which does not contain ghost points */ 312c4762a1bSJed Brown ierr = DMGlobalToLocalEnd (da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 313c4762a1bSJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); /* F is the right hand side function corresponds to center points */ 314c4762a1bSJed Brown ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 315c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 316c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 317c4762a1bSJed Brown ierr = PetscMalloc4(dof,&r,dof,&min,dof,&alpha,dof,&gamma);CHKERRQ(ierr); 318c4762a1bSJed Brown 319c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 320c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 321c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 322c4762a1bSJed Brown } 323c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 324c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 325c4762a1bSJed Brown } 326c4762a1bSJed Brown } 327c4762a1bSJed Brown 328c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 329c4762a1bSJed Brown PetscReal maxspeed; 330c4762a1bSJed Brown PetscScalar *u; 331c4762a1bSJed Brown if (i < sf) { 332c4762a1bSJed Brown u = &ctx->u[0]; 333c4762a1bSJed Brown alpha[0] = 1.0/6.0; 334c4762a1bSJed Brown gamma[0] = 1.0/3.0; 335c4762a1bSJed Brown for (j=0; j<dof; j++) { 336c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 337c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 338c4762a1bSJed 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]); 339c4762a1bSJed Brown } 340c4762a1bSJed Brown ierr = (*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); 341c4762a1bSJed Brown if (i > xs) { 342c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hs; 343c4762a1bSJed Brown } 344c4762a1bSJed Brown if (i < xs+xm) { 345c2a16db8SHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hs; 346c4762a1bSJed Brown islow++; 347c4762a1bSJed Brown } 348c4762a1bSJed Brown } else if (i == sf) { 349c4762a1bSJed Brown u = &ctx->u[0]; 350c4762a1bSJed Brown alpha[0] = hs*hf/(hs+hs)/(hs+hs+hf); 351c4762a1bSJed Brown gamma[0] = hs*(hs+hs)/(hs+hf)/(hs+hs+hf); 352c4762a1bSJed Brown for (j=0; j<dof; j++) { 353c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 354c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 355c4762a1bSJed 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]); 356c4762a1bSJed Brown } 357c4762a1bSJed Brown ierr = (*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); 358c2a16db8SHong Zhang if (i > xs) { 359c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hs; 360c4762a1bSJed Brown } 361c4762a1bSJed Brown } else if (i == fs) { 362c4762a1bSJed Brown u = &ctx->u[0]; 363c4762a1bSJed Brown alpha[0] = hf*hs/(hf+hf)/(hf+hf+hs); 364c4762a1bSJed Brown gamma[0] = hf*(hf+hf)/(hf+hs)/(hf+hf+hs); 365c4762a1bSJed Brown for (j=0; j<dof; j++) { 366c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 367c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 368c4762a1bSJed 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]); 369c4762a1bSJed Brown } 370c4762a1bSJed Brown ierr = (*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); 371c4762a1bSJed Brown if (i < xs+xm) { 372c2a16db8SHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hs; 373c4762a1bSJed Brown islow++; 374c4762a1bSJed Brown } 375c4762a1bSJed Brown } else if (i == fs+1) { 376c4762a1bSJed Brown u = &ctx->u[0]; 377c4762a1bSJed Brown alpha[0] = hs*hs/(hf+hs)/(hf+hs+hs); 378c4762a1bSJed Brown gamma[0] = hs*(hf+hs)/(hs+hs)/(hf+hs+hs); 379c4762a1bSJed Brown for (j=0; j<dof; j++) { 380c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 381c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 382c4762a1bSJed 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]); 383c4762a1bSJed Brown } 384c4762a1bSJed Brown ierr = (*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); 385c4762a1bSJed Brown if (i > xs) { 386c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hs; 387c4762a1bSJed Brown } 388c4762a1bSJed Brown if (i < xs+xm) { 389c2a16db8SHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hs; 390c2a16db8SHong Zhang islow++; 391c4762a1bSJed Brown } 392c4762a1bSJed Brown } else if (i > fs+1) { 393c4762a1bSJed Brown u = &ctx->u[0]; 394c4762a1bSJed Brown alpha[0] = 1.0/6.0; 395c4762a1bSJed Brown gamma[0] = 1.0/3.0; 396c4762a1bSJed Brown for (j=0; j<dof; j++) { 397c4762a1bSJed Brown r[j] = (x[i*dof+j] - x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 398c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 399c4762a1bSJed 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]); 400c4762a1bSJed Brown } 401c4762a1bSJed Brown ierr = (*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); 402c4762a1bSJed Brown if (i > xs) { 403c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hs; 404c4762a1bSJed Brown } 405c4762a1bSJed Brown if (i < xs+xm) { 406c2a16db8SHong Zhang for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hs; 407c2a16db8SHong Zhang islow++; 408c4762a1bSJed Brown } 409c4762a1bSJed Brown } 410c4762a1bSJed Brown } 411c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 412c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 413c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 414c4762a1bSJed Brown ierr = PetscFree4(r,min,alpha,gamma);CHKERRQ(ierr); 415c4762a1bSJed Brown PetscFunctionReturn(0); 416c4762a1bSJed Brown } 417c4762a1bSJed Brown 418c4762a1bSJed Brown static PetscErrorCode FVRHSFunctionfast(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 419c4762a1bSJed Brown { 420c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 421c4762a1bSJed Brown PetscErrorCode ierr; 422c4762a1bSJed Brown PetscInt i,j,Mx,dof,xs,xm,ifast = 0,sf = ctx->sf,fs = ctx->fs; 423c4762a1bSJed Brown PetscReal hf,hs; 424c4762a1bSJed Brown PetscScalar *x,*f,*r,*min,*alpha,*gamma; 425c4762a1bSJed Brown Vec Xloc; 426c4762a1bSJed Brown DM da; 427c4762a1bSJed Brown 428c4762a1bSJed Brown PetscFunctionBeginUser; 429c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 430c4762a1bSJed Brown ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); /* Xloc contains ghost points */ 431c4762a1bSJed Brown ierr = DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); /* Mx is the number of center points */ 432c4762a1bSJed Brown hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; 433c4762a1bSJed Brown hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 434c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); /* X is solution vector which does not contain ghost points */ 435c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 436c4762a1bSJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); /* F is the right hand side function corresponds to center points */ 437c4762a1bSJed Brown ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 438c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 439c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 440c4762a1bSJed Brown ierr = PetscMalloc4(dof,&r,dof,&min,dof,&alpha,dof,&gamma);CHKERRQ(ierr); 441c4762a1bSJed Brown 442c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 443c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 444c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 445c4762a1bSJed Brown } 446c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 447c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 448c4762a1bSJed Brown } 449c4762a1bSJed Brown } 450c4762a1bSJed Brown 451c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 452c4762a1bSJed Brown PetscReal maxspeed; 453c4762a1bSJed Brown PetscScalar *u; 454c4762a1bSJed Brown if (i == sf) { 455c4762a1bSJed Brown u = &ctx->u[0]; 456c4762a1bSJed Brown alpha[0] = hs*hf/(hs+hs)/(hs+hs+hf); 457c4762a1bSJed Brown gamma[0] = hs*(hs+hs)/(hs+hf)/(hs+hs+hf); 458c4762a1bSJed Brown for (j=0; j<dof; j++) { 459c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 460c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 461c4762a1bSJed 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]); 462c4762a1bSJed Brown } 463c4762a1bSJed Brown ierr = (*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); 464c4762a1bSJed Brown if (i < xs+xm) { 465c2a16db8SHong Zhang for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hf; 466c4762a1bSJed Brown ifast++; 467c4762a1bSJed Brown } 468c4762a1bSJed Brown } else if (i == sf+1) { 469c4762a1bSJed Brown u = &ctx->u[0]; 470c4762a1bSJed Brown alpha[0] = hf*hf/(hs+hf)/(hs+hf+hf); 471c4762a1bSJed Brown gamma[0] = hf*(hs+hf)/(hf+hf)/(hs+hf+hf); 472c4762a1bSJed Brown for (j=0; j<dof; j++) { 473c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 474c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 475c4762a1bSJed 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]); 476c4762a1bSJed Brown } 477c4762a1bSJed Brown ierr = (*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); 478c4762a1bSJed Brown if (i > xs) { 479c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hf; 480c4762a1bSJed Brown } 481c4762a1bSJed Brown if (i < xs+xm) { 482c2a16db8SHong Zhang for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hf; 483c4762a1bSJed Brown ifast++; 484c4762a1bSJed Brown } 485c4762a1bSJed Brown } else if (i > sf+1 && i < fs) { 486c4762a1bSJed Brown u = &ctx->u[0]; 487c4762a1bSJed Brown alpha[0] = 1.0/6.0; 488c4762a1bSJed Brown gamma[0] = 1.0/3.0; 489c4762a1bSJed Brown for (j=0; j<dof; j++) { 490c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 491c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 492c4762a1bSJed 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]); 493c4762a1bSJed Brown } 494c4762a1bSJed Brown ierr = (*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); 495c4762a1bSJed Brown if (i > xs) { 496c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hf; 497c4762a1bSJed Brown } 498c4762a1bSJed Brown if (i < xs+xm) { 499c2a16db8SHong Zhang for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hf; 500c4762a1bSJed Brown ifast++; 501c4762a1bSJed Brown } 502c4762a1bSJed Brown } else if (i == fs) { 503c4762a1bSJed Brown u = &ctx->u[0]; 504c4762a1bSJed Brown alpha[0] = hf*hs/(hf+hf)/(hf+hf+hs); 505c4762a1bSJed Brown gamma[0] = hf*(hf+hf)/(hf+hs)/(hf+hf+hs); 506c4762a1bSJed Brown for (j=0; j<dof; j++) { 507c4762a1bSJed Brown r[j] = (x[i*dof+j]-x[(i-1)*dof+j])/(x[(i-1)*dof+j]-x[(i-2)*dof+j]); 508c4762a1bSJed Brown min[j] = PetscMin(r[j],2.0); 509c4762a1bSJed 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]); 510c4762a1bSJed Brown } 511c4762a1bSJed Brown ierr = (*ctx->physics.flux)(ctx->physics.user,u,ctx->flux,&maxspeed);CHKERRQ(ierr); 512c4762a1bSJed Brown if (i > xs) { 513c2a16db8SHong Zhang for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hf; 514c4762a1bSJed Brown } 515c4762a1bSJed Brown } 516c4762a1bSJed Brown } 517c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 518c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 519c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 520c4762a1bSJed Brown ierr = PetscFree4(r,min,alpha,gamma);CHKERRQ(ierr); 521c4762a1bSJed Brown PetscFunctionReturn(0); 522c4762a1bSJed Brown } 523c4762a1bSJed Brown 524c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */ 525c4762a1bSJed Brown 526c4762a1bSJed Brown PetscErrorCode FVSample(FVCtx *ctx,DM da,PetscReal time,Vec U) 527c4762a1bSJed Brown { 528c4762a1bSJed Brown PetscErrorCode ierr; 529c4762a1bSJed Brown PetscScalar *u,*uj,xj,xi; 530c4762a1bSJed Brown PetscInt i,j,k,dof,xs,xm,Mx,count_slow,count_fast; 531c4762a1bSJed Brown const PetscInt N=200; 532c4762a1bSJed Brown 533c4762a1bSJed Brown PetscFunctionBeginUser; 534*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!ctx->physics.sample,PETSC_COMM_SELF,PETSC_ERR_SUP,"Physics has not provided a sampling function"); 535c4762a1bSJed Brown ierr = DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 536c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 537c4762a1bSJed Brown ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr); 538c4762a1bSJed Brown ierr = PetscMalloc1(dof,&uj);CHKERRQ(ierr); 539c4762a1bSJed Brown const PetscReal hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; 540c4762a1bSJed Brown const PetscReal hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 541c4762a1bSJed Brown count_slow = Mx/(1+ctx->hratio); 542c4762a1bSJed Brown count_fast = Mx-count_slow; 543c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 544c4762a1bSJed Brown if (i*hs+0.5*hs<(ctx->xmax-ctx->xmin)*0.25) { 545c4762a1bSJed Brown xi = ctx->xmin+0.5*hs+i*hs; 546c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 547c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 548c4762a1bSJed Brown for (j=0; j<N+1; j++) { 549c4762a1bSJed Brown xj = xi+hs*(j-N/2)/(PetscReal)N; 550c4762a1bSJed Brown ierr = (*ctx->physics.sample)(ctx->physics.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); 551c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 552c4762a1bSJed Brown } 553c4762a1bSJed Brown } else if ((ctx->xmax-ctx->xmin)*0.25+(i-count_slow/2)*hf+0.5*hf<(ctx->xmax-ctx->xmin)*0.75) { 554c4762a1bSJed Brown xi = ctx->xmin+(ctx->xmax-ctx->xmin)*0.25+0.5*hf+(i-count_slow/2)*hf; 555c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 556c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 557c4762a1bSJed Brown for (j=0; j<N+1; j++) { 558c4762a1bSJed Brown xj = xi+hf*(j-N/2)/(PetscReal)N; 559c4762a1bSJed Brown ierr = (*ctx->physics.sample)(ctx->physics.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); 560c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 561c4762a1bSJed Brown } 562c4762a1bSJed Brown } else { 563c4762a1bSJed Brown xi = ctx->xmin+(ctx->xmax-ctx->xmin)*0.75+0.5*hs+(i-count_slow/2-count_fast)*hs; 564c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 565c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 566c4762a1bSJed Brown for (j=0; j<N+1; j++) { 567c4762a1bSJed Brown xj = xi+hs*(j-N/2)/(PetscReal)N; 568c4762a1bSJed Brown ierr = (*ctx->physics.sample)(ctx->physics.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); 569c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 570c4762a1bSJed Brown } 571c4762a1bSJed Brown } 572c4762a1bSJed Brown } 573c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr); 574c4762a1bSJed Brown ierr = PetscFree(uj);CHKERRQ(ierr); 575c4762a1bSJed Brown PetscFunctionReturn(0); 576c4762a1bSJed Brown } 577c4762a1bSJed Brown 578c4762a1bSJed Brown static PetscErrorCode SolutionStatsView(DM da,Vec X,PetscViewer viewer) 579c4762a1bSJed Brown { 580c4762a1bSJed Brown PetscErrorCode ierr; 581c4762a1bSJed Brown PetscReal xmin,xmax; 582c4762a1bSJed Brown PetscScalar sum,tvsum,tvgsum; 583c4762a1bSJed Brown const PetscScalar *x; 584c4762a1bSJed Brown PetscInt imin,imax,Mx,i,j,xs,xm,dof; 585c4762a1bSJed Brown Vec Xloc; 586c4762a1bSJed Brown PetscBool iascii; 587c4762a1bSJed Brown 588c4762a1bSJed Brown PetscFunctionBeginUser; 589c4762a1bSJed Brown ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 590c4762a1bSJed Brown if (iascii) { 591c4762a1bSJed Brown /* PETSc lacks a function to compute total variation norm (difficult in multiple dimensions), we do it here */ 592c4762a1bSJed Brown ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 593c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 594c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 595c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,Xloc,(void*)&x);CHKERRQ(ierr); 596c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 597c4762a1bSJed Brown ierr = DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 598c4762a1bSJed Brown tvsum = 0; 599c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 600c4762a1bSJed Brown for (j=0; j<dof; j++) tvsum += PetscAbsScalar(x[i*dof+j]-x[(i-1)*dof+j]); 601c4762a1bSJed Brown } 602ffc4695bSBarry Smith ierr = MPI_Allreduce(&tvsum,&tvgsum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr); 603c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,Xloc,(void*)&x);CHKERRQ(ierr); 604c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 605c4762a1bSJed Brown 606c4762a1bSJed Brown ierr = VecMin(X,&imin,&xmin);CHKERRQ(ierr); 607c4762a1bSJed Brown ierr = VecMax(X,&imax,&xmax);CHKERRQ(ierr); 608c4762a1bSJed Brown ierr = VecSum(X,&sum);CHKERRQ(ierr); 609c4762a1bSJed Brown ierr = 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));CHKERRQ(ierr); 610c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type not supported"); 611c4762a1bSJed Brown PetscFunctionReturn(0); 612c4762a1bSJed Brown } 613c4762a1bSJed Brown 614c4762a1bSJed Brown static PetscErrorCode SolutionErrorNorms(FVCtx *ctx,DM da,PetscReal t,Vec X,PetscReal *nrm1) 615c4762a1bSJed Brown { 616c4762a1bSJed Brown PetscErrorCode ierr; 617c4762a1bSJed Brown Vec Y; 618c4762a1bSJed Brown PetscInt i,Mx,count_slow=0,count_fast=0; 619c4762a1bSJed Brown const PetscScalar *ptr_X,*ptr_Y; 620c4762a1bSJed Brown 621c4762a1bSJed Brown PetscFunctionBeginUser; 622c4762a1bSJed Brown ierr = VecGetSize(X,&Mx);CHKERRQ(ierr); 623c4762a1bSJed Brown ierr = VecDuplicate(X,&Y);CHKERRQ(ierr); 624c4762a1bSJed Brown ierr = FVSample(ctx,da,t,Y);CHKERRQ(ierr); 625c4762a1bSJed Brown const PetscReal hs = (ctx->xmax-ctx->xmin)/2.0*(ctx->hratio+1.0)/Mx; 626c4762a1bSJed Brown const PetscReal hf = (ctx->xmax-ctx->xmin)/2.0*(1.0+1.0/ctx->hratio)/Mx; 627c4762a1bSJed Brown count_slow = (PetscReal)Mx/(1.0+ctx->hratio); 628c4762a1bSJed Brown count_fast = Mx-count_slow; 629c4762a1bSJed Brown ierr = VecGetArrayRead(X,&ptr_X);CHKERRQ(ierr); 630c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&ptr_Y);CHKERRQ(ierr); 631c4762a1bSJed Brown for (i=0; i<Mx; i++) { 632c4762a1bSJed Brown if (i < count_slow/2 || i > count_slow/2+count_fast-1) *nrm1 += hs*PetscAbs(ptr_X[i]-ptr_Y[i]); 633c4762a1bSJed Brown else *nrm1 += hf*PetscAbs(ptr_X[i]-ptr_Y[i]); 634c4762a1bSJed Brown } 635c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&ptr_X);CHKERRQ(ierr); 636c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&ptr_Y);CHKERRQ(ierr); 637c4762a1bSJed Brown ierr = VecDestroy(&Y);CHKERRQ(ierr); 638c4762a1bSJed Brown PetscFunctionReturn(0); 639c4762a1bSJed Brown } 640c4762a1bSJed Brown 641c4762a1bSJed Brown int main(int argc,char *argv[]) 642c4762a1bSJed Brown { 643c4762a1bSJed Brown char physname[256] = "advect",final_fname[256] = "solution.m"; 644c4762a1bSJed Brown PetscFunctionList physics = 0; 645c4762a1bSJed Brown MPI_Comm comm; 646c4762a1bSJed Brown TS ts; 647c4762a1bSJed Brown DM da; 648c4762a1bSJed Brown Vec X,X0,R; 649c4762a1bSJed Brown FVCtx ctx; 650c4762a1bSJed Brown PetscInt i,k,dof,xs,xm,Mx,draw = 0,count_slow,count_fast,islow = 0,ifast = 0,*index_slow,*index_fast; 651c4762a1bSJed Brown PetscBool view_final = PETSC_FALSE; 652c4762a1bSJed Brown PetscReal ptime; 653c4762a1bSJed Brown PetscErrorCode ierr; 654c4762a1bSJed Brown 655c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr; 656c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 657c4762a1bSJed Brown ierr = PetscMemzero(&ctx,sizeof(ctx));CHKERRQ(ierr); 658c4762a1bSJed Brown 659c4762a1bSJed Brown /* Register physical models to be available on the command line */ 660c4762a1bSJed Brown ierr = PetscFunctionListAdd(&physics,"advect",PhysicsCreate_Advect);CHKERRQ(ierr); 661c4762a1bSJed Brown 662c4762a1bSJed Brown ctx.comm = comm; 663c4762a1bSJed Brown ctx.cfl = 0.9; 664c4762a1bSJed Brown ctx.bctype = FVBC_PERIODIC; 665c4762a1bSJed Brown ctx.xmin = -1.0; 666c4762a1bSJed Brown ctx.xmax = 1.0; 667c4762a1bSJed Brown ierr = PetscOptionsBegin(comm,NULL,"Finite Volume solver options","");CHKERRQ(ierr); 668c4762a1bSJed Brown ierr = PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL);CHKERRQ(ierr); 669c4762a1bSJed Brown ierr = PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL);CHKERRQ(ierr); 670c4762a1bSJed Brown ierr = PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL);CHKERRQ(ierr); 671c4762a1bSJed Brown ierr = PetscOptionsString("-view_final","Write final solution in ASCII MATLAB format to given file name","",final_fname,final_fname,sizeof(final_fname),&view_final);CHKERRQ(ierr); 672c4762a1bSJed Brown ierr = PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL);CHKERRQ(ierr); 673c4762a1bSJed Brown ierr = PetscOptionsBool("-exact","Compare errors with exact solution","",ctx.exact,&ctx.exact,NULL);CHKERRQ(ierr); 674c4762a1bSJed Brown ierr = PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL);CHKERRQ(ierr); 675c4762a1bSJed Brown ierr = PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL);CHKERRQ(ierr); 676c4762a1bSJed Brown ierr = PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL);CHKERRQ(ierr); 677c4762a1bSJed Brown ierr = PetscOptionsInt("-hratio","Spacing ratio","",ctx.hratio,&ctx.hratio,NULL);CHKERRQ(ierr); 678c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 679c4762a1bSJed Brown 680c4762a1bSJed Brown /* Choose the physics from the list of registered models */ 681c4762a1bSJed Brown { 682c4762a1bSJed Brown PetscErrorCode (*r)(FVCtx*); 683c4762a1bSJed Brown ierr = PetscFunctionListFind(physics,physname,&r);CHKERRQ(ierr); 684*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Physics '%s' not found",physname); 685c4762a1bSJed Brown /* Create the physics, will set the number of fields and their names */ 686c4762a1bSJed Brown ierr = (*r)(&ctx);CHKERRQ(ierr); 687c4762a1bSJed Brown } 688c4762a1bSJed Brown 689c4762a1bSJed Brown /* Create a DMDA to manage the parallel grid */ 690c4762a1bSJed Brown ierr = DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics.dof,2,NULL,&da);CHKERRQ(ierr); 691c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 692c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 693c4762a1bSJed Brown /* Inform the DMDA of the field names provided by the physics. */ 694c4762a1bSJed Brown /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ 695c4762a1bSJed Brown for (i=0; i<ctx.physics.dof; i++) { 696c4762a1bSJed Brown ierr = DMDASetFieldName(da,i,ctx.physics.fieldname[i]);CHKERRQ(ierr); 697c4762a1bSJed Brown } 698c4762a1bSJed Brown ierr = DMDAGetInfo(da,0,&Mx,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 699c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 700c4762a1bSJed Brown 701c4762a1bSJed Brown /* Set coordinates of cell centers */ 702c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(da,ctx.xmin+0.5*(ctx.xmax-ctx.xmin)/Mx,ctx.xmax+0.5*(ctx.xmax-ctx.xmin)/Mx,0,0,0,0);CHKERRQ(ierr); 703c4762a1bSJed Brown 704c4762a1bSJed Brown /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */ 705c4762a1bSJed Brown ierr = PetscMalloc3(dof,&ctx.u,dof,&ctx.flux,dof,&ctx.speeds);CHKERRQ(ierr); 706c4762a1bSJed Brown 707c4762a1bSJed Brown /* Create a vector to store the solution and to save the initial state */ 708c4762a1bSJed Brown ierr = DMCreateGlobalVector(da,&X);CHKERRQ(ierr); 709c4762a1bSJed Brown ierr = VecDuplicate(X,&X0);CHKERRQ(ierr); 710c4762a1bSJed Brown ierr = VecDuplicate(X,&R);CHKERRQ(ierr); 711c4762a1bSJed Brown 712c4762a1bSJed Brown /* create index for slow parts and fast parts*/ 713c4762a1bSJed Brown count_slow = Mx/(1+ctx.hratio); 714*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(count_slow%2,PETSC_COMM_WORLD,PETSC_ERR_USER,"Please adjust grid size Mx (-da_grid_x) and hratio (-hratio) so that Mx/(1+hartio) is even"); 715c4762a1bSJed Brown count_fast = Mx-count_slow; 716c4762a1bSJed Brown ctx.sf = count_slow/2; 717c4762a1bSJed Brown ctx.fs = ctx.sf + count_fast; 718c4762a1bSJed Brown ierr = PetscMalloc1(xm*dof,&index_slow);CHKERRQ(ierr); 719c4762a1bSJed Brown ierr = PetscMalloc1(xm*dof,&index_fast);CHKERRQ(ierr); 720c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 721c4762a1bSJed Brown if (i < count_slow/2 || i > count_slow/2+count_fast-1) 722c4762a1bSJed Brown for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k; 723c4762a1bSJed Brown else 724c4762a1bSJed Brown for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k; 725c4762a1bSJed Brown } 726c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss);CHKERRQ(ierr); 727c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf);CHKERRQ(ierr); 728c4762a1bSJed Brown 729c4762a1bSJed Brown /* Create a time-stepping object */ 730c4762a1bSJed Brown ierr = TSCreate(comm,&ts);CHKERRQ(ierr); 731c4762a1bSJed Brown ierr = TSSetDM(ts,da);CHKERRQ(ierr); 732c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,R,FVRHSFunction,&ctx);CHKERRQ(ierr); 733c4762a1bSJed Brown ierr = TSRHSSplitSetIS(ts,"slow",ctx.iss);CHKERRQ(ierr); 734c4762a1bSJed Brown ierr = TSRHSSplitSetIS(ts,"fast",ctx.isf);CHKERRQ(ierr); 735c4762a1bSJed Brown ierr = TSRHSSplitSetRHSFunction(ts,"slow",NULL,FVRHSFunctionslow,&ctx);CHKERRQ(ierr); 736c4762a1bSJed Brown ierr = TSRHSSplitSetRHSFunction(ts,"fast",NULL,FVRHSFunctionfast,&ctx);CHKERRQ(ierr); 737c4762a1bSJed Brown 738c4762a1bSJed Brown ierr = TSSetType(ts,TSMPRK);CHKERRQ(ierr); 739c4762a1bSJed Brown ierr = TSSetMaxTime(ts,10);CHKERRQ(ierr); 740c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 741c4762a1bSJed Brown 742c4762a1bSJed Brown /* Compute initial conditions and starting time step */ 743c4762a1bSJed Brown ierr = FVSample(&ctx,da,0,X0);CHKERRQ(ierr); 744c4762a1bSJed Brown ierr = FVRHSFunction(ts,0,X0,X,(void*)&ctx);CHKERRQ(ierr); /* Initial function evaluation, only used to determine max speed */ 745c4762a1bSJed Brown ierr = VecCopy(X0,X);CHKERRQ(ierr); /* The function value was not used so we set X=X0 again */ 746c4762a1bSJed Brown ierr = TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt);CHKERRQ(ierr); 747c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); /* Take runtime options */ 748c4762a1bSJed Brown ierr = SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 749c4762a1bSJed Brown { 750c4762a1bSJed Brown PetscInt steps; 751c4762a1bSJed Brown PetscScalar mass_initial,mass_final,mass_difference,mass_differenceg; 752c4762a1bSJed Brown const PetscScalar *ptr_X,*ptr_X0; 753c2a16db8SHong Zhang const PetscReal hs = (ctx.xmax-ctx.xmin)/2.0/count_slow; 754c2a16db8SHong Zhang const PetscReal hf = (ctx.xmax-ctx.xmin)/2.0/count_fast; 755c4762a1bSJed Brown ierr = TSSolve(ts,X);CHKERRQ(ierr); 756c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&ptime);CHKERRQ(ierr); 757c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 758c4762a1bSJed Brown /* calculate the total mass at initial time and final time */ 759c4762a1bSJed Brown mass_initial = 0.0; 760c4762a1bSJed Brown mass_final = 0.0; 761c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,X0,(void*)&ptr_X0);CHKERRQ(ierr); 762c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,X,(void*)&ptr_X);CHKERRQ(ierr); 763c2a16db8SHong Zhang for (i=xs; i<xs+xm; i++) { 764c2a16db8SHong Zhang if (i < ctx.sf || i > ctx.fs-1) { 765c2a16db8SHong Zhang for (k=0; k<dof; k++) { 766c2a16db8SHong Zhang mass_initial = mass_initial+hs*ptr_X0[i*dof+k]; 767c2a16db8SHong Zhang mass_final = mass_final+hs*ptr_X[i*dof+k]; 768c2a16db8SHong Zhang } 769c4762a1bSJed Brown } else { 770c2a16db8SHong Zhang for (k=0; k<dof; k++) { 771c2a16db8SHong Zhang mass_initial = mass_initial+hf*ptr_X0[i*dof+k]; 772c2a16db8SHong Zhang mass_final = mass_final+hf*ptr_X[i*dof+k]; 773c2a16db8SHong Zhang } 774c4762a1bSJed Brown } 775c4762a1bSJed Brown } 776c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,X0,(void*)&ptr_X0);CHKERRQ(ierr); 777c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,X,(void*)&ptr_X);CHKERRQ(ierr); 778c4762a1bSJed Brown mass_difference = mass_final-mass_initial; 779ffc4695bSBarry Smith ierr = MPI_Allreduce(&mass_difference,&mass_differenceg,1,MPIU_SCALAR,MPIU_SUM,comm);CHKERRMPI(ierr); 780c4762a1bSJed Brown ierr = PetscPrintf(comm,"Mass difference %g\n",(double)mass_differenceg);CHKERRQ(ierr); 781c4762a1bSJed Brown ierr = PetscPrintf(comm,"Final time %g, steps %D\n",(double)ptime,steps);CHKERRQ(ierr); 782c4762a1bSJed Brown if (ctx.exact) { 783c4762a1bSJed Brown PetscReal nrm1 = 0; 784c4762a1bSJed Brown ierr = SolutionErrorNorms(&ctx,da,ptime,X,&nrm1);CHKERRQ(ierr); 785c4762a1bSJed Brown ierr = PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1);CHKERRQ(ierr); 786c4762a1bSJed Brown } 787c4762a1bSJed Brown if (ctx.simulation) { 788c4762a1bSJed Brown PetscReal nrm1 = 0; 789c4762a1bSJed Brown PetscViewer fd; 790c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = "binaryoutput"; 791c4762a1bSJed Brown Vec XR; 792c4762a1bSJed Brown PetscBool flg; 793c4762a1bSJed Brown const PetscScalar *ptr_XR; 794c4762a1bSJed Brown ierr = PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg);CHKERRQ(ierr); 795*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 796c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&fd);CHKERRQ(ierr); 797c4762a1bSJed Brown ierr = VecDuplicate(X0,&XR);CHKERRQ(ierr); 798c4762a1bSJed Brown ierr = VecLoad(XR,fd);CHKERRQ(ierr); 799c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 800c4762a1bSJed Brown ierr = VecGetArrayRead(X,&ptr_X);CHKERRQ(ierr); 801c4762a1bSJed Brown ierr = VecGetArrayRead(XR,&ptr_XR);CHKERRQ(ierr); 802c4762a1bSJed Brown for (i=0; i<Mx; i++) { 803c4762a1bSJed Brown if (i < count_slow/2 || i > count_slow/2+count_fast-1) nrm1 = nrm1 + hs*PetscAbs(ptr_X[i]-ptr_XR[i]); 804c4762a1bSJed Brown else nrm1 = nrm1 + hf*PetscAbs(ptr_X[i]-ptr_XR[i]); 805c4762a1bSJed Brown } 806c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&ptr_X);CHKERRQ(ierr); 807c4762a1bSJed Brown ierr = VecRestoreArrayRead(XR,&ptr_XR);CHKERRQ(ierr); 808c4762a1bSJed Brown ierr = PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1);CHKERRQ(ierr); 809c4762a1bSJed Brown ierr = VecDestroy(&XR);CHKERRQ(ierr); 810c4762a1bSJed Brown } 811c4762a1bSJed Brown } 812c4762a1bSJed Brown 813c4762a1bSJed Brown ierr = SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 814c4762a1bSJed Brown if (draw & 0x1) { ierr = VecView(X0,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr); } 815c4762a1bSJed Brown if (draw & 0x2) { ierr = VecView(X,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr); } 816c4762a1bSJed Brown if (draw & 0x4) { 817c4762a1bSJed Brown Vec Y; 818c4762a1bSJed Brown ierr = VecDuplicate(X,&Y);CHKERRQ(ierr); 819c4762a1bSJed Brown ierr = FVSample(&ctx,da,ptime,Y);CHKERRQ(ierr); 820c4762a1bSJed Brown ierr = VecAYPX(Y,-1,X);CHKERRQ(ierr); 821c4762a1bSJed Brown ierr = VecView(Y,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr); 822c4762a1bSJed Brown ierr = VecDestroy(&Y);CHKERRQ(ierr); 823c4762a1bSJed Brown } 824c4762a1bSJed Brown 825c4762a1bSJed Brown if (view_final) { 826c4762a1bSJed Brown PetscViewer viewer; 827c4762a1bSJed Brown ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer);CHKERRQ(ierr); 828c4762a1bSJed Brown ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 829c4762a1bSJed Brown ierr = VecView(X,viewer);CHKERRQ(ierr); 830c4762a1bSJed Brown ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 831c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 832c4762a1bSJed Brown } 833c4762a1bSJed Brown 834c4762a1bSJed Brown /* Clean up */ 835c4762a1bSJed Brown ierr = (*ctx.physics.destroy)(ctx.physics.user);CHKERRQ(ierr); 836c4762a1bSJed Brown for (i=0; i<ctx.physics.dof; i++) {ierr = PetscFree(ctx.physics.fieldname[i]);CHKERRQ(ierr);} 837c4762a1bSJed Brown ierr = PetscFree3(ctx.u,ctx.flux,ctx.speeds);CHKERRQ(ierr); 838c4762a1bSJed Brown ierr = ISDestroy(&ctx.iss);CHKERRQ(ierr); 839c4762a1bSJed Brown ierr = ISDestroy(&ctx.isf);CHKERRQ(ierr); 840c4762a1bSJed Brown ierr = VecDestroy(&X);CHKERRQ(ierr); 841c4762a1bSJed Brown ierr = VecDestroy(&X0);CHKERRQ(ierr); 842c4762a1bSJed Brown ierr = VecDestroy(&R);CHKERRQ(ierr); 843c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 844c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 845c4762a1bSJed Brown ierr = PetscFree(index_slow);CHKERRQ(ierr); 846c4762a1bSJed Brown ierr = PetscFree(index_fast);CHKERRQ(ierr); 847c4762a1bSJed Brown ierr = PetscFunctionListDestroy(&physics);CHKERRQ(ierr); 848c4762a1bSJed Brown ierr = PetscFinalize(); 849c4762a1bSJed Brown return ierr; 850c4762a1bSJed Brown } 851c4762a1bSJed Brown 852c4762a1bSJed Brown /*TEST 853c4762a1bSJed Brown 854c4762a1bSJed Brown build: 855f56ea12dSJed Brown requires: !complex 856c4762a1bSJed Brown 857c4762a1bSJed Brown test: 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 rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 0 859c4762a1bSJed Brown 860c4762a1bSJed Brown test: 861c4762a1bSJed Brown suffix: 2 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 rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 1 863c4762a1bSJed Brown output_file: output/ex7_1.out 864c4762a1bSJed Brown 865c4762a1bSJed Brown test: 866c4762a1bSJed Brown suffix: 3 867c4762a1bSJed 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 868c4762a1bSJed Brown 869c4762a1bSJed Brown test: 870c4762a1bSJed Brown suffix: 4 871c4762a1bSJed 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 872c4762a1bSJed Brown output_file: output/ex7_3.out 873c2a16db8SHong Zhang 874c2a16db8SHong Zhang test: 875c2a16db8SHong Zhang suffix: 5 876c2a16db8SHong Zhang nsize: 2 877c2a16db8SHong 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 878c2a16db8SHong Zhang output_file: output/ex7_3.out 879c4762a1bSJed Brown TEST*/ 880