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