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