1c4762a1bSJed Brown static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n" 2c4762a1bSJed Brown " advection - Constant coefficient scalar advection\n" 3c4762a1bSJed Brown " u_t + (a*u)_x = 0\n" 4c4762a1bSJed Brown " for this toy problem, we choose different meshsizes for different sub-domains (slow-medium-fast-medium-slow), \n" 5c4762a1bSJed Brown " the meshsize ratio between two adjacient sub-domains is controlled with -hratio,\n" 6c4762a1bSJed Brown " exact - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n" 7c4762a1bSJed Brown " the states across shocks and rarefactions\n" 8c4762a1bSJed Brown " simulation - use reference solution which is generated by smaller time step size to be true solution,\n" 9c4762a1bSJed Brown " also the reference solution should be generated by user and stored in a binary file.\n" 10c4762a1bSJed Brown " characteristic - Limit the characteristic variables, this is usually preferred (default)\n" 11c4762a1bSJed Brown "Several initial conditions can be chosen with -initial N\n\n" 12c4762a1bSJed Brown "The problem size should be set with -da_grid_x M\n\n"; 13c4762a1bSJed Brown 14c4762a1bSJed Brown #include <petscts.h> 15c4762a1bSJed Brown #include <petscdm.h> 16c4762a1bSJed Brown #include <petscdmda.h> 17c4762a1bSJed Brown #include <petscdraw.h> 18c4762a1bSJed Brown #include "finitevolume1d.h" 19c4762a1bSJed Brown 209fbee547SJacob Faibussowitsch static inline PetscReal RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax) { PetscReal range = xmax-xmin; return xmin +PetscFmodReal(range+PetscFmodReal(a,range),range); } 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* --------------------------------- Advection ----------------------------------- */ 23c4762a1bSJed Brown typedef struct { 24c4762a1bSJed Brown PetscReal a; /* advective velocity */ 25c4762a1bSJed Brown } AdvectCtx; 26c4762a1bSJed Brown 27c4762a1bSJed Brown static PetscErrorCode PhysicsRiemann_Advect(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed) 28c4762a1bSJed Brown { 29c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx*)vctx; 30c4762a1bSJed Brown PetscReal speed; 31c4762a1bSJed Brown 32c4762a1bSJed Brown PetscFunctionBeginUser; 33c4762a1bSJed Brown speed = ctx->a; 34c4762a1bSJed Brown flux[0] = PetscMax(0,speed)*uL[0] + PetscMin(0,speed)*uR[0]; 35c4762a1bSJed Brown *maxspeed = speed; 36c4762a1bSJed Brown PetscFunctionReturn(0); 37c4762a1bSJed Brown } 38c4762a1bSJed Brown 39c4762a1bSJed Brown static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds) 40c4762a1bSJed Brown { 41c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx*)vctx; 42c4762a1bSJed Brown 43c4762a1bSJed Brown PetscFunctionBeginUser; 44c4762a1bSJed Brown X[0] = 1.; 45c4762a1bSJed Brown Xi[0] = 1.; 46c4762a1bSJed Brown speeds[0] = ctx->a; 47c4762a1bSJed Brown PetscFunctionReturn(0); 48c4762a1bSJed Brown } 49c4762a1bSJed Brown 50c4762a1bSJed Brown static PetscErrorCode PhysicsSample_Advect(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u) 51c4762a1bSJed Brown { 52c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx*)vctx; 53c4762a1bSJed Brown PetscReal a = ctx->a,x0; 54c4762a1bSJed Brown 55c4762a1bSJed Brown PetscFunctionBeginUser; 56c4762a1bSJed Brown switch (bctype) { 57c4762a1bSJed Brown case FVBC_OUTFLOW: x0 = x-a*t; break; 58c4762a1bSJed Brown case FVBC_PERIODIC: x0 = RangeMod(x-a*t,xmin,xmax); break; 59c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown BCType"); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown switch (initial) { 62c4762a1bSJed Brown case 0: u[0] = (x0 < 0) ? 1 : -1; break; 63c4762a1bSJed Brown case 1: u[0] = (x0 < 0) ? -1 : 1; break; 64c4762a1bSJed Brown case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break; 65c4762a1bSJed Brown case 3: u[0] = PetscSinReal(2*PETSC_PI*x0); break; 66c4762a1bSJed Brown case 4: u[0] = PetscAbs(x0); break; 67c4762a1bSJed Brown case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2*PETSC_PI*x0)); break; 68c4762a1bSJed Brown case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2-x0 : 0)); break; 69c4762a1bSJed Brown case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI*x0),10.0);break; 70c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition"); 71c4762a1bSJed Brown } 72c4762a1bSJed Brown PetscFunctionReturn(0); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 75c4762a1bSJed Brown static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx) 76c4762a1bSJed Brown { 77c4762a1bSJed Brown AdvectCtx *user; 78c4762a1bSJed Brown 79c4762a1bSJed Brown PetscFunctionBeginUser; 809566063dSJacob Faibussowitsch PetscCall(PetscNew(&user)); 81c4762a1bSJed Brown ctx->physics2.sample2 = PhysicsSample_Advect; 82c4762a1bSJed Brown ctx->physics2.riemann2 = PhysicsRiemann_Advect; 83c4762a1bSJed Brown ctx->physics2.characteristic2 = PhysicsCharacteristic_Advect; 84c4762a1bSJed Brown ctx->physics2.destroy = PhysicsDestroy_SimpleFree; 85c4762a1bSJed Brown ctx->physics2.user = user; 86c4762a1bSJed Brown ctx->physics2.dof = 1; 879566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("u",&ctx->physics2.fieldname[0])); 88c4762a1bSJed Brown user->a = 1; 89d0609cedSBarry Smith PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection",""); 90c4762a1bSJed Brown { 919566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-physics_advect_a","Speed","",user->a,&user->a,NULL)); 92c4762a1bSJed Brown } 93d0609cedSBarry Smith PetscOptionsEnd(); 94c4762a1bSJed Brown PetscFunctionReturn(0); 95c4762a1bSJed Brown } 96c4762a1bSJed Brown 97c4762a1bSJed Brown PetscErrorCode FVSample_3WaySplit(FVCtx *ctx,DM da,PetscReal time,Vec U) 98c4762a1bSJed Brown { 99c4762a1bSJed Brown PetscScalar *u,*uj,xj,xi; 100c4762a1bSJed Brown PetscInt i,j,k,dof,xs,xm,Mx; 101c4762a1bSJed Brown const PetscInt N = 200; 102c4762a1bSJed Brown PetscReal hs,hm,hf; 103c4762a1bSJed Brown 104c4762a1bSJed Brown PetscFunctionBeginUser; 1053c633725SBarry Smith PetscCheck(ctx->physics2.sample2,PETSC_COMM_SELF,PETSC_ERR_SUP,"Physics has not provided a sampling function"); 1069566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 1079566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 1089566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,U,&u)); 1099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof,&uj)); 110c4762a1bSJed Brown 111c4762a1bSJed Brown hs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 112c4762a1bSJed Brown hm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 113c4762a1bSJed Brown hf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 114c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 115c4762a1bSJed Brown if (i < ctx->sm) { 116c4762a1bSJed Brown xi = ctx->xmin+0.5*hs+i*hs; 117c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 118c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 119c4762a1bSJed Brown for (j=0; j<N+1; j++) { 120c4762a1bSJed Brown xj = xi+hs*(j-N/2)/(PetscReal)N; 1219566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj)); 122c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 123c4762a1bSJed Brown } 124c4762a1bSJed Brown } else if (i < ctx->mf) { 125c4762a1bSJed Brown xi = ctx->xmin+ctx->sm*hs+0.5*hm+(i-ctx->sm)*hm; 126c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 127c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 128c4762a1bSJed Brown for (j=0; j<N+1; j++) { 129c4762a1bSJed Brown xj = xi+hm*(j-N/2)/(PetscReal)N; 1309566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj)); 131c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 132c4762a1bSJed Brown } 133c4762a1bSJed Brown } else if (i < ctx->fm) { 134c4762a1bSJed Brown xi = ctx->xmin+ctx->sm*hs+(ctx->mf-ctx->sm)*hm+0.5*hf+(i-ctx->mf)*hf; 135c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 136c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 137c4762a1bSJed Brown for (j=0; j<N+1; j++) { 138c4762a1bSJed Brown xj = xi+hf*(j-N/2)/(PetscReal)N; 1399566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj)); 140c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 141c4762a1bSJed Brown } 142c4762a1bSJed Brown } else if (i < ctx->ms) { 143c4762a1bSJed Brown xi = ctx->xmin+ctx->sm*hs+(ctx->mf-ctx->sm)*hm+(ctx->fm-ctx->mf)*hf+0.5*hm+(i-ctx->fm)*hm; 144c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 145c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 146c4762a1bSJed Brown for (j=0; j<N+1; j++) { 147c4762a1bSJed Brown xj = xi+hm*(j-N/2)/(PetscReal)N; 1489566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj)); 149c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 150c4762a1bSJed Brown } 151c4762a1bSJed Brown } else { 152c4762a1bSJed Brown xi = ctx->xmin+ctx->sm*hs+(ctx->mf-ctx->sm)*hm+(ctx->fm-ctx->mf)*hf+(ctx->ms-ctx->fm)*hm+0.5*hs+(i-ctx->ms)*hs; 153c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 154c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] = 0; 155c4762a1bSJed Brown for (j=0; j<N+1; j++) { 156c4762a1bSJed Brown xj = xi+hs*(j-N/2)/(PetscReal)N; 1579566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj)); 158c4762a1bSJed Brown for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 159c4762a1bSJed Brown } 160c4762a1bSJed Brown } 161c4762a1bSJed Brown } 1629566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,U,&u)); 1639566063dSJacob Faibussowitsch PetscCall(PetscFree(uj)); 164c4762a1bSJed Brown PetscFunctionReturn(0); 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167c4762a1bSJed Brown static PetscErrorCode SolutionErrorNorms_3WaySplit(FVCtx *ctx,DM da,PetscReal t,Vec X,PetscReal *nrm1) 168c4762a1bSJed Brown { 169c4762a1bSJed Brown Vec Y; 170c4762a1bSJed Brown PetscInt i,Mx; 171c4762a1bSJed Brown const PetscScalar *ptr_X,*ptr_Y; 172c4762a1bSJed Brown PetscReal hs,hm,hf; 173c4762a1bSJed Brown 174c4762a1bSJed Brown PetscFunctionBeginUser; 1759566063dSJacob Faibussowitsch PetscCall(VecGetSize(X,&Mx)); 176c4762a1bSJed Brown hs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 177c4762a1bSJed Brown hm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 178c4762a1bSJed Brown hf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 1799566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&Y)); 1809566063dSJacob Faibussowitsch PetscCall(FVSample_3WaySplit(ctx,da,t,Y)); 1819566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&ptr_X)); 1829566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&ptr_Y)); 183c4762a1bSJed Brown for (i=0;i<Mx;i++) { 184c4762a1bSJed Brown if (i < ctx->sm || i > ctx->ms - 1) *nrm1 += hs*PetscAbs(ptr_X[i]-ptr_Y[i]); 185c4762a1bSJed Brown else if (i < ctx->mf || i > ctx->fm - 1) *nrm1 += hm*PetscAbs(ptr_X[i]-ptr_Y[i]); 186c4762a1bSJed Brown else *nrm1 += hf*PetscAbs(ptr_X[i]-ptr_Y[i]); 187c4762a1bSJed Brown } 1889566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&ptr_X)); 1899566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&ptr_Y)); 1909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 191c4762a1bSJed Brown PetscFunctionReturn(0); 192c4762a1bSJed Brown } 193c4762a1bSJed Brown 194c4762a1bSJed Brown PetscErrorCode FVRHSFunction_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 195c4762a1bSJed Brown { 196c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 197c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,sm = ctx->sm,mf = ctx->mf,fm = ctx->fm,ms = ctx->ms; 198c4762a1bSJed Brown PetscReal hxf,hxm,hxs,cfl_idt = 0; 199c4762a1bSJed Brown PetscScalar *x,*f,*slope; 200c4762a1bSJed Brown Vec Xloc; 201c4762a1bSJed Brown DM da; 202c4762a1bSJed Brown 203c4762a1bSJed Brown PetscFunctionBeginUser; 2049566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 2059566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); /* Xloc contains ghost points */ 2069566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); /* Mx is the number of center points */ 207c4762a1bSJed Brown hxs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 208c4762a1bSJed Brown hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 209c4762a1bSJed Brown hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 2109566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); /* X is solution vector which does not contain ghost points */ 2119566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 212c4762a1bSJed Brown 2139566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); /* F is the right hand side function corresponds to center points */ 214c4762a1bSJed Brown 2159566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 2169566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,F,&f)); 2179566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da,PETSC_TRUE,&slope)); /* contains ghost points */ 218c4762a1bSJed Brown 2199566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 220c4762a1bSJed Brown 221c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 222c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 223c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 224c4762a1bSJed Brown } 225c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 226c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 227c4762a1bSJed Brown } 228c4762a1bSJed Brown } 229c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 230c4762a1bSJed Brown struct _LimitInfo info; 231c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 232c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 2339566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds)); 234c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 2359566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR,2*dof)); 236c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 237c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 238c4762a1bSJed Brown for (j=0; j<dof; j++) { 239c4762a1bSJed Brown PetscScalar jmpL,jmpR; 240c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 241c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 242c4762a1bSJed Brown for (k=0; k<dof; k++) { 243c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 244c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 245c4762a1bSJed Brown } 246c4762a1bSJed Brown } 247c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 248c4762a1bSJed Brown info.m = dof; 249c4762a1bSJed Brown info.hxs = hxs; 250c4762a1bSJed Brown info.hxm = hxm; 251c4762a1bSJed Brown info.hxf = hxf; 252c4762a1bSJed Brown (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope); 253c4762a1bSJed Brown for (j=0; j<dof; j++) { 254c4762a1bSJed Brown PetscScalar tmp = 0; 255c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 256c4762a1bSJed Brown slope[i*dof+j] = tmp; 257c4762a1bSJed Brown } 258c4762a1bSJed Brown } 259c4762a1bSJed Brown 260c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 261c4762a1bSJed Brown PetscReal maxspeed; 262c4762a1bSJed Brown PetscScalar *uL,*uR; 263c4762a1bSJed Brown uL = &ctx->uLR[0]; 264c4762a1bSJed Brown uR = &ctx->uLR[dof]; 265c4762a1bSJed Brown if (i < sm || i > ms) { /* slow region */ 266c4762a1bSJed Brown for (j=0; j<dof; j++) { 267c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 268c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 269c4762a1bSJed Brown } 2709566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 271c4762a1bSJed Brown if (i > xs) { 272c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs; 273c4762a1bSJed Brown } 274c4762a1bSJed Brown if (i < xs+xm) { 275c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs; 276c4762a1bSJed Brown } 277c4762a1bSJed Brown } else if (i == sm) { /* interface between slow and medium component */ 278c4762a1bSJed Brown for (j=0; j<dof; j++) { 279c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 280c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 281c4762a1bSJed Brown } 2829566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 283c4762a1bSJed Brown if (i > xs) { 284c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs; 285c4762a1bSJed Brown } 286c4762a1bSJed Brown if (i < xs+xm) { 287c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxm; 288c4762a1bSJed Brown } 289c4762a1bSJed Brown } else if (i == ms) { /* interface between medium and slow regions */ 290c4762a1bSJed Brown for (j=0; j<dof; j++) { 291c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 292c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 293c4762a1bSJed Brown } 2949566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 295c4762a1bSJed Brown if (i > xs) { 296c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxm; 297c4762a1bSJed Brown } 298c4762a1bSJed Brown if (i < xs+xm) { 299c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs; 300c4762a1bSJed Brown } 301c4762a1bSJed Brown } else if (i < mf || i > fm) { /* medium region */ 302c4762a1bSJed Brown for (j=0; j<dof; j++) { 303c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 304c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 305c4762a1bSJed Brown } 3069566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 307c4762a1bSJed Brown if (i > xs) { 308c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxm; 309c4762a1bSJed Brown } 310c4762a1bSJed Brown if (i < xs+xm) { 311c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxm; 312c4762a1bSJed Brown } 313c4762a1bSJed Brown } else if (i == mf) { /* interface between medium and fast regions */ 314c4762a1bSJed Brown for (j=0; j<dof; j++) { 315c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 316c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 317c4762a1bSJed Brown } 3189566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 319c4762a1bSJed Brown if (i > xs) { 320c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxm; 321c4762a1bSJed Brown } 322c4762a1bSJed Brown if (i < xs+xm) { 323c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxf; 324c4762a1bSJed Brown } 325c4762a1bSJed Brown } else if (i == fm) { /* interface between fast and medium regions */ 326c4762a1bSJed Brown for (j=0; j<dof; j++) { 327c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 328c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 329c4762a1bSJed Brown } 3309566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 331c4762a1bSJed Brown if (i > xs) { 332c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf; 333c4762a1bSJed Brown } 334c4762a1bSJed Brown if (i < xs+xm) { 335c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxm; 336c4762a1bSJed Brown } 337c4762a1bSJed Brown } else { /* fast region */ 338c4762a1bSJed Brown for (j=0; j<dof; j++) { 339c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 340c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 341c4762a1bSJed Brown } 3429566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 343c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */ 344c4762a1bSJed Brown if (i > xs) { 345c4762a1bSJed Brown for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf; 346c4762a1bSJed Brown } 347c4762a1bSJed Brown if (i < xs+xm) { 348c4762a1bSJed Brown for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxf; 349c4762a1bSJed Brown } 350c4762a1bSJed Brown } 351c4762a1bSJed Brown } 3529566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 3539566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,F,&f)); 3549566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da,PETSC_TRUE,&slope)); 3559566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 3569566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da))); 357c4762a1bSJed Brown if (0) { 358c4762a1bSJed Brown /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */ 359c4762a1bSJed Brown PetscReal dt,tnow; 3609566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&dt)); 3619566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts,&tnow)); 362c4762a1bSJed Brown if (dt > 0.5/ctx->cfl_idt) { 3639566063dSJacob 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))); 364c4762a1bSJed Brown } 365c4762a1bSJed Brown } 366c4762a1bSJed Brown PetscFunctionReturn(0); 367c4762a1bSJed Brown } 368c4762a1bSJed Brown 369c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */ 370c4762a1bSJed Brown PetscErrorCode FVRHSFunctionslow_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 371c4762a1bSJed Brown { 372c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 373c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,islow = 0,sm = ctx->sm,ms = ctx->ms,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth; 374c4762a1bSJed Brown PetscReal hxs,hxm,hxf,cfl_idt = 0; 375c4762a1bSJed Brown PetscScalar *x,*f,*slope; 376c4762a1bSJed Brown Vec Xloc; 377c4762a1bSJed Brown DM da; 378c4762a1bSJed Brown 379c4762a1bSJed Brown PetscFunctionBeginUser; 3809566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 3819566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); 3829566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 383c4762a1bSJed Brown hxs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 384c4762a1bSJed Brown hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 385c4762a1bSJed Brown hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 3869566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 3879566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 3889566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 3899566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 3909566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 3919566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da,PETSC_TRUE,&slope)); 3929566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 393c4762a1bSJed Brown 394c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 395c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 396c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 397c4762a1bSJed Brown } 398c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 399c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 400c4762a1bSJed Brown } 401c4762a1bSJed Brown } 402c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 403c4762a1bSJed Brown struct _LimitInfo info; 404c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 405c4762a1bSJed Brown if (i < sm-lsbwidth+1 || i > ms+rsbwidth-2) { /* slow components and the first and last fast components */ 406c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 4079566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds)); 408c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 4099566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR,2*dof)); 410c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 411c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 412c4762a1bSJed Brown for (j=0; j<dof; j++) { 413c4762a1bSJed Brown PetscScalar jmpL,jmpR; 414c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 415c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 416c4762a1bSJed Brown for (k=0; k<dof; k++) { 417c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 418c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 419c4762a1bSJed Brown } 420c4762a1bSJed Brown } 421c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 422c4762a1bSJed Brown info.m = dof; 423c4762a1bSJed Brown info.hxs = hxs; 424c4762a1bSJed Brown info.hxm = hxm; 425c4762a1bSJed Brown info.hxf = hxf; 426c4762a1bSJed Brown (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope); 427c4762a1bSJed Brown for (j=0; j<dof; j++) { 428c4762a1bSJed Brown PetscScalar tmp = 0; 429c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 430c4762a1bSJed Brown slope[i*dof+j] = tmp; 431c4762a1bSJed Brown } 432c4762a1bSJed Brown } 433c4762a1bSJed Brown } 434c4762a1bSJed Brown 435c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 436c4762a1bSJed Brown PetscReal maxspeed; 437c4762a1bSJed Brown PetscScalar *uL,*uR; 438c4762a1bSJed Brown uL = &ctx->uLR[0]; 439c4762a1bSJed Brown uR = &ctx->uLR[dof]; 440c4762a1bSJed Brown if (i < sm-lsbwidth) { /* 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 } 4459566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 446c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */ 447c4762a1bSJed Brown if (i > xs) { 448c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 449c4762a1bSJed Brown } 450c4762a1bSJed Brown if (i < xs+xm) { 451c4762a1bSJed Brown for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 452c4762a1bSJed Brown islow++; 453c4762a1bSJed Brown } 454c4762a1bSJed Brown } 455c4762a1bSJed Brown if (i == sm-lsbwidth) { /* interface between slow and medium regions */ 456c4762a1bSJed Brown for (j=0; j<dof; j++) { 457c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 458c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 459c4762a1bSJed Brown } 4609566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 461c4762a1bSJed Brown if (i > xs) { 462c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 463c4762a1bSJed Brown } 464c4762a1bSJed Brown } 465c4762a1bSJed Brown if (i == ms+rsbwidth) { /* interface between medium and slow regions */ 466c4762a1bSJed Brown for (j=0; j<dof; j++) { 467c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 468c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 469c4762a1bSJed Brown } 4709566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 471c4762a1bSJed Brown if (i < xs+xm) { 472c4762a1bSJed Brown for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 473c4762a1bSJed Brown islow++; 474c4762a1bSJed Brown } 475c4762a1bSJed Brown } 476c4762a1bSJed Brown if (i > ms+rsbwidth) { /* slow region */ 477c4762a1bSJed Brown for (j=0; j<dof; j++) { 478c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 479c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 480c4762a1bSJed Brown } 4819566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 482c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */ 483c4762a1bSJed Brown if (i > xs) { 484c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 485c4762a1bSJed Brown } 486c4762a1bSJed Brown if (i < xs+xm) { 487c4762a1bSJed Brown for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 488c4762a1bSJed Brown islow++; 489c4762a1bSJed Brown } 490c4762a1bSJed Brown } 491c4762a1bSJed Brown } 4929566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 4939566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 4949566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da,PETSC_TRUE,&slope)); 4959566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 4969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da))); 497c4762a1bSJed Brown PetscFunctionReturn(0); 498c4762a1bSJed Brown } 499c4762a1bSJed Brown 500c4762a1bSJed Brown PetscErrorCode FVRHSFunctionslowbuffer_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 501c4762a1bSJed Brown { 502c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 503c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,islowbuffer = 0,sm = ctx->sm,ms = ctx->ms,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth; 504c4762a1bSJed Brown PetscReal hxs,hxm,hxf; 505c4762a1bSJed Brown PetscScalar *x,*f,*slope; 506c4762a1bSJed Brown Vec Xloc; 507c4762a1bSJed Brown DM da; 508c4762a1bSJed Brown 509c4762a1bSJed Brown PetscFunctionBeginUser; 5109566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 5119566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); 5129566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 513c4762a1bSJed Brown hxs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 514c4762a1bSJed Brown hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 515c4762a1bSJed Brown hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 5169566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 5179566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 5189566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 5199566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 5209566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 5219566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da,PETSC_TRUE,&slope)); 5229566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 523c4762a1bSJed Brown 524c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 525c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 526c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 527c4762a1bSJed Brown } 528c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 529c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 530c4762a1bSJed Brown } 531c4762a1bSJed Brown } 532c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 533c4762a1bSJed Brown struct _LimitInfo info; 534c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 535c4762a1bSJed Brown if ((i > sm-lsbwidth-2 && i < sm+1) || (i > ms-2 && i < ms+rsbwidth+1)) { 536c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 5379566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds)); 538c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 5399566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR,2*dof)); 540c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 541c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 542c4762a1bSJed Brown for (j=0; j<dof; j++) { 543c4762a1bSJed Brown PetscScalar jmpL,jmpR; 544c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 545c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 546c4762a1bSJed Brown for (k=0; k<dof; k++) { 547c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 548c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 549c4762a1bSJed Brown } 550c4762a1bSJed Brown } 551c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 552c4762a1bSJed Brown info.m = dof; 553c4762a1bSJed Brown info.hxs = hxs; 554c4762a1bSJed Brown info.hxm = hxm; 555c4762a1bSJed Brown info.hxf = hxf; 556c4762a1bSJed Brown (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope); 557c4762a1bSJed Brown for (j=0; j<dof; j++) { 558c4762a1bSJed Brown PetscScalar tmp = 0; 559c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 560c4762a1bSJed Brown slope[i*dof+j] = tmp; 561c4762a1bSJed Brown } 562c4762a1bSJed Brown } 563c4762a1bSJed Brown } 564c4762a1bSJed Brown 565c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 566c4762a1bSJed Brown PetscReal maxspeed; 567c4762a1bSJed Brown PetscScalar *uL,*uR; 568c4762a1bSJed Brown uL = &ctx->uLR[0]; 569c4762a1bSJed Brown uR = &ctx->uLR[dof]; 570c4762a1bSJed Brown if (i == sm-lsbwidth) { 571c4762a1bSJed Brown for (j=0; j<dof; j++) { 572c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 573c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 574c4762a1bSJed Brown } 5759566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 576c4762a1bSJed Brown if (i < xs+xm) { 577c4762a1bSJed Brown for (j=0; j<dof; j++) f[islowbuffer*dof+j] += ctx->flux[j]/hxs; 578c4762a1bSJed Brown islowbuffer++; 579c4762a1bSJed Brown } 580c4762a1bSJed Brown } 581c4762a1bSJed Brown if (i > sm-lsbwidth && i < sm) { 582c4762a1bSJed Brown for (j=0; j<dof; j++) { 583c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 584c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 585c4762a1bSJed Brown } 5869566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 587c4762a1bSJed Brown if (i > xs) { 588c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islowbuffer-1)*dof+j] -= ctx->flux[j]/hxs; 589c4762a1bSJed Brown } 590c4762a1bSJed Brown if (i < xs+xm) { 591c4762a1bSJed Brown for (j=0; j<dof; j++) f[islowbuffer*dof+j] += ctx->flux[j]/hxs; 592c4762a1bSJed Brown islowbuffer++; 593c4762a1bSJed Brown } 594c4762a1bSJed Brown } 595c4762a1bSJed Brown if (i == sm) { /* interface between the slow region and the medium region */ 596c4762a1bSJed Brown for (j=0; j<dof; j++) { 597c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 598c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 599c4762a1bSJed Brown } 6009566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 601c4762a1bSJed Brown if (i > xs) { 602c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islowbuffer-1)*dof+j] -= ctx->flux[j]/hxs; 603c4762a1bSJed Brown } 604c4762a1bSJed Brown } 605c4762a1bSJed Brown if (i == ms) { /* interface between the medium region and the slow region */ 606c4762a1bSJed Brown for (j=0; j<dof; j++) { 607c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 608c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 609c4762a1bSJed Brown } 6109566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 611c4762a1bSJed Brown if (i < xs+xm) { 612c4762a1bSJed Brown for (j=0; j<dof; j++) f[islowbuffer*dof+j] += ctx->flux[j]/hxs; 613c4762a1bSJed Brown islowbuffer++; 614c4762a1bSJed Brown } 615c4762a1bSJed Brown } 616c4762a1bSJed Brown if (i > ms && i < ms+rsbwidth) { 617c4762a1bSJed Brown for (j=0; j<dof; j++) { 618c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 619c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 620c4762a1bSJed Brown } 6219566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 622c4762a1bSJed Brown if (i > xs) { 623c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islowbuffer-1)*dof+j] -= ctx->flux[j]/hxs; 624c4762a1bSJed Brown } 625c4762a1bSJed Brown if (i < xs+xm) { 626c4762a1bSJed Brown for (j=0; j<dof; j++) f[islowbuffer*dof+j] += ctx->flux[j]/hxs; 627c4762a1bSJed Brown islowbuffer++; 628c4762a1bSJed Brown } 629c4762a1bSJed Brown } 630c4762a1bSJed Brown if (i == ms+rsbwidth) { 631c4762a1bSJed Brown for (j=0; j<dof; j++) { 632c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 633c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 634c4762a1bSJed Brown } 6359566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 636c4762a1bSJed Brown if (i > xs) { 637c4762a1bSJed Brown for (j=0; j<dof; j++) f[(islowbuffer-1)*dof+j] -= ctx->flux[j]/hxs; 638c4762a1bSJed Brown } 639c4762a1bSJed Brown } 640c4762a1bSJed Brown } 6419566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 6429566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 6439566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da,PETSC_TRUE,&slope)); 6449566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 645c4762a1bSJed Brown PetscFunctionReturn(0); 646c4762a1bSJed Brown } 647c4762a1bSJed Brown 648c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for medium components ----------------------------------- */ 649c4762a1bSJed Brown PetscErrorCode FVRHSFunctionmedium_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 650c4762a1bSJed Brown { 651c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 652c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,imedium = 0,sm = ctx->sm,mf = ctx->mf,fm = ctx->fm,ms = ctx->ms,lmbwidth = ctx->lmbwidth,rmbwidth = ctx->rmbwidth; 653c4762a1bSJed Brown PetscReal hxs,hxm,hxf; 654c4762a1bSJed Brown PetscScalar *x,*f,*slope; 655c4762a1bSJed Brown Vec Xloc; 656c4762a1bSJed Brown DM da; 657c4762a1bSJed Brown 658c4762a1bSJed Brown PetscFunctionBeginUser; 6599566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 6609566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); 6619566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 662c4762a1bSJed Brown hxs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 663c4762a1bSJed Brown hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 664c4762a1bSJed Brown hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 6659566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 6669566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 6679566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 6689566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 6699566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 6709566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da,PETSC_TRUE,&slope)); 6719566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 672c4762a1bSJed Brown 673c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 674c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 675c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 676c4762a1bSJed Brown } 677c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 678c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 679c4762a1bSJed Brown } 680c4762a1bSJed Brown } 681c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 682c4762a1bSJed Brown struct _LimitInfo info; 683c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 684c4762a1bSJed Brown if ((i > sm-2 && i < mf-lmbwidth+1) || (i > fm+rmbwidth-2 && i < ms+1)) { /* slow components and the first and last fast components */ 685c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 6869566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds)); 687c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 6889566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR,2*dof)); 689c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 690c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 691c4762a1bSJed Brown for (j=0; j<dof; j++) { 692c4762a1bSJed Brown PetscScalar jmpL,jmpR; 693c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 694c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 695c4762a1bSJed Brown for (k=0; k<dof; k++) { 696c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 697c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 698c4762a1bSJed Brown } 699c4762a1bSJed Brown } 700c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 701c4762a1bSJed Brown info.m = dof; 702c4762a1bSJed Brown info.hxs = hxs; 703c4762a1bSJed Brown info.hxm = hxm; 704c4762a1bSJed Brown info.hxf = hxf; 705c4762a1bSJed Brown (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope); 706c4762a1bSJed Brown for (j=0; j<dof; j++) { 707c4762a1bSJed Brown PetscScalar tmp = 0; 708c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 709c4762a1bSJed Brown slope[i*dof+j] = tmp; 710c4762a1bSJed Brown } 711c4762a1bSJed Brown } 712c4762a1bSJed Brown } 713c4762a1bSJed Brown 714c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 715c4762a1bSJed Brown PetscReal maxspeed; 716c4762a1bSJed Brown PetscScalar *uL,*uR; 717c4762a1bSJed Brown uL = &ctx->uLR[0]; 718c4762a1bSJed Brown uR = &ctx->uLR[dof]; 719c4762a1bSJed Brown if (i == sm) { /* interface between slow and medium regions */ 720c4762a1bSJed Brown for (j=0; j<dof; j++) { 721c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 722c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 723c4762a1bSJed Brown } 7249566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 725c4762a1bSJed Brown if (i < xs+xm) { 726c4762a1bSJed Brown for (j=0; j<dof; j++) f[imedium*dof+j] += ctx->flux[j]/hxm; 727c4762a1bSJed Brown imedium++; 728c4762a1bSJed Brown } 729c4762a1bSJed Brown } 730c4762a1bSJed Brown if (i > sm && i < mf-lmbwidth) { /* medium region */ 731c4762a1bSJed Brown for (j=0; j<dof; j++) { 732c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 733c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 734c4762a1bSJed Brown } 7359566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 736c4762a1bSJed Brown if (i > xs) { 737c4762a1bSJed Brown for (j=0; j<dof; j++) f[(imedium-1)*dof+j] -= ctx->flux[j]/hxm; 738c4762a1bSJed Brown } 739c4762a1bSJed Brown if (i < xs+xm) { 740c4762a1bSJed Brown for (j=0; j<dof; j++) f[imedium*dof+j] += ctx->flux[j]/hxm; 741c4762a1bSJed Brown imedium++; 742c4762a1bSJed Brown } 743c4762a1bSJed Brown } 744c4762a1bSJed Brown if (i == mf-lmbwidth) { /* interface between medium and fast regions */ 745c4762a1bSJed Brown for (j=0; j<dof; j++) { 746c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 747c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 748c4762a1bSJed Brown } 7499566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 750c4762a1bSJed Brown if (i > xs) { 751c4762a1bSJed Brown for (j=0; j<dof; j++) f[(imedium-1)*dof+j] -= ctx->flux[j]/hxm; 752c4762a1bSJed Brown } 753c4762a1bSJed Brown } 754c4762a1bSJed Brown if (i == fm+rmbwidth) { /* interface between fast and medium regions */ 755c4762a1bSJed Brown for (j=0; j<dof; j++) { 756c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 757c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 758c4762a1bSJed Brown } 7599566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 760c4762a1bSJed Brown if (i < xs+xm) { 761c4762a1bSJed Brown for (j=0; j<dof; j++) f[imedium*dof+j] += ctx->flux[j]/hxm; 762c4762a1bSJed Brown imedium++; 763c4762a1bSJed Brown } 764c4762a1bSJed Brown } 765c4762a1bSJed Brown if (i > fm+rmbwidth && i < ms) { /* medium region */ 766c4762a1bSJed Brown for (j=0; j<dof; j++) { 767c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 768c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 769c4762a1bSJed Brown } 7709566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 771c4762a1bSJed Brown if (i > xs) { 772c4762a1bSJed Brown for (j=0; j<dof; j++) f[(imedium-1)*dof+j] -= ctx->flux[j]/hxm; 773c4762a1bSJed Brown } 774c4762a1bSJed Brown if (i < xs+xm) { 775c4762a1bSJed Brown for (j=0; j<dof; j++) f[imedium*dof+j] += ctx->flux[j]/hxm; 776c4762a1bSJed Brown imedium++; 777c4762a1bSJed Brown } 778c4762a1bSJed Brown } 779c4762a1bSJed Brown if (i == ms) { /* interface between medium and slow regions */ 780c4762a1bSJed Brown for (j=0; j<dof; j++) { 781c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 782c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 783c4762a1bSJed Brown } 7849566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 785c4762a1bSJed Brown if (i > xs) { 786c4762a1bSJed Brown for (j=0; j<dof; j++) f[(imedium-1)*dof+j] -= ctx->flux[j]/hxm; 787c4762a1bSJed Brown } 788c4762a1bSJed Brown } 789c4762a1bSJed Brown } 7909566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 7919566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 7929566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da,PETSC_TRUE,&slope)); 7939566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 794c4762a1bSJed Brown PetscFunctionReturn(0); 795c4762a1bSJed Brown } 796c4762a1bSJed Brown 797c4762a1bSJed Brown PetscErrorCode FVRHSFunctionmediumbuffer_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 798c4762a1bSJed Brown { 799c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 800c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,imediumbuffer = 0,mf = ctx->mf,fm = ctx->fm,lmbwidth = ctx->lmbwidth,rmbwidth = ctx->rmbwidth; 801c4762a1bSJed Brown PetscReal hxs,hxm,hxf; 802c4762a1bSJed Brown PetscScalar *x,*f,*slope; 803c4762a1bSJed Brown Vec Xloc; 804c4762a1bSJed Brown DM da; 805c4762a1bSJed Brown 806c4762a1bSJed Brown PetscFunctionBeginUser; 8079566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 8089566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); 8099566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 810c4762a1bSJed Brown hxs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 811c4762a1bSJed Brown hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 812c4762a1bSJed Brown hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 8139566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 8149566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 8159566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 8169566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 8179566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 8189566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da,PETSC_TRUE,&slope)); 8199566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 820c4762a1bSJed Brown 821c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 822c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 823c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 824c4762a1bSJed Brown } 825c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 826c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 827c4762a1bSJed Brown } 828c4762a1bSJed Brown } 829c4762a1bSJed Brown for (i=xs-1; i<xs+xm+1; i++) { 830c4762a1bSJed Brown struct _LimitInfo info; 831c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 832c4762a1bSJed Brown if ((i > mf-lmbwidth-2 && i < mf+1) || (i > fm-2 && i < fm+rmbwidth+1)) { 833c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 8349566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds)); 835c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 8369566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR,2*dof)); 837c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 838c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 839c4762a1bSJed Brown for (j=0; j<dof; j++) { 840c4762a1bSJed Brown PetscScalar jmpL,jmpR; 841c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 842c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 843c4762a1bSJed Brown for (k=0; k<dof; k++) { 844c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 845c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 846c4762a1bSJed Brown } 847c4762a1bSJed Brown } 848c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 849c4762a1bSJed Brown info.m = dof; 850c4762a1bSJed Brown info.hxs = hxs; 851c4762a1bSJed Brown info.hxm = hxm; 852c4762a1bSJed Brown info.hxf = hxf; 853c4762a1bSJed Brown (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope); 854c4762a1bSJed Brown for (j=0; j<dof; j++) { 855c4762a1bSJed Brown PetscScalar tmp = 0; 856c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 857c4762a1bSJed Brown slope[i*dof+j] = tmp; 858c4762a1bSJed Brown } 859c4762a1bSJed Brown } 860c4762a1bSJed Brown } 861c4762a1bSJed Brown 862c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 863c4762a1bSJed Brown PetscReal maxspeed; 864c4762a1bSJed Brown PetscScalar *uL,*uR; 865c4762a1bSJed Brown uL = &ctx->uLR[0]; 866c4762a1bSJed Brown uR = &ctx->uLR[dof]; 867c4762a1bSJed Brown if (i == mf-lmbwidth) { 868c4762a1bSJed Brown for (j=0; j<dof; j++) { 869c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 870c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 871c4762a1bSJed Brown } 8729566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 873c4762a1bSJed Brown if (i < xs+xm) { 874c4762a1bSJed Brown for (j=0; j<dof; j++) f[imediumbuffer*dof+j] += ctx->flux[j]/hxm; 875c4762a1bSJed Brown imediumbuffer++; 876c4762a1bSJed Brown } 877c4762a1bSJed Brown } 878c4762a1bSJed Brown if (i > mf-lmbwidth && i < mf) { 879c4762a1bSJed Brown for (j=0; j<dof; j++) { 880c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 881c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 882c4762a1bSJed Brown } 8839566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 884c4762a1bSJed Brown if (i > xs) { 885c4762a1bSJed Brown for (j=0; j<dof; j++) f[(imediumbuffer-1)*dof+j] -= ctx->flux[j]/hxm; 886c4762a1bSJed Brown } 887c4762a1bSJed Brown if (i < xs+xm) { 888c4762a1bSJed Brown for (j=0; j<dof; j++) f[imediumbuffer*dof+j] += ctx->flux[j]/hxm; 889c4762a1bSJed Brown imediumbuffer++; 890c4762a1bSJed Brown } 891c4762a1bSJed Brown } 892c4762a1bSJed Brown if (i == mf) { /* interface between the medium region and the fast region */ 893c4762a1bSJed Brown for (j=0; j<dof; j++) { 894c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 895c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 896c4762a1bSJed Brown } 8979566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 898c4762a1bSJed Brown if (i > xs) { 899c4762a1bSJed Brown for (j=0; j<dof; j++) f[(imediumbuffer-1)*dof+j] -= ctx->flux[j]/hxm; 900c4762a1bSJed Brown } 901c4762a1bSJed Brown } 902c4762a1bSJed Brown if (i == fm) { /* interface between the fast region and the medium region */ 903c4762a1bSJed Brown for (j=0; j<dof; j++) { 904c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 905c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 906c4762a1bSJed Brown } 9079566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 908c4762a1bSJed Brown if (i < xs+xm) { 909c4762a1bSJed Brown for (j=0; j<dof; j++) f[imediumbuffer*dof+j] += ctx->flux[j]/hxm; 910c4762a1bSJed Brown imediumbuffer++; 911c4762a1bSJed Brown } 912c4762a1bSJed Brown } 913c4762a1bSJed Brown if (i > fm && i < fm+rmbwidth) { 914c4762a1bSJed Brown for (j=0; j<dof; j++) { 915c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 916c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 917c4762a1bSJed Brown } 9189566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 919c4762a1bSJed Brown if (i > xs) { 920c4762a1bSJed Brown for (j=0; j<dof; j++) f[(imediumbuffer-1)*dof+j] -= ctx->flux[j]/hxm; 921c4762a1bSJed Brown } 922c4762a1bSJed Brown if (i < xs+xm) { 923c4762a1bSJed Brown for (j=0; j<dof; j++) f[imediumbuffer*dof+j] += ctx->flux[j]/hxm; 924c4762a1bSJed Brown imediumbuffer++; 925c4762a1bSJed Brown } 926c4762a1bSJed Brown } 927c4762a1bSJed Brown if (i == fm+rmbwidth) { 928c4762a1bSJed Brown for (j=0; j<dof; j++) { 929c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 930c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 931c4762a1bSJed Brown } 9329566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 933c4762a1bSJed Brown if (i > xs) { 934c4762a1bSJed Brown for (j=0; j<dof; j++) f[(imediumbuffer-1)*dof+j] -= ctx->flux[j]/hxm; 935c4762a1bSJed Brown } 936c4762a1bSJed Brown } 937c4762a1bSJed Brown } 9389566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 9399566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 9409566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da,PETSC_TRUE,&slope)); 9419566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 942c4762a1bSJed Brown PetscFunctionReturn(0); 943c4762a1bSJed Brown } 944c4762a1bSJed Brown 945c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for fast parts ----------------------------------- */ 946c4762a1bSJed Brown PetscErrorCode FVRHSFunctionfast_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 947c4762a1bSJed Brown { 948c4762a1bSJed Brown FVCtx *ctx = (FVCtx*)vctx; 949c4762a1bSJed Brown PetscInt i,j,k,Mx,dof,xs,xm,ifast = 0,mf = ctx->mf,fm = ctx->fm; 950c4762a1bSJed Brown PetscReal hxs,hxm,hxf; 951c4762a1bSJed Brown PetscScalar *x,*f,*slope; 952c4762a1bSJed Brown Vec Xloc; 953c4762a1bSJed Brown DM da; 954c4762a1bSJed Brown 955c4762a1bSJed Brown PetscFunctionBeginUser; 9569566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 9579566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xloc)); 9589566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 959c4762a1bSJed Brown hxs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 960c4762a1bSJed Brown hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 961c4762a1bSJed Brown hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 9629566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 9639566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 9649566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 9659566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xloc,&x)); 9669566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 9679566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da,PETSC_TRUE,&slope)); 9689566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 969c4762a1bSJed Brown 970c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 971c4762a1bSJed Brown for (i=xs-2; i<0; i++) { 972c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 973c4762a1bSJed Brown } 974c4762a1bSJed Brown for (i=Mx; i<xs+xm+2; i++) { 975c4762a1bSJed Brown for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 976c4762a1bSJed Brown } 977c4762a1bSJed Brown } 9786aad120cSJose E. Roman for (i=xs-1; i<xs+xm+1; i++) { /* fast components and the last slow components before fast components and the first slow component after fast components */ 979c4762a1bSJed Brown struct _LimitInfo info; 980c4762a1bSJed Brown PetscScalar *cjmpL,*cjmpR; 981c4762a1bSJed Brown if (i > mf-2 && i < fm+1) { 9829566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds)); 9839566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR,2*dof)); 984c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 985c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 986c4762a1bSJed Brown for (j=0; j<dof; j++) { 987c4762a1bSJed Brown PetscScalar jmpL,jmpR; 988c4762a1bSJed Brown jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 989c4762a1bSJed Brown jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 990c4762a1bSJed Brown for (k=0; k<dof; k++) { 991c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 992c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 993c4762a1bSJed Brown } 994c4762a1bSJed Brown } 995c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 996c4762a1bSJed Brown info.m = dof; 997c4762a1bSJed Brown info.hxs = hxs; 998c4762a1bSJed Brown info.hxm = hxm; 999c4762a1bSJed Brown info.hxf = hxf; 1000c4762a1bSJed Brown (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope); 1001c4762a1bSJed Brown for (j=0; j<dof; j++) { 1002c4762a1bSJed Brown PetscScalar tmp = 0; 1003c4762a1bSJed Brown for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 1004c4762a1bSJed Brown slope[i*dof+j] = tmp; 1005c4762a1bSJed Brown } 1006c4762a1bSJed Brown } 1007c4762a1bSJed Brown } 1008c4762a1bSJed Brown 1009c4762a1bSJed Brown for (i=xs; i<xs+xm+1; i++) { 1010c4762a1bSJed Brown PetscReal maxspeed; 1011c4762a1bSJed Brown PetscScalar *uL,*uR; 1012c4762a1bSJed Brown uL = &ctx->uLR[0]; 1013c4762a1bSJed Brown uR = &ctx->uLR[dof]; 1014c4762a1bSJed Brown if (i == mf) { /* interface between medium and fast regions */ 1015c4762a1bSJed Brown for (j=0; j<dof; j++) { 1016c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 1017c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 1018c4762a1bSJed Brown } 10199566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 1020c4762a1bSJed Brown if (i < xs+xm) { 1021c4762a1bSJed Brown for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf; 1022c4762a1bSJed Brown ifast++; 1023c4762a1bSJed Brown } 1024c4762a1bSJed Brown } 1025c4762a1bSJed Brown if (i > mf && i < fm) { /* fast region */ 1026c4762a1bSJed Brown for (j=0; j<dof; j++) { 1027c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 1028c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 1029c4762a1bSJed Brown } 10309566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 1031c4762a1bSJed Brown if (i > xs) { 1032c4762a1bSJed Brown for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf; 1033c4762a1bSJed Brown } 1034c4762a1bSJed Brown if (i < xs+xm) { 1035c4762a1bSJed Brown for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf; 1036c4762a1bSJed Brown ifast++; 1037c4762a1bSJed Brown } 1038c4762a1bSJed Brown } 1039c4762a1bSJed Brown if (i == fm) { /* interface between fast and medium regions */ 1040c4762a1bSJed Brown for (j=0; j<dof; j++) { 1041c4762a1bSJed Brown uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 1042c4762a1bSJed Brown uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 1043c4762a1bSJed Brown } 10449566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed)); 1045c4762a1bSJed Brown if (i > xs) { 1046c4762a1bSJed Brown for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf; 1047c4762a1bSJed Brown } 1048c4762a1bSJed Brown } 1049c4762a1bSJed Brown } 10509566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xloc,&x)); 10519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 10529566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da,PETSC_TRUE,&slope)); 10539566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xloc)); 1054c4762a1bSJed Brown PetscFunctionReturn(0); 1055c4762a1bSJed Brown } 1056c4762a1bSJed Brown 1057c4762a1bSJed Brown int main(int argc,char *argv[]) 1058c4762a1bSJed Brown { 1059c4762a1bSJed Brown char lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m"; 1060c4762a1bSJed Brown PetscFunctionList limiters = 0,physics = 0; 1061c4762a1bSJed Brown MPI_Comm comm; 1062c4762a1bSJed Brown TS ts; 1063c4762a1bSJed Brown DM da; 1064c4762a1bSJed Brown Vec X,X0,R; 1065c4762a1bSJed Brown FVCtx ctx; 1066c4762a1bSJed Brown PetscInt i,k,dof,xs,xm,Mx,draw = 0,count_slow,count_medium,count_fast,islow = 0,imedium = 0,ifast = 0,*index_slow,*index_medium,*index_fast,islowbuffer = 0,*index_slowbuffer,imediumbuffer = 0,*index_mediumbuffer; 1067c4762a1bSJed Brown PetscBool view_final = PETSC_FALSE; 1068c4762a1bSJed Brown PetscReal ptime; 1069c4762a1bSJed Brown 1070*327415f7SBarry Smith PetscFunctionBeginUser; 10719566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,0,help)); 1072c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 10739566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&ctx,sizeof(ctx))); 1074c4762a1bSJed Brown 1075c4762a1bSJed Brown /* Register limiters to be available on the command line */ 10769566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"upwind" ,Limit3_Upwind)); 10779566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"lax-wendroff" ,Limit3_LaxWendroff)); 10789566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"beam-warming" ,Limit3_BeamWarming)); 10799566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"fromm" ,Limit3_Fromm)); 10809566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"minmod" ,Limit3_Minmod)); 10819566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"superbee" ,Limit3_Superbee)); 10829566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"mc" ,Limit3_MC)); 10839566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters,"koren3" ,Limit3_Koren3)); 1084c4762a1bSJed Brown 1085c4762a1bSJed Brown /* Register physical models to be available on the command line */ 10869566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&physics,"advect" ,PhysicsCreate_Advect)); 1087c4762a1bSJed Brown 1088c4762a1bSJed Brown ctx.comm = comm; 1089c4762a1bSJed Brown ctx.cfl = 0.9; 1090c4762a1bSJed Brown ctx.bctype = FVBC_PERIODIC; 1091c4762a1bSJed Brown ctx.xmin = -1.0; 1092c4762a1bSJed Brown ctx.xmax = 1.0; 1093d0609cedSBarry Smith PetscOptionsBegin(comm,NULL,"Finite Volume solver options",""); 10949566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL)); 10959566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL)); 10969566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-limit","Name of flux imiter to use","",limiters,lname,lname,sizeof(lname),NULL)); 10979566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL)); 10989566063dSJacob 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)); 10999566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL)); 11009566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-exact","Compare errors with exact solution","",ctx.exact,&ctx.exact,NULL)); 11019566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL)); 11029566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL)); 11039566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL)); 11049566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-hratio","Spacing ratio","",ctx.hratio,&ctx.hratio,NULL)); 1105d0609cedSBarry Smith PetscOptionsEnd(); 1106c4762a1bSJed Brown 1107c4762a1bSJed Brown /* Choose the limiter from the list of registered limiters */ 11089566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(limiters,lname,&ctx.limit3)); 11093c633725SBarry Smith PetscCheck(ctx.limit3,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Limiter '%s' not found",lname); 1110c4762a1bSJed Brown 1111c4762a1bSJed Brown /* Choose the physics from the list of registered models */ 1112c4762a1bSJed Brown { 1113c4762a1bSJed Brown PetscErrorCode (*r)(FVCtx*); 11149566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(physics,physname,&r)); 11153c633725SBarry Smith PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Physics '%s' not found",physname); 1116c4762a1bSJed Brown /* Create the physics, will set the number of fields and their names */ 11179566063dSJacob Faibussowitsch PetscCall((*r)(&ctx)); 1118c4762a1bSJed Brown } 1119c4762a1bSJed Brown 1120c4762a1bSJed Brown /* Create a DMDA to manage the parallel grid */ 11219566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics2.dof,2,NULL,&da)); 11229566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 11239566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1124c4762a1bSJed Brown /* Inform the DMDA of the field names provided by the physics. */ 1125c4762a1bSJed Brown /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ 1126c4762a1bSJed Brown for (i=0; i<ctx.physics2.dof; i++) { 11279566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da,i,ctx.physics2.fieldname[i])); 1128c4762a1bSJed Brown } 11299566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0)); 11309566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0)); 1131c4762a1bSJed Brown 1132c4762a1bSJed Brown /* Set coordinates of cell centers */ 11339566063dSJacob 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)); 1134c4762a1bSJed Brown 1135c4762a1bSJed Brown /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */ 11369566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dof*dof,&ctx.R,dof*dof,&ctx.Rinv,2*dof,&ctx.cjmpLR,1*dof,&ctx.cslope)); 11379566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(2*dof,&ctx.uLR,dof,&ctx.flux,dof,&ctx.speeds)); 1138c4762a1bSJed Brown 1139c4762a1bSJed Brown /* Create a vector to store the solution and to save the initial state */ 11409566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da,&X)); 11419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&X0)); 11429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&R)); 1143c4762a1bSJed Brown 1144c4762a1bSJed Brown /* create index for slow parts and fast parts, 1145c4762a1bSJed Brown count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */ 1146c4762a1bSJed Brown count_slow = Mx/(1+ctx.hratio)/(1+ctx.hratio); 1147c4762a1bSJed Brown count_medium = 2*ctx.hratio*count_slow; 1148cad9d221SBarry Smith PetscCheck((count_slow%2) == 0 && (count_medium%2) == 0,PETSC_COMM_WORLD,PETSC_ERR_USER,"Please adjust grid size Mx (-da_grid_x) and hratio (-hratio) so that Mx/(1+hartio)^2 and Mx*2*hratio/(1+hratio)^2 is even"); 1149c4762a1bSJed Brown count_fast = ctx.hratio*ctx.hratio*count_slow; 1150c4762a1bSJed Brown ctx.sm = count_slow/2; 1151c4762a1bSJed Brown ctx.mf = ctx.sm + count_medium/2; 1152c4762a1bSJed Brown ctx.fm = ctx.mf + count_fast; 1153c4762a1bSJed Brown ctx.ms = ctx.fm + count_medium/2; 11549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm*dof,&index_slow)); 11559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm*dof,&index_medium)); 11569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm*dof,&index_fast)); 11579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(6*dof,&index_slowbuffer)); 11589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(6*dof,&index_mediumbuffer)); 1159c4762a1bSJed Brown if (((AdvectCtx*)ctx.physics2.user)->a > 0) { 1160c4762a1bSJed Brown ctx.lsbwidth = 2; 1161c4762a1bSJed Brown ctx.rsbwidth = 4; 1162c4762a1bSJed Brown ctx.lmbwidth = 2; 1163c4762a1bSJed Brown ctx.rmbwidth = 4; 1164c4762a1bSJed Brown } else { 1165c4762a1bSJed Brown ctx.lsbwidth = 4; 1166c4762a1bSJed Brown ctx.rsbwidth = 2; 1167c4762a1bSJed Brown ctx.lmbwidth = 4; 1168c4762a1bSJed Brown ctx.rmbwidth = 2; 1169c4762a1bSJed Brown } 1170c4762a1bSJed Brown 1171c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 1172c4762a1bSJed Brown if (i < ctx.sm-ctx.lsbwidth || i > ctx.ms+ctx.rsbwidth-1) 1173c4762a1bSJed Brown for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k; 1174c4762a1bSJed Brown else if ((i >= ctx.sm-ctx.lsbwidth && i < ctx.sm) || (i > ctx.ms-1 && i <= ctx.ms+ctx.rsbwidth-1)) 1175c4762a1bSJed Brown for (k=0; k<dof; k++) index_slowbuffer[islowbuffer++] = i*dof+k; 1176c4762a1bSJed Brown else if (i < ctx.mf-ctx.lmbwidth || i > ctx.fm+ctx.rmbwidth-1) 1177c4762a1bSJed Brown for (k=0; k<dof; k++) index_medium[imedium++] = i*dof+k; 1178c4762a1bSJed Brown else if ((i >= ctx.mf-ctx.lmbwidth && i < ctx.mf) || (i > ctx.fm-1 && i <= ctx.fm+ctx.rmbwidth-1)) 1179c4762a1bSJed Brown for (k=0; k<dof; k++) index_mediumbuffer[imediumbuffer++] = i*dof+k; 1180c4762a1bSJed Brown else 1181c4762a1bSJed Brown for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k; 1182c4762a1bSJed Brown } 11839566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss)); 11849566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,imedium,index_medium,PETSC_COPY_VALUES,&ctx.ism)); 11859566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf)); 11869566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,islowbuffer,index_slowbuffer,PETSC_COPY_VALUES,&ctx.issb)); 11879566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,imediumbuffer,index_mediumbuffer,PETSC_COPY_VALUES,&ctx.ismb)); 1188c4762a1bSJed Brown 1189c4762a1bSJed Brown /* Create a time-stepping object */ 11909566063dSJacob Faibussowitsch PetscCall(TSCreate(comm,&ts)); 11919566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts,da)); 11929566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,R,FVRHSFunction_3WaySplit,&ctx)); 11939566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"slow",ctx.iss)); 11949566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"medium",ctx.ism)); 11959566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"fast",ctx.isf)); 11969566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"slowbuffer",ctx.issb)); 11979566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"mediumbuffer",ctx.ismb)); 11989566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"slow",NULL,FVRHSFunctionslow_3WaySplit,&ctx)); 11999566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"medium",NULL,FVRHSFunctionmedium_3WaySplit,&ctx)); 12009566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"fast",NULL,FVRHSFunctionfast_3WaySplit,&ctx)); 12019566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"slowbuffer",NULL,FVRHSFunctionslowbuffer_3WaySplit,&ctx)); 12029566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"mediumbuffer",NULL,FVRHSFunctionmediumbuffer_3WaySplit,&ctx)); 1203c4762a1bSJed Brown 12049566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSSSP)); 12059566063dSJacob Faibussowitsch /*PetscCall(TSSetType(ts,TSMPRK));*/ 12069566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,10)); 12079566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 1208c4762a1bSJed Brown 1209c4762a1bSJed Brown /* Compute initial conditions and starting time step */ 12109566063dSJacob Faibussowitsch PetscCall(FVSample_3WaySplit(&ctx,da,0,X0)); 12119566063dSJacob Faibussowitsch PetscCall(FVRHSFunction_3WaySplit(ts,0,X0,X,(void*)&ctx)); /* Initial function evaluation, only used to determine max speed */ 12129566063dSJacob Faibussowitsch PetscCall(VecCopy(X0,X)); /* The function value was not used so we set X=X0 again */ 12139566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt)); 12149566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); /* Take runtime options */ 12159566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD)); 1216c4762a1bSJed Brown { 1217c4762a1bSJed Brown PetscInt steps; 1218c4762a1bSJed Brown PetscScalar mass_initial,mass_final,mass_difference,mass_differenceg; 1219c4762a1bSJed Brown const PetscScalar *ptr_X,*ptr_X0; 1220c4762a1bSJed Brown const PetscReal hs = (ctx.xmax-ctx.xmin)/4.0/count_slow; 1221c4762a1bSJed Brown const PetscReal hm = (ctx.xmax-ctx.xmin)/2.0/count_medium; 1222c4762a1bSJed Brown const PetscReal hf = (ctx.xmax-ctx.xmin)/4.0/count_fast; 1223c4762a1bSJed Brown 12249566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,X)); 12259566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&ptime)); 12269566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&steps)); 1227c4762a1bSJed Brown /* calculate the total mass at initial time and final time */ 1228c4762a1bSJed Brown mass_initial = 0.0; 1229c4762a1bSJed Brown mass_final = 0.0; 12309566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,X0,(void*)&ptr_X0)); 12319566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,X,(void*)&ptr_X)); 1232c4762a1bSJed Brown for (i=xs;i<xs+xm;i++) { 1233c4762a1bSJed Brown if (i < ctx.sm || i > ctx.ms-1) 1234c4762a1bSJed Brown for (k=0; k<dof; k++) { 1235c4762a1bSJed Brown mass_initial = mass_initial + hs*ptr_X0[i*dof+k]; 1236c4762a1bSJed Brown mass_final = mass_final + hs*ptr_X[i*dof+k]; 1237c4762a1bSJed Brown } 1238c4762a1bSJed Brown else if (i < ctx.mf || i > ctx.fm-1) 1239c4762a1bSJed Brown for (k=0; k<dof; k++) { 1240c4762a1bSJed Brown mass_initial = mass_initial + hm*ptr_X0[i*dof+k]; 1241c4762a1bSJed Brown mass_final = mass_final + hm*ptr_X[i*dof+k]; 1242c4762a1bSJed Brown } 1243c4762a1bSJed Brown else { 1244c4762a1bSJed Brown for (k=0; k<dof; k++) { 1245c4762a1bSJed Brown mass_initial = mass_initial + hf*ptr_X0[i*dof+k]; 1246c4762a1bSJed Brown mass_final = mass_final + hf*ptr_X[i*dof+k]; 1247c4762a1bSJed Brown } 1248c4762a1bSJed Brown } 1249c4762a1bSJed Brown } 12509566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,X0,(void*)&ptr_X0)); 12519566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,X,(void*)&ptr_X)); 1252c4762a1bSJed Brown mass_difference = mass_final - mass_initial; 12539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&mass_difference,&mass_differenceg,1,MPIU_SCALAR,MPIU_SUM,comm)); 12549566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Mass difference %g\n",(double)mass_differenceg)); 125563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Final time %g, steps %" PetscInt_FMT "\n",(double)ptime,steps)); 12569566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Maximum allowable stepsize according to CFL %g\n",(double)(1.0/ctx.cfl_idt))); 1257c4762a1bSJed Brown if (ctx.exact) { 1258c4762a1bSJed Brown PetscReal nrm1=0; 12599566063dSJacob Faibussowitsch PetscCall(SolutionErrorNorms_3WaySplit(&ctx,da,ptime,X,&nrm1)); 12609566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1)); 1261c4762a1bSJed Brown } 1262c4762a1bSJed Brown if (ctx.simulation) { 1263c4762a1bSJed Brown PetscReal nrm1=0; 1264c4762a1bSJed Brown PetscViewer fd; 1265c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = "binaryoutput"; 1266c4762a1bSJed Brown Vec XR; 1267c4762a1bSJed Brown PetscBool flg; 1268c4762a1bSJed Brown const PetscScalar *ptr_XR; 12699566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg)); 12703c633725SBarry Smith PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 12719566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&fd)); 12729566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X0,&XR)); 12739566063dSJacob Faibussowitsch PetscCall(VecLoad(XR,fd)); 12749566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 12759566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&ptr_X)); 12769566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XR,&ptr_XR)); 1277c4762a1bSJed Brown for (i=xs;i<xs+xm;i++) { 1278c4762a1bSJed Brown if (i < ctx.sm || i > ctx.ms-1) 1279c4762a1bSJed Brown for (k=0; k<dof; k++) nrm1 = nrm1 + hs*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]); 1280c4762a1bSJed Brown else if (i < ctx.mf || i < ctx.fm-1) 1281c4762a1bSJed Brown for (k=0; k<dof; k++) nrm1 = nrm1 + hm*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]); 1282c4762a1bSJed Brown else 1283c4762a1bSJed Brown for (k=0; k<dof; k++) nrm1 = nrm1 + hf*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]); 1284c4762a1bSJed Brown } 12859566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&ptr_X)); 12869566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XR,&ptr_XR)); 12879566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1)); 12889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XR)); 1289c4762a1bSJed Brown } 1290c4762a1bSJed Brown } 1291c4762a1bSJed Brown 12929566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD)); 12939566063dSJacob Faibussowitsch if (draw & 0x1) PetscCall(VecView(X0,PETSC_VIEWER_DRAW_WORLD)); 12949566063dSJacob Faibussowitsch if (draw & 0x2) PetscCall(VecView(X,PETSC_VIEWER_DRAW_WORLD)); 1295c4762a1bSJed Brown if (draw & 0x4) { 1296c4762a1bSJed Brown Vec Y; 12979566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&Y)); 12989566063dSJacob Faibussowitsch PetscCall(FVSample_3WaySplit(&ctx,da,ptime,Y)); 12999566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y,-1,X)); 13009566063dSJacob Faibussowitsch PetscCall(VecView(Y,PETSC_VIEWER_DRAW_WORLD)); 13019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 1302c4762a1bSJed Brown } 1303c4762a1bSJed Brown 1304c4762a1bSJed Brown if (view_final) { 1305c4762a1bSJed Brown PetscViewer viewer; 13069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer)); 13079566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB)); 13089566063dSJacob Faibussowitsch PetscCall(VecView(X,viewer)); 13099566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 13109566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1311c4762a1bSJed Brown } 1312c4762a1bSJed Brown 1313c4762a1bSJed Brown /* Clean up */ 13149566063dSJacob Faibussowitsch PetscCall((*ctx.physics2.destroy)(ctx.physics2.user)); 13159566063dSJacob Faibussowitsch for (i=0; i<ctx.physics2.dof; i++) PetscCall(PetscFree(ctx.physics2.fieldname[i])); 13169566063dSJacob Faibussowitsch PetscCall(PetscFree4(ctx.R,ctx.Rinv,ctx.cjmpLR,ctx.cslope)); 13179566063dSJacob Faibussowitsch PetscCall(PetscFree3(ctx.uLR,ctx.flux,ctx.speeds)); 13189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 13199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X0)); 13209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&R)); 13219566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 13229566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 13239566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.iss)); 13249566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.ism)); 13259566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.isf)); 13269566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.issb)); 13279566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.ismb)); 13289566063dSJacob Faibussowitsch PetscCall(PetscFree(index_slow)); 13299566063dSJacob Faibussowitsch PetscCall(PetscFree(index_medium)); 13309566063dSJacob Faibussowitsch PetscCall(PetscFree(index_fast)); 13319566063dSJacob Faibussowitsch PetscCall(PetscFree(index_slowbuffer)); 13329566063dSJacob Faibussowitsch PetscCall(PetscFree(index_mediumbuffer)); 13339566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&limiters)); 13349566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&physics)); 13359566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 1336b122ec5aSJacob Faibussowitsch return 0; 1337c4762a1bSJed Brown } 1338c4762a1bSJed Brown 1339c4762a1bSJed Brown /*TEST 1340c4762a1bSJed Brown 1341c4762a1bSJed Brown build: 1342f56ea12dSJed Brown requires: !complex 1343c4762a1bSJed Brown depends: finitevolume1d.c 1344c4762a1bSJed Brown 1345c4762a1bSJed Brown test: 1346c4762a1bSJed Brown suffix: 1 1347c4762a1bSJed Brown args: -da_grid_x 90 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a23 -ts_use_splitrhsfunction 0 1348c4762a1bSJed Brown 1349c4762a1bSJed Brown test: 1350c4762a1bSJed Brown suffix: 2 1351c4762a1bSJed Brown args: -da_grid_x 90 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a23 -ts_use_splitrhsfunction 1 1352c4762a1bSJed Brown 1353c4762a1bSJed Brown TEST*/ 1354