1*c4762a1bSJed Brown /* 2*c4762a1bSJed Brown Note: 3*c4762a1bSJed Brown To performance a convergence study, a reference solution can be generated with a small stepsize and stored into a binary file, e.g. -ts_type ssp -ts_dt 1e-5 -ts_max_steps 30000 -ts_view_solution binary:reference.bin 4*c4762a1bSJed Brown Errors can be computed in the following runs with -simulation -f reference.bin 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown Multirate RK options: 7*c4762a1bSJed Brown -ts_rk_dtratio is the ratio between larger time step size and small time step size 8*c4762a1bSJed Brown -ts_rk_multirate_type has three choices: none (for single step RK) 9*c4762a1bSJed Brown combined (for multirate method and user just need to provide the same RHS with the single step RK) 10*c4762a1bSJed Brown partitioned (for multiraet method and user need to provide two RHS, one is for fast components and the orther is for slow components 11*c4762a1bSJed Brown */ 12*c4762a1bSJed Brown 13*c4762a1bSJed Brown static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n" 14*c4762a1bSJed Brown " advection - Variable coefficient scalar advection\n" 15*c4762a1bSJed Brown " u_t + (a*u)_x = 0\n" 16*c4762a1bSJed Brown " a is a piecewise function with two values say a0 and a1 (both a0 and a1 are constant).\n" 17*c4762a1bSJed Brown " in this toy problem we assume a=a0 when x<0 and a=a1 when x>0 with a0<a1 which has the same discontinuous point with initial condition.\n" 18*c4762a1bSJed Brown " we don't provide an exact solution, so you need to generate reference solution to do convergnece staudy,\n" 19*c4762a1bSJed Brown " more precisely, you need firstly generate a reference to a binary file say file.bin, then on commend line,\n" 20*c4762a1bSJed Brown " you should type -simulation -f file.bin.\n" 21*c4762a1bSJed Brown " you can choose the number of grids by -da_grid_x.\n" 22*c4762a1bSJed Brown " you can choose the value of a by -physics_advect_a1 and -physics_advect_a2.\n"; 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown #include <petscts.h> 25*c4762a1bSJed Brown #include <petscdm.h> 26*c4762a1bSJed Brown #include <petscdmda.h> 27*c4762a1bSJed Brown #include <petscdraw.h> 28*c4762a1bSJed Brown #include <petsc/private/tsimpl.h> 29*c4762a1bSJed Brown 30*c4762a1bSJed Brown #include <petsc/private/kernels/blockinvert.h> /* For the Kernel_*_gets_* stuff for BAIJ */ 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown #include "finitevolume1d.h" 33*c4762a1bSJed Brown 34*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); } 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown /* --------------------------------- Advection ----------------------------------- */ 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown typedef struct { 39*c4762a1bSJed Brown PetscReal a[2]; /* advective velocity */ 40*c4762a1bSJed Brown } AdvectCtx; 41*c4762a1bSJed Brown 42*c4762a1bSJed Brown static PetscErrorCode PhysicsRiemann_Advect(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed,PetscReal x,PetscReal xmin, PetscReal xmax) 43*c4762a1bSJed Brown { 44*c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx*)vctx; 45*c4762a1bSJed Brown PetscReal *speed; 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown PetscFunctionBeginUser; 48*c4762a1bSJed Brown speed = ctx->a; 49*c4762a1bSJed Brown if(x==0 || x == xmin || x == xmax) flux[0] = PetscMax(0,(speed[0]+speed[1])/2.0)*uL[0] + PetscMin(0,(speed[0]+speed[1])/2.0)*uR[0]; /* if condition need to be changed base on different problem, '0' is the discontinuous point of a */ 50*c4762a1bSJed Brown else if(x<0) flux[0] = PetscMax(0,speed[0])*uL[0] + PetscMin(0,speed[0])*uR[0]; /* else if condition need to be changed based on diferent problem, 'x = 0' is discontinuous point of a */ 51*c4762a1bSJed Brown else flux[0] = PetscMax(0,speed[1])*uL[0] + PetscMin(0,speed[1])*uR[0]; 52*c4762a1bSJed Brown *maxspeed = *speed; 53*c4762a1bSJed Brown PetscFunctionReturn(0); 54*c4762a1bSJed Brown } 55*c4762a1bSJed Brown 56*c4762a1bSJed Brown static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds,PetscReal x) 57*c4762a1bSJed Brown { 58*c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx*)vctx; 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown PetscFunctionBeginUser; 61*c4762a1bSJed Brown X[0] = 1.; 62*c4762a1bSJed Brown Xi[0] = 1.; 63*c4762a1bSJed Brown if(x<0) speeds[0] = ctx->a[0]; /* x is discontinuous point of a */ 64*c4762a1bSJed Brown else speeds[0] = ctx->a[1]; 65*c4762a1bSJed Brown PetscFunctionReturn(0); 66*c4762a1bSJed Brown } 67*c4762a1bSJed Brown 68*c4762a1bSJed Brown static PetscErrorCode PhysicsSample_Advect(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u) 69*c4762a1bSJed Brown { 70*c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx*)vctx; 71*c4762a1bSJed Brown PetscReal *a = ctx->a,x0; 72*c4762a1bSJed Brown 73*c4762a1bSJed Brown PetscFunctionBeginUser; 74*c4762a1bSJed Brown if(x<0){ /* x is cell center */ 75*c4762a1bSJed Brown switch (bctype) { 76*c4762a1bSJed Brown case FVBC_OUTFLOW: x0 = x-a[0]*t; break; 77*c4762a1bSJed Brown case FVBC_PERIODIC: x0 = RangeMod(x-a[0]*t,xmin,xmax); break; 78*c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown BCType"); 79*c4762a1bSJed Brown } 80*c4762a1bSJed Brown switch (initial) { 81*c4762a1bSJed Brown case 0: u[0] = (x0 < 0) ? 1 : -1; break; 82*c4762a1bSJed Brown case 1: u[0] = (x0 < 0) ? -1 : 1; break; 83*c4762a1bSJed Brown case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break; 84*c4762a1bSJed Brown case 3: u[0] = PetscSinReal(2*PETSC_PI*x0); break; 85*c4762a1bSJed Brown case 4: u[0] = PetscAbs(x0); break; 86*c4762a1bSJed Brown case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2*PETSC_PI*x0)); break; 87*c4762a1bSJed Brown case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2-x0 : 0)); break; 88*c4762a1bSJed Brown case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI*x0),10.0);break; 89*c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition"); 90*c4762a1bSJed Brown } 91*c4762a1bSJed Brown } 92*c4762a1bSJed Brown else{ 93*c4762a1bSJed Brown switch (bctype) { 94*c4762a1bSJed Brown case FVBC_OUTFLOW: x0 = x-a[1]*t; break; 95*c4762a1bSJed Brown case FVBC_PERIODIC: x0 = RangeMod(x-a[1]*t,xmin,xmax); break; 96*c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown BCType"); 97*c4762a1bSJed Brown } 98*c4762a1bSJed Brown switch (initial) { 99*c4762a1bSJed Brown case 0: u[0] = (x0 < 0) ? 1 : -1; break; 100*c4762a1bSJed Brown case 1: u[0] = (x0 < 0) ? -1 : 1; break; 101*c4762a1bSJed Brown case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break; 102*c4762a1bSJed Brown case 3: u[0] = PetscSinReal(2*PETSC_PI*x0); break; 103*c4762a1bSJed Brown case 4: u[0] = PetscAbs(x0); break; 104*c4762a1bSJed Brown case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2*PETSC_PI*x0)); break; 105*c4762a1bSJed Brown case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2-x0 : 0)); break; 106*c4762a1bSJed Brown case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI*x0),10.0);break; 107*c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition"); 108*c4762a1bSJed Brown } 109*c4762a1bSJed Brown } 110*c4762a1bSJed Brown PetscFunctionReturn(0); 111*c4762a1bSJed Brown } 112*c4762a1bSJed Brown 113*c4762a1bSJed Brown static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx) 114*c4762a1bSJed Brown { 115*c4762a1bSJed Brown PetscErrorCode ierr; 116*c4762a1bSJed Brown AdvectCtx *user; 117*c4762a1bSJed Brown 118*c4762a1bSJed Brown PetscFunctionBeginUser; 119*c4762a1bSJed Brown ierr = PetscNew(&user);CHKERRQ(ierr); 120*c4762a1bSJed Brown ctx->physics.sample = PhysicsSample_Advect; 121*c4762a1bSJed Brown ctx->physics.riemann = PhysicsRiemann_Advect; 122*c4762a1bSJed Brown ctx->physics.characteristic = PhysicsCharacteristic_Advect; 123*c4762a1bSJed Brown ctx->physics.destroy = PhysicsDestroy_SimpleFree; 124*c4762a1bSJed Brown ctx->physics.user = user; 125*c4762a1bSJed Brown ctx->physics.dof = 1; 126*c4762a1bSJed Brown ierr = PetscStrallocpy("u",&ctx->physics.fieldname[0]);CHKERRQ(ierr); 127*c4762a1bSJed Brown user->a[0] = 1; 128*c4762a1bSJed Brown user->a[1] = 1; 129*c4762a1bSJed Brown ierr = PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");CHKERRQ(ierr); 130*c4762a1bSJed Brown { 131*c4762a1bSJed Brown ierr = PetscOptionsReal("-physics_advect_a1","Speed1","",user->a[0],&user->a[0],NULL);CHKERRQ(ierr); 132*c4762a1bSJed Brown ierr = PetscOptionsReal("-physics_advect_a2","Speed2","",user->a[1],&user->a[1],NULL);CHKERRQ(ierr); 133*c4762a1bSJed Brown } 134*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 135*c4762a1bSJed Brown PetscFunctionReturn(0); 136*c4762a1bSJed Brown } 137*c4762a1bSJed Brown 138*c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow parts ----------------------------------- */ 139*c4762a1bSJed Brown 140*c4762a1bSJed Brown PetscErrorCode FVRHSFunctionslow(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 141*c4762a1bSJed Brown { 142*c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 143*c4762a1bSJed Brown PetscErrorCode ierr; 144*c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,len_slow; 145*c4762a1bSJed Brown PetscReal hx,cfl_idt = 0; 146*c4762a1bSJed Brown PetscScalar *x,*f,*slope; 147*c4762a1bSJed Brown Vec Xloc; 148*c4762a1bSJed Brown DM da; 149*c4762a1bSJed Brown 150*c4762a1bSJed Brown PetscFunctionBeginUser; 151*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 152*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 153*c4762a1bSJed Brown ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 154*c4762a1bSJed Brown hx = (ctx->xmax-ctx->xmin)/Mx; 155*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 157*c4762a1bSJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); 158*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 159*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 160*c4762a1bSJed Brown ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 161*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 162*c4762a1bSJed Brown ierr = ISGetSize(ctx->iss,&len_slow);CHKERRQ(ierr); 163*c4762a1bSJed Brown 164*c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 165*c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 166*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 167*c4762a1bSJed Brown } 168*c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 169*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 170*c4762a1bSJed Brown } 171*c4762a1bSJed Brown } 172*c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 173*c4762a1bSJed Brown struct _LimitInfo info; 174*c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 175*c4762a1bSJed Brown if (i < len_slow+1) { 176*c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 177*c4762a1bSJed Brown ierr = (*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds,ctx->xmin+hx*i);CHKERRQ(ierr); 178*c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 179*c4762a1bSJed Brown ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 180*c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 181*c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 182*c4762a1bSJed Brown for (j=0; j<dof; j++) { 183*c4762a1bSJed Brown PetscScalar jmpL,jmpR; 184*c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 185*c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 186*c4762a1bSJed Brown for (k=0; k<dof; k++) { 187*c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 188*c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 189*c4762a1bSJed Brown } 190*c4762a1bSJed Brown } 191*c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 192*c4762a1bSJed Brown info.m = dof; 193*c4762a1bSJed Brown info.hx = hx; 194*c4762a1bSJed Brown (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope); 195*c4762a1bSJed Brown for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 196*c4762a1bSJed Brown for (j=0; j<dof; j++) { 197*c4762a1bSJed Brown PetscScalar tmp = 0; 198*c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 199*c4762a1bSJed Brown slope[i*dof+j] = tmp; 200*c4762a1bSJed Brown } 201*c4762a1bSJed Brown } 202*c4762a1bSJed Brown } 203*c4762a1bSJed Brown 204*c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 205*c4762a1bSJed Brown PetscReal maxspeed; 206*c4762a1bSJed Brown PetscScalar *uL,*uR; 207*c4762a1bSJed Brown if (i < len_slow) { /* slow parts can be changed based on a */ 208*c4762a1bSJed Brown uL = &ctx->uLR[0]; 209*c4762a1bSJed Brown uR = &ctx->uLR[dof]; 210*c4762a1bSJed Brown for (j=0; j<dof; j++) { 211*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 212*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 213*c4762a1bSJed Brown } 214*c4762a1bSJed Brown ierr = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr); 215*c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 216*c4762a1bSJed Brown if (i > xs) { 217*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hx; 218*c4762a1bSJed Brown } 219*c4762a1bSJed Brown if (i < xs+xm) { 220*c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hx; 221*c4762a1bSJed Brown } 222*c4762a1bSJed Brown } 223*c4762a1bSJed Brown if (i == len_slow) { /* slow parts can be changed based on a */ 224*c4762a1bSJed Brown uL = &ctx->uLR[0]; 225*c4762a1bSJed Brown uR = &ctx->uLR[dof]; 226*c4762a1bSJed Brown for (j=0; j<dof; j++) { 227*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 228*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 229*c4762a1bSJed Brown } 230*c4762a1bSJed Brown ierr = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr); 231*c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 232*c4762a1bSJed Brown if (i > xs) { 233*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hx; 234*c4762a1bSJed Brown } 235*c4762a1bSJed Brown } 236*c4762a1bSJed Brown } 237*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 238*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 239*c4762a1bSJed Brown ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 240*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 241*c4762a1bSJed Brown PetscFunctionReturn(0); 242*c4762a1bSJed Brown } 243*c4762a1bSJed Brown 244*c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for fast parts ----------------------------------- */ 245*c4762a1bSJed Brown 246*c4762a1bSJed Brown PetscErrorCode FVRHSFunctionfast(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 247*c4762a1bSJed Brown { 248*c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 249*c4762a1bSJed Brown PetscErrorCode ierr; 250*c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,len_slow; 251*c4762a1bSJed Brown PetscReal hx,cfl_idt = 0; 252*c4762a1bSJed Brown PetscScalar *x,*f,*slope; 253*c4762a1bSJed Brown Vec Xloc; 254*c4762a1bSJed Brown DM da; 255*c4762a1bSJed Brown 256*c4762a1bSJed Brown PetscFunctionBeginUser; 257*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 258*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 259*c4762a1bSJed Brown ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 260*c4762a1bSJed Brown hx = (ctx->xmax-ctx->xmin)/Mx; 261*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 262*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 263*c4762a1bSJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); 264*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 265*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 266*c4762a1bSJed Brown ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 267*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 268*c4762a1bSJed Brown ierr = ISGetSize(ctx->iss,&len_slow);CHKERRQ(ierr); 269*c4762a1bSJed Brown 270*c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 271*c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 272*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 273*c4762a1bSJed Brown } 274*c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 275*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 276*c4762a1bSJed Brown } 277*c4762a1bSJed Brown } 278*c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 279*c4762a1bSJed Brown struct _LimitInfo info; 280*c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 281*c4762a1bSJed Brown if (i > len_slow-2) { 282*c4762a1bSJed Brown ierr = (*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds,ctx->xmin+hx*i);CHKERRQ(ierr); 283*c4762a1bSJed Brown ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 284*c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 285*c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 286*c4762a1bSJed Brown for (j=0; j<dof; j++) { 287*c4762a1bSJed Brown PetscScalar jmpL,jmpR; 288*c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 289*c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 290*c4762a1bSJed Brown for (k=0; k<dof; k++) { 291*c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 292*c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 293*c4762a1bSJed Brown } 294*c4762a1bSJed Brown } 295*c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 296*c4762a1bSJed Brown info.m = dof; 297*c4762a1bSJed Brown info.hx = hx; 298*c4762a1bSJed Brown (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope); 299*c4762a1bSJed Brown for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 300*c4762a1bSJed Brown for (j=0; j<dof; j++) { 301*c4762a1bSJed Brown PetscScalar tmp = 0; 302*c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 303*c4762a1bSJed Brown slope[i*dof+j] = tmp; 304*c4762a1bSJed Brown } 305*c4762a1bSJed Brown } 306*c4762a1bSJed Brown } 307*c4762a1bSJed Brown 308*c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 309*c4762a1bSJed Brown PetscReal maxspeed; 310*c4762a1bSJed Brown PetscScalar *uL,*uR; 311*c4762a1bSJed Brown if (i > len_slow) { 312*c4762a1bSJed Brown uL = &ctx->uLR[0]; 313*c4762a1bSJed Brown uR = &ctx->uLR[dof]; 314*c4762a1bSJed Brown for (j=0; j<dof; j++) { 315*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 316*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 317*c4762a1bSJed Brown } 318*c4762a1bSJed Brown ierr = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr); 319*c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 320*c4762a1bSJed Brown if (i > xs) { 321*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow-1)*dof+j] -= ctx->flux[j]/hx; 322*c4762a1bSJed Brown } 323*c4762a1bSJed Brown if (i < xs+xm) { 324*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow)*dof+j] += ctx->flux[j]/hx; 325*c4762a1bSJed Brown } 326*c4762a1bSJed Brown } 327*c4762a1bSJed Brown if (i == len_slow) { 328*c4762a1bSJed Brown uL = &ctx->uLR[0]; 329*c4762a1bSJed Brown uR = &ctx->uLR[dof]; 330*c4762a1bSJed Brown for (j=0; j<dof; j++) { 331*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 332*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 333*c4762a1bSJed Brown } 334*c4762a1bSJed Brown ierr = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr); 335*c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 336*c4762a1bSJed Brown if (i < xs+xm) { 337*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow)*dof+j] += ctx->flux[j]/hx; 338*c4762a1bSJed Brown } 339*c4762a1bSJed Brown } 340*c4762a1bSJed Brown } 341*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 342*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 343*c4762a1bSJed Brown ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 344*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 345*c4762a1bSJed Brown PetscFunctionReturn(0); 346*c4762a1bSJed Brown } 347*c4762a1bSJed Brown 348*c4762a1bSJed Brown PetscErrorCode FVRHSFunctionslow2(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 349*c4762a1bSJed Brown { 350*c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 351*c4762a1bSJed Brown PetscErrorCode ierr; 352*c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,len_slow1,len_slow2; 353*c4762a1bSJed Brown PetscReal hx,cfl_idt = 0; 354*c4762a1bSJed Brown PetscScalar *x,*f,*slope; 355*c4762a1bSJed Brown Vec Xloc; 356*c4762a1bSJed Brown DM da; 357*c4762a1bSJed Brown 358*c4762a1bSJed Brown PetscFunctionBeginUser; 359*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 360*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 361*c4762a1bSJed Brown ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 362*c4762a1bSJed Brown hx = (ctx->xmax-ctx->xmin)/Mx; 363*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 364*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 365*c4762a1bSJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); 366*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 367*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 368*c4762a1bSJed Brown ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 369*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 370*c4762a1bSJed Brown ierr = ISGetSize(ctx->iss,&len_slow1);CHKERRQ(ierr); 371*c4762a1bSJed Brown ierr = ISGetSize(ctx->iss2,&len_slow2);CHKERRQ(ierr); 372*c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 373*c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 374*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 375*c4762a1bSJed Brown } 376*c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 377*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 378*c4762a1bSJed Brown } 379*c4762a1bSJed Brown } 380*c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 381*c4762a1bSJed Brown struct _LimitInfo info; 382*c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 383*c4762a1bSJed Brown if (i < len_slow1+len_slow2+1 && i > len_slow1-2) { 384*c4762a1bSJed Brown ierr = (*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds,ctx->xmin+hx*i);CHKERRQ(ierr); 385*c4762a1bSJed Brown ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 386*c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 387*c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 388*c4762a1bSJed Brown for (j=0; j<dof; j++) { 389*c4762a1bSJed Brown PetscScalar jmpL,jmpR; 390*c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 391*c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 392*c4762a1bSJed Brown for (k=0; k<dof; k++) { 393*c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 394*c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 395*c4762a1bSJed Brown } 396*c4762a1bSJed Brown } 397*c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 398*c4762a1bSJed Brown info.m = dof; 399*c4762a1bSJed Brown info.hx = hx; 400*c4762a1bSJed Brown (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope); 401*c4762a1bSJed Brown for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 402*c4762a1bSJed Brown for (j=0; j<dof; j++) { 403*c4762a1bSJed Brown PetscScalar tmp = 0; 404*c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 405*c4762a1bSJed Brown slope[i*dof+j] = tmp; 406*c4762a1bSJed Brown } 407*c4762a1bSJed Brown } 408*c4762a1bSJed Brown } 409*c4762a1bSJed Brown 410*c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 411*c4762a1bSJed Brown PetscReal maxspeed; 412*c4762a1bSJed Brown PetscScalar *uL,*uR; 413*c4762a1bSJed Brown if (i < len_slow1+len_slow2 && i > len_slow1) { /* slow parts can be changed based on a */ 414*c4762a1bSJed Brown uL = &ctx->uLR[0]; 415*c4762a1bSJed Brown uR = &ctx->uLR[dof]; 416*c4762a1bSJed Brown for (j=0; j<dof; j++) { 417*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 418*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 419*c4762a1bSJed Brown } 420*c4762a1bSJed Brown ierr = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr); 421*c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 422*c4762a1bSJed Brown if (i > xs) { 423*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow1-1)*dof+j] -= ctx->flux[j]/hx; 424*c4762a1bSJed Brown } 425*c4762a1bSJed Brown if (i < xs+xm) { 426*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow1)*dof+j] += ctx->flux[j]/hx; 427*c4762a1bSJed Brown } 428*c4762a1bSJed Brown } 429*c4762a1bSJed Brown if (i == len_slow1+len_slow2) { /* slow parts can be changed based on a */ 430*c4762a1bSJed Brown uL = &ctx->uLR[0]; 431*c4762a1bSJed Brown uR = &ctx->uLR[dof]; 432*c4762a1bSJed Brown for (j=0; j<dof; j++) { 433*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 434*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 435*c4762a1bSJed Brown } 436*c4762a1bSJed Brown ierr = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr); 437*c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 438*c4762a1bSJed Brown if (i > xs) { 439*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow1-1)*dof+j] -= ctx->flux[j]/hx; 440*c4762a1bSJed Brown } 441*c4762a1bSJed Brown } 442*c4762a1bSJed Brown if (i == len_slow1) { /* slow parts can be changed based on a */ 443*c4762a1bSJed Brown uL = &ctx->uLR[0]; 444*c4762a1bSJed Brown uR = &ctx->uLR[dof]; 445*c4762a1bSJed Brown for (j=0; j<dof; j++) { 446*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 447*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 448*c4762a1bSJed Brown } 449*c4762a1bSJed Brown ierr = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr); 450*c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 451*c4762a1bSJed Brown if (i < xs+xm) { 452*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow1)*dof+j] += ctx->flux[j]/hx; 453*c4762a1bSJed Brown } 454*c4762a1bSJed Brown } 455*c4762a1bSJed Brown } 456*c4762a1bSJed Brown 457*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 458*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 459*c4762a1bSJed Brown ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 460*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 461*c4762a1bSJed Brown PetscFunctionReturn(0); 462*c4762a1bSJed Brown } 463*c4762a1bSJed Brown 464*c4762a1bSJed Brown PetscErrorCode FVRHSFunctionfast2(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 465*c4762a1bSJed Brown { 466*c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 467*c4762a1bSJed Brown PetscErrorCode ierr; 468*c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,len_slow1,len_slow2; 469*c4762a1bSJed Brown PetscReal hx,cfl_idt = 0; 470*c4762a1bSJed Brown PetscScalar *x,*f,*slope; 471*c4762a1bSJed Brown Vec Xloc; 472*c4762a1bSJed Brown DM da; 473*c4762a1bSJed Brown 474*c4762a1bSJed Brown PetscFunctionBeginUser; 475*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 476*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 477*c4762a1bSJed Brown ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 478*c4762a1bSJed Brown hx = (ctx->xmax-ctx->xmin)/Mx; 479*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 480*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 481*c4762a1bSJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); 482*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 483*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 484*c4762a1bSJed Brown ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 485*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 486*c4762a1bSJed Brown ierr = ISGetSize(ctx->iss,&len_slow1);CHKERRQ(ierr); 487*c4762a1bSJed Brown ierr = ISGetSize(ctx->iss2,&len_slow2);CHKERRQ(ierr); 488*c4762a1bSJed Brown 489*c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 490*c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 491*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 492*c4762a1bSJed Brown } 493*c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 494*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 495*c4762a1bSJed Brown } 496*c4762a1bSJed Brown } 497*c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 498*c4762a1bSJed Brown struct _LimitInfo info; 499*c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 500*c4762a1bSJed Brown if (i > len_slow1+len_slow2-2) { 501*c4762a1bSJed Brown ierr = (*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds,ctx->xmin+hx*i);CHKERRQ(ierr); 502*c4762a1bSJed Brown ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 503*c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 504*c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 505*c4762a1bSJed Brown for (j=0; j<dof; j++) { 506*c4762a1bSJed Brown PetscScalar jmpL,jmpR; 507*c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 508*c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 509*c4762a1bSJed Brown for (k=0; k<dof; k++) { 510*c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 511*c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 512*c4762a1bSJed Brown } 513*c4762a1bSJed Brown } 514*c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 515*c4762a1bSJed Brown info.m = dof; 516*c4762a1bSJed Brown info.hx = hx; 517*c4762a1bSJed Brown (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope); 518*c4762a1bSJed Brown for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 519*c4762a1bSJed Brown for (j=0; j<dof; j++) { 520*c4762a1bSJed Brown PetscScalar tmp = 0; 521*c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 522*c4762a1bSJed Brown slope[i*dof+j] = tmp; 523*c4762a1bSJed Brown } 524*c4762a1bSJed Brown } 525*c4762a1bSJed Brown } 526*c4762a1bSJed Brown 527*c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 528*c4762a1bSJed Brown PetscReal maxspeed; 529*c4762a1bSJed Brown PetscScalar *uL,*uR; 530*c4762a1bSJed Brown if (i > len_slow1+len_slow2) { 531*c4762a1bSJed Brown uL = &ctx->uLR[0]; 532*c4762a1bSJed Brown uR = &ctx->uLR[dof]; 533*c4762a1bSJed Brown for (j=0; j<dof; j++) { 534*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 535*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 536*c4762a1bSJed Brown } 537*c4762a1bSJed Brown ierr = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr); 538*c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 539*c4762a1bSJed Brown if (i > xs) { 540*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow1-len_slow2-1)*dof+j] -= ctx->flux[j]/hx; 541*c4762a1bSJed Brown } 542*c4762a1bSJed Brown if (i < xs+xm) { 543*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow1-len_slow2)*dof+j] += ctx->flux[j]/hx; 544*c4762a1bSJed Brown } 545*c4762a1bSJed Brown } 546*c4762a1bSJed Brown if (i == len_slow1+len_slow2) { 547*c4762a1bSJed Brown uL = &ctx->uLR[0]; 548*c4762a1bSJed Brown uR = &ctx->uLR[dof]; 549*c4762a1bSJed Brown for (j=0; j<dof; j++) { 550*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hx/2; 551*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hx/2; 552*c4762a1bSJed Brown } 553*c4762a1bSJed Brown ierr = (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed,ctx->xmin+hx*i,ctx->xmin,ctx->xmax);CHKERRQ(ierr); 554*c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */ 555*c4762a1bSJed Brown if (i < xs+xm) { 556*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-len_slow1-len_slow2)*dof+j] += ctx->flux[j]/hx; 557*c4762a1bSJed Brown } 558*c4762a1bSJed Brown } 559*c4762a1bSJed Brown } 560*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 561*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 562*c4762a1bSJed Brown ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 563*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 564*c4762a1bSJed Brown PetscFunctionReturn(0); 565*c4762a1bSJed Brown } 566*c4762a1bSJed Brown 567*c4762a1bSJed Brown int main(int argc,char *argv[]) 568*c4762a1bSJed Brown { 569*c4762a1bSJed Brown char lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m"; 570*c4762a1bSJed Brown PetscFunctionList limiters = 0,physics = 0; 571*c4762a1bSJed Brown MPI_Comm comm; 572*c4762a1bSJed Brown TS ts; 573*c4762a1bSJed Brown DM da; 574*c4762a1bSJed Brown Vec X,X0,R; 575*c4762a1bSJed Brown FVCtx ctx; 576*c4762a1bSJed Brown PetscInt i,k,dof,xs,xm,Mx,draw = 0,*index_slow,*index_fast,islow = 0,ifast = 0; 577*c4762a1bSJed Brown PetscBool view_final = PETSC_FALSE; 578*c4762a1bSJed Brown PetscReal ptime; 579*c4762a1bSJed Brown PetscErrorCode ierr; 580*c4762a1bSJed Brown 581*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr; 582*c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 583*c4762a1bSJed Brown ierr = PetscMemzero(&ctx,sizeof(ctx));CHKERRQ(ierr); 584*c4762a1bSJed Brown 585*c4762a1bSJed Brown /* Register limiters to be available on the command line */ 586*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"upwind" ,Limit_Upwind);CHKERRQ(ierr); 587*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"lax-wendroff" ,Limit_LaxWendroff);CHKERRQ(ierr); 588*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"beam-warming" ,Limit_BeamWarming);CHKERRQ(ierr); 589*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"fromm" ,Limit_Fromm);CHKERRQ(ierr); 590*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"minmod" ,Limit_Minmod);CHKERRQ(ierr); 591*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"superbee" ,Limit_Superbee);CHKERRQ(ierr); 592*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"mc" ,Limit_MC);CHKERRQ(ierr); 593*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"vanleer" ,Limit_VanLeer);CHKERRQ(ierr); 594*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"vanalbada" ,Limit_VanAlbada);CHKERRQ(ierr); 595*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"vanalbadatvd" ,Limit_VanAlbadaTVD);CHKERRQ(ierr); 596*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"koren" ,Limit_Koren);CHKERRQ(ierr); 597*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"korensym" ,Limit_KorenSym);CHKERRQ(ierr); 598*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"koren3" ,Limit_Koren3);CHKERRQ(ierr); 599*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"cada-torrilhon2" ,Limit_CadaTorrilhon2);CHKERRQ(ierr); 600*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"cada-torrilhon3-r0p1",Limit_CadaTorrilhon3R0p1);CHKERRQ(ierr); 601*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"cada-torrilhon3-r1" ,Limit_CadaTorrilhon3R1);CHKERRQ(ierr); 602*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"cada-torrilhon3-r10" ,Limit_CadaTorrilhon3R10);CHKERRQ(ierr); 603*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"cada-torrilhon3-r100",Limit_CadaTorrilhon3R100);CHKERRQ(ierr); 604*c4762a1bSJed Brown 605*c4762a1bSJed Brown /* Register physical models to be available on the command line */ 606*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&physics,"advect" ,PhysicsCreate_Advect);CHKERRQ(ierr); 607*c4762a1bSJed Brown 608*c4762a1bSJed Brown ctx.comm = comm; 609*c4762a1bSJed Brown ctx.cfl = 0.9; 610*c4762a1bSJed Brown ctx.bctype = FVBC_PERIODIC; 611*c4762a1bSJed Brown ctx.xmin = -1.0; 612*c4762a1bSJed Brown ctx.xmax = 1.0; 613*c4762a1bSJed Brown ierr = PetscOptionsBegin(comm,NULL,"Finite Volume solver options","");CHKERRQ(ierr); 614*c4762a1bSJed Brown ierr = PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL);CHKERRQ(ierr); 615*c4762a1bSJed Brown ierr = PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL);CHKERRQ(ierr); 616*c4762a1bSJed Brown ierr = PetscOptionsFList("-limit","Name of flux limiter to use","",limiters,lname,lname,sizeof(lname),NULL);CHKERRQ(ierr); 617*c4762a1bSJed Brown ierr = PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL);CHKERRQ(ierr); 618*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); 619*c4762a1bSJed Brown ierr = PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL);CHKERRQ(ierr); 620*c4762a1bSJed Brown ierr = PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL);CHKERRQ(ierr); 621*c4762a1bSJed Brown ierr = PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL);CHKERRQ(ierr); 622*c4762a1bSJed Brown ierr = PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL);CHKERRQ(ierr); 623*c4762a1bSJed Brown ierr = PetscOptionsBool("-recursive_split","Split the domain recursively","",ctx.recursive,&ctx.recursive,NULL);CHKERRQ(ierr); 624*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 625*c4762a1bSJed Brown 626*c4762a1bSJed Brown /* Choose the limiter from the list of registered limiters */ 627*c4762a1bSJed Brown ierr = PetscFunctionListFind(limiters,lname,&ctx.limit);CHKERRQ(ierr); 628*c4762a1bSJed Brown if (!ctx.limit) SETERRQ1(PETSC_COMM_SELF,1,"Limiter '%s' not found",lname); 629*c4762a1bSJed Brown 630*c4762a1bSJed Brown /* Choose the physics from the list of registered models */ 631*c4762a1bSJed Brown { 632*c4762a1bSJed Brown PetscErrorCode (*r)(FVCtx*); 633*c4762a1bSJed Brown ierr = PetscFunctionListFind(physics,physname,&r);CHKERRQ(ierr); 634*c4762a1bSJed Brown if (!r) SETERRQ1(PETSC_COMM_SELF,1,"Physics '%s' not found",physname); 635*c4762a1bSJed Brown /* Create the physics, will set the number of fields and their names */ 636*c4762a1bSJed Brown ierr = (*r)(&ctx);CHKERRQ(ierr); 637*c4762a1bSJed Brown } 638*c4762a1bSJed Brown 639*c4762a1bSJed Brown /* Create a DMDA to manage the parallel grid */ 640*c4762a1bSJed Brown ierr = DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics.dof,2,NULL,&da);CHKERRQ(ierr); 641*c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 642*c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 643*c4762a1bSJed Brown /* Inform the DMDA of the field names provided by the physics. */ 644*c4762a1bSJed Brown /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ 645*c4762a1bSJed Brown for (i=0; i<ctx.physics.dof; i++) { 646*c4762a1bSJed Brown ierr = DMDASetFieldName(da,i,ctx.physics.fieldname[i]);CHKERRQ(ierr); 647*c4762a1bSJed Brown } 648*c4762a1bSJed Brown ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 649*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 650*c4762a1bSJed Brown 651*c4762a1bSJed Brown /* Set coordinates of cell centers */ 652*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); 653*c4762a1bSJed Brown 654*c4762a1bSJed Brown /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */ 655*c4762a1bSJed Brown ierr = PetscMalloc4(dof*dof,&ctx.R,dof*dof,&ctx.Rinv,2*dof,&ctx.cjmpLR,1*dof,&ctx.cslope);CHKERRQ(ierr); 656*c4762a1bSJed Brown ierr = PetscMalloc3(2*dof,&ctx.uLR,dof,&ctx.flux,dof,&ctx.speeds);CHKERRQ(ierr); 657*c4762a1bSJed Brown 658*c4762a1bSJed Brown /* Create a vector to store the solution and to save the initial state */ 659*c4762a1bSJed Brown ierr = DMCreateGlobalVector(da,&X);CHKERRQ(ierr); 660*c4762a1bSJed Brown ierr = VecDuplicate(X,&X0);CHKERRQ(ierr); 661*c4762a1bSJed Brown ierr = VecDuplicate(X,&R);CHKERRQ(ierr); 662*c4762a1bSJed Brown 663*c4762a1bSJed Brown /*-------------------------------- create index for slow parts and fast parts ----------------------------------------*/ 664*c4762a1bSJed Brown ierr = PetscMalloc1(xm*dof,&index_slow);CHKERRQ(ierr); 665*c4762a1bSJed Brown ierr = PetscMalloc1(xm*dof,&index_fast);CHKERRQ(ierr); 666*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 667*c4762a1bSJed Brown if (ctx.xmin+i*(ctx.xmax-ctx.xmin)/(PetscReal)Mx+0.5*(ctx.xmax-ctx.xmin)/(PetscReal)Mx < 0) 668*c4762a1bSJed Brown for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k; 669*c4762a1bSJed Brown else 670*c4762a1bSJed Brown for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k; 671*c4762a1bSJed Brown } /* this step need to be changed based on discontinuous point of a */ 672*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss);CHKERRQ(ierr); 673*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf);CHKERRQ(ierr); 674*c4762a1bSJed Brown 675*c4762a1bSJed Brown /* Create a time-stepping object */ 676*c4762a1bSJed Brown ierr = TSCreate(comm,&ts);CHKERRQ(ierr); 677*c4762a1bSJed Brown ierr = TSSetDM(ts,da);CHKERRQ(ierr); 678*c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,R,FVRHSFunction,&ctx);CHKERRQ(ierr); 679*c4762a1bSJed Brown ierr = TSRHSSplitSetIS(ts,"slow",ctx.iss);CHKERRQ(ierr); 680*c4762a1bSJed Brown ierr = TSRHSSplitSetIS(ts,"fast",ctx.isf);CHKERRQ(ierr); 681*c4762a1bSJed Brown ierr = TSRHSSplitSetRHSFunction(ts,"slow",NULL,FVRHSFunctionslow,&ctx);CHKERRQ(ierr); 682*c4762a1bSJed Brown ierr = TSRHSSplitSetRHSFunction(ts,"fast",NULL,FVRHSFunctionfast,&ctx);CHKERRQ(ierr); 683*c4762a1bSJed Brown 684*c4762a1bSJed Brown if (ctx.recursive) { 685*c4762a1bSJed Brown TS subts; 686*c4762a1bSJed Brown islow = 0; 687*c4762a1bSJed Brown ifast = 0; 688*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 689*c4762a1bSJed Brown PetscReal coord = ctx.xmin+i*(ctx.xmax-ctx.xmin)/(PetscReal)Mx+0.5*(ctx.xmax-ctx.xmin)/(PetscReal)Mx; 690*c4762a1bSJed Brown if (coord >= 0 && coord < ctx.xmin+(ctx.xmax-ctx.xmin)*3/4.) 691*c4762a1bSJed Brown for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k; 692*c4762a1bSJed Brown if (coord >= ctx.xmin+(ctx.xmax-ctx.xmin)*3/4.) 693*c4762a1bSJed Brown for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k; 694*c4762a1bSJed Brown } /* this step need to be changed based on discontinuous point of a */ 695*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss2);CHKERRQ(ierr); 696*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf2);CHKERRQ(ierr); 697*c4762a1bSJed Brown 698*c4762a1bSJed Brown ierr = TSRHSSplitGetSubTS(ts,"fast",&subts);CHKERRQ(ierr); 699*c4762a1bSJed Brown ierr = TSRHSSplitSetIS(subts,"slow",ctx.iss2);CHKERRQ(ierr); 700*c4762a1bSJed Brown ierr = TSRHSSplitSetIS(subts,"fast",ctx.isf2);CHKERRQ(ierr); 701*c4762a1bSJed Brown ierr = TSRHSSplitSetRHSFunction(subts,"slow",NULL,FVRHSFunctionslow2,&ctx);CHKERRQ(ierr); 702*c4762a1bSJed Brown ierr = TSRHSSplitSetRHSFunction(subts,"fast",NULL,FVRHSFunctionfast2,&ctx);CHKERRQ(ierr); 703*c4762a1bSJed Brown } 704*c4762a1bSJed Brown 705*c4762a1bSJed Brown /*ierr = TSSetType(ts,TSSSP);CHKERRQ(ierr);*/ 706*c4762a1bSJed Brown ierr = TSSetType(ts,TSMPRK);CHKERRQ(ierr); 707*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,10);CHKERRQ(ierr); 708*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 709*c4762a1bSJed Brown 710*c4762a1bSJed Brown /* Compute initial conditions and starting time step */ 711*c4762a1bSJed Brown ierr = FVSample(&ctx,da,0,X0);CHKERRQ(ierr); 712*c4762a1bSJed Brown ierr = FVRHSFunction(ts,0,X0,X,(void*)&ctx);CHKERRQ(ierr); /* Initial function evaluation, only used to determine max speed */ 713*c4762a1bSJed Brown ierr = VecCopy(X0,X);CHKERRQ(ierr); /* The function value was not used so we set X=X0 again */ 714*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt);CHKERRQ(ierr); 715*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); /* Take runtime options */ 716*c4762a1bSJed Brown ierr = SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 717*c4762a1bSJed Brown { 718*c4762a1bSJed Brown PetscInt steps; 719*c4762a1bSJed Brown PetscScalar mass_initial,mass_final,mass_difference; 720*c4762a1bSJed Brown 721*c4762a1bSJed Brown ierr = TSSolve(ts,X);CHKERRQ(ierr); 722*c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&ptime);CHKERRQ(ierr); 723*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 724*c4762a1bSJed Brown ierr = PetscPrintf(comm,"Final time %g, steps %D\n",(double)ptime,steps);CHKERRQ(ierr); 725*c4762a1bSJed Brown /* calculate the total mass at initial time and final time */ 726*c4762a1bSJed Brown mass_initial = 0.0; 727*c4762a1bSJed Brown mass_final = 0.0; 728*c4762a1bSJed Brown ierr = VecSum(X0,&mass_initial);CHKERRQ(ierr); 729*c4762a1bSJed Brown ierr = VecSum(X,&mass_final);CHKERRQ(ierr); 730*c4762a1bSJed Brown mass_difference = (ctx.xmax-ctx.xmin)/(PetscScalar)Mx*(mass_final - mass_initial); 731*c4762a1bSJed Brown ierr = PetscPrintf(comm,"Mass difference %g\n",(double)mass_difference);CHKERRQ(ierr); 732*c4762a1bSJed Brown if (ctx.simulation) { 733*c4762a1bSJed Brown PetscViewer fd; 734*c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = "binaryoutput"; 735*c4762a1bSJed Brown Vec XR; 736*c4762a1bSJed Brown PetscReal nrm1,nrmsup; 737*c4762a1bSJed Brown PetscBool flg; 738*c4762a1bSJed Brown ierr = PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg);CHKERRQ(ierr); 739*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 740*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&fd);CHKERRQ(ierr); 741*c4762a1bSJed Brown ierr = VecDuplicate(X0,&XR);CHKERRQ(ierr); 742*c4762a1bSJed Brown ierr = VecLoad(XR,fd);CHKERRQ(ierr); 743*c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 744*c4762a1bSJed Brown ierr = VecAYPX(XR,-1,X);CHKERRQ(ierr); 745*c4762a1bSJed Brown ierr = VecNorm(XR,NORM_1,&nrm1);CHKERRQ(ierr); 746*c4762a1bSJed Brown ierr = VecNorm(XR,NORM_INFINITY,&nrmsup);CHKERRQ(ierr); 747*c4762a1bSJed Brown nrm1 /= Mx; 748*c4762a1bSJed Brown ierr = PetscPrintf(comm,"Error ||x-x_e||_1 %g ||x-x_e||_sup %g\n",(double)nrm1,(double)nrmsup);CHKERRQ(ierr); 749*c4762a1bSJed Brown ierr = VecDestroy(&XR);CHKERRQ(ierr); 750*c4762a1bSJed Brown } 751*c4762a1bSJed Brown } 752*c4762a1bSJed Brown 753*c4762a1bSJed Brown ierr = SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 754*c4762a1bSJed Brown if (draw & 0x1) {ierr = VecView(X0,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);} 755*c4762a1bSJed Brown if (draw & 0x2) {ierr = VecView(X,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);} 756*c4762a1bSJed Brown if (draw & 0x4) { 757*c4762a1bSJed Brown Vec Y; 758*c4762a1bSJed Brown ierr = VecDuplicate(X,&Y);CHKERRQ(ierr); 759*c4762a1bSJed Brown ierr = FVSample(&ctx,da,ptime,Y);CHKERRQ(ierr); 760*c4762a1bSJed Brown ierr = VecAYPX(Y,-1,X);CHKERRQ(ierr); 761*c4762a1bSJed Brown ierr = VecView(Y,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr); 762*c4762a1bSJed Brown ierr = VecDestroy(&Y);CHKERRQ(ierr); 763*c4762a1bSJed Brown } 764*c4762a1bSJed Brown 765*c4762a1bSJed Brown if (view_final) { 766*c4762a1bSJed Brown PetscViewer viewer; 767*c4762a1bSJed Brown ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer);CHKERRQ(ierr); 768*c4762a1bSJed Brown ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 769*c4762a1bSJed Brown ierr = VecView(X,viewer);CHKERRQ(ierr); 770*c4762a1bSJed Brown ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 771*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 772*c4762a1bSJed Brown } 773*c4762a1bSJed Brown 774*c4762a1bSJed Brown /* Clean up */ 775*c4762a1bSJed Brown ierr = (*ctx.physics.destroy)(ctx.physics.user);CHKERRQ(ierr); 776*c4762a1bSJed Brown for (i=0; i<ctx.physics.dof; i++) {ierr = PetscFree(ctx.physics.fieldname[i]);CHKERRQ(ierr);} 777*c4762a1bSJed Brown ierr = PetscFree4(ctx.R,ctx.Rinv,ctx.cjmpLR,ctx.cslope);CHKERRQ(ierr); 778*c4762a1bSJed Brown ierr = PetscFree3(ctx.uLR,ctx.flux,ctx.speeds);CHKERRQ(ierr); 779*c4762a1bSJed Brown ierr = ISDestroy(&ctx.iss);CHKERRQ(ierr); 780*c4762a1bSJed Brown ierr = ISDestroy(&ctx.isf);CHKERRQ(ierr); 781*c4762a1bSJed Brown ierr = ISDestroy(&ctx.iss2);CHKERRQ(ierr); 782*c4762a1bSJed Brown ierr = ISDestroy(&ctx.isf2);CHKERRQ(ierr); 783*c4762a1bSJed Brown ierr = VecDestroy(&X);CHKERRQ(ierr); 784*c4762a1bSJed Brown ierr = VecDestroy(&X0);CHKERRQ(ierr); 785*c4762a1bSJed Brown ierr = VecDestroy(&R);CHKERRQ(ierr); 786*c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 787*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 788*c4762a1bSJed Brown ierr = PetscFree(index_slow);CHKERRQ(ierr); 789*c4762a1bSJed Brown ierr = PetscFree(index_fast);CHKERRQ(ierr); 790*c4762a1bSJed Brown ierr = PetscFunctionListDestroy(&limiters);CHKERRQ(ierr); 791*c4762a1bSJed Brown ierr = PetscFunctionListDestroy(&physics);CHKERRQ(ierr); 792*c4762a1bSJed Brown ierr = PetscFinalize(); 793*c4762a1bSJed Brown return ierr; 794*c4762a1bSJed Brown } 795*c4762a1bSJed Brown 796*c4762a1bSJed Brown /*TEST 797*c4762a1bSJed Brown build: 798*c4762a1bSJed Brown requires: !complex c99 799*c4762a1bSJed Brown depends: finitevolume1d.c 800*c4762a1bSJed Brown 801*c4762a1bSJed Brown test: 802*c4762a1bSJed Brown args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 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 803*c4762a1bSJed Brown 804*c4762a1bSJed Brown test: 805*c4762a1bSJed Brown suffix: 2 806*c4762a1bSJed Brown args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 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 807*c4762a1bSJed Brown output_file: output/ex5_1.out 808*c4762a1bSJed Brown 809*c4762a1bSJed Brown test: 810*c4762a1bSJed Brown suffix: 3 811*c4762a1bSJed Brown args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 0 812*c4762a1bSJed Brown 813*c4762a1bSJed Brown test: 814*c4762a1bSJed Brown suffix: 4 815*c4762a1bSJed Brown args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 1 816*c4762a1bSJed Brown output_file: output/ex5_3.out 817*c4762a1bSJed Brown 818*c4762a1bSJed Brown test: 819*c4762a1bSJed Brown suffix: 5 820*c4762a1bSJed Brown args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 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 -recursive_split 821*c4762a1bSJed Brown 822*c4762a1bSJed Brown test: 823*c4762a1bSJed Brown suffix: 6 824*c4762a1bSJed Brown args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 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 -recursive_split 825*c4762a1bSJed Brown TEST*/ 826