1 static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n" 2 " advection - Constant coefficient scalar advection\n" 3 " u_t + (a*u)_x = 0\n" 4 " for this toy problem, we choose different meshsizes for different sub-domains (slow-medium-fast-medium-slow), \n" 5 " the meshsize ratio between two adjacient sub-domains is controlled with -hratio,\n" 6 " exact - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n" 7 " the states across shocks and rarefactions\n" 8 " simulation - use reference solution which is generated by smaller time step size to be true solution,\n" 9 " also the reference solution should be generated by user and stored in a binary file.\n" 10 " characteristic - Limit the characteristic variables, this is usually preferred (default)\n" 11 "Several initial conditions can be chosen with -initial N\n\n" 12 "The problem size should be set with -da_grid_x M\n\n"; 13 14 #include <petscts.h> 15 #include <petscdm.h> 16 #include <petscdmda.h> 17 #include <petscdraw.h> 18 #include "finitevolume1d.h" 19 20 PETSC_STATIC_INLINE PetscReal RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax) { PetscReal range = xmax-xmin; return xmin +PetscFmodReal(range+PetscFmodReal(a,range),range); } 21 22 /* --------------------------------- Advection ----------------------------------- */ 23 typedef struct { 24 PetscReal a; /* advective velocity */ 25 } AdvectCtx; 26 27 static PetscErrorCode PhysicsRiemann_Advect(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed) 28 { 29 AdvectCtx *ctx = (AdvectCtx*)vctx; 30 PetscReal speed; 31 32 PetscFunctionBeginUser; 33 speed = ctx->a; 34 flux[0] = PetscMax(0,speed)*uL[0] + PetscMin(0,speed)*uR[0]; 35 *maxspeed = speed; 36 PetscFunctionReturn(0); 37 } 38 39 static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds) 40 { 41 AdvectCtx *ctx = (AdvectCtx*)vctx; 42 43 PetscFunctionBeginUser; 44 X[0] = 1.; 45 Xi[0] = 1.; 46 speeds[0] = ctx->a; 47 PetscFunctionReturn(0); 48 } 49 50 static PetscErrorCode PhysicsSample_Advect(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u) 51 { 52 AdvectCtx *ctx = (AdvectCtx*)vctx; 53 PetscReal a = ctx->a,x0; 54 55 PetscFunctionBeginUser; 56 switch (bctype) { 57 case FVBC_OUTFLOW: x0 = x-a*t; break; 58 case FVBC_PERIODIC: x0 = RangeMod(x-a*t,xmin,xmax); break; 59 default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown BCType"); 60 } 61 switch (initial) { 62 case 0: u[0] = (x0 < 0) ? 1 : -1; break; 63 case 1: u[0] = (x0 < 0) ? -1 : 1; break; 64 case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break; 65 case 3: u[0] = PetscSinReal(2*PETSC_PI*x0); break; 66 case 4: u[0] = PetscAbs(x0); break; 67 case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2*PETSC_PI*x0)); break; 68 case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2-x0 : 0)); break; 69 case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI*x0),10.0);break; 70 default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition"); 71 } 72 PetscFunctionReturn(0); 73 } 74 75 static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx) 76 { 77 PetscErrorCode ierr; 78 AdvectCtx *user; 79 80 PetscFunctionBeginUser; 81 ierr = PetscNew(&user);CHKERRQ(ierr); 82 ctx->physics2.sample2 = PhysicsSample_Advect; 83 ctx->physics2.riemann2 = PhysicsRiemann_Advect; 84 ctx->physics2.characteristic2 = PhysicsCharacteristic_Advect; 85 ctx->physics2.destroy = PhysicsDestroy_SimpleFree; 86 ctx->physics2.user = user; 87 ctx->physics2.dof = 1; 88 ierr = PetscStrallocpy("u",&ctx->physics2.fieldname[0]);CHKERRQ(ierr); 89 user->a = 1; 90 ierr = PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");CHKERRQ(ierr); 91 { 92 ierr = PetscOptionsReal("-physics_advect_a","Speed","",user->a,&user->a,NULL);CHKERRQ(ierr); 93 } 94 ierr = PetscOptionsEnd();CHKERRQ(ierr); 95 PetscFunctionReturn(0); 96 } 97 98 PetscErrorCode FVSample_3WaySplit(FVCtx *ctx,DM da,PetscReal time,Vec U) 99 { 100 PetscErrorCode ierr; 101 PetscScalar *u,*uj,xj,xi; 102 PetscInt i,j,k,dof,xs,xm,Mx; 103 const PetscInt N = 200; 104 PetscReal hs,hm,hf; 105 106 PetscFunctionBeginUser; 107 if (!ctx->physics2.sample2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Physics has not provided a sampling function"); 108 ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 109 ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 110 ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr); 111 ierr = PetscMalloc1(dof,&uj);CHKERRQ(ierr); 112 113 hs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 114 hm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 115 hf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 116 for (i=xs; i<xs+xm; i++) { 117 if (i < ctx->sm) { 118 xi = ctx->xmin+0.5*hs+i*hs; 119 /* Integrate over cell i using trapezoid rule with N points. */ 120 for (k=0; k<dof; k++) u[i*dof+k] = 0; 121 for (j=0; j<N+1; j++) { 122 xj = xi+hs*(j-N/2)/(PetscReal)N; 123 ierr = (*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); 124 for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 125 } 126 } else if (i < ctx->mf) { 127 xi = ctx->xmin+ctx->sm*hs+0.5*hm+(i-ctx->sm)*hm; 128 /* Integrate over cell i using trapezoid rule with N points. */ 129 for (k=0; k<dof; k++) u[i*dof+k] = 0; 130 for (j=0; j<N+1; j++) { 131 xj = xi+hm*(j-N/2)/(PetscReal)N; 132 ierr = (*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); 133 for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 134 } 135 } else if (i < ctx->fm) { 136 xi = ctx->xmin+ctx->sm*hs+(ctx->mf-ctx->sm)*hm+0.5*hf+(i-ctx->mf)*hf; 137 /* Integrate over cell i using trapezoid rule with N points. */ 138 for (k=0; k<dof; k++) u[i*dof+k] = 0; 139 for (j=0; j<N+1; j++) { 140 xj = xi+hf*(j-N/2)/(PetscReal)N; 141 ierr = (*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); 142 for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 143 } 144 } else if (i < ctx->ms) { 145 xi = ctx->xmin+ctx->sm*hs+(ctx->mf-ctx->sm)*hm+(ctx->fm-ctx->mf)*hf+0.5*hm+(i-ctx->fm)*hm; 146 /* Integrate over cell i using trapezoid rule with N points. */ 147 for (k=0; k<dof; k++) u[i*dof+k] = 0; 148 for (j=0; j<N+1; j++) { 149 xj = xi+hm*(j-N/2)/(PetscReal)N; 150 ierr = (*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); 151 for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 152 } 153 } else { 154 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; 155 /* Integrate over cell i using trapezoid rule with N points. */ 156 for (k=0; k<dof; k++) u[i*dof+k] = 0; 157 for (j=0; j<N+1; j++) { 158 xj = xi+hs*(j-N/2)/(PetscReal)N; 159 ierr = (*ctx->physics2.sample2)(ctx->physics2.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);CHKERRQ(ierr); 160 for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N; 161 } 162 } 163 } 164 ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr); 165 ierr = PetscFree(uj);CHKERRQ(ierr); 166 PetscFunctionReturn(0); 167 } 168 169 static PetscErrorCode SolutionErrorNorms_3WaySplit(FVCtx *ctx,DM da,PetscReal t,Vec X,PetscReal *nrm1) 170 { 171 PetscErrorCode ierr; 172 Vec Y; 173 PetscInt i,Mx; 174 const PetscScalar *ptr_X,*ptr_Y; 175 PetscReal hs,hm,hf; 176 177 PetscFunctionBeginUser; 178 ierr = VecGetSize(X,&Mx);CHKERRQ(ierr); 179 hs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 180 hm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 181 hf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 182 ierr = VecDuplicate(X,&Y);CHKERRQ(ierr); 183 ierr = FVSample_3WaySplit(ctx,da,t,Y);CHKERRQ(ierr); 184 ierr = VecGetArrayRead(X,&ptr_X);CHKERRQ(ierr); 185 ierr = VecGetArrayRead(Y,&ptr_Y);CHKERRQ(ierr); 186 for (i=0;i<Mx;i++) { 187 if (i < ctx->sm || i > ctx->ms - 1) *nrm1 += hs*PetscAbs(ptr_X[i]-ptr_Y[i]); 188 else if (i < ctx->mf || i > ctx->fm - 1) *nrm1 += hm*PetscAbs(ptr_X[i]-ptr_Y[i]); 189 else *nrm1 += hf*PetscAbs(ptr_X[i]-ptr_Y[i]); 190 } 191 ierr = VecRestoreArrayRead(X,&ptr_X);CHKERRQ(ierr); 192 ierr = VecRestoreArrayRead(Y,&ptr_Y);CHKERRQ(ierr); 193 ierr = VecDestroy(&Y);CHKERRQ(ierr); 194 PetscFunctionReturn(0); 195 } 196 197 PetscErrorCode FVRHSFunction_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 198 { 199 FVCtx *ctx = (FVCtx*)vctx; 200 PetscErrorCode ierr; 201 PetscInt i,j,k,Mx,dof,xs,xm,sm = ctx->sm,mf = ctx->mf,fm = ctx->fm,ms = ctx->ms; 202 PetscReal hxf,hxm,hxs,cfl_idt = 0; 203 PetscScalar *x,*f,*slope; 204 Vec Xloc; 205 DM da; 206 207 PetscFunctionBeginUser; 208 ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 209 ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); /* Xloc contains ghost points */ 210 ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); /* Mx is the number of center points */ 211 hxs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 212 hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 213 hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 214 ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); /* X is solution vector which does not contain ghost points */ 215 ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 216 217 ierr = VecZeroEntries(F);CHKERRQ(ierr); /* F is the right hand side function corresponds to center points */ 218 219 ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 220 ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr); 221 ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); /* contains ghost points */ 222 223 ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 224 225 if (ctx->bctype == FVBC_OUTFLOW) { 226 for (i=xs-2; i<0; i++) { 227 for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 228 } 229 for (i=Mx; i<xs+xm+2; i++) { 230 for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 231 } 232 } 233 for (i=xs-1; i<xs+xm+1; i++) { 234 struct _LimitInfo info; 235 PetscScalar *cjmpL,*cjmpR; 236 /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 237 ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); 238 /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 239 ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 240 cjmpL = &ctx->cjmpLR[0]; 241 cjmpR = &ctx->cjmpLR[dof]; 242 for (j=0; j<dof; j++) { 243 PetscScalar jmpL,jmpR; 244 jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 245 jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 246 for (k=0; k<dof; k++) { 247 cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 248 cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 249 } 250 } 251 /* Apply limiter to the left and right characteristic jumps */ 252 info.m = dof; 253 info.hxs = hxs; 254 info.hxm = hxm; 255 info.hxf = hxf; 256 (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope); 257 for (j=0; j<dof; j++) { 258 PetscScalar tmp = 0; 259 for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 260 slope[i*dof+j] = tmp; 261 } 262 } 263 264 for (i=xs; i<xs+xm+1; i++) { 265 PetscReal maxspeed; 266 PetscScalar *uL,*uR; 267 uL = &ctx->uLR[0]; 268 uR = &ctx->uLR[dof]; 269 if (i < sm || i > ms) { /* slow region */ 270 for (j=0; j<dof; j++) { 271 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 272 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 273 } 274 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 275 if (i > xs) { 276 for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs; 277 } 278 if (i < xs+xm) { 279 for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs; 280 } 281 } else if (i == sm) { /* interface between slow and medium component */ 282 for (j=0; j<dof; j++) { 283 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 284 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 285 } 286 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 287 if (i > xs) { 288 for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxs; 289 } 290 if (i < xs+xm) { 291 for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxm; 292 } 293 } else if (i == ms) { /* interface between medium and slow regions */ 294 for (j=0; j<dof; j++) { 295 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 296 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 297 } 298 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 299 if (i > xs) { 300 for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxm; 301 } 302 if (i < xs+xm) { 303 for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxs; 304 } 305 } else if (i < mf || i > fm) { /* medium region */ 306 for (j=0; j<dof; j++) { 307 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 308 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 309 } 310 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 311 if (i > xs) { 312 for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxm; 313 } 314 if (i < xs+xm) { 315 for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxm; 316 } 317 } else if (i == mf) { /* interface between medium and fast regions */ 318 for (j=0; j<dof; j++) { 319 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 320 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 321 } 322 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 323 if (i > xs) { 324 for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxm; 325 } 326 if (i < xs+xm) { 327 for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxf; 328 } 329 } else if (i == fm) { /* interface between fast and medium regions */ 330 for (j=0; j<dof; j++) { 331 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 332 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 333 } 334 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 335 if (i > xs) { 336 for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf; 337 } 338 if (i < xs+xm) { 339 for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxm; 340 } 341 } else { /* fast region */ 342 for (j=0; j<dof; j++) { 343 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 344 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 345 } 346 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 347 cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */ 348 if (i > xs) { 349 for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hxf; 350 } 351 if (i < xs+xm) { 352 for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hxf; 353 } 354 } 355 } 356 ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 357 ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr); 358 ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 359 ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 360 ierr = MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr); 361 if (0) { 362 /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */ 363 PetscReal dt,tnow; 364 ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 365 ierr = TSGetTime(ts,&tnow);CHKERRQ(ierr); 366 if (dt > 0.5/ctx->cfl_idt) { 367 ierr = PetscPrintf(ctx->comm,"Stability constraint exceeded at t=%g, dt %g > %g\n",(double)tnow,(double)dt,(double)(0.5/ctx->cfl_idt));CHKERRQ(ierr); 368 } 369 } 370 PetscFunctionReturn(0); 371 } 372 373 /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */ 374 PetscErrorCode FVRHSFunctionslow_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 375 { 376 FVCtx *ctx = (FVCtx*)vctx; 377 PetscErrorCode ierr; 378 PetscInt i,j,k,Mx,dof,xs,xm,islow = 0,sm = ctx->sm,ms = ctx->ms,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth; 379 PetscReal hxs,hxm,hxf,cfl_idt = 0; 380 PetscScalar *x,*f,*slope; 381 Vec Xloc; 382 DM da; 383 384 PetscFunctionBeginUser; 385 ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 386 ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 387 ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 388 hxs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 389 hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 390 hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 391 ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 392 ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 393 ierr = VecZeroEntries(F);CHKERRQ(ierr); 394 ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 395 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 396 ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 397 ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 398 399 if (ctx->bctype == FVBC_OUTFLOW) { 400 for (i=xs-2; i<0; i++) { 401 for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 402 } 403 for (i=Mx; i<xs+xm+2; i++) { 404 for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 405 } 406 } 407 for (i=xs-1; i<xs+xm+1; i++) { 408 struct _LimitInfo info; 409 PetscScalar *cjmpL,*cjmpR; 410 if (i < sm-lsbwidth+1 || i > ms+rsbwidth-2) { /* slow components and the first and last fast components */ 411 /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 412 ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); 413 /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 414 ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 415 cjmpL = &ctx->cjmpLR[0]; 416 cjmpR = &ctx->cjmpLR[dof]; 417 for (j=0; j<dof; j++) { 418 PetscScalar jmpL,jmpR; 419 jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 420 jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 421 for (k=0; k<dof; k++) { 422 cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 423 cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 424 } 425 } 426 /* Apply limiter to the left and right characteristic jumps */ 427 info.m = dof; 428 info.hxs = hxs; 429 info.hxm = hxm; 430 info.hxf = hxf; 431 (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope); 432 for (j=0; j<dof; j++) { 433 PetscScalar tmp = 0; 434 for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 435 slope[i*dof+j] = tmp; 436 } 437 } 438 } 439 440 for (i=xs; i<xs+xm+1; i++) { 441 PetscReal maxspeed; 442 PetscScalar *uL,*uR; 443 uL = &ctx->uLR[0]; 444 uR = &ctx->uLR[dof]; 445 if (i < sm-lsbwidth) { /* slow region */ 446 for (j=0; j<dof; j++) { 447 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 448 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 449 } 450 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 451 cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */ 452 if (i > xs) { 453 for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 454 } 455 if (i < xs+xm) { 456 for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 457 islow++; 458 } 459 } 460 if (i == sm-lsbwidth) { /* interface between slow and medium regions */ 461 for (j=0; j<dof; j++) { 462 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 463 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 464 } 465 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 466 if (i > xs) { 467 for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 468 } 469 } 470 if (i == ms+rsbwidth) { /* interface between medium and slow regions */ 471 for (j=0; j<dof; j++) { 472 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 473 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 474 } 475 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 476 if (i < xs+xm) { 477 for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 478 islow++; 479 } 480 } 481 if (i > ms+rsbwidth) { /* slow region */ 482 for (j=0; j<dof; j++) { 483 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 484 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 485 } 486 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 487 cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */ 488 if (i > xs) { 489 for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 490 } 491 if (i < xs+xm) { 492 for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 493 islow++; 494 } 495 } 496 } 497 ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 498 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 499 ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 500 ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 501 ierr = MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr); 502 PetscFunctionReturn(0); 503 } 504 505 PetscErrorCode FVRHSFunctionslowbuffer_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 506 { 507 FVCtx *ctx = (FVCtx*)vctx; 508 PetscErrorCode ierr; 509 PetscInt i,j,k,Mx,dof,xs,xm,islowbuffer = 0,sm = ctx->sm,ms = ctx->ms,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth; 510 PetscReal hxs,hxm,hxf; 511 PetscScalar *x,*f,*slope; 512 Vec Xloc; 513 DM da; 514 515 PetscFunctionBeginUser; 516 ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 517 ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 518 ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 519 hxs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 520 hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 521 hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 522 ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 523 ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 524 ierr = VecZeroEntries(F);CHKERRQ(ierr); 525 ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 526 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 527 ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 528 ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 529 530 if (ctx->bctype == FVBC_OUTFLOW) { 531 for (i=xs-2; i<0; i++) { 532 for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 533 } 534 for (i=Mx; i<xs+xm+2; i++) { 535 for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 536 } 537 } 538 for (i=xs-1; i<xs+xm+1; i++) { 539 struct _LimitInfo info; 540 PetscScalar *cjmpL,*cjmpR; 541 if ((i > sm-lsbwidth-2 && i < sm+1) || (i > ms-2 && i < ms+rsbwidth+1)) { 542 /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 543 ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); 544 /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 545 ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 546 cjmpL = &ctx->cjmpLR[0]; 547 cjmpR = &ctx->cjmpLR[dof]; 548 for (j=0; j<dof; j++) { 549 PetscScalar jmpL,jmpR; 550 jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 551 jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 552 for (k=0; k<dof; k++) { 553 cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 554 cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 555 } 556 } 557 /* Apply limiter to the left and right characteristic jumps */ 558 info.m = dof; 559 info.hxs = hxs; 560 info.hxm = hxm; 561 info.hxf = hxf; 562 (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope); 563 for (j=0; j<dof; j++) { 564 PetscScalar tmp = 0; 565 for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 566 slope[i*dof+j] = tmp; 567 } 568 } 569 } 570 571 for (i=xs; i<xs+xm+1; i++) { 572 PetscReal maxspeed; 573 PetscScalar *uL,*uR; 574 uL = &ctx->uLR[0]; 575 uR = &ctx->uLR[dof]; 576 if (i == sm-lsbwidth) { 577 for (j=0; j<dof; j++) { 578 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 579 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 580 } 581 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 582 if (i < xs+xm) { 583 for (j=0; j<dof; j++) f[islowbuffer*dof+j] += ctx->flux[j]/hxs; 584 islowbuffer++; 585 } 586 } 587 if (i > sm-lsbwidth && i < sm) { 588 for (j=0; j<dof; j++) { 589 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 590 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 591 } 592 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 593 if (i > xs) { 594 for (j=0; j<dof; j++) f[(islowbuffer-1)*dof+j] -= ctx->flux[j]/hxs; 595 } 596 if (i < xs+xm) { 597 for (j=0; j<dof; j++) f[islowbuffer*dof+j] += ctx->flux[j]/hxs; 598 islowbuffer++; 599 } 600 } 601 if (i == sm) { /* interface between the slow region and the medium region */ 602 for (j=0; j<dof; j++) { 603 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 604 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 605 } 606 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 607 if (i > xs) { 608 for (j=0; j<dof; j++) f[(islowbuffer-1)*dof+j] -= ctx->flux[j]/hxs; 609 } 610 } 611 if (i == ms) { /* interface between the medium region and the slow region */ 612 for (j=0; j<dof; j++) { 613 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 614 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 615 } 616 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 617 if (i < xs+xm) { 618 for (j=0; j<dof; j++) f[islowbuffer*dof+j] += ctx->flux[j]/hxs; 619 islowbuffer++; 620 } 621 } 622 if (i > ms && i < ms+rsbwidth) { 623 for (j=0; j<dof; j++) { 624 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 625 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 626 } 627 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 628 if (i > xs) { 629 for (j=0; j<dof; j++) f[(islowbuffer-1)*dof+j] -= ctx->flux[j]/hxs; 630 } 631 if (i < xs+xm) { 632 for (j=0; j<dof; j++) f[islowbuffer*dof+j] += ctx->flux[j]/hxs; 633 islowbuffer++; 634 } 635 } 636 if (i == ms+rsbwidth) { 637 for (j=0; j<dof; j++) { 638 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 639 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 640 } 641 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 642 if (i > xs) { 643 for (j=0; j<dof; j++) f[(islowbuffer-1)*dof+j] -= ctx->flux[j]/hxs; 644 } 645 } 646 } 647 ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 648 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 649 ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 650 ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 651 PetscFunctionReturn(0); 652 } 653 654 /* --------------------------------- Finite Volume Solver for medium components ----------------------------------- */ 655 PetscErrorCode FVRHSFunctionmedium_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 656 { 657 FVCtx *ctx = (FVCtx*)vctx; 658 PetscErrorCode ierr; 659 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; 660 PetscReal hxs,hxm,hxf; 661 PetscScalar *x,*f,*slope; 662 Vec Xloc; 663 DM da; 664 665 PetscFunctionBeginUser; 666 ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 667 ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 668 ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 669 hxs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 670 hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 671 hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 672 ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 673 ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 674 ierr = VecZeroEntries(F);CHKERRQ(ierr); 675 ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 676 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 677 ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 678 ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 679 680 if (ctx->bctype == FVBC_OUTFLOW) { 681 for (i=xs-2; i<0; i++) { 682 for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 683 } 684 for (i=Mx; i<xs+xm+2; i++) { 685 for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 686 } 687 } 688 for (i=xs-1; i<xs+xm+1; i++) { 689 struct _LimitInfo info; 690 PetscScalar *cjmpL,*cjmpR; 691 if ((i > sm-2 && i < mf-lmbwidth+1) || (i > fm+rmbwidth-2 && i < ms+1)) { /* slow components and the first and last fast components */ 692 /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 693 ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); 694 /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 695 ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 696 cjmpL = &ctx->cjmpLR[0]; 697 cjmpR = &ctx->cjmpLR[dof]; 698 for (j=0; j<dof; j++) { 699 PetscScalar jmpL,jmpR; 700 jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 701 jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 702 for (k=0; k<dof; k++) { 703 cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 704 cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 705 } 706 } 707 /* Apply limiter to the left and right characteristic jumps */ 708 info.m = dof; 709 info.hxs = hxs; 710 info.hxm = hxm; 711 info.hxf = hxf; 712 (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope); 713 for (j=0; j<dof; j++) { 714 PetscScalar tmp = 0; 715 for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 716 slope[i*dof+j] = tmp; 717 } 718 } 719 } 720 721 for (i=xs; i<xs+xm+1; i++) { 722 PetscReal maxspeed; 723 PetscScalar *uL,*uR; 724 uL = &ctx->uLR[0]; 725 uR = &ctx->uLR[dof]; 726 if (i == sm) { /* interface between slow and medium regions */ 727 for (j=0; j<dof; j++) { 728 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 729 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 730 } 731 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 732 if (i < xs+xm) { 733 for (j=0; j<dof; j++) f[imedium*dof+j] += ctx->flux[j]/hxm; 734 imedium++; 735 } 736 } 737 if (i > sm && i < mf-lmbwidth) { /* medium region */ 738 for (j=0; j<dof; j++) { 739 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 740 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 741 } 742 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 743 if (i > xs) { 744 for (j=0; j<dof; j++) f[(imedium-1)*dof+j] -= ctx->flux[j]/hxm; 745 } 746 if (i < xs+xm) { 747 for (j=0; j<dof; j++) f[imedium*dof+j] += ctx->flux[j]/hxm; 748 imedium++; 749 } 750 } 751 if (i == mf-lmbwidth) { /* interface between medium and fast regions */ 752 for (j=0; j<dof; j++) { 753 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 754 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 755 } 756 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 757 if (i > xs) { 758 for (j=0; j<dof; j++) f[(imedium-1)*dof+j] -= ctx->flux[j]/hxm; 759 } 760 } 761 if (i == fm+rmbwidth) { /* interface between fast and medium regions */ 762 for (j=0; j<dof; j++) { 763 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 764 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 765 } 766 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 767 if (i < xs+xm) { 768 for (j=0; j<dof; j++) f[imedium*dof+j] += ctx->flux[j]/hxm; 769 imedium++; 770 } 771 } 772 if (i > fm+rmbwidth && i < ms) { /* medium region */ 773 for (j=0; j<dof; j++) { 774 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 775 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 776 } 777 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 778 if (i > xs) { 779 for (j=0; j<dof; j++) f[(imedium-1)*dof+j] -= ctx->flux[j]/hxm; 780 } 781 if (i < xs+xm) { 782 for (j=0; j<dof; j++) f[imedium*dof+j] += ctx->flux[j]/hxm; 783 imedium++; 784 } 785 } 786 if (i == ms) { /* interface between medium and slow regions */ 787 for (j=0; j<dof; j++) { 788 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 789 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 790 } 791 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 792 if (i > xs) { 793 for (j=0; j<dof; j++) f[(imedium-1)*dof+j] -= ctx->flux[j]/hxm; 794 } 795 } 796 } 797 ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 798 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 799 ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 800 ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 801 PetscFunctionReturn(0); 802 } 803 804 PetscErrorCode FVRHSFunctionmediumbuffer_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 805 { 806 FVCtx *ctx = (FVCtx*)vctx; 807 PetscErrorCode ierr; 808 PetscInt i,j,k,Mx,dof,xs,xm,imediumbuffer = 0,mf = ctx->mf,fm = ctx->fm,lmbwidth = ctx->lmbwidth,rmbwidth = ctx->rmbwidth; 809 PetscReal hxs,hxm,hxf; 810 PetscScalar *x,*f,*slope; 811 Vec Xloc; 812 DM da; 813 814 PetscFunctionBeginUser; 815 ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 816 ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 817 ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 818 hxs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 819 hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 820 hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 821 ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 822 ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 823 ierr = VecZeroEntries(F);CHKERRQ(ierr); 824 ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 825 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 826 ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 827 ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 828 829 if (ctx->bctype == FVBC_OUTFLOW) { 830 for (i=xs-2; i<0; i++) { 831 for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 832 } 833 for (i=Mx; i<xs+xm+2; i++) { 834 for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 835 } 836 } 837 for (i=xs-1; i<xs+xm+1; i++) { 838 struct _LimitInfo info; 839 PetscScalar *cjmpL,*cjmpR; 840 if ((i > mf-lmbwidth-2 && i < mf+1) || (i > fm-2 && i < fm+rmbwidth+1)) { 841 /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 842 ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); 843 /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 844 ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 845 cjmpL = &ctx->cjmpLR[0]; 846 cjmpR = &ctx->cjmpLR[dof]; 847 for (j=0; j<dof; j++) { 848 PetscScalar jmpL,jmpR; 849 jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 850 jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 851 for (k=0; k<dof; k++) { 852 cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 853 cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 854 } 855 } 856 /* Apply limiter to the left and right characteristic jumps */ 857 info.m = dof; 858 info.hxs = hxs; 859 info.hxm = hxm; 860 info.hxf = hxf; 861 (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope); 862 for (j=0; j<dof; j++) { 863 PetscScalar tmp = 0; 864 for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 865 slope[i*dof+j] = tmp; 866 } 867 } 868 } 869 870 for (i=xs; i<xs+xm+1; i++) { 871 PetscReal maxspeed; 872 PetscScalar *uL,*uR; 873 uL = &ctx->uLR[0]; 874 uR = &ctx->uLR[dof]; 875 if (i == mf-lmbwidth) { 876 for (j=0; j<dof; j++) { 877 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 878 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 879 } 880 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 881 if (i < xs+xm) { 882 for (j=0; j<dof; j++) f[imediumbuffer*dof+j] += ctx->flux[j]/hxm; 883 imediumbuffer++; 884 } 885 } 886 if (i > mf-lmbwidth && i < mf) { 887 for (j=0; j<dof; j++) { 888 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 889 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 890 } 891 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 892 if (i > xs) { 893 for (j=0; j<dof; j++) f[(imediumbuffer-1)*dof+j] -= ctx->flux[j]/hxm; 894 } 895 if (i < xs+xm) { 896 for (j=0; j<dof; j++) f[imediumbuffer*dof+j] += ctx->flux[j]/hxm; 897 imediumbuffer++; 898 } 899 } 900 if (i == mf) { /* interface between the medium region and the fast region */ 901 for (j=0; j<dof; j++) { 902 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 903 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 904 } 905 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 906 if (i > xs) { 907 for (j=0; j<dof; j++) f[(imediumbuffer-1)*dof+j] -= ctx->flux[j]/hxm; 908 } 909 } 910 if (i == fm) { /* interface between the fast region and the medium region */ 911 for (j=0; j<dof; j++) { 912 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 913 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 914 } 915 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 916 if (i < xs+xm) { 917 for (j=0; j<dof; j++) f[imediumbuffer*dof+j] += ctx->flux[j]/hxm; 918 imediumbuffer++; 919 } 920 } 921 if (i > fm && i < fm+rmbwidth) { 922 for (j=0; j<dof; j++) { 923 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 924 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 925 } 926 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 927 if (i > xs) { 928 for (j=0; j<dof; j++) f[(imediumbuffer-1)*dof+j] -= ctx->flux[j]/hxm; 929 } 930 if (i < xs+xm) { 931 for (j=0; j<dof; j++) f[imediumbuffer*dof+j] += ctx->flux[j]/hxm; 932 imediumbuffer++; 933 } 934 } 935 if (i == fm+rmbwidth) { 936 for (j=0; j<dof; j++) { 937 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 938 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 939 } 940 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 941 if (i > xs) { 942 for (j=0; j<dof; j++) f[(imediumbuffer-1)*dof+j] -= ctx->flux[j]/hxm; 943 } 944 } 945 } 946 ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 947 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 948 ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 949 ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 950 PetscFunctionReturn(0); 951 } 952 953 /* --------------------------------- Finite Volume Solver for fast parts ----------------------------------- */ 954 PetscErrorCode FVRHSFunctionfast_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 955 { 956 FVCtx *ctx = (FVCtx*)vctx; 957 PetscErrorCode ierr; 958 PetscInt i,j,k,Mx,dof,xs,xm,ifast = 0,mf = ctx->mf,fm = ctx->fm; 959 PetscReal hxs,hxm,hxf; 960 PetscScalar *x,*f,*slope; 961 Vec Xloc; 962 DM da; 963 964 PetscFunctionBeginUser; 965 ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 966 ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 967 ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 968 hxs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 969 hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 970 hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 971 ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 972 ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 973 ierr = VecZeroEntries(F);CHKERRQ(ierr); 974 ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 975 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 976 ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 977 ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 978 979 if (ctx->bctype == FVBC_OUTFLOW) { 980 for (i=xs-2; i<0; i++) { 981 for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 982 } 983 for (i=Mx; i<xs+xm+2; i++) { 984 for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 985 } 986 } 987 for (i=xs-1; i<xs+xm+1; i++) { /* fast components and the last slow componets before fast components and the first slow component after fast components */ 988 struct _LimitInfo info; 989 PetscScalar *cjmpL,*cjmpR; 990 if (i > mf-2 && i < fm+1) { 991 ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); 992 ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 993 cjmpL = &ctx->cjmpLR[0]; 994 cjmpR = &ctx->cjmpLR[dof]; 995 for (j=0; j<dof; j++) { 996 PetscScalar jmpL,jmpR; 997 jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 998 jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 999 for (k=0; k<dof; k++) { 1000 cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 1001 cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 1002 } 1003 } 1004 /* Apply limiter to the left and right characteristic jumps */ 1005 info.m = dof; 1006 info.hxs = hxs; 1007 info.hxm = hxm; 1008 info.hxf = hxf; 1009 (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope); 1010 for (j=0; j<dof; j++) { 1011 PetscScalar tmp = 0; 1012 for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 1013 slope[i*dof+j] = tmp; 1014 } 1015 } 1016 } 1017 1018 for (i=xs; i<xs+xm+1; i++) { 1019 PetscReal maxspeed; 1020 PetscScalar *uL,*uR; 1021 uL = &ctx->uLR[0]; 1022 uR = &ctx->uLR[dof]; 1023 if (i == mf) { /* interface between medium and fast regions */ 1024 for (j=0; j<dof; j++) { 1025 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 1026 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 1027 } 1028 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 1029 if (i < xs+xm) { 1030 for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf; 1031 ifast++; 1032 } 1033 } 1034 if (i > mf && i < fm) { /* fast region */ 1035 for (j=0; j<dof; j++) { 1036 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 1037 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 1038 } 1039 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 1040 if (i > xs) { 1041 for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf; 1042 } 1043 if (i < xs+xm) { 1044 for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf; 1045 ifast++; 1046 } 1047 } 1048 if (i == fm) { /* interface between fast and medium regions */ 1049 for (j=0; j<dof; j++) { 1050 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 1051 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 1052 } 1053 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 1054 if (i > xs) { 1055 for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf; 1056 } 1057 } 1058 } 1059 ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 1060 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 1061 ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 1062 ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 1063 PetscFunctionReturn(0); 1064 } 1065 1066 int main(int argc,char *argv[]) 1067 { 1068 char lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m"; 1069 PetscFunctionList limiters = 0,physics = 0; 1070 MPI_Comm comm; 1071 TS ts; 1072 DM da; 1073 Vec X,X0,R; 1074 FVCtx ctx; 1075 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; 1076 PetscBool view_final = PETSC_FALSE; 1077 PetscReal ptime; 1078 PetscErrorCode ierr; 1079 1080 ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr; 1081 comm = PETSC_COMM_WORLD; 1082 ierr = PetscMemzero(&ctx,sizeof(ctx));CHKERRQ(ierr); 1083 1084 /* Register limiters to be available on the command line */ 1085 ierr = PetscFunctionListAdd(&limiters,"upwind" ,Limit3_Upwind);CHKERRQ(ierr); 1086 ierr = PetscFunctionListAdd(&limiters,"lax-wendroff" ,Limit3_LaxWendroff);CHKERRQ(ierr); 1087 ierr = PetscFunctionListAdd(&limiters,"beam-warming" ,Limit3_BeamWarming);CHKERRQ(ierr); 1088 ierr = PetscFunctionListAdd(&limiters,"fromm" ,Limit3_Fromm);CHKERRQ(ierr); 1089 ierr = PetscFunctionListAdd(&limiters,"minmod" ,Limit3_Minmod);CHKERRQ(ierr); 1090 ierr = PetscFunctionListAdd(&limiters,"superbee" ,Limit3_Superbee);CHKERRQ(ierr); 1091 ierr = PetscFunctionListAdd(&limiters,"mc" ,Limit3_MC);CHKERRQ(ierr); 1092 ierr = PetscFunctionListAdd(&limiters,"koren3" ,Limit3_Koren3);CHKERRQ(ierr); 1093 1094 /* Register physical models to be available on the command line */ 1095 ierr = PetscFunctionListAdd(&physics,"advect" ,PhysicsCreate_Advect);CHKERRQ(ierr); 1096 1097 ctx.comm = comm; 1098 ctx.cfl = 0.9; 1099 ctx.bctype = FVBC_PERIODIC; 1100 ctx.xmin = -1.0; 1101 ctx.xmax = 1.0; 1102 ierr = PetscOptionsBegin(comm,NULL,"Finite Volume solver options","");CHKERRQ(ierr); 1103 ierr = PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL);CHKERRQ(ierr); 1104 ierr = PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL);CHKERRQ(ierr); 1105 ierr = PetscOptionsFList("-limit","Name of flux imiter to use","",limiters,lname,lname,sizeof(lname),NULL);CHKERRQ(ierr); 1106 ierr = PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL);CHKERRQ(ierr); 1107 ierr = PetscOptionsString("-view_final","Write final solution in ASCII MATLAB format to given file name","",final_fname,final_fname,sizeof(final_fname),&view_final);CHKERRQ(ierr); 1108 ierr = PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL);CHKERRQ(ierr); 1109 ierr = PetscOptionsBool("-exact","Compare errors with exact solution","",ctx.exact,&ctx.exact,NULL);CHKERRQ(ierr); 1110 ierr = PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL);CHKERRQ(ierr); 1111 ierr = PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL);CHKERRQ(ierr); 1112 ierr = PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL);CHKERRQ(ierr); 1113 ierr = PetscOptionsInt("-hratio","Spacing ratio","",ctx.hratio,&ctx.hratio,NULL);CHKERRQ(ierr); 1114 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1115 1116 /* Choose the limiter from the list of registered limiters */ 1117 ierr = PetscFunctionListFind(limiters,lname,&ctx.limit3);CHKERRQ(ierr); 1118 if (!ctx.limit3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Limiter '%s' not found",lname); 1119 1120 /* Choose the physics from the list of registered models */ 1121 { 1122 PetscErrorCode (*r)(FVCtx*); 1123 ierr = PetscFunctionListFind(physics,physname,&r);CHKERRQ(ierr); 1124 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Physics '%s' not found",physname); 1125 /* Create the physics, will set the number of fields and their names */ 1126 ierr = (*r)(&ctx);CHKERRQ(ierr); 1127 } 1128 1129 /* Create a DMDA to manage the parallel grid */ 1130 ierr = DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics2.dof,2,NULL,&da);CHKERRQ(ierr); 1131 ierr = DMSetFromOptions(da);CHKERRQ(ierr); 1132 ierr = DMSetUp(da);CHKERRQ(ierr); 1133 /* Inform the DMDA of the field names provided by the physics. */ 1134 /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ 1135 for (i=0; i<ctx.physics2.dof; i++) { 1136 ierr = DMDASetFieldName(da,i,ctx.physics2.fieldname[i]);CHKERRQ(ierr); 1137 } 1138 ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 1139 ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 1140 1141 /* Set coordinates of cell centers */ 1142 ierr = DMDASetUniformCoordinates(da,ctx.xmin+0.5*(ctx.xmax-ctx.xmin)/Mx,ctx.xmax+0.5*(ctx.xmax-ctx.xmin)/Mx,0,0,0,0);CHKERRQ(ierr); 1143 1144 /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */ 1145 ierr = PetscMalloc4(dof*dof,&ctx.R,dof*dof,&ctx.Rinv,2*dof,&ctx.cjmpLR,1*dof,&ctx.cslope);CHKERRQ(ierr); 1146 ierr = PetscMalloc3(2*dof,&ctx.uLR,dof,&ctx.flux,dof,&ctx.speeds);CHKERRQ(ierr); 1147 1148 /* Create a vector to store the solution and to save the initial state */ 1149 ierr = DMCreateGlobalVector(da,&X);CHKERRQ(ierr); 1150 ierr = VecDuplicate(X,&X0);CHKERRQ(ierr); 1151 ierr = VecDuplicate(X,&R);CHKERRQ(ierr); 1152 1153 /* create index for slow parts and fast parts, 1154 count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */ 1155 count_slow = Mx/(1+ctx.hratio)/(1+ctx.hratio); 1156 count_medium = 2*ctx.hratio*count_slow; 1157 if (count_slow%2 || count_medium%2) SETERRQ(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"); 1158 count_fast = ctx.hratio*ctx.hratio*count_slow; 1159 ctx.sm = count_slow/2; 1160 ctx.mf = ctx.sm + count_medium/2; 1161 ctx.fm = ctx.mf + count_fast; 1162 ctx.ms = ctx.fm + count_medium/2; 1163 ierr = PetscMalloc1(xm*dof,&index_slow);CHKERRQ(ierr); 1164 ierr = PetscMalloc1(xm*dof,&index_medium);CHKERRQ(ierr); 1165 ierr = PetscMalloc1(xm*dof,&index_fast);CHKERRQ(ierr); 1166 ierr = PetscMalloc1(6*dof,&index_slowbuffer);CHKERRQ(ierr); 1167 ierr = PetscMalloc1(6*dof,&index_mediumbuffer);CHKERRQ(ierr); 1168 if (((AdvectCtx*)ctx.physics2.user)->a > 0) { 1169 ctx.lsbwidth = 2; 1170 ctx.rsbwidth = 4; 1171 ctx.lmbwidth = 2; 1172 ctx.rmbwidth = 4; 1173 } else { 1174 ctx.lsbwidth = 4; 1175 ctx.rsbwidth = 2; 1176 ctx.lmbwidth = 4; 1177 ctx.rmbwidth = 2; 1178 } 1179 1180 for (i=xs; i<xs+xm; i++) { 1181 if (i < ctx.sm-ctx.lsbwidth || i > ctx.ms+ctx.rsbwidth-1) 1182 for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k; 1183 else if ((i >= ctx.sm-ctx.lsbwidth && i < ctx.sm) || (i > ctx.ms-1 && i <= ctx.ms+ctx.rsbwidth-1)) 1184 for (k=0; k<dof; k++) index_slowbuffer[islowbuffer++] = i*dof+k; 1185 else if (i < ctx.mf-ctx.lmbwidth || i > ctx.fm+ctx.rmbwidth-1) 1186 for (k=0; k<dof; k++) index_medium[imedium++] = i*dof+k; 1187 else if ((i >= ctx.mf-ctx.lmbwidth && i < ctx.mf) || (i > ctx.fm-1 && i <= ctx.fm+ctx.rmbwidth-1)) 1188 for (k=0; k<dof; k++) index_mediumbuffer[imediumbuffer++] = i*dof+k; 1189 else 1190 for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k; 1191 } 1192 ierr = ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss);CHKERRQ(ierr); 1193 ierr = ISCreateGeneral(PETSC_COMM_WORLD,imedium,index_medium,PETSC_COPY_VALUES,&ctx.ism);CHKERRQ(ierr); 1194 ierr = ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf);CHKERRQ(ierr); 1195 ierr = ISCreateGeneral(PETSC_COMM_WORLD,islowbuffer,index_slowbuffer,PETSC_COPY_VALUES,&ctx.issb);CHKERRQ(ierr); 1196 ierr = ISCreateGeneral(PETSC_COMM_WORLD,imediumbuffer,index_mediumbuffer,PETSC_COPY_VALUES,&ctx.ismb);CHKERRQ(ierr); 1197 1198 /* Create a time-stepping object */ 1199 ierr = TSCreate(comm,&ts);CHKERRQ(ierr); 1200 ierr = TSSetDM(ts,da);CHKERRQ(ierr); 1201 ierr = TSSetRHSFunction(ts,R,FVRHSFunction_3WaySplit,&ctx);CHKERRQ(ierr); 1202 ierr = TSRHSSplitSetIS(ts,"slow",ctx.iss);CHKERRQ(ierr); 1203 ierr = TSRHSSplitSetIS(ts,"medium",ctx.ism);CHKERRQ(ierr); 1204 ierr = TSRHSSplitSetIS(ts,"fast",ctx.isf);CHKERRQ(ierr); 1205 ierr = TSRHSSplitSetIS(ts,"slowbuffer",ctx.issb);CHKERRQ(ierr); 1206 ierr = TSRHSSplitSetIS(ts,"mediumbuffer",ctx.ismb);CHKERRQ(ierr); 1207 ierr = TSRHSSplitSetRHSFunction(ts,"slow",NULL,FVRHSFunctionslow_3WaySplit,&ctx);CHKERRQ(ierr); 1208 ierr = TSRHSSplitSetRHSFunction(ts,"medium",NULL,FVRHSFunctionmedium_3WaySplit,&ctx);CHKERRQ(ierr); 1209 ierr = TSRHSSplitSetRHSFunction(ts,"fast",NULL,FVRHSFunctionfast_3WaySplit,&ctx);CHKERRQ(ierr); 1210 ierr = TSRHSSplitSetRHSFunction(ts,"slowbuffer",NULL,FVRHSFunctionslowbuffer_3WaySplit,&ctx);CHKERRQ(ierr); 1211 ierr = TSRHSSplitSetRHSFunction(ts,"mediumbuffer",NULL,FVRHSFunctionmediumbuffer_3WaySplit,&ctx);CHKERRQ(ierr); 1212 1213 ierr = TSSetType(ts,TSSSP);CHKERRQ(ierr); 1214 /*ierr = TSSetType(ts,TSMPRK);CHKERRQ(ierr);*/ 1215 ierr = TSSetMaxTime(ts,10);CHKERRQ(ierr); 1216 ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 1217 1218 /* Compute initial conditions and starting time step */ 1219 ierr = FVSample_3WaySplit(&ctx,da,0,X0);CHKERRQ(ierr); 1220 ierr = FVRHSFunction_3WaySplit(ts,0,X0,X,(void*)&ctx);CHKERRQ(ierr); /* Initial function evaluation, only used to determine max speed */ 1221 ierr = VecCopy(X0,X);CHKERRQ(ierr); /* The function value was not used so we set X=X0 again */ 1222 ierr = TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt);CHKERRQ(ierr); 1223 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); /* Take runtime options */ 1224 ierr = SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1225 { 1226 PetscInt steps; 1227 PetscScalar mass_initial,mass_final,mass_difference,mass_differenceg; 1228 const PetscScalar *ptr_X,*ptr_X0; 1229 const PetscReal hs = (ctx.xmax-ctx.xmin)/4.0/count_slow; 1230 const PetscReal hm = (ctx.xmax-ctx.xmin)/2.0/count_medium; 1231 const PetscReal hf = (ctx.xmax-ctx.xmin)/4.0/count_fast; 1232 1233 ierr = TSSolve(ts,X);CHKERRQ(ierr); 1234 ierr = TSGetSolveTime(ts,&ptime);CHKERRQ(ierr); 1235 ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 1236 /* calculate the total mass at initial time and final time */ 1237 mass_initial = 0.0; 1238 mass_final = 0.0; 1239 ierr = DMDAVecGetArrayRead(da,X0,(void*)&ptr_X0);CHKERRQ(ierr); 1240 ierr = DMDAVecGetArrayRead(da,X,(void*)&ptr_X);CHKERRQ(ierr); 1241 for (i=xs;i<xs+xm;i++) { 1242 if (i < ctx.sm || i > ctx.ms-1) 1243 for (k=0; k<dof; k++) { 1244 mass_initial = mass_initial + hs*ptr_X0[i*dof+k]; 1245 mass_final = mass_final + hs*ptr_X[i*dof+k]; 1246 } 1247 else if (i < ctx.mf || i > ctx.fm-1) 1248 for (k=0; k<dof; k++) { 1249 mass_initial = mass_initial + hm*ptr_X0[i*dof+k]; 1250 mass_final = mass_final + hm*ptr_X[i*dof+k]; 1251 } 1252 else { 1253 for (k=0; k<dof; k++) { 1254 mass_initial = mass_initial + hf*ptr_X0[i*dof+k]; 1255 mass_final = mass_final + hf*ptr_X[i*dof+k]; 1256 } 1257 } 1258 } 1259 ierr = DMDAVecRestoreArrayRead(da,X0,(void*)&ptr_X0);CHKERRQ(ierr); 1260 ierr = DMDAVecRestoreArrayRead(da,X,(void*)&ptr_X);CHKERRQ(ierr); 1261 mass_difference = mass_final - mass_initial; 1262 ierr = MPI_Allreduce(&mass_difference,&mass_differenceg,1,MPIU_SCALAR,MPIU_SUM,comm);CHKERRMPI(ierr); 1263 ierr = PetscPrintf(comm,"Mass difference %g\n",(double)mass_differenceg);CHKERRQ(ierr); 1264 ierr = PetscPrintf(comm,"Final time %g, steps %D\n",(double)ptime,steps);CHKERRQ(ierr); 1265 ierr = PetscPrintf(comm,"Maximum allowable stepsize according to CFL %g\n",(double)(1.0/ctx.cfl_idt));CHKERRQ(ierr); 1266 if (ctx.exact) { 1267 PetscReal nrm1=0; 1268 ierr = SolutionErrorNorms_3WaySplit(&ctx,da,ptime,X,&nrm1);CHKERRQ(ierr); 1269 ierr = PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1);CHKERRQ(ierr); 1270 } 1271 if (ctx.simulation) { 1272 PetscReal nrm1=0; 1273 PetscViewer fd; 1274 char filename[PETSC_MAX_PATH_LEN] = "binaryoutput"; 1275 Vec XR; 1276 PetscBool flg; 1277 const PetscScalar *ptr_XR; 1278 ierr = PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg);CHKERRQ(ierr); 1279 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 1280 ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&fd);CHKERRQ(ierr); 1281 ierr = VecDuplicate(X0,&XR);CHKERRQ(ierr); 1282 ierr = VecLoad(XR,fd);CHKERRQ(ierr); 1283 ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 1284 ierr = VecGetArrayRead(X,&ptr_X);CHKERRQ(ierr); 1285 ierr = VecGetArrayRead(XR,&ptr_XR);CHKERRQ(ierr); 1286 for (i=xs;i<xs+xm;i++) { 1287 if (i < ctx.sm || i > ctx.ms-1) 1288 for (k=0; k<dof; k++) nrm1 = nrm1 + hs*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]); 1289 else if (i < ctx.mf || i < ctx.fm-1) 1290 for (k=0; k<dof; k++) nrm1 = nrm1 + hm*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]); 1291 else 1292 for (k=0; k<dof; k++) nrm1 = nrm1 + hf*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]); 1293 } 1294 ierr = VecRestoreArrayRead(X,&ptr_X);CHKERRQ(ierr); 1295 ierr = VecRestoreArrayRead(XR,&ptr_XR);CHKERRQ(ierr); 1296 ierr = PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1);CHKERRQ(ierr); 1297 ierr = VecDestroy(&XR);CHKERRQ(ierr); 1298 } 1299 } 1300 1301 ierr = SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1302 if (draw & 0x1) {ierr = VecView(X0,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);} 1303 if (draw & 0x2) {ierr = VecView(X,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);} 1304 if (draw & 0x4) { 1305 Vec Y; 1306 ierr = VecDuplicate(X,&Y);CHKERRQ(ierr); 1307 ierr = FVSample_3WaySplit(&ctx,da,ptime,Y);CHKERRQ(ierr); 1308 ierr = VecAYPX(Y,-1,X);CHKERRQ(ierr); 1309 ierr = VecView(Y,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr); 1310 ierr = VecDestroy(&Y);CHKERRQ(ierr); 1311 } 1312 1313 if (view_final) { 1314 PetscViewer viewer; 1315 ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer);CHKERRQ(ierr); 1316 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 1317 ierr = VecView(X,viewer);CHKERRQ(ierr); 1318 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1319 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1320 } 1321 1322 /* Clean up */ 1323 ierr = (*ctx.physics2.destroy)(ctx.physics2.user);CHKERRQ(ierr); 1324 for (i=0; i<ctx.physics2.dof; i++) {ierr = PetscFree(ctx.physics2.fieldname[i]);CHKERRQ(ierr);} 1325 ierr = PetscFree4(ctx.R,ctx.Rinv,ctx.cjmpLR,ctx.cslope);CHKERRQ(ierr); 1326 ierr = PetscFree3(ctx.uLR,ctx.flux,ctx.speeds);CHKERRQ(ierr); 1327 ierr = VecDestroy(&X);CHKERRQ(ierr); 1328 ierr = VecDestroy(&X0);CHKERRQ(ierr); 1329 ierr = VecDestroy(&R);CHKERRQ(ierr); 1330 ierr = DMDestroy(&da);CHKERRQ(ierr); 1331 ierr = TSDestroy(&ts);CHKERRQ(ierr); 1332 ierr = ISDestroy(&ctx.iss);CHKERRQ(ierr); 1333 ierr = ISDestroy(&ctx.ism);CHKERRQ(ierr); 1334 ierr = ISDestroy(&ctx.isf);CHKERRQ(ierr); 1335 ierr = ISDestroy(&ctx.issb);CHKERRQ(ierr); 1336 ierr = ISDestroy(&ctx.ismb);CHKERRQ(ierr); 1337 ierr = PetscFree(index_slow);CHKERRQ(ierr); 1338 ierr = PetscFree(index_medium);CHKERRQ(ierr); 1339 ierr = PetscFree(index_fast);CHKERRQ(ierr); 1340 ierr = PetscFree(index_slowbuffer);CHKERRQ(ierr); 1341 ierr = PetscFree(index_mediumbuffer);CHKERRQ(ierr); 1342 ierr = PetscFunctionListDestroy(&limiters);CHKERRQ(ierr); 1343 ierr = PetscFunctionListDestroy(&physics);CHKERRQ(ierr); 1344 ierr = PetscFinalize(); 1345 return ierr; 1346 } 1347 1348 /*TEST 1349 1350 build: 1351 requires: !complex 1352 depends: finitevolume1d.c 1353 1354 test: 1355 suffix: 1 1356 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 1357 1358 test: 1359 suffix: 2 1360 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 1361 1362 TEST*/ 1363