1*c4762a1bSJed Brown /* 2*c4762a1bSJed Brown Note: 3*c4762a1bSJed Brown -hratio is the ratio between mesh size of corse grids and fine grids 4*c4762a1bSJed Brown -ts_rk_dtratio is the ratio between the large stepsize and the small stepsize 5*c4762a1bSJed Brown */ 6*c4762a1bSJed Brown 7*c4762a1bSJed Brown static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n" 8*c4762a1bSJed Brown " advection - Constant coefficient scalar advection\n" 9*c4762a1bSJed Brown " u_t + (a*u)_x = 0\n" 10*c4762a1bSJed Brown " for this toy problem, we choose different meshsizes for different sub-domains (slow-fast-slow), say\n" 11*c4762a1bSJed Brown " hxs = hratio*hxf \n" 12*c4762a1bSJed Brown " where hxs and hxf are the grid spacings for coarse and fine grids respectively.\n" 13*c4762a1bSJed Brown " exact - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n" 14*c4762a1bSJed Brown " the states across shocks and rarefactions\n" 15*c4762a1bSJed Brown " simulation - use reference solution which is generated by smaller time step size to be true solution,\n" 16*c4762a1bSJed Brown " also the reference solution should be generated by user and stored in a binary file.\n" 17*c4762a1bSJed Brown " characteristic - Limit the characteristic variables, this is usually preferred (default)\n" 18*c4762a1bSJed Brown "Several initial conditions can be chosen with -initial N\n\n" 19*c4762a1bSJed Brown "The problem size should be set with -da_grid_x M\n\n"; 20*c4762a1bSJed Brown 21*c4762a1bSJed Brown #include <petscts.h> 22*c4762a1bSJed Brown #include <petscdm.h> 23*c4762a1bSJed Brown #include <petscdmda.h> 24*c4762a1bSJed Brown #include <petscdraw.h> 25*c4762a1bSJed Brown #include "finitevolume1d.h" 26*c4762a1bSJed Brown 27*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); } 28*c4762a1bSJed Brown 29*c4762a1bSJed Brown /* --------------------------------- Advection ----------------------------------- */ 30*c4762a1bSJed Brown typedef struct { 31*c4762a1bSJed Brown PetscReal a; /* advective velocity */ 32*c4762a1bSJed Brown } AdvectCtx; 33*c4762a1bSJed Brown 34*c4762a1bSJed Brown static PetscErrorCode PhysicsRiemann_Advect(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed) 35*c4762a1bSJed Brown { 36*c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx*)vctx; 37*c4762a1bSJed Brown PetscReal speed; 38*c4762a1bSJed Brown 39*c4762a1bSJed Brown PetscFunctionBeginUser; 40*c4762a1bSJed Brown speed = ctx->a; 41*c4762a1bSJed Brown flux[0] = PetscMax(0,speed)*uL[0] + PetscMin(0,speed)*uR[0]; 42*c4762a1bSJed Brown *maxspeed = speed; 43*c4762a1bSJed Brown PetscFunctionReturn(0); 44*c4762a1bSJed Brown } 45*c4762a1bSJed Brown 46*c4762a1bSJed Brown static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds) 47*c4762a1bSJed Brown { 48*c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx*)vctx; 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown PetscFunctionBeginUser; 51*c4762a1bSJed Brown X[0] = 1.; 52*c4762a1bSJed Brown Xi[0] = 1.; 53*c4762a1bSJed Brown speeds[0] = ctx->a; 54*c4762a1bSJed Brown PetscFunctionReturn(0); 55*c4762a1bSJed Brown } 56*c4762a1bSJed Brown 57*c4762a1bSJed Brown static PetscErrorCode PhysicsSample_Advect(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u) 58*c4762a1bSJed Brown { 59*c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx*)vctx; 60*c4762a1bSJed Brown PetscReal a = ctx->a,x0; 61*c4762a1bSJed Brown 62*c4762a1bSJed Brown PetscFunctionBeginUser; 63*c4762a1bSJed Brown switch (bctype) { 64*c4762a1bSJed Brown case FVBC_OUTFLOW: x0 = x-a*t; break; 65*c4762a1bSJed Brown case FVBC_PERIODIC: x0 = RangeMod(x-a*t,xmin,xmax); break; 66*c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown BCType"); 67*c4762a1bSJed Brown } 68*c4762a1bSJed Brown switch (initial) { 69*c4762a1bSJed Brown case 0: u[0] = (x0 < 0) ? 1 : -1; break; 70*c4762a1bSJed Brown case 1: u[0] = (x0 < 0) ? -1 : 1; break; 71*c4762a1bSJed Brown case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break; 72*c4762a1bSJed Brown case 3: u[0] = PetscSinReal(2*PETSC_PI*x0); break; 73*c4762a1bSJed Brown case 4: u[0] = PetscAbs(x0); break; 74*c4762a1bSJed Brown case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2*PETSC_PI*x0)); break; 75*c4762a1bSJed Brown case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2-x0 : 0)); break; 76*c4762a1bSJed Brown case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI*x0),10.0);break; 77*c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition"); 78*c4762a1bSJed Brown } 79*c4762a1bSJed Brown PetscFunctionReturn(0); 80*c4762a1bSJed Brown } 81*c4762a1bSJed Brown 82*c4762a1bSJed Brown static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx) 83*c4762a1bSJed Brown { 84*c4762a1bSJed Brown PetscErrorCode ierr; 85*c4762a1bSJed Brown AdvectCtx *user; 86*c4762a1bSJed Brown 87*c4762a1bSJed Brown PetscFunctionBeginUser; 88*c4762a1bSJed Brown ierr = PetscNew(&user);CHKERRQ(ierr); 89*c4762a1bSJed Brown ctx->physics2.sample2 = PhysicsSample_Advect; 90*c4762a1bSJed Brown ctx->physics2.riemann2 = PhysicsRiemann_Advect; 91*c4762a1bSJed Brown ctx->physics2.characteristic2 = PhysicsCharacteristic_Advect; 92*c4762a1bSJed Brown ctx->physics2.destroy = PhysicsDestroy_SimpleFree; 93*c4762a1bSJed Brown ctx->physics2.user = user; 94*c4762a1bSJed Brown ctx->physics2.dof = 1; 95*c4762a1bSJed Brown ierr = PetscStrallocpy("u",&ctx->physics2.fieldname[0]);CHKERRQ(ierr); 96*c4762a1bSJed Brown user->a = 1; 97*c4762a1bSJed Brown ierr = PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");CHKERRQ(ierr); 98*c4762a1bSJed Brown { 99*c4762a1bSJed Brown ierr = PetscOptionsReal("-physics_advect_a","Speed","",user->a,&user->a,NULL);CHKERRQ(ierr); 100*c4762a1bSJed Brown } 101*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 102*c4762a1bSJed Brown PetscFunctionReturn(0); 103*c4762a1bSJed Brown } 104*c4762a1bSJed Brown 105*c4762a1bSJed Brown PetscErrorCode FVSample_2WaySplit(FVCtx *ctx,DM da,PetscReal time,Vec U) 106*c4762a1bSJed Brown { 107*c4762a1bSJed Brown PetscErrorCode ierr; 108*c4762a1bSJed Brown PetscScalar *u,*uj,xj,xi; 109*c4762a1bSJed Brown PetscInt i,j,k,dof,xs,xm,Mx; 110*c4762a1bSJed Brown const PetscInt N = 200; 111*c4762a1bSJed Brown PetscReal hs,hf; 112*c4762a1bSJed Brown 113*c4762a1bSJed Brown PetscFunctionBeginUser; 114*c4762a1bSJed Brown if (!ctx->physics2.sample2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Physics has not provided a sampling function"); 115*c4762a1bSJed Brown ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 116*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 117*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr); 118*c4762a1bSJed Brown ierr = PetscMalloc1(dof,&uj);CHKERRQ(ierr); 119*c4762a1bSJed Brown hs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 120*c4762a1bSJed Brown hf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 121*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 122*c4762a1bSJed Brown if (i < ctx->sf) { 123*c4762a1bSJed Brown xi = ctx->xmin+0.5*hs+i*hs; 124*c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 125*c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 126*c4762a1bSJed Brown for (j=0; j<N+1; j++) { 127*c4762a1bSJed Brown xj = xi+hs*(j-N/2)/(PetscReal)N; 128*c4762a1bSJed Brown ierr = (*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); 129*c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 130*c4762a1bSJed Brown } 131*c4762a1bSJed Brown } else if (i < ctx->fs) { 132*c4762a1bSJed Brown xi = ctx->xmin+ctx->sf*hs+0.5*hf+(i-ctx->sf)*hf; 133*c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 134*c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 135*c4762a1bSJed Brown for (j=0; j<N+1; j++) { 136*c4762a1bSJed Brown xj = xi+hf*(j-N/2)/(PetscReal)N; 137*c4762a1bSJed Brown ierr = (*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); 138*c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 139*c4762a1bSJed Brown } 140*c4762a1bSJed Brown } else { 141*c4762a1bSJed Brown xi = ctx->xmin+ctx->sf*hs+(ctx->fs-ctx->sf)*hf+0.5*hs+(i-ctx->fs)*hs; 142*c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 143*c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 144*c4762a1bSJed Brown for (j=0; j<N+1; j++) { 145*c4762a1bSJed Brown xj = xi+hs*(j-N/2)/(PetscReal)N; 146*c4762a1bSJed Brown ierr = (*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); 147*c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 148*c4762a1bSJed Brown } 149*c4762a1bSJed Brown } 150*c4762a1bSJed Brown } 151*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr); 152*c4762a1bSJed Brown ierr = PetscFree(uj);CHKERRQ(ierr); 153*c4762a1bSJed Brown PetscFunctionReturn(0); 154*c4762a1bSJed Brown } 155*c4762a1bSJed Brown 156*c4762a1bSJed Brown static PetscErrorCode SolutionErrorNorms_2WaySplit(FVCtx *ctx,DM da,PetscReal t,Vec X,PetscReal *nrm1) 157*c4762a1bSJed Brown { 158*c4762a1bSJed Brown PetscErrorCode ierr; 159*c4762a1bSJed Brown Vec Y; 160*c4762a1bSJed Brown PetscInt i,Mx; 161*c4762a1bSJed Brown const PetscScalar *ptr_X,*ptr_Y; 162*c4762a1bSJed Brown PetscReal hs,hf; 163*c4762a1bSJed Brown 164*c4762a1bSJed Brown PetscFunctionBeginUser; 165*c4762a1bSJed Brown ierr = VecGetSize(X,&Mx);CHKERRQ(ierr); 166*c4762a1bSJed Brown ierr = VecDuplicate(X,&Y);CHKERRQ(ierr); 167*c4762a1bSJed Brown ierr = FVSample_2WaySplit(ctx,da,t,Y);CHKERRQ(ierr); 168*c4762a1bSJed Brown hs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 169*c4762a1bSJed Brown hf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 170*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&ptr_X);CHKERRQ(ierr); 171*c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&ptr_Y);CHKERRQ(ierr); 172*c4762a1bSJed Brown for (i=0; i<Mx; i++) { 173*c4762a1bSJed Brown if (i < ctx->sf || i > ctx->fs -1) *nrm1 += hs*PetscAbs(ptr_X[i]-ptr_Y[i]); 174*c4762a1bSJed Brown else *nrm1 += hf*PetscAbs(ptr_X[i]-ptr_Y[i]); 175*c4762a1bSJed Brown } 176*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&ptr_X);CHKERRQ(ierr); 177*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&ptr_Y);CHKERRQ(ierr); 178*c4762a1bSJed Brown ierr = VecDestroy(&Y);CHKERRQ(ierr); 179*c4762a1bSJed Brown PetscFunctionReturn(0); 180*c4762a1bSJed Brown } 181*c4762a1bSJed Brown 182*c4762a1bSJed Brown PetscErrorCode FVRHSFunction_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 183*c4762a1bSJed Brown { 184*c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 185*c4762a1bSJed Brown PetscErrorCode ierr; 186*c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,sf = ctx->sf,fs = ctx->fs; 187*c4762a1bSJed Brown PetscReal hxf,hxs,cfl_idt = 0; 188*c4762a1bSJed Brown PetscScalar *x,*f,*slope; 189*c4762a1bSJed Brown Vec Xloc; 190*c4762a1bSJed Brown DM da; 191*c4762a1bSJed Brown 192*c4762a1bSJed Brown PetscFunctionBeginUser; 193*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 194*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); /* Xloc contains ghost points */ 195*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 */ 196*c4762a1bSJed Brown hxs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 197*c4762a1bSJed Brown hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 198*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); /* X is solution vector which does not contain ghost points */ 199*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 200*c4762a1bSJed Brown 201*c4762a1bSJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); /* F is the right hand side function corresponds to center points */ 202*c4762a1bSJed Brown 203*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 204*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr); 205*c4762a1bSJed Brown ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); /* contains ghost points */ 206*c4762a1bSJed Brown 207*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 208*c4762a1bSJed Brown 209*c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 210*c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 211*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 212*c4762a1bSJed Brown } 213*c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 214*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 215*c4762a1bSJed Brown } 216*c4762a1bSJed Brown } 217*c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 218*c4762a1bSJed Brown struct _LimitInfo info; 219*c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 220*c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 221*c4762a1bSJed Brown ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); 222*c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 223*c4762a1bSJed Brown ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 224*c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 225*c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 226*c4762a1bSJed Brown for (j=0; j<dof; j++) { 227*c4762a1bSJed Brown PetscScalar jmpL,jmpR; 228*c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 229*c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 230*c4762a1bSJed Brown for (k=0; k<dof; k++) { 231*c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 232*c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 233*c4762a1bSJed Brown } 234*c4762a1bSJed Brown } 235*c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 236*c4762a1bSJed Brown info.m = dof; 237*c4762a1bSJed Brown info.hxs = hxs; 238*c4762a1bSJed Brown info.hxf = hxf; 239*c4762a1bSJed Brown (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope); 240*c4762a1bSJed Brown for (j=0; j<dof; j++) { 241*c4762a1bSJed Brown PetscScalar tmp = 0; 242*c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 243*c4762a1bSJed Brown slope[i*dof+j] = tmp; 244*c4762a1bSJed Brown } 245*c4762a1bSJed Brown } 246*c4762a1bSJed Brown 247*c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 248*c4762a1bSJed Brown PetscReal maxspeed; 249*c4762a1bSJed Brown PetscScalar *uL,*uR; 250*c4762a1bSJed Brown uL = &ctx->uLR[0]; 251*c4762a1bSJed Brown uR = &ctx->uLR[dof]; 252*c4762a1bSJed Brown if (i < sf) { /* slow region */ 253*c4762a1bSJed Brown for (j=0; j<dof; j++) { 254*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 255*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 256*c4762a1bSJed Brown } 257*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 258*c4762a1bSJed Brown if (i > xs) { 259*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs; 260*c4762a1bSJed Brown } 261*c4762a1bSJed Brown if (i < xs+xm) { 262*c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs; 263*c4762a1bSJed Brown } 264*c4762a1bSJed Brown } else if (i == sf) { /* interface between the slow region and the fast region */ 265*c4762a1bSJed Brown for (j=0; j<dof; j++) { 266*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 267*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 268*c4762a1bSJed Brown } 269*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,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]/hxs; 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]/hxf; 275*c4762a1bSJed Brown } 276*c4762a1bSJed Brown } else if (i > sf && i < fs) { /* fast region */ 277*c4762a1bSJed Brown for (j=0; j<dof; j++) { 278*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 279*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 280*c4762a1bSJed Brown } 281*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 282*c4762a1bSJed Brown if (i > xs) { 283*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf; 284*c4762a1bSJed Brown } 285*c4762a1bSJed Brown if (i < xs+xm) { 286*c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxf; 287*c4762a1bSJed Brown } 288*c4762a1bSJed Brown } else if (i == fs) { /* interface between the fast region and the slow region */ 289*c4762a1bSJed Brown for (j=0; j<dof; j++) { 290*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 291*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 292*c4762a1bSJed Brown } 293*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 294*c4762a1bSJed Brown if (i > xs) { 295*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf; 296*c4762a1bSJed Brown } 297*c4762a1bSJed Brown if (i < xs+xm) { 298*c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs; 299*c4762a1bSJed Brown } 300*c4762a1bSJed Brown } else { /* slow region */ 301*c4762a1bSJed Brown for (j=0; j<dof; j++) { 302*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 303*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 304*c4762a1bSJed Brown } 305*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 306*c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */ 307*c4762a1bSJed Brown if (i > xs) { 308*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs; 309*c4762a1bSJed Brown } 310*c4762a1bSJed Brown if (i < xs+xm) { 311*c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs; 312*c4762a1bSJed Brown } 313*c4762a1bSJed Brown } 314*c4762a1bSJed Brown } 315*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 316*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr); 317*c4762a1bSJed Brown ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 318*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 319*c4762a1bSJed Brown ierr = MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da));CHKERRQ(ierr); 320*c4762a1bSJed Brown if (0) { 321*c4762a1bSJed Brown /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */ 322*c4762a1bSJed Brown PetscReal dt,tnow; 323*c4762a1bSJed Brown ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 324*c4762a1bSJed Brown ierr = TSGetTime(ts,&tnow);CHKERRQ(ierr); 325*c4762a1bSJed Brown if (dt > 0.5/ctx->cfl_idt) { 326*c4762a1bSJed Brown if (1) { 327*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); 328*c4762a1bSJed Brown } else SETERRQ2(PETSC_COMM_SELF,1,"Stability constraint exceeded, %g > %g",(double)dt,(double)(ctx->cfl/ctx->cfl_idt)); 329*c4762a1bSJed Brown } 330*c4762a1bSJed Brown } 331*c4762a1bSJed Brown PetscFunctionReturn(0); 332*c4762a1bSJed Brown } 333*c4762a1bSJed Brown 334*c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */ 335*c4762a1bSJed Brown PetscErrorCode FVRHSFunctionslow_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 336*c4762a1bSJed Brown { 337*c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 338*c4762a1bSJed Brown PetscErrorCode ierr; 339*c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth; 340*c4762a1bSJed Brown PetscReal hxs,hxf,cfl_idt = 0; 341*c4762a1bSJed Brown PetscScalar *x,*f,*slope; 342*c4762a1bSJed Brown Vec Xloc; 343*c4762a1bSJed Brown DM da; 344*c4762a1bSJed Brown 345*c4762a1bSJed Brown PetscFunctionBeginUser; 346*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 347*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 348*c4762a1bSJed Brown ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 349*c4762a1bSJed Brown hxs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 350*c4762a1bSJed Brown hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 351*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 352*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 353*c4762a1bSJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); 354*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 355*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 356*c4762a1bSJed Brown ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 357*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 358*c4762a1bSJed Brown 359*c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 360*c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 361*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 362*c4762a1bSJed Brown } 363*c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 364*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 365*c4762a1bSJed Brown } 366*c4762a1bSJed Brown } 367*c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 368*c4762a1bSJed Brown struct _LimitInfo info; 369*c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 370*c4762a1bSJed Brown if (i < sf-lsbwidth+1 || i > fs+rsbwidth-2) { /* slow components and the first and last fast components */ 371*c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 372*c4762a1bSJed Brown ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); 373*c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 374*c4762a1bSJed Brown ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 375*c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 376*c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 377*c4762a1bSJed Brown for (j=0; j<dof; j++) { 378*c4762a1bSJed Brown PetscScalar jmpL,jmpR; 379*c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 380*c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 381*c4762a1bSJed Brown for (k=0; k<dof; k++) { 382*c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 383*c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 384*c4762a1bSJed Brown } 385*c4762a1bSJed Brown } 386*c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 387*c4762a1bSJed Brown info.m = dof; 388*c4762a1bSJed Brown info.hxs = hxs; 389*c4762a1bSJed Brown info.hxf = hxf; 390*c4762a1bSJed Brown (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope); 391*c4762a1bSJed Brown for (j=0; j<dof; j++) { 392*c4762a1bSJed Brown PetscScalar tmp = 0; 393*c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 394*c4762a1bSJed Brown slope[i*dof+j] = tmp; 395*c4762a1bSJed Brown } 396*c4762a1bSJed Brown } 397*c4762a1bSJed Brown } 398*c4762a1bSJed Brown 399*c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 400*c4762a1bSJed Brown PetscReal maxspeed; 401*c4762a1bSJed Brown PetscScalar *uL,*uR; 402*c4762a1bSJed Brown uL = &ctx->uLR[0]; 403*c4762a1bSJed Brown uR = &ctx->uLR[dof]; 404*c4762a1bSJed Brown if (i < sf-lsbwidth) { /* slow region */ 405*c4762a1bSJed Brown for (j=0; j<dof; j++) { 406*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 407*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 408*c4762a1bSJed Brown } 409*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 410*c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */ 411*c4762a1bSJed Brown if (i > xs) { 412*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 413*c4762a1bSJed Brown } 414*c4762a1bSJed Brown if (i < xs+xm) { 415*c4762a1bSJed Brown for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 416*c4762a1bSJed Brown islow++; 417*c4762a1bSJed Brown } 418*c4762a1bSJed Brown } 419*c4762a1bSJed Brown if (i == sf-lsbwidth) { /* interface between the slow region and the fast region */ 420*c4762a1bSJed Brown for (j=0; j<dof; j++) { 421*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 422*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 423*c4762a1bSJed Brown } 424*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 425*c4762a1bSJed Brown if (i > xs) { 426*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 427*c4762a1bSJed Brown } 428*c4762a1bSJed Brown } 429*c4762a1bSJed Brown if (i == fs+rsbwidth) { /* slow region */ 430*c4762a1bSJed Brown for (j=0; j<dof; j++) { 431*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 432*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 433*c4762a1bSJed Brown } 434*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 435*c4762a1bSJed Brown if (i < xs+xm) { 436*c4762a1bSJed Brown for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 437*c4762a1bSJed Brown islow++; 438*c4762a1bSJed Brown } 439*c4762a1bSJed Brown } 440*c4762a1bSJed Brown if (i > fs+rsbwidth) { /* slow region */ 441*c4762a1bSJed Brown for (j=0; j<dof; j++) { 442*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 443*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 444*c4762a1bSJed Brown } 445*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 446*c4762a1bSJed Brown if (i > xs) { 447*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 448*c4762a1bSJed Brown } 449*c4762a1bSJed Brown if (i < xs+xm) { 450*c4762a1bSJed Brown for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 451*c4762a1bSJed Brown islow++; 452*c4762a1bSJed Brown } 453*c4762a1bSJed Brown } 454*c4762a1bSJed Brown } 455*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 456*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 457*c4762a1bSJed Brown ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 458*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 459*c4762a1bSJed Brown ierr = MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da));CHKERRQ(ierr); 460*c4762a1bSJed Brown PetscFunctionReturn(0); 461*c4762a1bSJed Brown } 462*c4762a1bSJed Brown 463*c4762a1bSJed Brown PetscErrorCode FVRHSFunctionslowbuffer_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 464*c4762a1bSJed Brown { 465*c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 466*c4762a1bSJed Brown PetscErrorCode ierr; 467*c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth; 468*c4762a1bSJed Brown PetscReal hxs,hxf; 469*c4762a1bSJed Brown PetscScalar *x,*f,*slope; 470*c4762a1bSJed Brown Vec Xloc; 471*c4762a1bSJed Brown DM da; 472*c4762a1bSJed Brown 473*c4762a1bSJed Brown PetscFunctionBeginUser; 474*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 475*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 476*c4762a1bSJed Brown ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 477*c4762a1bSJed Brown hxs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 478*c4762a1bSJed Brown hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 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 487*c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 488*c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 489*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 490*c4762a1bSJed Brown } 491*c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 492*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 493*c4762a1bSJed Brown } 494*c4762a1bSJed Brown } 495*c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 496*c4762a1bSJed Brown struct _LimitInfo info; 497*c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 498*c4762a1bSJed Brown if ((i > sf-lsbwidth-2 && i < sf+1) || (i > fs-2 && i < fs+rsbwidth+1)) { 499*c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 500*c4762a1bSJed Brown ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); 501*c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 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.hxs = hxs; 517*c4762a1bSJed Brown info.hxf = hxf; 518*c4762a1bSJed Brown (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope); 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 uL = &ctx->uLR[0]; 531*c4762a1bSJed Brown uR = &ctx->uLR[dof]; 532*c4762a1bSJed Brown if (i == sf-lsbwidth) { 533*c4762a1bSJed Brown for (j=0; j<dof; j++) { 534*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 535*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 536*c4762a1bSJed Brown } 537*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 538*c4762a1bSJed Brown if (i < xs+xm) { 539*c4762a1bSJed Brown for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 540*c4762a1bSJed Brown islow++; 541*c4762a1bSJed Brown } 542*c4762a1bSJed Brown } 543*c4762a1bSJed Brown if (i > sf-lsbwidth && i < sf) { 544*c4762a1bSJed Brown for (j=0; j<dof; j++) { 545*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 546*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 547*c4762a1bSJed Brown } 548*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 549*c4762a1bSJed Brown if (i > xs) { 550*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 551*c4762a1bSJed Brown } 552*c4762a1bSJed Brown if (i < xs+xm) { 553*c4762a1bSJed Brown for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 554*c4762a1bSJed Brown islow++; 555*c4762a1bSJed Brown } 556*c4762a1bSJed Brown } 557*c4762a1bSJed Brown if (i == sf) { /* interface between the slow region and the fast region */ 558*c4762a1bSJed Brown for (j=0; j<dof; j++) { 559*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 560*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 561*c4762a1bSJed Brown } 562*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 563*c4762a1bSJed Brown if (i > xs) { 564*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 565*c4762a1bSJed Brown } 566*c4762a1bSJed Brown } 567*c4762a1bSJed Brown if (i == fs) { /* interface between the fast region and the slow region */ 568*c4762a1bSJed Brown for (j=0; j<dof; j++) { 569*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 570*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 571*c4762a1bSJed Brown } 572*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 573*c4762a1bSJed Brown if (i < xs+xm) { 574*c4762a1bSJed Brown for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 575*c4762a1bSJed Brown islow++; 576*c4762a1bSJed Brown } 577*c4762a1bSJed Brown } 578*c4762a1bSJed Brown if (i > fs && i < fs+rsbwidth) { 579*c4762a1bSJed Brown for (j=0; j<dof; j++) { 580*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 581*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 582*c4762a1bSJed Brown } 583*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 584*c4762a1bSJed Brown if (i > xs) { 585*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 586*c4762a1bSJed Brown } 587*c4762a1bSJed Brown if (i < xs+xm) { 588*c4762a1bSJed Brown for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 589*c4762a1bSJed Brown islow++; 590*c4762a1bSJed Brown } 591*c4762a1bSJed Brown } 592*c4762a1bSJed Brown if (i == fs+rsbwidth) { 593*c4762a1bSJed Brown for (j=0; j<dof; j++) { 594*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 595*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 596*c4762a1bSJed Brown } 597*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 598*c4762a1bSJed Brown if (i > xs) { 599*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 600*c4762a1bSJed Brown } 601*c4762a1bSJed Brown } 602*c4762a1bSJed Brown } 603*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 604*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 605*c4762a1bSJed Brown ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 606*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 607*c4762a1bSJed Brown PetscFunctionReturn(0); 608*c4762a1bSJed Brown } 609*c4762a1bSJed Brown 610*c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for fast parts ----------------------------------- */ 611*c4762a1bSJed Brown PetscErrorCode FVRHSFunctionfast_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 612*c4762a1bSJed Brown { 613*c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 614*c4762a1bSJed Brown PetscErrorCode ierr; 615*c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,ifast = 0,sf = ctx->sf,fs = ctx->fs; 616*c4762a1bSJed Brown PetscReal hxs,hxf; 617*c4762a1bSJed Brown PetscScalar *x,*f,*slope; 618*c4762a1bSJed Brown Vec Xloc; 619*c4762a1bSJed Brown DM da; 620*c4762a1bSJed Brown 621*c4762a1bSJed Brown PetscFunctionBeginUser; 622*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 623*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 624*c4762a1bSJed Brown ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 625*c4762a1bSJed Brown hxs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 626*c4762a1bSJed Brown hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 627*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 628*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 629*c4762a1bSJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); 630*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 631*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 632*c4762a1bSJed Brown ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 633*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 634*c4762a1bSJed Brown 635*c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 636*c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 637*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 638*c4762a1bSJed Brown } 639*c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 640*c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 641*c4762a1bSJed Brown } 642*c4762a1bSJed Brown } 643*c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 644*c4762a1bSJed Brown struct _LimitInfo info; 645*c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 646*c4762a1bSJed Brown if (i > sf-2 && i < fs+1) { 647*c4762a1bSJed Brown ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); 648*c4762a1bSJed Brown ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 649*c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 650*c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 651*c4762a1bSJed Brown for (j=0; j<dof; j++) { 652*c4762a1bSJed Brown PetscScalar jmpL,jmpR; 653*c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 654*c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 655*c4762a1bSJed Brown for (k=0; k<dof; k++) { 656*c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 657*c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 658*c4762a1bSJed Brown } 659*c4762a1bSJed Brown } 660*c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 661*c4762a1bSJed Brown info.m = dof; 662*c4762a1bSJed Brown info.hxs = hxs; 663*c4762a1bSJed Brown info.hxf = hxf; 664*c4762a1bSJed Brown (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope); 665*c4762a1bSJed Brown for (j=0; j<dof; j++) { 666*c4762a1bSJed Brown PetscScalar tmp = 0; 667*c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 668*c4762a1bSJed Brown slope[i*dof+j] = tmp; 669*c4762a1bSJed Brown } 670*c4762a1bSJed Brown } 671*c4762a1bSJed Brown } 672*c4762a1bSJed Brown 673*c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 674*c4762a1bSJed Brown PetscReal maxspeed; 675*c4762a1bSJed Brown PetscScalar *uL,*uR; 676*c4762a1bSJed Brown uL = &ctx->uLR[0]; 677*c4762a1bSJed Brown uR = &ctx->uLR[dof]; 678*c4762a1bSJed Brown if (i == sf) { /* interface between the slow region and the fast region */ 679*c4762a1bSJed Brown for (j=0; j<dof; j++) { 680*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 681*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 682*c4762a1bSJed Brown } 683*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 684*c4762a1bSJed Brown if (i < xs+xm) { 685*c4762a1bSJed Brown for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf; 686*c4762a1bSJed Brown ifast++; 687*c4762a1bSJed Brown } 688*c4762a1bSJed Brown } 689*c4762a1bSJed Brown if (i > sf && i < fs) { /* fast region */ 690*c4762a1bSJed Brown for (j=0; j<dof; j++) { 691*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 692*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 693*c4762a1bSJed Brown } 694*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 695*c4762a1bSJed Brown if (i > xs) { 696*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf; 697*c4762a1bSJed Brown } 698*c4762a1bSJed Brown if (i < xs+xm) { 699*c4762a1bSJed Brown for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf; 700*c4762a1bSJed Brown ifast++; 701*c4762a1bSJed Brown } 702*c4762a1bSJed Brown } 703*c4762a1bSJed Brown if (i == fs) { /* interface between the fast region and the slow region */ 704*c4762a1bSJed Brown for (j=0; j<dof; j++) { 705*c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 706*c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 707*c4762a1bSJed Brown } 708*c4762a1bSJed Brown ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 709*c4762a1bSJed Brown if (i > xs) { 710*c4762a1bSJed Brown for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf; 711*c4762a1bSJed Brown } 712*c4762a1bSJed Brown } 713*c4762a1bSJed Brown } 714*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 715*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 716*c4762a1bSJed Brown ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 717*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 718*c4762a1bSJed Brown PetscFunctionReturn(0); 719*c4762a1bSJed Brown } 720*c4762a1bSJed Brown 721*c4762a1bSJed Brown int main(int argc,char *argv[]) 722*c4762a1bSJed Brown { 723*c4762a1bSJed Brown char lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m"; 724*c4762a1bSJed Brown PetscFunctionList limiters = 0,physics = 0; 725*c4762a1bSJed Brown MPI_Comm comm; 726*c4762a1bSJed Brown TS ts; 727*c4762a1bSJed Brown DM da; 728*c4762a1bSJed Brown Vec X,X0,R; 729*c4762a1bSJed Brown FVCtx ctx; 730*c4762a1bSJed Brown PetscInt i,k,dof,xs,xm,Mx,draw = 0,count_slow,count_fast,islow = 0,ifast =0,islowbuffer = 0,*index_slow,*index_fast,*index_slowbuffer; 731*c4762a1bSJed Brown PetscBool view_final = PETSC_FALSE; 732*c4762a1bSJed Brown PetscReal ptime; 733*c4762a1bSJed Brown PetscErrorCode ierr; 734*c4762a1bSJed Brown 735*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr; 736*c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 737*c4762a1bSJed Brown ierr = PetscMemzero(&ctx,sizeof(ctx));CHKERRQ(ierr); 738*c4762a1bSJed Brown 739*c4762a1bSJed Brown /* Register limiters to be available on the command line */ 740*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"upwind" ,Limit2_Upwind);CHKERRQ(ierr); 741*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"lax-wendroff" ,Limit2_LaxWendroff);CHKERRQ(ierr); 742*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"beam-warming" ,Limit2_BeamWarming);CHKERRQ(ierr); 743*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"fromm" ,Limit2_Fromm);CHKERRQ(ierr); 744*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"minmod" ,Limit2_Minmod);CHKERRQ(ierr); 745*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"superbee" ,Limit2_Superbee);CHKERRQ(ierr); 746*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"mc" ,Limit2_MC);CHKERRQ(ierr); 747*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&limiters,"koren3" ,Limit2_Koren3);CHKERRQ(ierr); 748*c4762a1bSJed Brown 749*c4762a1bSJed Brown /* Register physical models to be available on the command line */ 750*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&physics,"advect" ,PhysicsCreate_Advect);CHKERRQ(ierr); 751*c4762a1bSJed Brown 752*c4762a1bSJed Brown ctx.comm = comm; 753*c4762a1bSJed Brown ctx.cfl = 0.9; 754*c4762a1bSJed Brown ctx.bctype = FVBC_PERIODIC; 755*c4762a1bSJed Brown ctx.xmin = -1.0; 756*c4762a1bSJed Brown ctx.xmax = 1.0; 757*c4762a1bSJed Brown ierr = PetscOptionsBegin(comm,NULL,"Finite Volume solver options","");CHKERRQ(ierr); 758*c4762a1bSJed Brown ierr = PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL);CHKERRQ(ierr); 759*c4762a1bSJed Brown ierr = PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL);CHKERRQ(ierr); 760*c4762a1bSJed Brown ierr = PetscOptionsFList("-limit","Name of flux imiter to use","",limiters,lname,lname,sizeof(lname),NULL);CHKERRQ(ierr); 761*c4762a1bSJed Brown ierr = PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL);CHKERRQ(ierr); 762*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); 763*c4762a1bSJed Brown ierr = PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL);CHKERRQ(ierr); 764*c4762a1bSJed Brown ierr = PetscOptionsBool("-exact","Compare errors with exact solution","",ctx.exact,&ctx.exact,NULL);CHKERRQ(ierr); 765*c4762a1bSJed Brown ierr = PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL);CHKERRQ(ierr); 766*c4762a1bSJed Brown ierr = PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL);CHKERRQ(ierr); 767*c4762a1bSJed Brown ierr = PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL);CHKERRQ(ierr); 768*c4762a1bSJed Brown ierr = PetscOptionsInt("-hratio","Spacing ratio","",ctx.hratio,&ctx.hratio,NULL);CHKERRQ(ierr); 769*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 770*c4762a1bSJed Brown 771*c4762a1bSJed Brown /* Choose the limiter from the list of registered limiters */ 772*c4762a1bSJed Brown ierr = PetscFunctionListFind(limiters,lname,&ctx.limit2);CHKERRQ(ierr); 773*c4762a1bSJed Brown if (!ctx.limit2) SETERRQ1(PETSC_COMM_SELF,1,"Limiter '%s' not found",lname); 774*c4762a1bSJed Brown 775*c4762a1bSJed Brown /* Choose the physics from the list of registered models */ 776*c4762a1bSJed Brown { 777*c4762a1bSJed Brown PetscErrorCode (*r)(FVCtx*); 778*c4762a1bSJed Brown ierr = PetscFunctionListFind(physics,physname,&r);CHKERRQ(ierr); 779*c4762a1bSJed Brown if (!r) SETERRQ1(PETSC_COMM_SELF,1,"Physics '%s' not found",physname); 780*c4762a1bSJed Brown /* Create the physics, will set the number of fields and their names */ 781*c4762a1bSJed Brown ierr = (*r)(&ctx);CHKERRQ(ierr); 782*c4762a1bSJed Brown } 783*c4762a1bSJed Brown 784*c4762a1bSJed Brown /* Create a DMDA to manage the parallel grid */ 785*c4762a1bSJed Brown ierr = DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics2.dof,2,NULL,&da);CHKERRQ(ierr); 786*c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 787*c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 788*c4762a1bSJed Brown /* Inform the DMDA of the field names provided by the physics. */ 789*c4762a1bSJed Brown /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ 790*c4762a1bSJed Brown for (i=0; i<ctx.physics2.dof; i++) { 791*c4762a1bSJed Brown ierr = DMDASetFieldName(da,i,ctx.physics2.fieldname[i]);CHKERRQ(ierr); 792*c4762a1bSJed Brown } 793*c4762a1bSJed Brown ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 794*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 795*c4762a1bSJed Brown 796*c4762a1bSJed Brown /* Set coordinates of cell centers */ 797*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); 798*c4762a1bSJed Brown 799*c4762a1bSJed Brown /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */ 800*c4762a1bSJed Brown ierr = PetscMalloc4(dof*dof,&ctx.R,dof*dof,&ctx.Rinv,2*dof,&ctx.cjmpLR,1*dof,&ctx.cslope);CHKERRQ(ierr); 801*c4762a1bSJed Brown ierr = PetscMalloc3(2*dof,&ctx.uLR,dof,&ctx.flux,dof,&ctx.speeds);CHKERRQ(ierr); 802*c4762a1bSJed Brown 803*c4762a1bSJed Brown /* Create a vector to store the solution and to save the initial state */ 804*c4762a1bSJed Brown ierr = DMCreateGlobalVector(da,&X);CHKERRQ(ierr); 805*c4762a1bSJed Brown ierr = VecDuplicate(X,&X0);CHKERRQ(ierr); 806*c4762a1bSJed Brown ierr = VecDuplicate(X,&R);CHKERRQ(ierr); 807*c4762a1bSJed Brown 808*c4762a1bSJed Brown /* create index for slow parts and fast parts, 809*c4762a1bSJed Brown count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */ 810*c4762a1bSJed Brown count_slow = Mx/(1.0+ctx.hratio/3.0); 811*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/3) is even"); 812*c4762a1bSJed Brown count_fast = Mx-count_slow; 813*c4762a1bSJed Brown ctx.sf = count_slow/2; 814*c4762a1bSJed Brown ctx.fs = ctx.sf+count_fast; 815*c4762a1bSJed Brown ierr = PetscMalloc1(xm*dof,&index_slow);CHKERRQ(ierr); 816*c4762a1bSJed Brown ierr = PetscMalloc1(xm*dof,&index_fast);CHKERRQ(ierr); 817*c4762a1bSJed Brown ierr = PetscMalloc1(6*dof,&index_slowbuffer);CHKERRQ(ierr); 818*c4762a1bSJed Brown if (((AdvectCtx*)ctx.physics2.user)->a > 0) { 819*c4762a1bSJed Brown ctx.lsbwidth = 2; 820*c4762a1bSJed Brown ctx.rsbwidth = 4; 821*c4762a1bSJed Brown } else { 822*c4762a1bSJed Brown ctx.lsbwidth = 4; 823*c4762a1bSJed Brown ctx.rsbwidth = 2; 824*c4762a1bSJed Brown } 825*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 826*c4762a1bSJed Brown if (i < ctx.sf-ctx.lsbwidth || i > ctx.fs+ctx.rsbwidth-1) 827*c4762a1bSJed Brown for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k; 828*c4762a1bSJed Brown else if ((i >= ctx.sf-ctx.lsbwidth && i < ctx.sf) || (i > ctx.fs-1 && i <= ctx.fs+ctx.rsbwidth-1)) 829*c4762a1bSJed Brown for (k=0; k<dof; k++) index_slowbuffer[islowbuffer++] = i*dof+k; 830*c4762a1bSJed Brown else 831*c4762a1bSJed Brown for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k; 832*c4762a1bSJed Brown } 833*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss);CHKERRQ(ierr); 834*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf);CHKERRQ(ierr); 835*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,islowbuffer,index_slowbuffer,PETSC_COPY_VALUES,&ctx.issb);CHKERRQ(ierr); 836*c4762a1bSJed Brown 837*c4762a1bSJed Brown /* Create a time-stepping object */ 838*c4762a1bSJed Brown ierr = TSCreate(comm,&ts);CHKERRQ(ierr); 839*c4762a1bSJed Brown ierr = TSSetDM(ts,da);CHKERRQ(ierr); 840*c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,R,FVRHSFunction_2WaySplit,&ctx);CHKERRQ(ierr); 841*c4762a1bSJed Brown ierr = TSRHSSplitSetIS(ts,"slow",ctx.iss);CHKERRQ(ierr); 842*c4762a1bSJed Brown ierr = TSRHSSplitSetIS(ts,"slowbuffer",ctx.issb);CHKERRQ(ierr); 843*c4762a1bSJed Brown ierr = TSRHSSplitSetIS(ts,"fast",ctx.isf);CHKERRQ(ierr); 844*c4762a1bSJed Brown ierr = TSRHSSplitSetRHSFunction(ts,"slow",NULL,FVRHSFunctionslow_2WaySplit,&ctx);CHKERRQ(ierr); 845*c4762a1bSJed Brown ierr = TSRHSSplitSetRHSFunction(ts,"fast",NULL,FVRHSFunctionfast_2WaySplit,&ctx);CHKERRQ(ierr); 846*c4762a1bSJed Brown ierr = TSRHSSplitSetRHSFunction(ts,"slowbuffer",NULL,FVRHSFunctionslowbuffer_2WaySplit,&ctx);CHKERRQ(ierr); 847*c4762a1bSJed Brown 848*c4762a1bSJed Brown ierr = TSSetType(ts,TSSSP);CHKERRQ(ierr); 849*c4762a1bSJed Brown /*ierr = TSSetType(ts,TSMPRK);CHKERRQ(ierr);*/ 850*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,10);CHKERRQ(ierr); 851*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 852*c4762a1bSJed Brown 853*c4762a1bSJed Brown /* Compute initial conditions and starting time step */ 854*c4762a1bSJed Brown ierr = FVSample_2WaySplit(&ctx,da,0,X0);CHKERRQ(ierr); 855*c4762a1bSJed Brown ierr = FVRHSFunction_2WaySplit(ts,0,X0,X,(void*)&ctx);CHKERRQ(ierr); /* Initial function evaluation, only used to determine max speed */ 856*c4762a1bSJed Brown ierr = VecCopy(X0,X);CHKERRQ(ierr); /* The function value was not used so we set X=X0 again */ 857*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt);CHKERRQ(ierr); 858*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); /* Take runtime options */ 859*c4762a1bSJed Brown ierr = SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 860*c4762a1bSJed Brown { 861*c4762a1bSJed Brown PetscInt steps; 862*c4762a1bSJed Brown PetscScalar mass_initial,mass_final,mass_difference,mass_differenceg; 863*c4762a1bSJed Brown const PetscScalar *ptr_X,*ptr_X0; 864*c4762a1bSJed Brown const PetscReal hs = (ctx.xmax-ctx.xmin)*3.0/4.0/count_slow; 865*c4762a1bSJed Brown const PetscReal hf = (ctx.xmax-ctx.xmin)/4.0/count_fast; 866*c4762a1bSJed Brown 867*c4762a1bSJed Brown ierr = TSSolve(ts,X);CHKERRQ(ierr); 868*c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&ptime);CHKERRQ(ierr); 869*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 870*c4762a1bSJed Brown /* calculate the total mass at initial time and final time */ 871*c4762a1bSJed Brown mass_initial = 0.0; 872*c4762a1bSJed Brown mass_final = 0.0; 873*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,X0,(void*)&ptr_X0);CHKERRQ(ierr); 874*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,X,(void*)&ptr_X);CHKERRQ(ierr); 875*c4762a1bSJed Brown for (i=xs;i<xs+xm;i++) { 876*c4762a1bSJed Brown if (i < ctx.sf || i > ctx.fs-1) { 877*c4762a1bSJed Brown for (k=0; k<dof; k++) { 878*c4762a1bSJed Brown mass_initial = mass_initial + hs*ptr_X0[i*dof+k]; 879*c4762a1bSJed Brown mass_final = mass_final + hs*ptr_X[i*dof+k]; 880*c4762a1bSJed Brown } 881*c4762a1bSJed Brown } else { 882*c4762a1bSJed Brown for (k=0; k<dof; k++) { 883*c4762a1bSJed Brown mass_initial = mass_initial + hf*ptr_X0[i*dof+k]; 884*c4762a1bSJed Brown mass_final = mass_final + hf*ptr_X[i*dof+k]; 885*c4762a1bSJed Brown } 886*c4762a1bSJed Brown } 887*c4762a1bSJed Brown } 888*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,X0,(void*)&ptr_X0);CHKERRQ(ierr); 889*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,X,(void*)&ptr_X);CHKERRQ(ierr); 890*c4762a1bSJed Brown mass_difference = mass_final - mass_initial; 891*c4762a1bSJed Brown ierr = MPI_Allreduce(&mass_difference,&mass_differenceg,1,MPIU_SCALAR,MPIU_SUM,comm);CHKERRQ(ierr); 892*c4762a1bSJed Brown ierr = PetscPrintf(comm,"Mass difference %g\n",(double)mass_differenceg);CHKERRQ(ierr); 893*c4762a1bSJed Brown ierr = PetscPrintf(comm,"Final time %g, steps %D\n",(double)ptime,steps);CHKERRQ(ierr); 894*c4762a1bSJed Brown ierr = PetscPrintf(comm,"Maximum allowable stepsize according to CFL %g\n",(double)(1.0/ctx.cfl_idt));CHKERRQ(ierr); 895*c4762a1bSJed Brown if (ctx.exact) { 896*c4762a1bSJed Brown PetscReal nrm1=0; 897*c4762a1bSJed Brown ierr = SolutionErrorNorms_2WaySplit(&ctx,da,ptime,X,&nrm1);CHKERRQ(ierr); 898*c4762a1bSJed Brown ierr = PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1);CHKERRQ(ierr); 899*c4762a1bSJed Brown } 900*c4762a1bSJed Brown if (ctx.simulation) { 901*c4762a1bSJed Brown PetscReal nrm1=0; 902*c4762a1bSJed Brown PetscViewer fd; 903*c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = "binaryoutput"; 904*c4762a1bSJed Brown Vec XR; 905*c4762a1bSJed Brown PetscBool flg; 906*c4762a1bSJed Brown const PetscScalar *ptr_XR; 907*c4762a1bSJed Brown ierr = PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg);CHKERRQ(ierr); 908*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 909*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&fd);CHKERRQ(ierr); 910*c4762a1bSJed Brown ierr = VecDuplicate(X0,&XR);CHKERRQ(ierr); 911*c4762a1bSJed Brown ierr = VecLoad(XR,fd);CHKERRQ(ierr); 912*c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 913*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&ptr_X);CHKERRQ(ierr); 914*c4762a1bSJed Brown ierr = VecGetArrayRead(XR,&ptr_XR);CHKERRQ(ierr); 915*c4762a1bSJed Brown for(i=xs;i<xs+xm;i++) { 916*c4762a1bSJed Brown if(i < ctx.sf || i > ctx.fs-1) 917*c4762a1bSJed Brown for (k=0; k<dof; k++) nrm1 = nrm1 + hs*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]); 918*c4762a1bSJed Brown else 919*c4762a1bSJed Brown for (k=0; k<dof; k++) nrm1 = nrm1 + hf*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]); 920*c4762a1bSJed Brown } 921*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&ptr_X);CHKERRQ(ierr); 922*c4762a1bSJed Brown ierr = VecRestoreArrayRead(XR,&ptr_XR);CHKERRQ(ierr); 923*c4762a1bSJed Brown ierr = PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1);CHKERRQ(ierr); 924*c4762a1bSJed Brown ierr = VecDestroy(&XR);CHKERRQ(ierr); 925*c4762a1bSJed Brown } 926*c4762a1bSJed Brown } 927*c4762a1bSJed Brown 928*c4762a1bSJed Brown ierr = SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 929*c4762a1bSJed Brown if (draw & 0x1) {ierr = VecView(X0,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);} 930*c4762a1bSJed Brown if (draw & 0x2) {ierr = VecView(X,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);} 931*c4762a1bSJed Brown if (draw & 0x4) { 932*c4762a1bSJed Brown Vec Y; 933*c4762a1bSJed Brown ierr = VecDuplicate(X,&Y);CHKERRQ(ierr); 934*c4762a1bSJed Brown ierr = FVSample_2WaySplit(&ctx,da,ptime,Y);CHKERRQ(ierr); 935*c4762a1bSJed Brown ierr = VecAYPX(Y,-1,X);CHKERRQ(ierr); 936*c4762a1bSJed Brown ierr = VecView(Y,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr); 937*c4762a1bSJed Brown ierr = VecDestroy(&Y);CHKERRQ(ierr); 938*c4762a1bSJed Brown } 939*c4762a1bSJed Brown 940*c4762a1bSJed Brown if (view_final) { 941*c4762a1bSJed Brown PetscViewer viewer; 942*c4762a1bSJed Brown ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer);CHKERRQ(ierr); 943*c4762a1bSJed Brown ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 944*c4762a1bSJed Brown ierr = VecView(X,viewer);CHKERRQ(ierr); 945*c4762a1bSJed Brown ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 946*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 947*c4762a1bSJed Brown } 948*c4762a1bSJed Brown 949*c4762a1bSJed Brown /* Clean up */ 950*c4762a1bSJed Brown ierr = (*ctx.physics2.destroy)(ctx.physics2.user);CHKERRQ(ierr); 951*c4762a1bSJed Brown for (i=0; i<ctx.physics2.dof; i++) {ierr = PetscFree(ctx.physics2.fieldname[i]);CHKERRQ(ierr);} 952*c4762a1bSJed Brown ierr = PetscFree4(ctx.R,ctx.Rinv,ctx.cjmpLR,ctx.cslope);CHKERRQ(ierr); 953*c4762a1bSJed Brown ierr = PetscFree3(ctx.uLR,ctx.flux,ctx.speeds);CHKERRQ(ierr); 954*c4762a1bSJed Brown ierr = VecDestroy(&X);CHKERRQ(ierr); 955*c4762a1bSJed Brown ierr = VecDestroy(&X0);CHKERRQ(ierr); 956*c4762a1bSJed Brown ierr = VecDestroy(&R);CHKERRQ(ierr); 957*c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 958*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 959*c4762a1bSJed Brown ierr = ISDestroy(&ctx.iss);CHKERRQ(ierr); 960*c4762a1bSJed Brown ierr = ISDestroy(&ctx.isf);CHKERRQ(ierr); 961*c4762a1bSJed Brown ierr = ISDestroy(&ctx.issb);CHKERRQ(ierr); 962*c4762a1bSJed Brown ierr = PetscFree(index_slow);CHKERRQ(ierr); 963*c4762a1bSJed Brown ierr = PetscFree(index_fast);CHKERRQ(ierr); 964*c4762a1bSJed Brown ierr = PetscFree(index_slowbuffer);CHKERRQ(ierr); 965*c4762a1bSJed Brown ierr = PetscFunctionListDestroy(&limiters);CHKERRQ(ierr); 966*c4762a1bSJed Brown ierr = PetscFunctionListDestroy(&physics);CHKERRQ(ierr); 967*c4762a1bSJed Brown ierr = PetscFinalize(); 968*c4762a1bSJed Brown return ierr; 969*c4762a1bSJed Brown } 970*c4762a1bSJed Brown 971*c4762a1bSJed Brown /*TEST 972*c4762a1bSJed Brown 973*c4762a1bSJed Brown build: 974*c4762a1bSJed Brown requires: !complex c99 975*c4762a1bSJed Brown depends: finitevolume1d.c 976*c4762a1bSJed Brown 977*c4762a1bSJed Brown test: 978*c4762a1bSJed Brown suffix: 1 979*c4762a1bSJed Brown args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 980*c4762a1bSJed Brown 981*c4762a1bSJed Brown test: 982*c4762a1bSJed Brown suffix: 2 983*c4762a1bSJed Brown args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 0 984*c4762a1bSJed Brown output_file: output/ex6_1.out 985*c4762a1bSJed Brown 986*c4762a1bSJed Brown TEST*/ 987