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));CHKERRQ(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 if (1) { 368 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); 369 } else SETERRQ2(PETSC_COMM_SELF,1,"Stability constraint exceeded, %g > %g",(double)dt,(double)(ctx->cfl/ctx->cfl_idt)); 370 } 371 } 372 PetscFunctionReturn(0); 373 } 374 375 /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */ 376 PetscErrorCode FVRHSFunctionslow_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 377 { 378 FVCtx *ctx = (FVCtx*)vctx; 379 PetscErrorCode ierr; 380 PetscInt i,j,k,Mx,dof,xs,xm,islow = 0,sm = ctx->sm,ms = ctx->ms,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth; 381 PetscReal hxs,hxm,hxf,cfl_idt = 0; 382 PetscScalar *x,*f,*slope; 383 Vec Xloc; 384 DM da; 385 386 PetscFunctionBeginUser; 387 ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 388 ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 389 ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 390 hxs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 391 hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 392 hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 393 ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 394 ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 395 ierr = VecZeroEntries(F);CHKERRQ(ierr); 396 ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 397 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 398 ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 399 ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 400 401 if (ctx->bctype == FVBC_OUTFLOW) { 402 for (i=xs-2; i<0; i++) { 403 for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 404 } 405 for (i=Mx; i<xs+xm+2; i++) { 406 for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 407 } 408 } 409 for (i=xs-1; i<xs+xm+1; i++) { 410 struct _LimitInfo info; 411 PetscScalar *cjmpL,*cjmpR; 412 if (i < sm-lsbwidth+1 || i > ms+rsbwidth-2) { /* slow components and the first and last fast components */ 413 /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 414 ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); 415 /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 416 ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 417 cjmpL = &ctx->cjmpLR[0]; 418 cjmpR = &ctx->cjmpLR[dof]; 419 for (j=0; j<dof; j++) { 420 PetscScalar jmpL,jmpR; 421 jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 422 jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 423 for (k=0; k<dof; k++) { 424 cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 425 cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 426 } 427 } 428 /* Apply limiter to the left and right characteristic jumps */ 429 info.m = dof; 430 info.hxs = hxs; 431 info.hxm = hxm; 432 info.hxf = hxf; 433 (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope); 434 for (j=0; j<dof; j++) { 435 PetscScalar tmp = 0; 436 for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 437 slope[i*dof+j] = tmp; 438 } 439 } 440 } 441 442 for (i=xs; i<xs+xm+1; i++) { 443 PetscReal maxspeed; 444 PetscScalar *uL,*uR; 445 uL = &ctx->uLR[0]; 446 uR = &ctx->uLR[dof]; 447 if (i < sm-lsbwidth) { /* slow region */ 448 for (j=0; j<dof; j++) { 449 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 450 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 451 } 452 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 453 cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */ 454 if (i > xs) { 455 for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 456 } 457 if (i < xs+xm) { 458 for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 459 islow++; 460 } 461 } 462 if (i == sm-lsbwidth) { /* interface between slow and medium regions */ 463 for (j=0; j<dof; j++) { 464 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 465 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 466 } 467 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 468 if (i > xs) { 469 for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 470 } 471 } 472 if (i == ms+rsbwidth) { /* interface between medium and slow regions */ 473 for (j=0; j<dof; j++) { 474 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 475 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 476 } 477 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 478 if (i < xs+xm) { 479 for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 480 islow++; 481 } 482 } 483 if (i > ms+rsbwidth) { /* slow region */ 484 for (j=0; j<dof; j++) { 485 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 486 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 487 } 488 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 489 cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hxs)); /* Max allowable value of 1/Delta t */ 490 if (i > xs) { 491 for (j=0; j<dof; j++) f[(islow-1)*dof+j] -= ctx->flux[j]/hxs; 492 } 493 if (i < xs+xm) { 494 for (j=0; j<dof; j++) f[islow*dof+j] += ctx->flux[j]/hxs; 495 islow++; 496 } 497 } 498 } 499 ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 500 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 501 ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 502 ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 503 ierr = MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_SCALAR,MPIU_MAX,PetscObjectComm((PetscObject)da));CHKERRQ(ierr); 504 PetscFunctionReturn(0); 505 } 506 507 PetscErrorCode FVRHSFunctionslowbuffer_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 508 { 509 FVCtx *ctx = (FVCtx*)vctx; 510 PetscErrorCode ierr; 511 PetscInt i,j,k,Mx,dof,xs,xm,islowbuffer = 0,sm = ctx->sm,ms = ctx->ms,lsbwidth = ctx->lsbwidth,rsbwidth = ctx->rsbwidth; 512 PetscReal hxs,hxm,hxf; 513 PetscScalar *x,*f,*slope; 514 Vec Xloc; 515 DM da; 516 517 PetscFunctionBeginUser; 518 ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 519 ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 520 ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 521 hxs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 522 hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 523 hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 524 ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 525 ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 526 ierr = VecZeroEntries(F);CHKERRQ(ierr); 527 ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 528 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 529 ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 530 ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 531 532 if (ctx->bctype == FVBC_OUTFLOW) { 533 for (i=xs-2; i<0; i++) { 534 for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 535 } 536 for (i=Mx; i<xs+xm+2; i++) { 537 for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 538 } 539 } 540 for (i=xs-1; i<xs+xm+1; i++) { 541 struct _LimitInfo info; 542 PetscScalar *cjmpL,*cjmpR; 543 if ((i > sm-lsbwidth-2 && i < sm+1) || (i > ms-2 && i < ms+rsbwidth+1)) { 544 /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 545 ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); 546 /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 547 ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 548 cjmpL = &ctx->cjmpLR[0]; 549 cjmpR = &ctx->cjmpLR[dof]; 550 for (j=0; j<dof; j++) { 551 PetscScalar jmpL,jmpR; 552 jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 553 jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 554 for (k=0; k<dof; k++) { 555 cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 556 cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 557 } 558 } 559 /* Apply limiter to the left and right characteristic jumps */ 560 info.m = dof; 561 info.hxs = hxs; 562 info.hxm = hxm; 563 info.hxf = hxf; 564 (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope); 565 for (j=0; j<dof; j++) { 566 PetscScalar tmp = 0; 567 for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 568 slope[i*dof+j] = tmp; 569 } 570 } 571 } 572 573 for (i=xs; i<xs+xm+1; i++) { 574 PetscReal maxspeed; 575 PetscScalar *uL,*uR; 576 uL = &ctx->uLR[0]; 577 uR = &ctx->uLR[dof]; 578 if (i == sm-lsbwidth) { 579 for (j=0; j<dof; j++) { 580 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 581 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 582 } 583 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 584 if (i < xs+xm) { 585 for (j=0; j<dof; j++) f[islowbuffer*dof+j] += ctx->flux[j]/hxs; 586 islowbuffer++; 587 } 588 } 589 if (i > sm-lsbwidth && i < sm) { 590 for (j=0; j<dof; j++) { 591 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 592 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 593 } 594 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 595 if (i > xs) { 596 for (j=0; j<dof; j++) f[(islowbuffer-1)*dof+j] -= ctx->flux[j]/hxs; 597 } 598 if (i < xs+xm) { 599 for (j=0; j<dof; j++) f[islowbuffer*dof+j] += ctx->flux[j]/hxs; 600 islowbuffer++; 601 } 602 } 603 if (i == sm) { /* interface between the slow region and the medium region */ 604 for (j=0; j<dof; j++) { 605 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 606 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 607 } 608 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 609 if (i > xs) { 610 for (j=0; j<dof; j++) f[(islowbuffer-1)*dof+j] -= ctx->flux[j]/hxs; 611 } 612 } 613 if (i == ms) { /* interface between the medium region and the slow region */ 614 for (j=0; j<dof; j++) { 615 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 616 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 617 } 618 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 619 if (i < xs+xm) { 620 for (j=0; j<dof; j++) f[islowbuffer*dof+j] += ctx->flux[j]/hxs; 621 islowbuffer++; 622 } 623 } 624 if (i > ms && i < ms+rsbwidth) { 625 for (j=0; j<dof; j++) { 626 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 627 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 628 } 629 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 630 if (i > xs) { 631 for (j=0; j<dof; j++) f[(islowbuffer-1)*dof+j] -= ctx->flux[j]/hxs; 632 } 633 if (i < xs+xm) { 634 for (j=0; j<dof; j++) f[islowbuffer*dof+j] += ctx->flux[j]/hxs; 635 islowbuffer++; 636 } 637 } 638 if (i == ms+rsbwidth) { 639 for (j=0; j<dof; j++) { 640 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 641 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 642 } 643 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 644 if (i > xs) { 645 for (j=0; j<dof; j++) f[(islowbuffer-1)*dof+j] -= ctx->flux[j]/hxs; 646 } 647 } 648 } 649 ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 650 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 651 ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 652 ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 653 PetscFunctionReturn(0); 654 } 655 656 /* --------------------------------- Finite Volume Solver for medium components ----------------------------------- */ 657 PetscErrorCode FVRHSFunctionmedium_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 658 { 659 FVCtx *ctx = (FVCtx*)vctx; 660 PetscErrorCode ierr; 661 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; 662 PetscReal hxs,hxm,hxf; 663 PetscScalar *x,*f,*slope; 664 Vec Xloc; 665 DM da; 666 667 PetscFunctionBeginUser; 668 ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 669 ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 670 ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 671 hxs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 672 hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 673 hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 674 ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 675 ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 676 ierr = VecZeroEntries(F);CHKERRQ(ierr); 677 ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 678 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 679 ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 680 ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 681 682 if (ctx->bctype == FVBC_OUTFLOW) { 683 for (i=xs-2; i<0; i++) { 684 for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 685 } 686 for (i=Mx; i<xs+xm+2; i++) { 687 for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 688 } 689 } 690 for (i=xs-1; i<xs+xm+1; i++) { 691 struct _LimitInfo info; 692 PetscScalar *cjmpL,*cjmpR; 693 if ((i > sm-2 && i < mf-lmbwidth+1) || (i > fm+rmbwidth-2 && i < ms+1)) { /* slow components and the first and last fast components */ 694 /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 695 ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); 696 /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 697 ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 698 cjmpL = &ctx->cjmpLR[0]; 699 cjmpR = &ctx->cjmpLR[dof]; 700 for (j=0; j<dof; j++) { 701 PetscScalar jmpL,jmpR; 702 jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 703 jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 704 for (k=0; k<dof; k++) { 705 cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 706 cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 707 } 708 } 709 /* Apply limiter to the left and right characteristic jumps */ 710 info.m = dof; 711 info.hxs = hxs; 712 info.hxm = hxm; 713 info.hxf = hxf; 714 (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope); 715 for (j=0; j<dof; j++) { 716 PetscScalar tmp = 0; 717 for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 718 slope[i*dof+j] = tmp; 719 } 720 } 721 } 722 723 for (i=xs; i<xs+xm+1; i++) { 724 PetscReal maxspeed; 725 PetscScalar *uL,*uR; 726 uL = &ctx->uLR[0]; 727 uR = &ctx->uLR[dof]; 728 if (i == sm) { /* interface between slow and medium regions */ 729 for (j=0; j<dof; j++) { 730 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxs/2; 731 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 732 } 733 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 734 if (i < xs+xm) { 735 for (j=0; j<dof; j++) f[imedium*dof+j] += ctx->flux[j]/hxm; 736 imedium++; 737 } 738 } 739 if (i > sm && i < mf-lmbwidth) { /* medium region */ 740 for (j=0; j<dof; j++) { 741 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 742 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 743 } 744 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 745 if (i > xs) { 746 for (j=0; j<dof; j++) f[(imedium-1)*dof+j] -= ctx->flux[j]/hxm; 747 } 748 if (i < xs+xm) { 749 for (j=0; j<dof; j++) f[imedium*dof+j] += ctx->flux[j]/hxm; 750 imedium++; 751 } 752 } 753 if (i == mf-lmbwidth) { /* interface between medium and fast regions */ 754 for (j=0; j<dof; j++) { 755 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 756 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 757 } 758 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 759 if (i > xs) { 760 for (j=0; j<dof; j++) f[(imedium-1)*dof+j] -= ctx->flux[j]/hxm; 761 } 762 } 763 if (i == fm+rmbwidth) { /* interface between fast and medium regions */ 764 for (j=0; j<dof; j++) { 765 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 766 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 767 } 768 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 769 if (i < xs+xm) { 770 for (j=0; j<dof; j++) f[imedium*dof+j] += ctx->flux[j]/hxm; 771 imedium++; 772 } 773 } 774 if (i > fm+rmbwidth && i < ms) { /* medium region */ 775 for (j=0; j<dof; j++) { 776 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 777 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 778 } 779 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 780 if (i > xs) { 781 for (j=0; j<dof; j++) f[(imedium-1)*dof+j] -= ctx->flux[j]/hxm; 782 } 783 if (i < xs+xm) { 784 for (j=0; j<dof; j++) f[imedium*dof+j] += ctx->flux[j]/hxm; 785 imedium++; 786 } 787 } 788 if (i == ms) { /* interface between medium and slow regions */ 789 for (j=0; j<dof; j++) { 790 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 791 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxs/2; 792 } 793 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 794 if (i > xs) { 795 for (j=0; j<dof; j++) f[(imedium-1)*dof+j] -= ctx->flux[j]/hxm; 796 } 797 } 798 } 799 ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 800 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 801 ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 802 ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 803 PetscFunctionReturn(0); 804 } 805 806 PetscErrorCode FVRHSFunctionmediumbuffer_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 807 { 808 FVCtx *ctx = (FVCtx*)vctx; 809 PetscErrorCode ierr; 810 PetscInt i,j,k,Mx,dof,xs,xm,imediumbuffer = 0,mf = ctx->mf,fm = ctx->fm,lmbwidth = ctx->lmbwidth,rmbwidth = ctx->rmbwidth; 811 PetscReal hxs,hxm,hxf; 812 PetscScalar *x,*f,*slope; 813 Vec Xloc; 814 DM da; 815 816 PetscFunctionBeginUser; 817 ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 818 ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 819 ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 820 hxs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 821 hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 822 hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 823 ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 824 ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 825 ierr = VecZeroEntries(F);CHKERRQ(ierr); 826 ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 827 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 828 ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 829 ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 830 831 if (ctx->bctype == FVBC_OUTFLOW) { 832 for (i=xs-2; i<0; i++) { 833 for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 834 } 835 for (i=Mx; i<xs+xm+2; i++) { 836 for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 837 } 838 } 839 for (i=xs-1; i<xs+xm+1; i++) { 840 struct _LimitInfo info; 841 PetscScalar *cjmpL,*cjmpR; 842 if ((i > mf-lmbwidth-2 && i < mf+1) || (i > fm-2 && i < fm+rmbwidth+1)) { 843 /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 844 ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); 845 /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 846 ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 847 cjmpL = &ctx->cjmpLR[0]; 848 cjmpR = &ctx->cjmpLR[dof]; 849 for (j=0; j<dof; j++) { 850 PetscScalar jmpL,jmpR; 851 jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 852 jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 853 for (k=0; k<dof; k++) { 854 cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 855 cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 856 } 857 } 858 /* Apply limiter to the left and right characteristic jumps */ 859 info.m = dof; 860 info.hxs = hxs; 861 info.hxm = hxm; 862 info.hxf = hxf; 863 (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope); 864 for (j=0; j<dof; j++) { 865 PetscScalar tmp = 0; 866 for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 867 slope[i*dof+j] = tmp; 868 } 869 } 870 } 871 872 for (i=xs; i<xs+xm+1; i++) { 873 PetscReal maxspeed; 874 PetscScalar *uL,*uR; 875 uL = &ctx->uLR[0]; 876 uR = &ctx->uLR[dof]; 877 if (i == mf-lmbwidth) { 878 for (j=0; j<dof; j++) { 879 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 880 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 881 } 882 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 883 if (i < xs+xm) { 884 for (j=0; j<dof; j++) f[imediumbuffer*dof+j] += ctx->flux[j]/hxm; 885 imediumbuffer++; 886 } 887 } 888 if (i > mf-lmbwidth && i < mf) { 889 for (j=0; j<dof; j++) { 890 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 891 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 892 } 893 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 894 if (i > xs) { 895 for (j=0; j<dof; j++) f[(imediumbuffer-1)*dof+j] -= ctx->flux[j]/hxm; 896 } 897 if (i < xs+xm) { 898 for (j=0; j<dof; j++) f[imediumbuffer*dof+j] += ctx->flux[j]/hxm; 899 imediumbuffer++; 900 } 901 } 902 if (i == mf) { /* interface between the medium region and the fast region */ 903 for (j=0; j<dof; j++) { 904 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 905 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 906 } 907 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 908 if (i > xs) { 909 for (j=0; j<dof; j++) f[(imediumbuffer-1)*dof+j] -= ctx->flux[j]/hxm; 910 } 911 } 912 if (i == fm) { /* interface between the fast region and the medium region */ 913 for (j=0; j<dof; j++) { 914 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 915 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 916 } 917 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 918 if (i < xs+xm) { 919 for (j=0; j<dof; j++) f[imediumbuffer*dof+j] += ctx->flux[j]/hxm; 920 imediumbuffer++; 921 } 922 } 923 if (i > fm && i < fm+rmbwidth) { 924 for (j=0; j<dof; j++) { 925 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 926 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 927 } 928 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 929 if (i > xs) { 930 for (j=0; j<dof; j++) f[(imediumbuffer-1)*dof+j] -= ctx->flux[j]/hxm; 931 } 932 if (i < xs+xm) { 933 for (j=0; j<dof; j++) f[imediumbuffer*dof+j] += ctx->flux[j]/hxm; 934 imediumbuffer++; 935 } 936 } 937 if (i == fm+rmbwidth) { 938 for (j=0; j<dof; j++) { 939 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 940 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 941 } 942 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 943 if (i > xs) { 944 for (j=0; j<dof; j++) f[(imediumbuffer-1)*dof+j] -= ctx->flux[j]/hxm; 945 } 946 } 947 } 948 ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 949 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 950 ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 951 ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 952 PetscFunctionReturn(0); 953 } 954 955 /* --------------------------------- Finite Volume Solver for fast parts ----------------------------------- */ 956 PetscErrorCode FVRHSFunctionfast_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void *vctx) 957 { 958 FVCtx *ctx = (FVCtx*)vctx; 959 PetscErrorCode ierr; 960 PetscInt i,j,k,Mx,dof,xs,xm,ifast = 0,mf = ctx->mf,fm = ctx->fm; 961 PetscReal hxs,hxm,hxf; 962 PetscScalar *x,*f,*slope; 963 Vec Xloc; 964 DM da; 965 966 PetscFunctionBeginUser; 967 ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 968 ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 969 ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 970 hxs = (ctx->xmax-ctx->xmin)/8.0/ctx->sm; 971 hxm = (ctx->xmax-ctx->xmin)/4.0/(ctx->mf-ctx->sm); 972 hxf = (ctx->xmax-ctx->xmin)/4.0/(ctx->fm-ctx->mf); 973 ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 974 ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 975 ierr = VecZeroEntries(F);CHKERRQ(ierr); 976 ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr); 977 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 978 ierr = DMDAGetArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 979 ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 980 981 if (ctx->bctype == FVBC_OUTFLOW) { 982 for (i=xs-2; i<0; i++) { 983 for (j=0; j<dof; j++) x[i*dof+j] = x[j]; 984 } 985 for (i=Mx; i<xs+xm+2; i++) { 986 for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j]; 987 } 988 } 989 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 */ 990 struct _LimitInfo info; 991 PetscScalar *cjmpL,*cjmpR; 992 if (i > mf-2 && i < fm+1) { 993 ierr = (*ctx->physics2.characteristic2)(ctx->physics2.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);CHKERRQ(ierr); 994 ierr = PetscArrayzero(ctx->cjmpLR,2*dof);CHKERRQ(ierr); 995 cjmpL = &ctx->cjmpLR[0]; 996 cjmpR = &ctx->cjmpLR[dof]; 997 for (j=0; j<dof; j++) { 998 PetscScalar jmpL,jmpR; 999 jmpL = x[(i+0)*dof+j]-x[(i-1)*dof+j]; 1000 jmpR = x[(i+1)*dof+j]-x[(i+0)*dof+j]; 1001 for (k=0; k<dof; k++) { 1002 cjmpL[k] += ctx->Rinv[k+j*dof]*jmpL; 1003 cjmpR[k] += ctx->Rinv[k+j*dof]*jmpR; 1004 } 1005 } 1006 /* Apply limiter to the left and right characteristic jumps */ 1007 info.m = dof; 1008 info.hxs = hxs; 1009 info.hxm = hxm; 1010 info.hxf = hxf; 1011 (*ctx->limit3)(&info,cjmpL,cjmpR,ctx->sm,ctx->mf,ctx->fm,ctx->ms,i,ctx->cslope); 1012 for (j=0; j<dof; j++) { 1013 PetscScalar tmp = 0; 1014 for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof]*ctx->cslope[k]; 1015 slope[i*dof+j] = tmp; 1016 } 1017 } 1018 } 1019 1020 for (i=xs; i<xs+xm+1; i++) { 1021 PetscReal maxspeed; 1022 PetscScalar *uL,*uR; 1023 uL = &ctx->uLR[0]; 1024 uR = &ctx->uLR[dof]; 1025 if (i == mf) { /* interface between medium and fast regions */ 1026 for (j=0; j<dof; j++) { 1027 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxm/2; 1028 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 1029 } 1030 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 1031 if (i < xs+xm) { 1032 for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf; 1033 ifast++; 1034 } 1035 } 1036 if (i > mf && i < fm) { /* fast region */ 1037 for (j=0; j<dof; j++) { 1038 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 1039 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxf/2; 1040 } 1041 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 1042 if (i > xs) { 1043 for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf; 1044 } 1045 if (i < xs+xm) { 1046 for (j=0; j<dof; j++) f[ifast*dof+j] += ctx->flux[j]/hxf; 1047 ifast++; 1048 } 1049 } 1050 if (i == fm) { /* interface between fast and medium regions */ 1051 for (j=0; j<dof; j++) { 1052 uL[j] = x[(i-1)*dof+j]+slope[(i-1)*dof+j]*hxf/2; 1053 uR[j] = x[(i-0)*dof+j]-slope[(i-0)*dof+j]*hxm/2; 1054 } 1055 ierr = (*ctx->physics2.riemann2)(ctx->physics2.user,dof,uL,uR,ctx->flux,&maxspeed);CHKERRQ(ierr); 1056 if (i > xs) { 1057 for (j=0; j<dof; j++) f[(ifast-1)*dof+j] -= ctx->flux[j]/hxf; 1058 } 1059 } 1060 } 1061 ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr); 1062 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 1063 ierr = DMDARestoreArray(da,PETSC_TRUE,&slope);CHKERRQ(ierr); 1064 ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 1065 PetscFunctionReturn(0); 1066 } 1067 1068 int main(int argc,char *argv[]) 1069 { 1070 char lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m"; 1071 PetscFunctionList limiters = 0,physics = 0; 1072 MPI_Comm comm; 1073 TS ts; 1074 DM da; 1075 Vec X,X0,R; 1076 FVCtx ctx; 1077 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; 1078 PetscBool view_final = PETSC_FALSE; 1079 PetscReal ptime; 1080 PetscErrorCode ierr; 1081 1082 ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr; 1083 comm = PETSC_COMM_WORLD; 1084 ierr = PetscMemzero(&ctx,sizeof(ctx));CHKERRQ(ierr); 1085 1086 /* Register limiters to be available on the command line */ 1087 ierr = PetscFunctionListAdd(&limiters,"upwind" ,Limit3_Upwind);CHKERRQ(ierr); 1088 ierr = PetscFunctionListAdd(&limiters,"lax-wendroff" ,Limit3_LaxWendroff);CHKERRQ(ierr); 1089 ierr = PetscFunctionListAdd(&limiters,"beam-warming" ,Limit3_BeamWarming);CHKERRQ(ierr); 1090 ierr = PetscFunctionListAdd(&limiters,"fromm" ,Limit3_Fromm);CHKERRQ(ierr); 1091 ierr = PetscFunctionListAdd(&limiters,"minmod" ,Limit3_Minmod);CHKERRQ(ierr); 1092 ierr = PetscFunctionListAdd(&limiters,"superbee" ,Limit3_Superbee);CHKERRQ(ierr); 1093 ierr = PetscFunctionListAdd(&limiters,"mc" ,Limit3_MC);CHKERRQ(ierr); 1094 ierr = PetscFunctionListAdd(&limiters,"koren3" ,Limit3_Koren3);CHKERRQ(ierr); 1095 1096 /* Register physical models to be available on the command line */ 1097 ierr = PetscFunctionListAdd(&physics,"advect" ,PhysicsCreate_Advect);CHKERRQ(ierr); 1098 1099 ctx.comm = comm; 1100 ctx.cfl = 0.9; 1101 ctx.bctype = FVBC_PERIODIC; 1102 ctx.xmin = -1.0; 1103 ctx.xmax = 1.0; 1104 ierr = PetscOptionsBegin(comm,NULL,"Finite Volume solver options","");CHKERRQ(ierr); 1105 ierr = PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL);CHKERRQ(ierr); 1106 ierr = PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL);CHKERRQ(ierr); 1107 ierr = PetscOptionsFList("-limit","Name of flux imiter to use","",limiters,lname,lname,sizeof(lname),NULL);CHKERRQ(ierr); 1108 ierr = PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL);CHKERRQ(ierr); 1109 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); 1110 ierr = PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL);CHKERRQ(ierr); 1111 ierr = PetscOptionsBool("-exact","Compare errors with exact solution","",ctx.exact,&ctx.exact,NULL);CHKERRQ(ierr); 1112 ierr = PetscOptionsBool("-simulation","Compare errors with reference solution","",ctx.simulation,&ctx.simulation,NULL);CHKERRQ(ierr); 1113 ierr = PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL);CHKERRQ(ierr); 1114 ierr = PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL);CHKERRQ(ierr); 1115 ierr = PetscOptionsInt("-hratio","Spacing ratio","",ctx.hratio,&ctx.hratio,NULL);CHKERRQ(ierr); 1116 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1117 1118 /* Choose the limiter from the list of registered limiters */ 1119 ierr = PetscFunctionListFind(limiters,lname,&ctx.limit3);CHKERRQ(ierr); 1120 if (!ctx.limit3) SETERRQ1(PETSC_COMM_SELF,1,"Limiter '%s' not found",lname); 1121 1122 /* Choose the physics from the list of registered models */ 1123 { 1124 PetscErrorCode (*r)(FVCtx*); 1125 ierr = PetscFunctionListFind(physics,physname,&r);CHKERRQ(ierr); 1126 if (!r) SETERRQ1(PETSC_COMM_SELF,1,"Physics '%s' not found",physname); 1127 /* Create the physics, will set the number of fields and their names */ 1128 ierr = (*r)(&ctx);CHKERRQ(ierr); 1129 } 1130 1131 /* Create a DMDA to manage the parallel grid */ 1132 ierr = DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics2.dof,2,NULL,&da);CHKERRQ(ierr); 1133 ierr = DMSetFromOptions(da);CHKERRQ(ierr); 1134 ierr = DMSetUp(da);CHKERRQ(ierr); 1135 /* Inform the DMDA of the field names provided by the physics. */ 1136 /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ 1137 for (i=0; i<ctx.physics2.dof; i++) { 1138 ierr = DMDASetFieldName(da,i,ctx.physics2.fieldname[i]);CHKERRQ(ierr); 1139 } 1140 ierr = DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);CHKERRQ(ierr); 1141 ierr = DMDAGetCorners(da,&xs,0,0,&xm,0,0);CHKERRQ(ierr); 1142 1143 /* Set coordinates of cell centers */ 1144 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); 1145 1146 /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */ 1147 ierr = PetscMalloc4(dof*dof,&ctx.R,dof*dof,&ctx.Rinv,2*dof,&ctx.cjmpLR,1*dof,&ctx.cslope);CHKERRQ(ierr); 1148 ierr = PetscMalloc3(2*dof,&ctx.uLR,dof,&ctx.flux,dof,&ctx.speeds);CHKERRQ(ierr); 1149 1150 /* Create a vector to store the solution and to save the initial state */ 1151 ierr = DMCreateGlobalVector(da,&X);CHKERRQ(ierr); 1152 ierr = VecDuplicate(X,&X0);CHKERRQ(ierr); 1153 ierr = VecDuplicate(X,&R);CHKERRQ(ierr); 1154 1155 /* create index for slow parts and fast parts, 1156 count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */ 1157 count_slow = Mx/(1+ctx.hratio)/(1+ctx.hratio); 1158 count_medium = 2*ctx.hratio*count_slow; 1159 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"); 1160 count_fast = ctx.hratio*ctx.hratio*count_slow; 1161 ctx.sm = count_slow/2; 1162 ctx.mf = ctx.sm + count_medium/2; 1163 ctx.fm = ctx.mf + count_fast; 1164 ctx.ms = ctx.fm + count_medium/2; 1165 ierr = PetscMalloc1(xm*dof,&index_slow);CHKERRQ(ierr); 1166 ierr = PetscMalloc1(xm*dof,&index_medium);CHKERRQ(ierr); 1167 ierr = PetscMalloc1(xm*dof,&index_fast);CHKERRQ(ierr); 1168 ierr = PetscMalloc1(6*dof,&index_slowbuffer);CHKERRQ(ierr); 1169 ierr = PetscMalloc1(6*dof,&index_mediumbuffer);CHKERRQ(ierr); 1170 if (((AdvectCtx*)ctx.physics2.user)->a > 0) { 1171 ctx.lsbwidth = 2; 1172 ctx.rsbwidth = 4; 1173 ctx.lmbwidth = 2; 1174 ctx.rmbwidth = 4; 1175 } else { 1176 ctx.lsbwidth = 4; 1177 ctx.rsbwidth = 2; 1178 ctx.lmbwidth = 4; 1179 ctx.rmbwidth = 2; 1180 } 1181 1182 for (i=xs; i<xs+xm; i++) { 1183 if (i < ctx.sm-ctx.lsbwidth || i > ctx.ms+ctx.rsbwidth-1) 1184 for (k=0; k<dof; k++) index_slow[islow++] = i*dof+k; 1185 else if ((i >= ctx.sm-ctx.lsbwidth && i < ctx.sm) || (i > ctx.ms-1 && i <= ctx.ms+ctx.rsbwidth-1)) 1186 for (k=0; k<dof; k++) index_slowbuffer[islowbuffer++] = i*dof+k; 1187 else if (i < ctx.mf-ctx.lmbwidth || i > ctx.fm+ctx.rmbwidth-1) 1188 for (k=0; k<dof; k++) index_medium[imedium++] = i*dof+k; 1189 else if ((i >= ctx.mf-ctx.lmbwidth && i < ctx.mf) || (i > ctx.fm-1 && i <= ctx.fm+ctx.rmbwidth-1)) 1190 for (k=0; k<dof; k++) index_mediumbuffer[imediumbuffer++] = i*dof+k; 1191 else 1192 for (k=0; k<dof; k++) index_fast[ifast++] = i*dof+k; 1193 } 1194 ierr = ISCreateGeneral(PETSC_COMM_WORLD,islow,index_slow,PETSC_COPY_VALUES,&ctx.iss);CHKERRQ(ierr); 1195 ierr = ISCreateGeneral(PETSC_COMM_WORLD,imedium,index_medium,PETSC_COPY_VALUES,&ctx.ism);CHKERRQ(ierr); 1196 ierr = ISCreateGeneral(PETSC_COMM_WORLD,ifast,index_fast,PETSC_COPY_VALUES,&ctx.isf);CHKERRQ(ierr); 1197 ierr = ISCreateGeneral(PETSC_COMM_WORLD,islowbuffer,index_slowbuffer,PETSC_COPY_VALUES,&ctx.issb);CHKERRQ(ierr); 1198 ierr = ISCreateGeneral(PETSC_COMM_WORLD,imediumbuffer,index_mediumbuffer,PETSC_COPY_VALUES,&ctx.ismb);CHKERRQ(ierr); 1199 1200 /* Create a time-stepping object */ 1201 ierr = TSCreate(comm,&ts);CHKERRQ(ierr); 1202 ierr = TSSetDM(ts,da);CHKERRQ(ierr); 1203 ierr = TSSetRHSFunction(ts,R,FVRHSFunction_3WaySplit,&ctx);CHKERRQ(ierr); 1204 ierr = TSRHSSplitSetIS(ts,"slow",ctx.iss);CHKERRQ(ierr); 1205 ierr = TSRHSSplitSetIS(ts,"medium",ctx.ism);CHKERRQ(ierr); 1206 ierr = TSRHSSplitSetIS(ts,"fast",ctx.isf);CHKERRQ(ierr); 1207 ierr = TSRHSSplitSetIS(ts,"slowbuffer",ctx.issb);CHKERRQ(ierr); 1208 ierr = TSRHSSplitSetIS(ts,"mediumbuffer",ctx.ismb);CHKERRQ(ierr); 1209 ierr = TSRHSSplitSetRHSFunction(ts,"slow",NULL,FVRHSFunctionslow_3WaySplit,&ctx);CHKERRQ(ierr); 1210 ierr = TSRHSSplitSetRHSFunction(ts,"medium",NULL,FVRHSFunctionmedium_3WaySplit,&ctx);CHKERRQ(ierr); 1211 ierr = TSRHSSplitSetRHSFunction(ts,"fast",NULL,FVRHSFunctionfast_3WaySplit,&ctx);CHKERRQ(ierr); 1212 ierr = TSRHSSplitSetRHSFunction(ts,"slowbuffer",NULL,FVRHSFunctionslowbuffer_3WaySplit,&ctx);CHKERRQ(ierr); 1213 ierr = TSRHSSplitSetRHSFunction(ts,"mediumbuffer",NULL,FVRHSFunctionmediumbuffer_3WaySplit,&ctx);CHKERRQ(ierr); 1214 1215 ierr = TSSetType(ts,TSSSP);CHKERRQ(ierr); 1216 /*ierr = TSSetType(ts,TSMPRK);CHKERRQ(ierr);*/ 1217 ierr = TSSetMaxTime(ts,10);CHKERRQ(ierr); 1218 ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 1219 1220 /* Compute initial conditions and starting time step */ 1221 ierr = FVSample_3WaySplit(&ctx,da,0,X0);CHKERRQ(ierr); 1222 ierr = FVRHSFunction_3WaySplit(ts,0,X0,X,(void*)&ctx);CHKERRQ(ierr); /* Initial function evaluation, only used to determine max speed */ 1223 ierr = VecCopy(X0,X);CHKERRQ(ierr); /* The function value was not used so we set X=X0 again */ 1224 ierr = TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt);CHKERRQ(ierr); 1225 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); /* Take runtime options */ 1226 ierr = SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1227 { 1228 PetscInt steps; 1229 PetscScalar mass_initial,mass_final,mass_difference,mass_differenceg; 1230 const PetscScalar *ptr_X,*ptr_X0; 1231 const PetscReal hs = (ctx.xmax-ctx.xmin)/4.0/count_slow; 1232 const PetscReal hm = (ctx.xmax-ctx.xmin)/2.0/count_medium; 1233 const PetscReal hf = (ctx.xmax-ctx.xmin)/4.0/count_fast; 1234 1235 ierr = TSSolve(ts,X);CHKERRQ(ierr); 1236 ierr = TSGetSolveTime(ts,&ptime);CHKERRQ(ierr); 1237 ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 1238 /* calculate the total mass at initial time and final time */ 1239 mass_initial = 0.0; 1240 mass_final = 0.0; 1241 ierr = DMDAVecGetArrayRead(da,X0,(void*)&ptr_X0);CHKERRQ(ierr); 1242 ierr = DMDAVecGetArrayRead(da,X,(void*)&ptr_X);CHKERRQ(ierr); 1243 for (i=xs;i<xs+xm;i++) { 1244 if (i < ctx.sm || i > ctx.ms-1) 1245 for (k=0; k<dof; k++) { 1246 mass_initial = mass_initial + hs*ptr_X0[i*dof+k]; 1247 mass_final = mass_final + hs*ptr_X[i*dof+k]; 1248 } 1249 else if (i < ctx.mf || i > ctx.fm-1) 1250 for (k=0; k<dof; k++) { 1251 mass_initial = mass_initial + hm*ptr_X0[i*dof+k]; 1252 mass_final = mass_final + hm*ptr_X[i*dof+k]; 1253 } 1254 else { 1255 for (k=0; k<dof; k++) { 1256 mass_initial = mass_initial + hf*ptr_X0[i*dof+k]; 1257 mass_final = mass_final + hf*ptr_X[i*dof+k]; 1258 } 1259 } 1260 } 1261 ierr = DMDAVecRestoreArrayRead(da,X0,(void*)&ptr_X0);CHKERRQ(ierr); 1262 ierr = DMDAVecRestoreArrayRead(da,X,(void*)&ptr_X);CHKERRQ(ierr); 1263 mass_difference = mass_final - mass_initial; 1264 ierr = MPI_Allreduce(&mass_difference,&mass_differenceg,1,MPIU_SCALAR,MPIU_SUM,comm);CHKERRQ(ierr); 1265 ierr = PetscPrintf(comm,"Mass difference %g\n",(double)mass_differenceg);CHKERRQ(ierr); 1266 ierr = PetscPrintf(comm,"Final time %g, steps %D\n",(double)ptime,steps);CHKERRQ(ierr); 1267 ierr = PetscPrintf(comm,"Maximum allowable stepsize according to CFL %g\n",(double)(1.0/ctx.cfl_idt));CHKERRQ(ierr); 1268 if (ctx.exact) { 1269 PetscReal nrm1=0; 1270 ierr = SolutionErrorNorms_3WaySplit(&ctx,da,ptime,X,&nrm1);CHKERRQ(ierr); 1271 ierr = PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1);CHKERRQ(ierr); 1272 } 1273 if (ctx.simulation) { 1274 PetscReal nrm1=0; 1275 PetscViewer fd; 1276 char filename[PETSC_MAX_PATH_LEN] = "binaryoutput"; 1277 Vec XR; 1278 PetscBool flg; 1279 const PetscScalar *ptr_XR; 1280 ierr = PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg);CHKERRQ(ierr); 1281 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 1282 ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&fd);CHKERRQ(ierr); 1283 ierr = VecDuplicate(X0,&XR);CHKERRQ(ierr); 1284 ierr = VecLoad(XR,fd);CHKERRQ(ierr); 1285 ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 1286 ierr = VecGetArrayRead(X,&ptr_X);CHKERRQ(ierr); 1287 ierr = VecGetArrayRead(XR,&ptr_XR);CHKERRQ(ierr); 1288 for (i=xs;i<xs+xm;i++) { 1289 if (i < ctx.sm || i > ctx.ms-1) 1290 for (k=0; k<dof; k++) nrm1 = nrm1 + hs*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]); 1291 else if (i < ctx.mf || i < ctx.fm-1) 1292 for (k=0; k<dof; k++) nrm1 = nrm1 + hm*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]); 1293 else 1294 for (k=0; k<dof; k++) nrm1 = nrm1 + hf*PetscAbs(ptr_X[i*dof+k]-ptr_XR[i*dof+k]); 1295 } 1296 ierr = VecRestoreArrayRead(X,&ptr_X);CHKERRQ(ierr); 1297 ierr = VecRestoreArrayRead(XR,&ptr_XR);CHKERRQ(ierr); 1298 ierr = PetscPrintf(comm,"Error ||x-x_e||_1 %g\n",(double)nrm1);CHKERRQ(ierr); 1299 ierr = VecDestroy(&XR);CHKERRQ(ierr); 1300 } 1301 } 1302 1303 ierr = SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1304 if (draw & 0x1) {ierr = VecView(X0,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);} 1305 if (draw & 0x2) {ierr = VecView(X,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);} 1306 if (draw & 0x4) { 1307 Vec Y; 1308 ierr = VecDuplicate(X,&Y);CHKERRQ(ierr); 1309 ierr = FVSample_3WaySplit(&ctx,da,ptime,Y);CHKERRQ(ierr); 1310 ierr = VecAYPX(Y,-1,X);CHKERRQ(ierr); 1311 ierr = VecView(Y,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr); 1312 ierr = VecDestroy(&Y);CHKERRQ(ierr); 1313 } 1314 1315 if (view_final) { 1316 PetscViewer viewer; 1317 ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer);CHKERRQ(ierr); 1318 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 1319 ierr = VecView(X,viewer);CHKERRQ(ierr); 1320 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1321 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1322 } 1323 1324 /* Clean up */ 1325 ierr = (*ctx.physics2.destroy)(ctx.physics2.user);CHKERRQ(ierr); 1326 for (i=0; i<ctx.physics2.dof; i++) {ierr = PetscFree(ctx.physics2.fieldname[i]);CHKERRQ(ierr);} 1327 ierr = PetscFree4(ctx.R,ctx.Rinv,ctx.cjmpLR,ctx.cslope);CHKERRQ(ierr); 1328 ierr = PetscFree3(ctx.uLR,ctx.flux,ctx.speeds);CHKERRQ(ierr); 1329 ierr = VecDestroy(&X);CHKERRQ(ierr); 1330 ierr = VecDestroy(&X0);CHKERRQ(ierr); 1331 ierr = VecDestroy(&R);CHKERRQ(ierr); 1332 ierr = DMDestroy(&da);CHKERRQ(ierr); 1333 ierr = TSDestroy(&ts);CHKERRQ(ierr); 1334 ierr = ISDestroy(&ctx.iss);CHKERRQ(ierr); 1335 ierr = ISDestroy(&ctx.ism);CHKERRQ(ierr); 1336 ierr = ISDestroy(&ctx.isf);CHKERRQ(ierr); 1337 ierr = ISDestroy(&ctx.issb);CHKERRQ(ierr); 1338 ierr = ISDestroy(&ctx.ismb);CHKERRQ(ierr); 1339 ierr = PetscFree(index_slow);CHKERRQ(ierr); 1340 ierr = PetscFree(index_medium);CHKERRQ(ierr); 1341 ierr = PetscFree(index_fast);CHKERRQ(ierr); 1342 ierr = PetscFree(index_slowbuffer);CHKERRQ(ierr); 1343 ierr = PetscFree(index_mediumbuffer);CHKERRQ(ierr); 1344 ierr = PetscFunctionListDestroy(&limiters);CHKERRQ(ierr); 1345 ierr = PetscFunctionListDestroy(&physics);CHKERRQ(ierr); 1346 ierr = PetscFinalize(); 1347 return ierr; 1348 } 1349 1350 /*TEST 1351 1352 build: 1353 requires: !complex c99 1354 depends: finitevolume1d.c 1355 1356 test: 1357 suffix: 1 1358 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 1359 1360 test: 1361 suffix: 2 1362 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 1363 1364 TEST*/ 1365