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 279fbee547SJacob Faibussowitsch 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 AdvectCtx *user; 85c4762a1bSJed Brown 86c4762a1bSJed Brown PetscFunctionBeginUser; 879566063dSJacob Faibussowitsch PetscCall(PetscNew(&user)); 88c4762a1bSJed Brown ctx->physics2.sample2 = PhysicsSample_Advect; 89c4762a1bSJed Brown ctx->physics2.riemann2 = PhysicsRiemann_Advect; 90c4762a1bSJed Brown ctx->physics2.characteristic2 = PhysicsCharacteristic_Advect; 91c4762a1bSJed Brown ctx->physics2.destroy = PhysicsDestroy_SimpleFree; 92c4762a1bSJed Brown ctx->physics2.user = user; 93c4762a1bSJed Brown ctx->physics2.dof = 1; 949566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("u",&ctx->physics2.fieldname[0])); 95c4762a1bSJed Brown user->a = 1; 96*d0609cedSBarry Smith PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection",""); 97c4762a1bSJed Brown { 989566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-physics_advect_a","Speed","",user->a,&user->a,NULL)); 99c4762a1bSJed Brown } 100*d0609cedSBarry Smith PetscOptionsEnd(); 101c4762a1bSJed Brown PetscFunctionReturn(0); 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown PetscErrorCode FVSample_2WaySplit(FVCtx *ctx,DM da,PetscReal time,Vec U) 105c4762a1bSJed Brown { 106c4762a1bSJed Brown PetscScalar *u,*uj,xj,xi; 107c4762a1bSJed Brown PetscInt i,j,k,dof,xs,xm,Mx; 108c4762a1bSJed Brown const PetscInt N = 200; 109c4762a1bSJed Brown PetscReal hs,hf; 110c4762a1bSJed Brown 111c4762a1bSJed Brown PetscFunctionBeginUser; 1123c633725SBarry Smith PetscCheck(ctx->physics2.sample2,PETSC_COMM_SELF,PETSC_ERR_SUP,"Physics has not provided a sampling function"); 1139566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 1149566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 1159566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,U,&u)); 1169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof,&uj)); 117c4762a1bSJed Brown hs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 118c4762a1bSJed Brown hf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 119c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 120c4762a1bSJed Brown if (i < ctx->sf) { 121c4762a1bSJed Brown xi = ctx->xmin+0.5*hs+i*hs; 122c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 123c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 124c4762a1bSJed Brown for (j=0; j<N+1; j++) { 125c4762a1bSJed Brown xj = xi+hs*(j-N/2)/(PetscReal)N; 1269566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj)); 127c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 128c4762a1bSJed Brown } 129c4762a1bSJed Brown } else if (i < ctx->fs) { 130c4762a1bSJed Brown xi = ctx->xmin+ctx->sf*hs+0.5*hf+(i-ctx->sf)*hf; 131c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 132c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 133c4762a1bSJed Brown for (j=0; j<N+1; j++) { 134c4762a1bSJed Brown xj = xi+hf*(j-N/2)/(PetscReal)N; 1359566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj)); 136c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 137c4762a1bSJed Brown } 138c4762a1bSJed Brown } else { 139c4762a1bSJed Brown xi = ctx->xmin+ctx->sf*hs+(ctx->fs-ctx->sf)*hf+0.5*hs+(i-ctx->fs)*hs; 140c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 141c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 142c4762a1bSJed Brown for (j=0; j<N+1; j++) { 143c4762a1bSJed Brown xj = xi+hs*(j-N/2)/(PetscReal)N; 1449566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj)); 145c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 146c4762a1bSJed Brown } 147c4762a1bSJed Brown } 148c4762a1bSJed Brown } 1499566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,U,&u)); 1509566063dSJacob Faibussowitsch PetscCall(PetscFree(uj)); 151c4762a1bSJed Brown PetscFunctionReturn(0); 152c4762a1bSJed Brown } 153c4762a1bSJed Brown 154c4762a1bSJed Brown static PetscErrorCode SolutionErrorNorms_2WaySplit(FVCtx *ctx,DM da,PetscReal t,Vec X,PetscReal *nrm1) 155c4762a1bSJed Brown { 156c4762a1bSJed Brown Vec Y; 157c4762a1bSJed Brown PetscInt i,Mx; 158c4762a1bSJed Brown const PetscScalar *ptr_X,*ptr_Y; 159c4762a1bSJed Brown PetscReal hs,hf; 160c4762a1bSJed Brown 161c4762a1bSJed Brown PetscFunctionBeginUser; 1629566063dSJacob Faibussowitsch PetscCall(VecGetSize(X,&Mx)); 1639566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&Y)); 1649566063dSJacob Faibussowitsch PetscCall(FVSample_2WaySplit(ctx,da,t,Y)); 165c4762a1bSJed Brown hs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 166c4762a1bSJed Brown hf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 1679566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&ptr_X)); 1689566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&ptr_Y)); 169c4762a1bSJed Brown for (i=0; i<Mx; i++) { 170c4762a1bSJed Brown if (i < ctx->sf || i > ctx->fs -1) *nrm1 += hs*PetscAbs(ptr_X[i]-ptr_Y[i]); 171c4762a1bSJed Brown else *nrm1 += hf*PetscAbs(ptr_X[i]-ptr_Y[i]); 172c4762a1bSJed Brown } 1739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&ptr_X)); 1749566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&ptr_Y)); 1759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 176c4762a1bSJed Brown PetscFunctionReturn(0); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179c4762a1bSJed Brown PetscErrorCode FVRHSFunction_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 180c4762a1bSJed Brown { 181c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 182c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,sf = ctx->sf,fs = ctx->fs; 183c4762a1bSJed Brown PetscReal hxf,hxs,cfl_idt = 0; 184c4762a1bSJed Brown PetscScalar *x,*f,*slope; 185c4762a1bSJed Brown Vec Xloc; 186c4762a1bSJed Brown DM da; 187c4762a1bSJed Brown 188c4762a1bSJed Brown PetscFunctionBeginUser; 1899566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 1909566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); /* Xloc contains ghost points */ 1919566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); /* Mx is the number of center points */ 192c4762a1bSJed Brown hxs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 193c4762a1bSJed Brown hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 1949566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); /* X is solution vector which does not contain ghost points */ 1959566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 196c4762a1bSJed Brown 1979566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); /* F is the right hand side function corresponds to center points */ 198c4762a1bSJed Brown 1999566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 2009566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,F,&f)); 2019566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da,PETSC_TRUE,&slope)); /* contains ghost points */ 202c4762a1bSJed Brown 2039566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 204c4762a1bSJed Brown 205c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 206c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 207c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 208c4762a1bSJed Brown } 209c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 210c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 211c4762a1bSJed Brown } 212c4762a1bSJed Brown } 213c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 214c4762a1bSJed Brown struct _LimitInfo info; 215c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 216c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 2179566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds)); 218c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 2199566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR,2*dof)); 220c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 221c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 222c4762a1bSJed Brown for (j=0; j<dof; j++) { 223c4762a1bSJed Brown PetscScalar jmpL,jmpR; 224c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 225c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 226c4762a1bSJed Brown for (k=0; k<dof; k++) { 227c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 228c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 229c4762a1bSJed Brown } 230c4762a1bSJed Brown } 231c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 232c4762a1bSJed Brown info.m = dof; 233c4762a1bSJed Brown info.hxs = hxs; 234c4762a1bSJed Brown info.hxf = hxf; 235c4762a1bSJed Brown (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope); 236c4762a1bSJed Brown for (j=0; j<dof; j++) { 237c4762a1bSJed Brown PetscScalar tmp = 0; 238c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 239c4762a1bSJed Brown slope[i*dof+j] = tmp; 240c4762a1bSJed Brown } 241c4762a1bSJed Brown } 242c4762a1bSJed Brown 243c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 244c4762a1bSJed Brown PetscReal maxspeed; 245c4762a1bSJed Brown PetscScalar *uL,*uR; 246c4762a1bSJed Brown uL = &ctx->uLR[0]; 247c4762a1bSJed Brown uR = &ctx->uLR[dof]; 248c4762a1bSJed Brown if (i < sf) { /* slow region */ 249c4762a1bSJed Brown for (j=0; j<dof; j++) { 250c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 251c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 252c4762a1bSJed Brown } 2539566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 254c4762a1bSJed Brown if (i > xs) { 255c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs; 256c4762a1bSJed Brown } 257c4762a1bSJed Brown if (i < xs+xm) { 258c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs; 259c4762a1bSJed Brown } 260c4762a1bSJed Brown } else if (i == sf) { /* interface between the slow region and the fast region */ 261c4762a1bSJed Brown for (j=0; j<dof; j++) { 262c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 263c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 264c4762a1bSJed Brown } 2659566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 266c4762a1bSJed Brown if (i > xs) { 267c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs; 268c4762a1bSJed Brown } 269c4762a1bSJed Brown if (i < xs+xm) { 270c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxf; 271c4762a1bSJed Brown } 272c4762a1bSJed Brown } else if (i > sf && i < fs) { /* fast region */ 273c4762a1bSJed Brown for (j=0; j<dof; j++) { 274c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 275c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 276c4762a1bSJed Brown } 2779566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 278c4762a1bSJed Brown if (i > xs) { 279c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf; 280c4762a1bSJed Brown } 281c4762a1bSJed Brown if (i < xs+xm) { 282c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxf; 283c4762a1bSJed Brown } 284c4762a1bSJed Brown } else if (i == fs) { /* interface between the fast region and the slow region */ 285c4762a1bSJed Brown for (j=0; j<dof; j++) { 286c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 287c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 288c4762a1bSJed Brown } 2899566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 290c4762a1bSJed Brown if (i > xs) { 291c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf; 292c4762a1bSJed Brown } 293c4762a1bSJed Brown if (i < xs+xm) { 294c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs; 295c4762a1bSJed Brown } 296c4762a1bSJed Brown } else { /* slow region */ 297c4762a1bSJed Brown for (j=0; j<dof; j++) { 298c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 299c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 300c4762a1bSJed Brown } 3019566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 302c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */ 303c4762a1bSJed Brown if (i > xs) { 304c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs; 305c4762a1bSJed Brown } 306c4762a1bSJed Brown if (i < xs+xm) { 307c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs; 308c4762a1bSJed Brown } 309c4762a1bSJed Brown } 310c4762a1bSJed Brown } 3119566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 3129566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,F,&f)); 3139566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da,PETSC_TRUE,&slope)); 3149566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 3159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da))); 316c4762a1bSJed Brown if (0) { 317c4762a1bSJed Brown /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */ 318c4762a1bSJed Brown PetscReal dt,tnow; 3199566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&dt)); 3209566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts,&tnow)); 321c4762a1bSJed Brown if (dt > 0.5/ctx->cfl_idt) { 322c4762a1bSJed Brown if (1) { 3239566063dSJacob Faibussowitsch PetscCall(PetscPrintf(ctx->comm,"Stability constraint exceeded at t=%g, dt %g > %g\n",(double)tnow,(double)dt,(double)(0.5/ctx->cfl_idt))); 32498921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Stability constraint exceeded, %g > %g",(double)dt,(double)(ctx->cfl/ctx->cfl_idt)); 325c4762a1bSJed Brown } 326c4762a1bSJed Brown } 327c4762a1bSJed Brown PetscFunctionReturn(0); 328c4762a1bSJed Brown } 329c4762a1bSJed Brown 330c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */ 331c4762a1bSJed Brown PetscErrorCode FVRHSFunctionslow_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 332c4762a1bSJed Brown { 333c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 334c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth; 335c4762a1bSJed Brown PetscReal hxs,hxf,cfl_idt = 0; 336c4762a1bSJed Brown PetscScalar *x,*f,*slope; 337c4762a1bSJed Brown Vec Xloc; 338c4762a1bSJed Brown DM da; 339c4762a1bSJed Brown 340c4762a1bSJed Brown PetscFunctionBeginUser; 3419566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 3429566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); 3439566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 344c4762a1bSJed Brown hxs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 345c4762a1bSJed Brown hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 3469566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 3479566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 3489566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 3499566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 3509566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 3519566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da,PETSC_TRUE,&slope)); 3529566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 353c4762a1bSJed Brown 354c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 355c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 356c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 357c4762a1bSJed Brown } 358c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 359c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 360c4762a1bSJed Brown } 361c4762a1bSJed Brown } 362c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 363c4762a1bSJed Brown struct _LimitInfo info; 364c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 365c4762a1bSJed Brown if (i < sf-lsbwidth+1 || i > fs+rsbwidth-2) { /* slow components and the first and last fast components */ 366c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 3679566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds)); 368c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 3699566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR,2*dof)); 370c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 371c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 372c4762a1bSJed Brown for (j=0; j<dof; j++) { 373c4762a1bSJed Brown PetscScalar jmpL,jmpR; 374c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 375c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 376c4762a1bSJed Brown for (k=0; k<dof; k++) { 377c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 378c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 379c4762a1bSJed Brown } 380c4762a1bSJed Brown } 381c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 382c4762a1bSJed Brown info.m = dof; 383c4762a1bSJed Brown info.hxs = hxs; 384c4762a1bSJed Brown info.hxf = hxf; 385c4762a1bSJed Brown (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope); 386c4762a1bSJed Brown for (j=0; j<dof; j++) { 387c4762a1bSJed Brown PetscScalar tmp = 0; 388c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 389c4762a1bSJed Brown slope[i*dof+j] = tmp; 390c4762a1bSJed Brown } 391c4762a1bSJed Brown } 392c4762a1bSJed Brown } 393c4762a1bSJed Brown 394c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 395c4762a1bSJed Brown PetscReal maxspeed; 396c4762a1bSJed Brown PetscScalar *uL,*uR; 397c4762a1bSJed Brown uL = &ctx->uLR[0]; 398c4762a1bSJed Brown uR = &ctx->uLR[dof]; 399c4762a1bSJed Brown if (i < sf-lsbwidth) { /* slow region */ 400c4762a1bSJed Brown for (j=0; j<dof; j++) { 401c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 402c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 403c4762a1bSJed Brown } 4049566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 405c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */ 406c4762a1bSJed Brown if (i > xs) { 407c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 408c4762a1bSJed Brown } 409c4762a1bSJed Brown if (i < xs+xm) { 410c4762a1bSJed Brown for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 411c4762a1bSJed Brown islow++; 412c4762a1bSJed Brown } 413c4762a1bSJed Brown } 414c4762a1bSJed Brown if (i == sf-lsbwidth) { /* interface between the slow region and the fast region */ 415c4762a1bSJed Brown for (j=0; j<dof; j++) { 416c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 417c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 418c4762a1bSJed Brown } 4199566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 420c4762a1bSJed Brown if (i > xs) { 421c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 422c4762a1bSJed Brown } 423c4762a1bSJed Brown } 424c4762a1bSJed Brown if (i == fs+rsbwidth) { /* slow region */ 425c4762a1bSJed Brown for (j=0; j<dof; j++) { 426c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 427c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 428c4762a1bSJed Brown } 4299566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 430c4762a1bSJed Brown if (i < xs+xm) { 431c4762a1bSJed Brown for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 432c4762a1bSJed Brown islow++; 433c4762a1bSJed Brown } 434c4762a1bSJed Brown } 435c4762a1bSJed Brown if (i > fs+rsbwidth) { /* slow region */ 436c4762a1bSJed Brown for (j=0; j<dof; j++) { 437c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 438c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 439c4762a1bSJed Brown } 4409566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 441c4762a1bSJed Brown if (i > xs) { 442c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 443c4762a1bSJed Brown } 444c4762a1bSJed Brown if (i < xs+xm) { 445c4762a1bSJed Brown for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 446c4762a1bSJed Brown islow++; 447c4762a1bSJed Brown } 448c4762a1bSJed Brown } 449c4762a1bSJed Brown } 4509566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 4519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 4529566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da,PETSC_TRUE,&slope)); 4539566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 4549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da))); 455c4762a1bSJed Brown PetscFunctionReturn(0); 456c4762a1bSJed Brown } 457c4762a1bSJed Brown 458c4762a1bSJed Brown PetscErrorCode FVRHSFunctionslowbuffer_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 459c4762a1bSJed Brown { 460c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 461c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,islow = 0,sf = ctx->sf,fs = ctx->fs,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth; 462c4762a1bSJed Brown PetscReal hxs,hxf; 463c4762a1bSJed Brown PetscScalar *x,*f,*slope; 464c4762a1bSJed Brown Vec Xloc; 465c4762a1bSJed Brown DM da; 466c4762a1bSJed Brown 467c4762a1bSJed Brown PetscFunctionBeginUser; 4689566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 4699566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); 4709566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 471c4762a1bSJed Brown hxs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 472c4762a1bSJed Brown hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 4739566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 4749566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 4759566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 4769566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 4779566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 4789566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da,PETSC_TRUE,&slope)); 4799566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 480c4762a1bSJed Brown 481c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 482c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 483c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 484c4762a1bSJed Brown } 485c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 486c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 487c4762a1bSJed Brown } 488c4762a1bSJed Brown } 489c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 490c4762a1bSJed Brown struct _LimitInfo info; 491c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 492c4762a1bSJed Brown if ((i > sf-lsbwidth-2 && i < sf+1) || (i > fs-2 && i < fs+rsbwidth+1)) { 493c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 4949566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds)); 495c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 4969566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR,2*dof)); 497c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 498c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 499c4762a1bSJed Brown for (j=0; j<dof; j++) { 500c4762a1bSJed Brown PetscScalar jmpL,jmpR; 501c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 502c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 503c4762a1bSJed Brown for (k=0; k<dof; k++) { 504c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 505c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 506c4762a1bSJed Brown } 507c4762a1bSJed Brown } 508c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 509c4762a1bSJed Brown info.m = dof; 510c4762a1bSJed Brown info.hxs = hxs; 511c4762a1bSJed Brown info.hxf = hxf; 512c4762a1bSJed Brown (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope); 513c4762a1bSJed Brown for (j=0; j<dof; j++) { 514c4762a1bSJed Brown PetscScalar tmp = 0; 515c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 516c4762a1bSJed Brown slope[i*dof+j] = tmp; 517c4762a1bSJed Brown } 518c4762a1bSJed Brown } 519c4762a1bSJed Brown } 520c4762a1bSJed Brown 521c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 522c4762a1bSJed Brown PetscReal maxspeed; 523c4762a1bSJed Brown PetscScalar *uL,*uR; 524c4762a1bSJed Brown uL = &ctx->uLR[0]; 525c4762a1bSJed Brown uR = &ctx->uLR[dof]; 526c4762a1bSJed Brown if (i == sf-lsbwidth) { 527c4762a1bSJed Brown for (j=0; j<dof; j++) { 528c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 529c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 530c4762a1bSJed Brown } 5319566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 532c4762a1bSJed Brown if (i < xs+xm) { 533c4762a1bSJed Brown for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 534c4762a1bSJed Brown islow++; 535c4762a1bSJed Brown } 536c4762a1bSJed Brown } 537c4762a1bSJed Brown if (i > sf-lsbwidth && i < sf) { 538c4762a1bSJed Brown for (j=0; j<dof; j++) { 539c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 540c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 541c4762a1bSJed Brown } 5429566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 543c4762a1bSJed Brown if (i > xs) { 544c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 545c4762a1bSJed Brown } 546c4762a1bSJed Brown if (i < xs+xm) { 547c4762a1bSJed Brown for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 548c4762a1bSJed Brown islow++; 549c4762a1bSJed Brown } 550c4762a1bSJed Brown } 551c4762a1bSJed Brown if (i == sf) { /* interface between the slow region and the fast region */ 552c4762a1bSJed Brown for (j=0; j<dof; j++) { 553c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 554c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 555c4762a1bSJed Brown } 5569566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 557c4762a1bSJed Brown if (i > xs) { 558c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 559c4762a1bSJed Brown } 560c4762a1bSJed Brown } 561c4762a1bSJed Brown if (i == fs) { /* interface between the fast region and the slow region */ 562c4762a1bSJed Brown for (j=0; j<dof; j++) { 563c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 564c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 565c4762a1bSJed Brown } 5669566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 567c4762a1bSJed Brown if (i < xs+xm) { 568c4762a1bSJed Brown for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 569c4762a1bSJed Brown islow++; 570c4762a1bSJed Brown } 571c4762a1bSJed Brown } 572c4762a1bSJed Brown if (i > fs && i < fs+rsbwidth) { 573c4762a1bSJed Brown for (j=0; j<dof; j++) { 574c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 575c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 576c4762a1bSJed Brown } 5779566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 578c4762a1bSJed Brown if (i > xs) { 579c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 580c4762a1bSJed Brown } 581c4762a1bSJed Brown if (i < xs+xm) { 582c4762a1bSJed Brown for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 583c4762a1bSJed Brown islow++; 584c4762a1bSJed Brown } 585c4762a1bSJed Brown } 586c4762a1bSJed Brown if (i == fs+rsbwidth) { 587c4762a1bSJed Brown for (j=0; j<dof; j++) { 588c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 589c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 590c4762a1bSJed Brown } 5919566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 592c4762a1bSJed Brown if (i > xs) { 593c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 594c4762a1bSJed Brown } 595c4762a1bSJed Brown } 596c4762a1bSJed Brown } 5979566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 5989566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 5999566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da,PETSC_TRUE,&slope)); 6009566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 601c4762a1bSJed Brown PetscFunctionReturn(0); 602c4762a1bSJed Brown } 603c4762a1bSJed Brown 604c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for fast parts ----------------------------------- */ 605c4762a1bSJed Brown PetscErrorCode FVRHSFunctionfast_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 606c4762a1bSJed Brown { 607c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 608c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,ifast = 0,sf = ctx->sf,fs = ctx->fs; 609c4762a1bSJed Brown PetscReal hxs,hxf; 610c4762a1bSJed Brown PetscScalar *x,*f,*slope; 611c4762a1bSJed Brown Vec Xloc; 612c4762a1bSJed Brown DM da; 613c4762a1bSJed Brown 614c4762a1bSJed Brown PetscFunctionBeginUser; 6159566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 6169566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); 6179566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 618c4762a1bSJed Brown hxs = (ctx->xmax-ctx->xmin)*3.0/8.0/ctx->sf; 619c4762a1bSJed Brown hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fs-ctx->sf); 6209566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 6219566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 6229566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 6239566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 6249566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 6259566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da,PETSC_TRUE,&slope)); 6269566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 627c4762a1bSJed Brown 628c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 629c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 630c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 631c4762a1bSJed Brown } 632c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 633c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 634c4762a1bSJed Brown } 635c4762a1bSJed Brown } 636c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 637c4762a1bSJed Brown struct _LimitInfo info; 638c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 639c4762a1bSJed Brown if (i > sf-2 && i < fs+1) { 6409566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds)); 6419566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR,2*dof)); 642c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 643c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 644c4762a1bSJed Brown for (j=0; j<dof; j++) { 645c4762a1bSJed Brown PetscScalar jmpL,jmpR; 646c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 647c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 648c4762a1bSJed Brown for (k=0; k<dof; k++) { 649c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 650c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 651c4762a1bSJed Brown } 652c4762a1bSJed Brown } 653c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 654c4762a1bSJed Brown info.m = dof; 655c4762a1bSJed Brown info.hxs = hxs; 656c4762a1bSJed Brown info.hxf = hxf; 657c4762a1bSJed Brown (*ctx->limit2)(&info,cjmpL,cjmpR,ctx->sf,ctx->fs,i,ctx->cslope); 658c4762a1bSJed Brown for (j=0; j<dof; j++) { 659c4762a1bSJed Brown PetscScalar tmp = 0; 660c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 661c4762a1bSJed Brown slope[i*dof+j] = tmp; 662c4762a1bSJed Brown } 663c4762a1bSJed Brown } 664c4762a1bSJed Brown } 665c4762a1bSJed Brown 666c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 667c4762a1bSJed Brown PetscReal maxspeed; 668c4762a1bSJed Brown PetscScalar *uL,*uR; 669c4762a1bSJed Brown uL = &ctx->uLR[0]; 670c4762a1bSJed Brown uR = &ctx->uLR[dof]; 671c4762a1bSJed Brown if (i == sf) { /* interface between the slow region and the fast region */ 672c4762a1bSJed Brown for (j=0; j<dof; j++) { 673c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 674c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 675c4762a1bSJed Brown } 6769566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 677c4762a1bSJed Brown if (i < xs+xm) { 678c4762a1bSJed Brown for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf; 679c4762a1bSJed Brown ifast++; 680c4762a1bSJed Brown } 681c4762a1bSJed Brown } 682c4762a1bSJed Brown if (i > sf && i < fs) { /* fast region */ 683c4762a1bSJed Brown for (j=0; j<dof; j++) { 684c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 685c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 686c4762a1bSJed Brown } 6879566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 688c4762a1bSJed Brown if (i > xs) { 689c4762a1bSJed Brown for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf; 690c4762a1bSJed Brown } 691c4762a1bSJed Brown if (i < xs+xm) { 692c4762a1bSJed Brown for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf; 693c4762a1bSJed Brown ifast++; 694c4762a1bSJed Brown } 695c4762a1bSJed Brown } 696c4762a1bSJed Brown if (i == fs) { /* interface between the fast region and the slow region */ 697c4762a1bSJed Brown for (j=0; j<dof; j++) { 698c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 699c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 700c4762a1bSJed Brown } 7019566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 702c4762a1bSJed Brown if (i > xs) { 703c4762a1bSJed Brown for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf; 704c4762a1bSJed Brown } 705c4762a1bSJed Brown } 706c4762a1bSJed Brown } 7079566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 7089566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 7099566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da,PETSC_TRUE,&slope)); 7109566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 711c4762a1bSJed Brown PetscFunctionReturn(0); 712c4762a1bSJed Brown } 713c4762a1bSJed Brown 714c4762a1bSJed Brown int main(int argc,char *argv[]) 715c4762a1bSJed Brown { 716c4762a1bSJed Brown char lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m"; 717c4762a1bSJed Brown PetscFunctionList limiters = 0,physics = 0; 718c4762a1bSJed Brown MPI_Comm comm; 719c4762a1bSJed Brown TS ts; 720c4762a1bSJed Brown DM da; 721c4762a1bSJed Brown Vec X,X0,R; 722c4762a1bSJed Brown FVCtx ctx; 723c4762a1bSJed 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; 724c4762a1bSJed Brown PetscBool view_final = PETSC_FALSE; 725c4762a1bSJed Brown PetscReal ptime; 726c4762a1bSJed Brown 7279566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,0,help)); 728c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 7299566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&ctx,sizeof(ctx))); 730c4762a1bSJed Brown 731c4762a1bSJed Brown /* Register limiters to be available on the command line */ 7329566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"upwind" ,Limit2_Upwind)); 7339566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"lax-wendroff" ,Limit2_LaxWendroff)); 7349566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"beam-warming" ,Limit2_BeamWarming)); 7359566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"fromm" ,Limit2_Fromm)); 7369566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"minmod" ,Limit2_Minmod)); 7379566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"superbee" ,Limit2_Superbee)); 7389566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"mc" ,Limit2_MC)); 7399566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"koren3" ,Limit2_Koren3)); 740c4762a1bSJed Brown 741c4762a1bSJed Brown /* Register physical models to be available on the command line */ 7429566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&physics,"advect" ,PhysicsCreate_Advect)); 743c4762a1bSJed Brown 744c4762a1bSJed Brown ctx.comm = comm; 745c4762a1bSJed Brown ctx.cfl = 0.9; 746c4762a1bSJed Brown ctx.bctype = FVBC_PERIODIC; 747c4762a1bSJed Brown ctx.xmin = -1.0; 748c4762a1bSJed Brown ctx.xmax = 1.0; 749*d0609cedSBarry Smith PetscOptionsBegin(comm,NULL,"Finite Volume solver options",""); 7509566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL)); 7519566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL)); 7529566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-limit","Name of flux imiter to use","",limiters,lname,lname,sizeof(lname),NULL)); 7539566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL)); 7549566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-view_final","Write final solution in ASCII MATLAB format to given file name","",final_fname,final_fname,sizeof(final_fname),&view_final)); 7559566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL)); 7569566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-exact","Compare errors with exact solution","",ctx.exact,&ctx.exact,NULL)); 7579566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL)); 7589566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL)); 7599566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL)); 7609566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-hratio","Spacing ratio","",ctx.hratio,&ctx.hratio,NULL)); 761*d0609cedSBarry Smith PetscOptionsEnd(); 762c4762a1bSJed Brown 763c4762a1bSJed Brown /* Choose the limiter from the list of registered limiters */ 7649566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(limiters,lname,&ctx.limit2)); 7653c633725SBarry Smith PetscCheck(ctx.limit2,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Limiter '%s' not found",lname); 766c4762a1bSJed Brown 767c4762a1bSJed Brown /* Choose the physics from the list of registered models */ 768c4762a1bSJed Brown { 769c4762a1bSJed Brown PetscErrorCode (*r)(FVCtx*); 7709566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(physics,physname,&r)); 7713c633725SBarry Smith PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Physics '%s' not found",physname); 772c4762a1bSJed Brown /* Create the physics, will set the number of fields and their names */ 7739566063dSJacob Faibussowitsch PetscCall((*r)(&ctx)); 774c4762a1bSJed Brown } 775c4762a1bSJed Brown 776c4762a1bSJed Brown /* Create a DMDA to manage the parallel grid */ 7779566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics2.dof,2,NULL,&da)); 7789566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 7799566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 780c4762a1bSJed Brown /* Inform the DMDA of the field names provided by the physics. */ 781c4762a1bSJed Brown /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ 782c4762a1bSJed Brown for (i=0; i<ctx.physics2.dof; i++) { 7839566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da,i,ctx.physics2.fieldname[i])); 784c4762a1bSJed Brown } 7859566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 7869566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 787c4762a1bSJed Brown 788c4762a1bSJed Brown /* Set coordinates of cell centers */ 7899566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da,ctx.xmin+0.5*(ctx.xmax-ctx.xmin)/Mx,ctx.xmax+0.5*(ctx.xmax-ctx.xmin)/Mx,0,0,0,0)); 790c4762a1bSJed Brown 791c4762a1bSJed Brown /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */ 7929566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dof*dof,&ctx.R,dof*dof,&ctx.Rinv,2*dof,&ctx.cjmpLR,1*dof,&ctx.cslope)); 7939566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(2*dof,&ctx.uLR,dof,&ctx.flux,dof,&ctx.speeds)); 794c4762a1bSJed Brown 795c4762a1bSJed Brown /* Create a vector to store the solution and to save the initial state */ 7969566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da,&X)); 7979566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&X0)); 7989566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&R)); 799c4762a1bSJed Brown 800c4762a1bSJed Brown /* create index for slow parts and fast parts, 801c4762a1bSJed Brown count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */ 802c4762a1bSJed Brown count_slow = Mx/(1.0+ctx.hratio/3.0); 80308401ef6SPierre Jolivet PetscCheck(count_slow%2 == 0,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"); 804c4762a1bSJed Brown count_fast = Mx-count_slow; 805c4762a1bSJed Brown ctx.sf = count_slow/2; 806c4762a1bSJed Brown ctx.fs = ctx.sf+count_fast; 8079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm*dof,&index_slow)); 8089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm*dof,&index_fast)); 8099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(6*dof,&index_slowbuffer)); 810c4762a1bSJed Brown if (((AdvectCtx*)ctx.physics2.user)->a > 0) { 811c4762a1bSJed Brown ctx.lsbwidth = 2; 812c4762a1bSJed Brown ctx.rsbwidth = 4; 813c4762a1bSJed Brown } else { 814c4762a1bSJed Brown ctx.lsbwidth = 4; 815c4762a1bSJed Brown ctx.rsbwidth = 2; 816c4762a1bSJed Brown } 817c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 818c4762a1bSJed Brown if (i < ctx.sf-ctx.lsbwidth || i > ctx.fs+ctx.rsbwidth-1) 819c4762a1bSJed Brown for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k; 820c4762a1bSJed Brown else if ((i >= ctx.sf-ctx.lsbwidth && i < ctx.sf) || (i > ctx.fs-1 && i <= ctx.fs+ctx.rsbwidth-1)) 821c4762a1bSJed Brown for (k=0; k<dof; k++) index_slowbuffer[islowbuffer++] = i*dof+k; 822c4762a1bSJed Brown else 823c4762a1bSJed Brown for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k; 824c4762a1bSJed Brown } 8259566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss)); 8269566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf)); 8279566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,islowbuffer,index_slowbuffer,PETSC_COPY_VALUES,&ctx.issb)); 828c4762a1bSJed Brown 829c4762a1bSJed Brown /* Create a time-stepping object */ 8309566063dSJacob Faibussowitsch PetscCall(TSCreate(comm,&ts)); 8319566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts,da)); 8329566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,R,FVRHSFunction_2WaySplit,&ctx)); 8339566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"slow",ctx.iss)); 8349566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"slowbuffer",ctx.issb)); 8359566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"fast",ctx.isf)); 8369566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"slow",NULL,FVRHSFunctionslow_2WaySplit,&ctx)); 8379566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"fast",NULL,FVRHSFunctionfast_2WaySplit,&ctx)); 8389566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"slowbuffer",NULL,FVRHSFunctionslowbuffer_2WaySplit,&ctx)); 839c4762a1bSJed Brown 8409566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSSSP)); 8419566063dSJacob Faibussowitsch /*PetscCall(TSSetType(ts,TSMPRK));*/ 8429566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,10)); 8439566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 844c4762a1bSJed Brown 845c4762a1bSJed Brown /* Compute initial conditions and starting time step */ 8469566063dSJacob Faibussowitsch PetscCall(FVSample_2WaySplit(&ctx,da,0,X0)); 8479566063dSJacob Faibussowitsch PetscCall(FVRHSFunction_2WaySplit(ts,0,X0,X,(void*)&ctx)); /* Initial function evaluation, only used to determine max speed */ 8489566063dSJacob Faibussowitsch PetscCall(VecCopy(X0,X)); /* The function value was not used so we set X=X0 again */ 8499566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt)); 8509566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); /* Take runtime options */ 8519566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD)); 852c4762a1bSJed Brown { 853c4762a1bSJed Brown PetscInt steps; 854c4762a1bSJed Brown PetscScalar mass_initial,mass_final,mass_difference,mass_differenceg; 855c4762a1bSJed Brown const PetscScalar *ptr_X,*ptr_X0; 856c4762a1bSJed Brown const PetscReal hs = (ctx.xmax-ctx.xmin)*3.0/4.0/count_slow; 857c4762a1bSJed Brown const PetscReal hf = (ctx.xmax-ctx.xmin)/4.0/count_fast; 858c4762a1bSJed Brown 8599566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,X)); 8609566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&ptime)); 8619566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&steps)); 862c4762a1bSJed Brown /* calculate the total mass at initial time and final time */ 863c4762a1bSJed Brown mass_initial = 0.0; 864c4762a1bSJed Brown mass_final = 0.0; 8659566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,X0,(void*)&ptr_X0)); 8669566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,X,(void*)&ptr_X)); 867c4762a1bSJed Brown for (i=xs;i<xs+xm;i++) { 868c4762a1bSJed Brown if (i < ctx.sf || i > ctx.fs-1) { 869c4762a1bSJed Brown for (k=0; k<dof; k++) { 870c4762a1bSJed Brown mass_initial = mass_initial + hs*ptr_X0[i*dof+k]; 871c4762a1bSJed Brown mass_final = mass_final + hs*ptr_X[i*dof+k]; 872c4762a1bSJed Brown } 873c4762a1bSJed Brown } else { 874c4762a1bSJed Brown for (k=0; k<dof; k++) { 875c4762a1bSJed Brown mass_initial = mass_initial + hf*ptr_X0[i*dof+k]; 876c4762a1bSJed Brown mass_final = mass_final + hf*ptr_X[i*dof+k]; 877c4762a1bSJed Brown } 878c4762a1bSJed Brown } 879c4762a1bSJed Brown } 8809566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,X0,(void*)&ptr_X0)); 8819566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,X,(void*)&ptr_X)); 882c4762a1bSJed Brown mass_difference = mass_final - mass_initial; 8839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&mass_difference,&mass_differenceg,1,MPIU_SCALAR,MPIU_SUM,comm)); 8849566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Mass difference %g\n",(double)mass_differenceg)); 8859566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Final time %g, steps %D\n",(double)ptime,steps)); 8869566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Maximum allowable stepsize according to CFL %g\n",(double)(1.0/ctx.cfl_idt))); 887c4762a1bSJed Brown if (ctx.exact) { 888c4762a1bSJed Brown PetscReal nrm1=0; 8899566063dSJacob Faibussowitsch PetscCall(SolutionErrorNorms_2WaySplit(&ctx,da,ptime,X,&nrm1)); 8909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1)); 891c4762a1bSJed Brown } 892c4762a1bSJed Brown if (ctx.simulation) { 893c4762a1bSJed Brown PetscReal nrm1=0; 894c4762a1bSJed Brown PetscViewer fd; 895c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = "binaryoutput"; 896c4762a1bSJed Brown Vec XR; 897c4762a1bSJed Brown PetscBool flg; 898c4762a1bSJed Brown const PetscScalar *ptr_XR; 8999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg)); 9003c633725SBarry Smith PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 9019566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&fd)); 9029566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X0,&XR)); 9039566063dSJacob Faibussowitsch PetscCall(VecLoad(XR,fd)); 9049566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 9059566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&ptr_X)); 9069566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XR,&ptr_XR)); 907c4762a1bSJed Brown for (i=xs;i<xs+xm;i++) { 908c4762a1bSJed Brown if (i < ctx.sf || i > ctx.fs-1) 909c4762a1bSJed Brown for (k=0; k<dof; k++) nrm1 = nrm1 + hs*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]); 910c4762a1bSJed Brown else 911c4762a1bSJed Brown for (k=0; k<dof; k++) nrm1 = nrm1 + hf*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]); 912c4762a1bSJed Brown } 9139566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&ptr_X)); 9149566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XR,&ptr_XR)); 9159566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1)); 9169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XR)); 917c4762a1bSJed Brown } 918c4762a1bSJed Brown } 919c4762a1bSJed Brown 9209566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD)); 9219566063dSJacob Faibussowitsch if (draw & 0x1) PetscCall(VecView(X0,PETSC_VIEWER_DRAW_WORLD)); 9229566063dSJacob Faibussowitsch if (draw & 0x2) PetscCall(VecView(X,PETSC_VIEWER_DRAW_WORLD)); 923c4762a1bSJed Brown if (draw & 0x4) { 924c4762a1bSJed Brown Vec Y; 9259566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&Y)); 9269566063dSJacob Faibussowitsch PetscCall(FVSample_2WaySplit(&ctx,da,ptime,Y)); 9279566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y,-1,X)); 9289566063dSJacob Faibussowitsch PetscCall(VecView(Y,PETSC_VIEWER_DRAW_WORLD)); 9299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 930c4762a1bSJed Brown } 931c4762a1bSJed Brown 932c4762a1bSJed Brown if (view_final) { 933c4762a1bSJed Brown PetscViewer viewer; 9349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer)); 9359566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB)); 9369566063dSJacob Faibussowitsch PetscCall(VecView(X,viewer)); 9379566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 9389566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 939c4762a1bSJed Brown } 940c4762a1bSJed Brown 941c4762a1bSJed Brown /* Clean up */ 9429566063dSJacob Faibussowitsch PetscCall((*ctx.physics2.destroy)(ctx.physics2.user)); 9439566063dSJacob Faibussowitsch for (i=0; i<ctx.physics2.dof; i++) PetscCall(PetscFree(ctx.physics2.fieldname[i])); 9449566063dSJacob Faibussowitsch PetscCall(PetscFree4(ctx.R,ctx.Rinv,ctx.cjmpLR,ctx.cslope)); 9459566063dSJacob Faibussowitsch PetscCall(PetscFree3(ctx.uLR,ctx.flux,ctx.speeds)); 9469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 9479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X0)); 9489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&R)); 9499566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 9509566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 9519566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.iss)); 9529566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.isf)); 9539566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.issb)); 9549566063dSJacob Faibussowitsch PetscCall(PetscFree(index_slow)); 9559566063dSJacob Faibussowitsch PetscCall(PetscFree(index_fast)); 9569566063dSJacob Faibussowitsch PetscCall(PetscFree(index_slowbuffer)); 9579566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&limiters)); 9589566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&physics)); 9599566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 960b122ec5aSJacob Faibussowitsch return 0; 961c4762a1bSJed Brown } 962c4762a1bSJed Brown 963c4762a1bSJed Brown /*TEST 964c4762a1bSJed Brown 965c4762a1bSJed Brown build: 966f56ea12dSJed Brown requires: !complex 967c4762a1bSJed Brown depends: finitevolume1d.c 968c4762a1bSJed Brown 969c4762a1bSJed Brown test: 970c4762a1bSJed Brown suffix: 1 971c4762a1bSJed 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 972c4762a1bSJed Brown 973c4762a1bSJed Brown test: 974c4762a1bSJed Brown suffix: 2 975c4762a1bSJed 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 976c4762a1bSJed Brown output_file: output/ex6_1.out 977c4762a1bSJed Brown 978c4762a1bSJed Brown TEST*/ 979