1 static char help[] = "Second Order TVD Finite Volume Example.\n"; 2 /*F 3 4 We use a second order TVD finite volume method to evolve a system of PDEs. Our simple upwinded residual evaluation loops 5 over all mesh faces and uses a Riemann solver to produce the flux given the face geometry and cell values, 6 \begin{equation} 7 f_i = \mathrm{riemann}(\mathrm{phys}, p_\mathrm{centroid}, \hat n, x^L, x^R) 8 \end{equation} 9 and then update the cell values given the cell volume. 10 \begin{eqnarray} 11 f^L_i &-=& \frac{f_i}{vol^L} \\ 12 f^R_i &+=& \frac{f_i}{vol^R} 13 \end{eqnarray} 14 15 As an example, we can consider the shallow water wave equation, 16 \begin{eqnarray} 17 h_t + \nabla\cdot \left( uh \right) &=& 0 \\ 18 (uh)_t + \nabla\cdot \left( u\otimes uh + \frac{g h^2}{2} I \right) &=& 0 19 \end{eqnarray} 20 where $h$ is wave height, $u$ is wave velocity, and $g$ is the acceleration due to gravity. 21 22 A representative Riemann solver for the shallow water equations is given in the PhysicsRiemann_SW() function, 23 \begin{eqnarray} 24 f^{L,R}_h &=& uh^{L,R} \cdot \hat n \\ 25 f^{L,R}_{uh} &=& \frac{f^{L,R}_h}{h^{L,R}} uh^{L,R} + g (h^{L,R})^2 \hat n \\ 26 c^{L,R} &=& \sqrt{g h^{L,R}} \\ 27 s &=& \max\left( \left|\frac{uh^L \cdot \hat n}{h^L}\right| + c^L, \left|\frac{uh^R \cdot \hat n}{h^R}\right| + c^R \right) \\ 28 f_i &=& \frac{A_\mathrm{face}}{2} \left( f^L_i + f^R_i + s \left( x^L_i - x^R_i \right) \right) 29 \end{eqnarray} 30 where $c$ is the local gravity wave speed and $f_i$ is a Rusanov flux. 31 32 The more sophisticated residual evaluation in RHSFunctionLocal_LS() uses a least-squares fit to a quadratic polynomial 33 over a neighborhood of the given element. 34 35 The mesh is read in from an ExodusII file, usually generated by Cubit. 36 F*/ 37 #include <petscdmplex.h> 38 #include <petscdmforest.h> 39 #include <petscds.h> 40 #include <petscts.h> 41 42 #define DIM 2 /* Geometric dimension */ 43 #define ALEN(a) (sizeof(a)/sizeof((a)[0])) 44 45 static PetscFunctionList PhysicsList, PhysicsRiemannList_SW; 46 47 /* Represents continuum physical equations. */ 48 typedef struct _n_Physics *Physics; 49 50 /* Physical model includes boundary conditions, initial conditions, and functionals of interest. It is 51 * discretization-independent, but its members depend on the scenario being solved. */ 52 typedef struct _n_Model *Model; 53 54 /* 'User' implements a discretization of a continuous model. */ 55 typedef struct _n_User *User; 56 typedef PetscErrorCode (*SolutionFunction)(Model,PetscReal,const PetscReal*,PetscScalar*,void*); 57 typedef PetscErrorCode (*SetUpBCFunction)(DM,PetscDS,Physics); 58 typedef PetscErrorCode (*FunctionalFunction)(Model,PetscReal,const PetscReal*,const PetscScalar*,PetscReal*,void*); 59 typedef PetscErrorCode (*SetupFields)(Physics,PetscSection); 60 static PetscErrorCode ModelSolutionSetDefault(Model,SolutionFunction,void*); 61 static PetscErrorCode ModelFunctionalRegister(Model,const char*,PetscInt*,FunctionalFunction,void*); 62 static PetscErrorCode OutputVTK(DM,const char*,PetscViewer*); 63 64 struct FieldDescription { 65 const char *name; 66 PetscInt dof; 67 }; 68 69 typedef struct _n_FunctionalLink *FunctionalLink; 70 struct _n_FunctionalLink { 71 char *name; 72 FunctionalFunction func; 73 void *ctx; 74 PetscInt offset; 75 FunctionalLink next; 76 }; 77 78 struct _n_Physics { 79 PetscRiemannFunc riemann; 80 PetscInt dof; /* number of degrees of freedom per cell */ 81 PetscReal maxspeed; /* kludge to pick initial time step, need to add monitoring and step control */ 82 void *data; 83 PetscInt nfields; 84 const struct FieldDescription *field_desc; 85 }; 86 87 struct _n_Model { 88 MPI_Comm comm; /* Does not do collective communicaton, but some error conditions can be collective */ 89 Physics physics; 90 FunctionalLink functionalRegistry; 91 PetscInt maxComputed; 92 PetscInt numMonitored; 93 FunctionalLink *functionalMonitored; 94 PetscInt numCall; 95 FunctionalLink *functionalCall; 96 SolutionFunction solution; 97 SetUpBCFunction setupbc; 98 void *solutionctx; 99 PetscReal maxspeed; /* estimate of global maximum speed (for CFL calculation) */ 100 PetscReal bounds[2*DIM]; 101 PetscErrorCode (*errorIndicator)(PetscInt, PetscReal, PetscInt, const PetscScalar[], const PetscScalar[], PetscReal *, void *); 102 void *errorCtx; 103 }; 104 105 struct _n_User { 106 PetscInt vtkInterval; /* For monitor */ 107 char outputBasename[PETSC_MAX_PATH_LEN]; /* Basename for output files */ 108 PetscInt monitorStepOffset; 109 Model model; 110 PetscBool vtkmon; 111 }; 112 113 PETSC_STATIC_INLINE PetscReal DotDIMReal(const PetscReal *x,const PetscReal *y) 114 { 115 PetscInt i; 116 PetscReal prod=0.0; 117 118 for (i=0; i<DIM; i++) prod += x[i]*y[i]; 119 return prod; 120 } 121 PETSC_STATIC_INLINE PetscReal NormDIM(const PetscReal *x) { return PetscSqrtReal(PetscAbsReal(DotDIMReal(x,x))); } 122 123 PETSC_STATIC_INLINE PetscReal Dot2Real(const PetscReal *x,const PetscReal *y) { return x[0]*y[0] + x[1]*y[1];} 124 PETSC_STATIC_INLINE PetscReal Norm2Real(const PetscReal *x) { return PetscSqrtReal(PetscAbsReal(Dot2Real(x,x)));} 125 PETSC_STATIC_INLINE void Normalize2Real(PetscReal *x) { PetscReal a = 1./Norm2Real(x); x[0] *= a; x[1] *= a; } 126 PETSC_STATIC_INLINE void Waxpy2Real(PetscReal a,const PetscReal *x,const PetscReal *y,PetscReal *w) { w[0] = a*x[0] + y[0]; w[1] = a*x[1] + y[1]; } 127 PETSC_STATIC_INLINE void Scale2Real(PetscReal a,const PetscReal *x,PetscReal *y) { y[0] = a*x[0]; y[1] = a*x[1]; } 128 129 /******************* Advect ********************/ 130 typedef enum {ADVECT_SOL_TILTED,ADVECT_SOL_BUMP,ADVECT_SOL_BUMP_CAVITY} AdvectSolType; 131 static const char *const AdvectSolTypes[] = {"TILTED","BUMP","BUMP_CAVITY","AdvectSolType","ADVECT_SOL_",0}; 132 typedef enum {ADVECT_SOL_BUMP_CONE,ADVECT_SOL_BUMP_COS} AdvectSolBumpType; 133 static const char *const AdvectSolBumpTypes[] = {"CONE","COS","AdvectSolBumpType","ADVECT_SOL_BUMP_",0}; 134 135 typedef struct { 136 PetscReal wind[DIM]; 137 } Physics_Advect_Tilted; 138 typedef struct { 139 PetscReal center[DIM]; 140 PetscReal radius; 141 AdvectSolBumpType type; 142 } Physics_Advect_Bump; 143 144 typedef struct { 145 PetscReal inflowState; 146 AdvectSolType soltype; 147 union { 148 Physics_Advect_Tilted tilted; 149 Physics_Advect_Bump bump; 150 } sol; 151 struct { 152 PetscInt Solution; 153 PetscInt Error; 154 } functional; 155 } Physics_Advect; 156 157 static const struct FieldDescription PhysicsFields_Advect[] = {{"U",1},{NULL,0}}; 158 159 static PetscErrorCode PhysicsBoundary_Advect_Inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 160 { 161 Physics phys = (Physics)ctx; 162 Physics_Advect *advect = (Physics_Advect*)phys->data; 163 164 PetscFunctionBeginUser; 165 xG[0] = advect->inflowState; 166 PetscFunctionReturn(0); 167 } 168 169 static PetscErrorCode PhysicsBoundary_Advect_Outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 170 { 171 PetscFunctionBeginUser; 172 xG[0] = xI[0]; 173 PetscFunctionReturn(0); 174 } 175 176 static void PhysicsRiemann_Advect(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys) 177 { 178 Physics_Advect *advect = (Physics_Advect*)phys->data; 179 PetscReal wind[DIM],wn; 180 181 switch (advect->soltype) { 182 case ADVECT_SOL_TILTED: { 183 Physics_Advect_Tilted *tilted = &advect->sol.tilted; 184 wind[0] = tilted->wind[0]; 185 wind[1] = tilted->wind[1]; 186 } break; 187 case ADVECT_SOL_BUMP: 188 wind[0] = -qp[1]; 189 wind[1] = qp[0]; 190 break; 191 case ADVECT_SOL_BUMP_CAVITY: 192 { 193 PetscInt i; 194 PetscReal comp2[3] = {0.,0.,0.}, rad2; 195 196 rad2 = 0.; 197 for (i = 0; i < dim; i++) { 198 comp2[i] = qp[i] * qp[i]; 199 rad2 += comp2[i]; 200 } 201 202 wind[0] = -qp[1]; 203 wind[1] = qp[0]; 204 if (rad2 > 1.) { 205 PetscInt maxI = 0; 206 PetscReal maxComp2 = comp2[0]; 207 208 for (i = 1; i < dim; i++) { 209 if (comp2[i] > maxComp2) { 210 maxI = i; 211 maxComp2 = comp2[i]; 212 } 213 } 214 wind[maxI] = 0.; 215 } 216 } 217 break; 218 default: 219 { 220 PetscInt i; 221 for (i = 0; i < DIM; ++i) wind[i] = 0.0; 222 } 223 /* default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for solution type %s",AdvectSolBumpTypes[advect->soltype]); */ 224 } 225 wn = Dot2Real(wind, n); 226 flux[0] = (wn > 0 ? xL[0] : xR[0]) * wn; 227 } 228 229 static PetscErrorCode PhysicsSolution_Advect(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx) 230 { 231 Physics phys = (Physics)ctx; 232 Physics_Advect *advect = (Physics_Advect*)phys->data; 233 234 PetscFunctionBeginUser; 235 switch (advect->soltype) { 236 case ADVECT_SOL_TILTED: { 237 PetscReal x0[DIM]; 238 Physics_Advect_Tilted *tilted = &advect->sol.tilted; 239 Waxpy2Real(-time,tilted->wind,x,x0); 240 if (x0[1] > 0) u[0] = 1.*x[0] + 3.*x[1]; 241 else u[0] = advect->inflowState; 242 } break; 243 case ADVECT_SOL_BUMP_CAVITY: 244 case ADVECT_SOL_BUMP: { 245 Physics_Advect_Bump *bump = &advect->sol.bump; 246 PetscReal x0[DIM],v[DIM],r,cost,sint; 247 cost = PetscCosReal(time); 248 sint = PetscSinReal(time); 249 x0[0] = cost*x[0] + sint*x[1]; 250 x0[1] = -sint*x[0] + cost*x[1]; 251 Waxpy2Real(-1,bump->center,x0,v); 252 r = Norm2Real(v); 253 switch (bump->type) { 254 case ADVECT_SOL_BUMP_CONE: 255 u[0] = PetscMax(1 - r/bump->radius,0); 256 break; 257 case ADVECT_SOL_BUMP_COS: 258 u[0] = 0.5 + 0.5*PetscCosReal(PetscMin(r/bump->radius,1)*PETSC_PI); 259 break; 260 } 261 } break; 262 default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown solution type"); 263 } 264 PetscFunctionReturn(0); 265 } 266 267 static PetscErrorCode PhysicsFunctional_Advect(Model mod,PetscReal time,const PetscReal *x,const PetscScalar *y,PetscReal *f,void *ctx) 268 { 269 Physics phys = (Physics)ctx; 270 Physics_Advect *advect = (Physics_Advect*)phys->data; 271 PetscScalar yexact[1] = {0.0}; 272 PetscErrorCode ierr; 273 274 PetscFunctionBeginUser; 275 ierr = PhysicsSolution_Advect(mod,time,x,yexact,phys);CHKERRQ(ierr); 276 f[advect->functional.Solution] = PetscRealPart(y[0]); 277 f[advect->functional.Error] = PetscAbsScalar(y[0]-yexact[0]); 278 PetscFunctionReturn(0); 279 } 280 281 static PetscErrorCode SetUpBC_Advect(DM dm, PetscDS prob, Physics phys) 282 { 283 PetscErrorCode ierr; 284 const PetscInt inflowids[] = {100,200,300},outflowids[] = {101}; 285 DMLabel label; 286 287 PetscFunctionBeginUser; 288 /* Register "canned" boundary conditions and defaults for where to apply. */ 289 ierr = DMGetLabel(dm, "Face Sets", &label);CHKERRQ(ierr); 290 ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "inflow", label, ALEN(inflowids), inflowids, 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Advect_Inflow, NULL, phys, NULL);CHKERRQ(ierr); 291 ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "outflow", label, ALEN(outflowids), outflowids, 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Advect_Outflow, NULL, phys, NULL);CHKERRQ(ierr); 292 PetscFunctionReturn(0); 293 } 294 295 static PetscErrorCode PhysicsCreate_Advect(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject) 296 { 297 Physics_Advect *advect; 298 PetscErrorCode ierr; 299 300 PetscFunctionBeginUser; 301 phys->field_desc = PhysicsFields_Advect; 302 phys->riemann = (PetscRiemannFunc)PhysicsRiemann_Advect; 303 ierr = PetscNew(&advect);CHKERRQ(ierr); 304 phys->data = advect; 305 mod->setupbc = SetUpBC_Advect; 306 307 ierr = PetscOptionsHead(PetscOptionsObject,"Advect options");CHKERRQ(ierr); 308 { 309 PetscInt two = 2,dof = 1; 310 advect->soltype = ADVECT_SOL_TILTED; 311 ierr = PetscOptionsEnum("-advect_sol_type","solution type","",AdvectSolTypes,(PetscEnum)advect->soltype,(PetscEnum*)&advect->soltype,NULL);CHKERRQ(ierr); 312 switch (advect->soltype) { 313 case ADVECT_SOL_TILTED: { 314 Physics_Advect_Tilted *tilted = &advect->sol.tilted; 315 two = 2; 316 tilted->wind[0] = 0.0; 317 tilted->wind[1] = 1.0; 318 ierr = PetscOptionsRealArray("-advect_tilted_wind","background wind vx,vy","",tilted->wind,&two,NULL);CHKERRQ(ierr); 319 advect->inflowState = -2.0; 320 ierr = PetscOptionsRealArray("-advect_tilted_inflow","Inflow state","",&advect->inflowState,&dof,NULL);CHKERRQ(ierr); 321 phys->maxspeed = Norm2Real(tilted->wind); 322 } break; 323 case ADVECT_SOL_BUMP_CAVITY: 324 case ADVECT_SOL_BUMP: { 325 Physics_Advect_Bump *bump = &advect->sol.bump; 326 two = 2; 327 bump->center[0] = 2.; 328 bump->center[1] = 0.; 329 ierr = PetscOptionsRealArray("-advect_bump_center","location of center of bump x,y","",bump->center,&two,NULL);CHKERRQ(ierr); 330 bump->radius = 0.9; 331 ierr = PetscOptionsReal("-advect_bump_radius","radius of bump","",bump->radius,&bump->radius,NULL);CHKERRQ(ierr); 332 bump->type = ADVECT_SOL_BUMP_CONE; 333 ierr = PetscOptionsEnum("-advect_bump_type","type of bump","",AdvectSolBumpTypes,(PetscEnum)bump->type,(PetscEnum*)&bump->type,NULL);CHKERRQ(ierr); 334 phys->maxspeed = 3.; /* radius of mesh, kludge */ 335 } break; 336 } 337 } 338 ierr = PetscOptionsTail();CHKERRQ(ierr); 339 /* Initial/transient solution with default boundary conditions */ 340 ierr = ModelSolutionSetDefault(mod,PhysicsSolution_Advect,phys);CHKERRQ(ierr); 341 /* Register "canned" functionals */ 342 ierr = ModelFunctionalRegister(mod,"Solution",&advect->functional.Solution,PhysicsFunctional_Advect,phys);CHKERRQ(ierr); 343 ierr = ModelFunctionalRegister(mod,"Error",&advect->functional.Error,PhysicsFunctional_Advect,phys);CHKERRQ(ierr); 344 PetscFunctionReturn(0); 345 } 346 347 /******************* Shallow Water ********************/ 348 typedef struct { 349 PetscReal gravity; 350 PetscReal boundaryHeight; 351 struct { 352 PetscInt Height; 353 PetscInt Speed; 354 PetscInt Energy; 355 } functional; 356 } Physics_SW; 357 typedef struct { 358 PetscReal h; 359 PetscReal uh[DIM]; 360 } SWNode; 361 typedef union { 362 SWNode swnode; 363 PetscReal vals[DIM+1]; 364 } SWNodeUnion; 365 366 static const struct FieldDescription PhysicsFields_SW[] = {{"Height",1},{"Momentum",DIM},{NULL,0}}; 367 368 /* 369 * h_t + div(uh) = 0 370 * (uh)_t + div (u\otimes uh + g h^2 / 2 I) = 0 371 * 372 * */ 373 static PetscErrorCode SWFlux(Physics phys,const PetscReal *n,const SWNode *x,SWNode *f) 374 { 375 Physics_SW *sw = (Physics_SW*)phys->data; 376 PetscReal uhn,u[DIM]; 377 PetscInt i; 378 379 PetscFunctionBeginUser; 380 Scale2Real(1./x->h,x->uh,u); 381 uhn = x->uh[0] * n[0] + x->uh[1] * n[1]; 382 f->h = uhn; 383 for (i=0; i<DIM; i++) f->uh[i] = u[i] * uhn + sw->gravity * PetscSqr(x->h) * n[i]; 384 PetscFunctionReturn(0); 385 } 386 387 static PetscErrorCode PhysicsBoundary_SW_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 388 { 389 PetscFunctionBeginUser; 390 xG[0] = xI[0]; 391 xG[1] = -xI[1]; 392 xG[2] = -xI[2]; 393 PetscFunctionReturn(0); 394 } 395 396 static void PhysicsRiemann_SW_HLL(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys) 397 { 398 Physics_SW *sw = (Physics_SW *) phys->data; 399 PetscReal aL, aR; 400 PetscReal nn[DIM]; 401 #if !defined(PETSC_USE_COMPLEX) 402 const SWNode *uL = (const SWNode *) xL, *uR = (const SWNode *) xR; 403 #else 404 SWNodeUnion uLreal, uRreal; 405 const SWNode *uL = &uLreal.swnode; 406 const SWNode *uR = &uRreal.swnode; 407 #endif 408 SWNodeUnion fL, fR; 409 PetscInt i; 410 PetscReal zero = 0.; 411 412 #if defined(PETSC_USE_COMPLEX) 413 uLreal.swnode.h = 0; uRreal.swnode.h = 0; 414 for (i = 0; i < 1+dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]); 415 for (i = 0; i < 1+dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]); 416 #endif 417 if (uL->h <= 0 || uR->h <= 0) { 418 for (i = 0; i < 1 + dim; i++) flux[i] = zero; 419 return; 420 } /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative"); */ 421 nn[0] = n[0]; 422 nn[1] = n[1]; 423 Normalize2Real(nn); 424 SWFlux(phys, nn, uL, &(fL.swnode)); 425 SWFlux(phys, nn, uR, &(fR.swnode)); 426 /* gravity wave speed */ 427 aL = PetscSqrtReal(sw->gravity * uL->h); 428 aR = PetscSqrtReal(sw->gravity * uR->h); 429 // Defining u_tilda and v_tilda as u and v 430 PetscReal u_L, u_R; 431 u_L = Dot2Real(uL->uh,nn)/uL->h; 432 u_R = Dot2Real(uR->uh,nn)/uR->h; 433 PetscReal sL, sR; 434 sL = PetscMin(u_L - aL, u_R - aR); 435 sR = PetscMax(u_L + aL, u_R + aR); 436 if (sL > zero) { 437 for (i = 0; i < dim + 1; i++) { 438 flux[i] = fL.vals[i] * Norm2Real(n); 439 } 440 } else if (sR < zero) { 441 for (i = 0; i < dim + 1; i++) { 442 flux[i] = fR.vals[i] * Norm2Real(n); 443 } 444 } else { 445 for (i = 0; i < dim + 1; i++) { 446 flux[i] = ((sR * fL.vals[i] - sL * fR.vals[i] + sR * sL * (xR[i] - xL[i])) / (sR - sL)) * Norm2Real(n); 447 } 448 } 449 } 450 451 static void PhysicsRiemann_SW_Rusanov(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys) 452 { 453 Physics_SW *sw = (Physics_SW*)phys->data; 454 PetscReal cL,cR,speed; 455 PetscReal nn[DIM]; 456 #if !defined(PETSC_USE_COMPLEX) 457 const SWNode *uL = (const SWNode*)xL,*uR = (const SWNode*)xR; 458 #else 459 SWNodeUnion uLreal, uRreal; 460 const SWNode *uL = &uLreal.swnode; 461 const SWNode *uR = &uRreal.swnode; 462 #endif 463 SWNodeUnion fL,fR; 464 PetscInt i; 465 PetscReal zero=0.; 466 467 #if defined(PETSC_USE_COMPLEX) 468 uLreal.swnode.h = 0; uRreal.swnode.h = 0; 469 for (i = 0; i < 1+dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]); 470 for (i = 0; i < 1+dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]); 471 #endif 472 if (uL->h < 0 || uR->h < 0) {for (i=0; i<1+dim; i++) flux[i] = zero/zero; return;} /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative"); */ 473 nn[0] = n[0]; 474 nn[1] = n[1]; 475 Normalize2Real(nn); 476 SWFlux(phys,nn,uL,&(fL.swnode)); 477 SWFlux(phys,nn,uR,&(fR.swnode)); 478 cL = PetscSqrtReal(sw->gravity*uL->h); 479 cR = PetscSqrtReal(sw->gravity*uR->h); /* gravity wave speed */ 480 speed = PetscMax(PetscAbsReal(Dot2Real(uL->uh,nn)/uL->h) + cL,PetscAbsReal(Dot2Real(uR->uh,nn)/uR->h) + cR); 481 for (i=0; i<1+dim; i++) flux[i] = (0.5*(fL.vals[i] + fR.vals[i]) + 0.5*speed*(xL[i] - xR[i])) * Norm2Real(n); 482 } 483 484 static PetscErrorCode PhysicsSolution_SW(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx) 485 { 486 PetscReal dx[2],r,sigma; 487 488 PetscFunctionBeginUser; 489 if (time != 0.0) SETERRQ1(mod->comm,PETSC_ERR_SUP,"No solution known for time %g",(double)time); 490 dx[0] = x[0] - 1.5; 491 dx[1] = x[1] - 1.0; 492 r = Norm2Real(dx); 493 sigma = 0.5; 494 u[0] = 1 + 2*PetscExpReal(-PetscSqr(r)/(2*PetscSqr(sigma))); 495 u[1] = 0.0; 496 u[2] = 0.0; 497 PetscFunctionReturn(0); 498 } 499 500 static PetscErrorCode PhysicsFunctional_SW(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx) 501 { 502 Physics phys = (Physics)ctx; 503 Physics_SW *sw = (Physics_SW*)phys->data; 504 const SWNode *x = (const SWNode*)xx; 505 PetscReal u[2]; 506 PetscReal h; 507 508 PetscFunctionBeginUser; 509 h = x->h; 510 Scale2Real(1./x->h,x->uh,u); 511 f[sw->functional.Height] = h; 512 f[sw->functional.Speed] = Norm2Real(u) + PetscSqrtReal(sw->gravity*h); 513 f[sw->functional.Energy] = 0.5*(Dot2Real(x->uh,u) + sw->gravity*PetscSqr(h)); 514 PetscFunctionReturn(0); 515 } 516 517 static PetscErrorCode SetUpBC_SW(DM dm, PetscDS prob,Physics phys) 518 { 519 PetscErrorCode ierr; 520 const PetscInt wallids[] = {100,101,200,300}; 521 DMLabel label; 522 523 PetscFunctionBeginUser; 524 ierr = DMGetLabel(dm, "Face Sets", &label);CHKERRQ(ierr); 525 ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, ALEN(wallids), wallids, 0, 0, NULL, (void (*)(void)) PhysicsBoundary_SW_Wall, NULL, phys, NULL);CHKERRQ(ierr); 526 PetscFunctionReturn(0); 527 } 528 529 static PetscErrorCode PhysicsCreate_SW(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject) 530 { 531 Physics_SW *sw; 532 char sw_riemann[64] = "rusanov"; 533 PetscErrorCode ierr; 534 535 PetscFunctionBeginUser; 536 phys->field_desc = PhysicsFields_SW; 537 ierr = PetscNew(&sw);CHKERRQ(ierr); 538 phys->data = sw; 539 mod->setupbc = SetUpBC_SW; 540 541 PetscFunctionListAdd(&PhysicsRiemannList_SW, "rusanov", PhysicsRiemann_SW_Rusanov); 542 PetscFunctionListAdd(&PhysicsRiemannList_SW, "hll", PhysicsRiemann_SW_HLL); 543 544 ierr = PetscOptionsHead(PetscOptionsObject,"SW options");CHKERRQ(ierr); 545 { 546 void (*PhysicsRiemann_SW)(PetscInt, PetscInt, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, PetscInt, const PetscScalar, PetscScalar *, Physics); 547 sw->gravity = 1.0; 548 ierr = PetscOptionsReal("-sw_gravity","Gravitational constant","",sw->gravity,&sw->gravity,NULL);CHKERRQ(ierr); 549 ierr = PetscOptionsFList("-sw_riemann","Riemann solver","",PhysicsRiemannList_SW,sw_riemann,sw_riemann,sizeof sw_riemann,NULL);CHKERRQ(ierr); 550 ierr = PetscFunctionListFind(PhysicsRiemannList_SW,sw_riemann,&PhysicsRiemann_SW);CHKERRQ(ierr); 551 phys->riemann = (PetscRiemannFunc) PhysicsRiemann_SW; 552 } 553 ierr = PetscOptionsTail();CHKERRQ(ierr); 554 phys->maxspeed = PetscSqrtReal(2.0*sw->gravity); /* Mach 1 for depth of 2 */ 555 556 ierr = ModelSolutionSetDefault(mod,PhysicsSolution_SW,phys);CHKERRQ(ierr); 557 ierr = ModelFunctionalRegister(mod,"Height",&sw->functional.Height,PhysicsFunctional_SW,phys);CHKERRQ(ierr); 558 ierr = ModelFunctionalRegister(mod,"Speed",&sw->functional.Speed,PhysicsFunctional_SW,phys);CHKERRQ(ierr); 559 ierr = ModelFunctionalRegister(mod,"Energy",&sw->functional.Energy,PhysicsFunctional_SW,phys);CHKERRQ(ierr); 560 561 PetscFunctionReturn(0); 562 } 563 564 /******************* Euler Density Shock (EULER_IV_SHOCK,EULER_SS_SHOCK) ********************/ 565 /* An initial-value and self-similar solutions of the compressible Euler equations */ 566 /* Ravi Samtaney and D. I. Pullin */ 567 /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */ 568 typedef enum {EULER_PAR_GAMMA,EULER_PAR_RHOR,EULER_PAR_AMACH,EULER_PAR_ITANA,EULER_PAR_SIZE} EulerParamIdx; 569 typedef enum {EULER_IV_SHOCK,EULER_SS_SHOCK,EULER_SHOCK_TUBE,EULER_LINEAR_WAVE} EulerType; 570 typedef struct { 571 PetscReal r; 572 PetscReal ru[DIM]; 573 PetscReal E; 574 } EulerNode; 575 typedef union { 576 EulerNode eulernode; 577 PetscReal vals[DIM+2]; 578 } EulerNodeUnion; 579 typedef PetscErrorCode (*EquationOfState)(const PetscReal*, const EulerNode*, PetscReal*); 580 typedef struct { 581 EulerType type; 582 PetscReal pars[EULER_PAR_SIZE]; 583 EquationOfState sound; 584 struct { 585 PetscInt Density; 586 PetscInt Momentum; 587 PetscInt Energy; 588 PetscInt Pressure; 589 PetscInt Speed; 590 } monitor; 591 } Physics_Euler; 592 593 static const struct FieldDescription PhysicsFields_Euler[] = {{"Density",1},{"Momentum",DIM},{"Energy",1},{NULL,0}}; 594 595 /* initial condition */ 596 int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx); 597 static PetscErrorCode PhysicsSolution_Euler(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx) 598 { 599 PetscInt i; 600 Physics phys = (Physics)ctx; 601 Physics_Euler *eu = (Physics_Euler*)phys->data; 602 EulerNode *uu = (EulerNode*)u; 603 PetscReal p0,gamma,c; 604 PetscFunctionBeginUser; 605 if (time != 0.0) SETERRQ1(mod->comm,PETSC_ERR_SUP,"No solution known for time %g",(double)time); 606 607 for (i=0; i<DIM; i++) uu->ru[i] = 0.0; /* zero out initial velocity */ 608 /* set E and rho */ 609 gamma = eu->pars[EULER_PAR_GAMMA]; 610 611 if (eu->type==EULER_IV_SHOCK || eu->type==EULER_SS_SHOCK) { 612 /******************* Euler Density Shock ********************/ 613 /* On initial-value and self-similar solutions of the compressible Euler equations */ 614 /* Ravi Samtaney and D. I. Pullin */ 615 /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */ 616 /* initial conditions 1: left of shock, 0: left of discontinuity 2: right of discontinuity, */ 617 p0 = 1.; 618 if (x[0] < 0.0 + x[1]*eu->pars[EULER_PAR_ITANA]) { 619 if (x[0] < mod->bounds[0]*0.5) { /* left of shock (1) */ 620 PetscReal amach,rho,press,gas1,p1; 621 amach = eu->pars[EULER_PAR_AMACH]; 622 rho = 1.; 623 press = p0; 624 p1 = press*(1.0+2.0*gamma/(gamma+1.0)*(amach*amach-1.0)); 625 gas1 = (gamma-1.0)/(gamma+1.0); 626 uu->r = rho*(p1/press+gas1)/(gas1*p1/press+1.0); 627 uu->ru[0] = ((uu->r - rho)*PetscSqrtReal(gamma*press/rho)*amach); 628 uu->E = p1/(gamma-1.0) + .5/uu->r*uu->ru[0]*uu->ru[0]; 629 } 630 else { /* left of discontinuity (0) */ 631 uu->r = 1.; /* rho = 1 */ 632 uu->E = p0/(gamma-1.0); 633 } 634 } 635 else { /* right of discontinuity (2) */ 636 uu->r = eu->pars[EULER_PAR_RHOR]; 637 uu->E = p0/(gamma-1.0); 638 } 639 } 640 else if (eu->type==EULER_SHOCK_TUBE) { 641 /* For (x<x0) set (rho,u,p)=(8,0,10) and for (x>x0) set (rho,u,p)=(1,0,1). Choose x0 to the midpoint of the domain in the x-direction. */ 642 if (x[0] < 0.0) { 643 uu->r = 8.; 644 uu->E = 10./(gamma-1.); 645 } 646 else { 647 uu->r = 1.; 648 uu->E = 1./(gamma-1.); 649 } 650 } 651 else if (eu->type==EULER_LINEAR_WAVE) { 652 initLinearWave( uu, gamma, x, mod->bounds[1] - mod->bounds[0]); 653 } 654 else SETERRQ1(mod->comm,PETSC_ERR_SUP,"Unknown type %d",eu->type); 655 656 /* set phys->maxspeed: (mod->maxspeed = phys->maxspeed) in main; */ 657 eu->sound(&gamma,uu,&c); 658 c = (uu->ru[0]/uu->r) + c; 659 if (c > phys->maxspeed) phys->maxspeed = c; 660 661 PetscFunctionReturn(0); 662 } 663 664 static PetscErrorCode Pressure_PG(const PetscReal gamma,const EulerNode *x,PetscReal *p) 665 { 666 PetscReal ru2; 667 668 PetscFunctionBeginUser; 669 ru2 = DotDIMReal(x->ru,x->ru); 670 (*p)=(x->E - 0.5*ru2/x->r)*(gamma - 1.0); /* (E - rho V^2/2)(gamma-1) = e rho (gamma-1) */ 671 PetscFunctionReturn(0); 672 } 673 674 static PetscErrorCode SpeedOfSound_PG(const PetscReal *gamma, const EulerNode *x, PetscReal *c) 675 { 676 PetscReal p; 677 678 PetscFunctionBeginUser; 679 Pressure_PG(*gamma,x,&p); 680 if (p<0.) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"negative pressure time %g -- NEED TO FIX!!!!!!",(double) p); 681 /* pars[EULER_PAR_GAMMA] = heat capacity ratio */ 682 (*c)=PetscSqrtReal(*gamma * p / x->r); 683 PetscFunctionReturn(0); 684 } 685 686 /* 687 * x = (rho,rho*(u_1),...,rho*e)^T 688 * x_t+div(f_1(x))+...+div(f_DIM(x)) = 0 689 * 690 * f_i(x) = u_i*x+(0,0,...,p,...,p*u_i)^T 691 * 692 */ 693 static PetscErrorCode EulerFlux(Physics phys,const PetscReal *n,const EulerNode *x,EulerNode *f) 694 { 695 Physics_Euler *eu = (Physics_Euler*)phys->data; 696 PetscReal nu,p; 697 PetscInt i; 698 699 PetscFunctionBeginUser; 700 Pressure_PG(eu->pars[EULER_PAR_GAMMA],x,&p); 701 nu = DotDIMReal(x->ru,n); 702 f->r = nu; /* A rho u */ 703 nu /= x->r; /* A u */ 704 for (i=0; i<DIM; i++) f->ru[i] = nu * x->ru[i] + n[i]*p; /* r u^2 + p */ 705 f->E = nu * (x->E + p); /* u(e+p) */ 706 PetscFunctionReturn(0); 707 } 708 709 /* PetscReal* => EulerNode* conversion */ 710 static PetscErrorCode PhysicsBoundary_Euler_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *a_xI, PetscScalar *a_xG, void *ctx) 711 { 712 PetscInt i; 713 const EulerNode *xI = (const EulerNode*)a_xI; 714 EulerNode *xG = (EulerNode*)a_xG; 715 Physics phys = (Physics)ctx; 716 Physics_Euler *eu = (Physics_Euler*)phys->data; 717 PetscFunctionBeginUser; 718 xG->r = xI->r; /* ghost cell density - same */ 719 xG->E = xI->E; /* ghost cell energy - same */ 720 if (n[1] != 0.) { /* top and bottom */ 721 xG->ru[0] = xI->ru[0]; /* copy tang to wall */ 722 xG->ru[1] = -xI->ru[1]; /* reflect perp to t/b wall */ 723 } 724 else { /* sides */ 725 for (i=0; i<DIM; i++) xG->ru[i] = xI->ru[i]; /* copy */ 726 } 727 if (eu->type == EULER_LINEAR_WAVE) { /* debug */ 728 #if 0 729 PetscPrintf(PETSC_COMM_WORLD,"%s coord=%g,%g\n",PETSC_FUNCTION_NAME,c[0],c[1]); 730 #endif 731 } 732 PetscFunctionReturn(0); 733 } 734 int godunovflux( const PetscScalar *ul, const PetscScalar *ur, PetscScalar *flux, const PetscReal *nn, const int *ndim, const PetscReal *gamma); 735 /* PetscReal* => EulerNode* conversion */ 736 static void PhysicsRiemann_Euler_Godunov( PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, 737 const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys) 738 { 739 Physics_Euler *eu = (Physics_Euler*)phys->data; 740 PetscReal cL,cR,speed,velL,velR,nn[DIM],s2; 741 PetscInt i; 742 PetscErrorCode ierr; 743 PetscFunctionBeginUser; 744 745 for (i=0,s2=0.; i<DIM; i++) { 746 nn[i] = n[i]; 747 s2 += nn[i]*nn[i]; 748 } 749 s2 = PetscSqrtReal(s2); /* |n|_2 = sum(n^2)^1/2 */ 750 for (i=0.; i<DIM; i++) nn[i] /= s2; 751 if (0) { /* Rusanov */ 752 const EulerNode *uL = (const EulerNode*)xL,*uR = (const EulerNode*)xR; 753 EulerNodeUnion fL,fR; 754 EulerFlux(phys,nn,uL,&(fL.eulernode)); 755 EulerFlux(phys,nn,uR,&(fR.eulernode)); 756 ierr = eu->sound(&eu->pars[EULER_PAR_GAMMA],uL,&cL);if (ierr) exit(13); 757 ierr = eu->sound(&eu->pars[EULER_PAR_GAMMA],uR,&cR);if (ierr) exit(14); 758 velL = DotDIMReal(uL->ru,nn)/uL->r; 759 velR = DotDIMReal(uR->ru,nn)/uR->r; 760 speed = PetscMax(velR + cR, velL + cL); 761 for (i=0; i<2+dim; i++) flux[i] = 0.5*((fL.vals[i]+fR.vals[i]) + speed*(xL[i] - xR[i]))*s2; 762 } 763 else { 764 int dim = DIM; 765 /* int iwave = */ 766 godunovflux(xL, xR, flux, nn, &dim, &eu->pars[EULER_PAR_GAMMA]); 767 for (i=0; i<2+dim; i++) flux[i] *= s2; 768 } 769 PetscFunctionReturnVoid(); 770 } 771 772 static PetscErrorCode PhysicsFunctional_Euler(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx) 773 { 774 Physics phys = (Physics)ctx; 775 Physics_Euler *eu = (Physics_Euler*)phys->data; 776 const EulerNode *x = (const EulerNode*)xx; 777 PetscReal p; 778 779 PetscFunctionBeginUser; 780 f[eu->monitor.Density] = x->r; 781 f[eu->monitor.Momentum] = NormDIM(x->ru); 782 f[eu->monitor.Energy] = x->E; 783 f[eu->monitor.Speed] = NormDIM(x->ru)/x->r; 784 Pressure_PG(eu->pars[EULER_PAR_GAMMA], x, &p); 785 f[eu->monitor.Pressure] = p; 786 PetscFunctionReturn(0); 787 } 788 789 static PetscErrorCode SetUpBC_Euler(DM dm, PetscDS prob,Physics phys) 790 { 791 PetscErrorCode ierr; 792 Physics_Euler *eu = (Physics_Euler *) phys->data; 793 DMLabel label; 794 795 PetscFunctionBeginUser; 796 ierr = DMGetLabel(dm, "Face Sets", &label);CHKERRQ(ierr); 797 if (eu->type == EULER_LINEAR_WAVE) { 798 const PetscInt wallids[] = {100,101}; 799 ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, ALEN(wallids), wallids, 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Euler_Wall, NULL, phys, NULL);CHKERRQ(ierr); 800 } 801 else { 802 const PetscInt wallids[] = {100,101,200,300}; 803 ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, ALEN(wallids), wallids, 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Euler_Wall, NULL, phys, NULL);CHKERRQ(ierr); 804 } 805 PetscFunctionReturn(0); 806 } 807 808 static PetscErrorCode PhysicsCreate_Euler(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject) 809 { 810 Physics_Euler *eu; 811 PetscErrorCode ierr; 812 813 PetscFunctionBeginUser; 814 phys->field_desc = PhysicsFields_Euler; 815 phys->riemann = (PetscRiemannFunc) PhysicsRiemann_Euler_Godunov; 816 ierr = PetscNew(&eu);CHKERRQ(ierr); 817 phys->data = eu; 818 mod->setupbc = SetUpBC_Euler; 819 ierr = PetscOptionsHead(PetscOptionsObject,"Euler options");CHKERRQ(ierr); 820 { 821 PetscReal alpha; 822 char type[64] = "linear_wave"; 823 PetscBool is; 824 eu->pars[EULER_PAR_GAMMA] = 1.4; 825 eu->pars[EULER_PAR_AMACH] = 2.02; 826 eu->pars[EULER_PAR_RHOR] = 3.0; 827 eu->pars[EULER_PAR_ITANA] = 0.57735026918963; /* angle of Euler self similar (SS) shock */ 828 ierr = PetscOptionsReal("-eu_gamma","Heat capacity ratio","",eu->pars[EULER_PAR_GAMMA],&eu->pars[EULER_PAR_GAMMA],NULL);CHKERRQ(ierr); 829 ierr = PetscOptionsReal("-eu_amach","Shock speed (Mach)","",eu->pars[EULER_PAR_AMACH],&eu->pars[EULER_PAR_AMACH],NULL);CHKERRQ(ierr); 830 ierr = PetscOptionsReal("-eu_rho2","Density right of discontinuity","",eu->pars[EULER_PAR_RHOR],&eu->pars[EULER_PAR_RHOR],NULL);CHKERRQ(ierr); 831 alpha = 60.; 832 ierr = PetscOptionsReal("-eu_alpha","Angle of discontinuity","",alpha,&alpha,NULL);CHKERRQ(ierr); 833 if (alpha<=0. || alpha>90.) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Alpha bust be > 0 and <= 90 (%g)",alpha); 834 eu->pars[EULER_PAR_ITANA] = 1./PetscTanReal( alpha * PETSC_PI / 180.0); 835 ierr = PetscOptionsString("-eu_type","Type of Euler test","",type,type,sizeof(type),NULL);CHKERRQ(ierr); 836 ierr = PetscStrcmp(type,"linear_wave", &is);CHKERRQ(ierr); 837 if (is) { 838 /* Remember this should be periodic */ 839 eu->type = EULER_LINEAR_WAVE; 840 ierr = PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"linear_wave");CHKERRQ(ierr); 841 } 842 else { 843 if (DIM != 2) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DIM must be 2 unless linear wave test %s",type); 844 ierr = PetscStrcmp(type,"iv_shock", &is);CHKERRQ(ierr); 845 if (is) { 846 eu->type = EULER_IV_SHOCK; 847 ierr = PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"iv_shock");CHKERRQ(ierr); 848 } 849 else { 850 ierr = PetscStrcmp(type,"ss_shock", &is);CHKERRQ(ierr); 851 if (is) { 852 eu->type = EULER_SS_SHOCK; 853 ierr = PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"ss_shock");CHKERRQ(ierr); 854 } 855 else { 856 ierr = PetscStrcmp(type,"shock_tube", &is);CHKERRQ(ierr); 857 if (is) eu->type = EULER_SHOCK_TUBE; 858 else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Unknown Euler type %s",type); 859 ierr = PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"shock_tube");CHKERRQ(ierr); 860 } 861 } 862 } 863 } 864 ierr = PetscOptionsTail();CHKERRQ(ierr); 865 eu->sound = SpeedOfSound_PG; 866 phys->maxspeed = 0.; /* will get set in solution */ 867 ierr = ModelSolutionSetDefault(mod,PhysicsSolution_Euler,phys);CHKERRQ(ierr); 868 ierr = ModelFunctionalRegister(mod,"Speed",&eu->monitor.Speed,PhysicsFunctional_Euler,phys);CHKERRQ(ierr); 869 ierr = ModelFunctionalRegister(mod,"Energy",&eu->monitor.Energy,PhysicsFunctional_Euler,phys);CHKERRQ(ierr); 870 ierr = ModelFunctionalRegister(mod,"Density",&eu->monitor.Density,PhysicsFunctional_Euler,phys);CHKERRQ(ierr); 871 ierr = ModelFunctionalRegister(mod,"Momentum",&eu->monitor.Momentum,PhysicsFunctional_Euler,phys);CHKERRQ(ierr); 872 ierr = ModelFunctionalRegister(mod,"Pressure",&eu->monitor.Pressure,PhysicsFunctional_Euler,phys);CHKERRQ(ierr); 873 874 PetscFunctionReturn(0); 875 } 876 877 static PetscErrorCode ErrorIndicator_Simple(PetscInt dim, PetscReal volume, PetscInt numComps, const PetscScalar u[], const PetscScalar grad[], PetscReal *error, void *ctx) 878 { 879 PetscReal err = 0.; 880 PetscInt i, j; 881 882 PetscFunctionBeginUser; 883 for (i = 0; i < numComps; i++) { 884 for (j = 0; j < dim; j++) { 885 err += PetscSqr(PetscRealPart(grad[i * dim + j])); 886 } 887 } 888 *error = volume * err; 889 PetscFunctionReturn(0); 890 } 891 892 PetscErrorCode CreatePartitionVec(DM dm, DM *dmCell, Vec *partition) 893 { 894 PetscSF sfPoint; 895 PetscSection coordSection; 896 Vec coordinates; 897 PetscSection sectionCell; 898 PetscScalar *part; 899 PetscInt cStart, cEnd, c; 900 PetscMPIInt rank; 901 PetscErrorCode ierr; 902 903 PetscFunctionBeginUser; 904 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 905 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 906 ierr = DMClone(dm, dmCell);CHKERRQ(ierr); 907 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 908 ierr = DMSetPointSF(*dmCell, sfPoint);CHKERRQ(ierr); 909 ierr = DMSetCoordinateSection(*dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 910 ierr = DMSetCoordinatesLocal(*dmCell, coordinates);CHKERRQ(ierr); 911 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRMPI(ierr); 912 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ionCell);CHKERRQ(ierr); 913 ierr = DMPlexGetHeightStratum(*dmCell, 0, &cStart, &cEnd);CHKERRQ(ierr); 914 ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 915 for (c = cStart; c < cEnd; ++c) { 916 ierr = PetscSectionSetDof(sectionCell, c, 1);CHKERRQ(ierr); 917 } 918 ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 919 ierr = DMSetLocalSection(*dmCell, sectionCell);CHKERRQ(ierr); 920 ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 921 ierr = DMCreateLocalVector(*dmCell, partition);CHKERRQ(ierr); 922 ierr = PetscObjectSetName((PetscObject)*partition, "partition");CHKERRQ(ierr); 923 ierr = VecGetArray(*partition, &part);CHKERRQ(ierr); 924 for (c = cStart; c < cEnd; ++c) { 925 PetscScalar *p; 926 927 ierr = DMPlexPointLocalRef(*dmCell, c, part, &p);CHKERRQ(ierr); 928 p[0] = rank; 929 } 930 ierr = VecRestoreArray(*partition, &part);CHKERRQ(ierr); 931 PetscFunctionReturn(0); 932 } 933 934 PetscErrorCode CreateMassMatrix(DM dm, Vec *massMatrix, User user) 935 { 936 DM plex, dmMass, dmFace, dmCell, dmCoord; 937 PetscSection coordSection; 938 Vec coordinates, facegeom, cellgeom; 939 PetscSection sectionMass; 940 PetscScalar *m; 941 const PetscScalar *fgeom, *cgeom, *coords; 942 PetscInt vStart, vEnd, v; 943 PetscErrorCode ierr; 944 945 PetscFunctionBeginUser; 946 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 947 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 948 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 949 ierr = DMClone(dm, &dmMass);CHKERRQ(ierr); 950 ierr = DMSetCoordinateSection(dmMass, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 951 ierr = DMSetCoordinatesLocal(dmMass, coordinates);CHKERRQ(ierr); 952 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ionMass);CHKERRQ(ierr); 953 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 954 ierr = PetscSectionSetChart(sectionMass, vStart, vEnd);CHKERRQ(ierr); 955 for (v = vStart; v < vEnd; ++v) { 956 PetscInt numFaces; 957 958 ierr = DMPlexGetSupportSize(dmMass, v, &numFaces);CHKERRQ(ierr); 959 ierr = PetscSectionSetDof(sectionMass, v, numFaces*numFaces);CHKERRQ(ierr); 960 } 961 ierr = PetscSectionSetUp(sectionMass);CHKERRQ(ierr); 962 ierr = DMSetLocalSection(dmMass, sectionMass);CHKERRQ(ierr); 963 ierr = PetscSectionDestroy(§ionMass);CHKERRQ(ierr); 964 ierr = DMGetLocalVector(dmMass, massMatrix);CHKERRQ(ierr); 965 ierr = VecGetArray(*massMatrix, &m);CHKERRQ(ierr); 966 ierr = DMPlexGetGeometryFVM(plex, &facegeom, &cellgeom, NULL);CHKERRQ(ierr); 967 ierr = VecGetDM(facegeom, &dmFace);CHKERRQ(ierr); 968 ierr = VecGetArrayRead(facegeom, &fgeom);CHKERRQ(ierr); 969 ierr = VecGetDM(cellgeom, &dmCell);CHKERRQ(ierr); 970 ierr = VecGetArrayRead(cellgeom, &cgeom);CHKERRQ(ierr); 971 ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr); 972 ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 973 for (v = vStart; v < vEnd; ++v) { 974 const PetscInt *faces; 975 PetscFVFaceGeom *fgA, *fgB, *cg; 976 PetscScalar *vertex; 977 PetscInt numFaces, sides[2], f, g; 978 979 ierr = DMPlexPointLocalRead(dmCoord, v, coords, &vertex);CHKERRQ(ierr); 980 ierr = DMPlexGetSupportSize(dmMass, v, &numFaces);CHKERRQ(ierr); 981 ierr = DMPlexGetSupport(dmMass, v, &faces);CHKERRQ(ierr); 982 for (f = 0; f < numFaces; ++f) { 983 sides[0] = faces[f]; 984 ierr = DMPlexPointLocalRead(dmFace, faces[f], fgeom, &fgA);CHKERRQ(ierr); 985 for (g = 0; g < numFaces; ++g) { 986 const PetscInt *cells = NULL; 987 PetscReal area = 0.0; 988 PetscInt numCells; 989 990 sides[1] = faces[g]; 991 ierr = DMPlexPointLocalRead(dmFace, faces[g], fgeom, &fgB);CHKERRQ(ierr); 992 ierr = DMPlexGetJoin(dmMass, 2, sides, &numCells, &cells);CHKERRQ(ierr); 993 if (numCells != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Invalid join for faces"); 994 ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cg);CHKERRQ(ierr); 995 area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgA->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgA->centroid[0] - cg->centroid[0])); 996 area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgB->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgB->centroid[0] - cg->centroid[0])); 997 m[f*numFaces+g] = Dot2Real(fgA->normal, fgB->normal)*area*0.5; 998 ierr = DMPlexRestoreJoin(dmMass, 2, sides, &numCells, &cells);CHKERRQ(ierr); 999 } 1000 } 1001 } 1002 ierr = VecRestoreArrayRead(facegeom, &fgeom);CHKERRQ(ierr); 1003 ierr = VecRestoreArrayRead(cellgeom, &cgeom);CHKERRQ(ierr); 1004 ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 1005 ierr = VecRestoreArray(*massMatrix, &m);CHKERRQ(ierr); 1006 ierr = DMDestroy(&dmMass);CHKERRQ(ierr); 1007 ierr = DMDestroy(&plex);CHKERRQ(ierr); 1008 PetscFunctionReturn(0); 1009 } 1010 1011 /* Behavior will be different for multi-physics or when using non-default boundary conditions */ 1012 static PetscErrorCode ModelSolutionSetDefault(Model mod,SolutionFunction func,void *ctx) 1013 { 1014 PetscFunctionBeginUser; 1015 mod->solution = func; 1016 mod->solutionctx = ctx; 1017 PetscFunctionReturn(0); 1018 } 1019 1020 static PetscErrorCode ModelFunctionalRegister(Model mod,const char *name,PetscInt *offset,FunctionalFunction func,void *ctx) 1021 { 1022 PetscErrorCode ierr; 1023 FunctionalLink link,*ptr; 1024 PetscInt lastoffset = -1; 1025 1026 PetscFunctionBeginUser; 1027 for (ptr=&mod->functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset; 1028 ierr = PetscNew(&link);CHKERRQ(ierr); 1029 ierr = PetscStrallocpy(name,&link->name);CHKERRQ(ierr); 1030 link->offset = lastoffset + 1; 1031 link->func = func; 1032 link->ctx = ctx; 1033 link->next = NULL; 1034 *ptr = link; 1035 *offset = link->offset; 1036 PetscFunctionReturn(0); 1037 } 1038 1039 static PetscErrorCode ModelFunctionalSetFromOptions(Model mod,PetscOptionItems *PetscOptionsObject) 1040 { 1041 PetscErrorCode ierr; 1042 PetscInt i,j; 1043 FunctionalLink link; 1044 char *names[256]; 1045 1046 PetscFunctionBeginUser; 1047 mod->numMonitored = ALEN(names); 1048 ierr = PetscOptionsStringArray("-monitor","list of functionals to monitor","",names,&mod->numMonitored,NULL);CHKERRQ(ierr); 1049 /* Create list of functionals that will be computed somehow */ 1050 ierr = PetscMalloc1(mod->numMonitored,&mod->functionalMonitored);CHKERRQ(ierr); 1051 /* Create index of calls that we will have to make to compute these functionals (over-allocation in general). */ 1052 ierr = PetscMalloc1(mod->numMonitored,&mod->functionalCall);CHKERRQ(ierr); 1053 mod->numCall = 0; 1054 for (i=0; i<mod->numMonitored; i++) { 1055 for (link=mod->functionalRegistry; link; link=link->next) { 1056 PetscBool match; 1057 ierr = PetscStrcasecmp(names[i],link->name,&match);CHKERRQ(ierr); 1058 if (match) break; 1059 } 1060 if (!link) SETERRQ1(mod->comm,PETSC_ERR_USER,"No known functional '%s'",names[i]); 1061 mod->functionalMonitored[i] = link; 1062 for (j=0; j<i; j++) { 1063 if (mod->functionalCall[j]->func == link->func && mod->functionalCall[j]->ctx == link->ctx) goto next_name; 1064 } 1065 mod->functionalCall[mod->numCall++] = link; /* Just points to the first link using the result. There may be more results. */ 1066 next_name: 1067 ierr = PetscFree(names[i]);CHKERRQ(ierr); 1068 } 1069 1070 /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */ 1071 mod->maxComputed = -1; 1072 for (link=mod->functionalRegistry; link; link=link->next) { 1073 for (i=0; i<mod->numCall; i++) { 1074 FunctionalLink call = mod->functionalCall[i]; 1075 if (link->func == call->func && link->ctx == call->ctx) { 1076 mod->maxComputed = PetscMax(mod->maxComputed,link->offset); 1077 } 1078 } 1079 } 1080 PetscFunctionReturn(0); 1081 } 1082 1083 static PetscErrorCode FunctionalLinkDestroy(FunctionalLink *link) 1084 { 1085 PetscErrorCode ierr; 1086 FunctionalLink l,next; 1087 1088 PetscFunctionBeginUser; 1089 if (!link) PetscFunctionReturn(0); 1090 l = *link; 1091 *link = NULL; 1092 for (; l; l=next) { 1093 next = l->next; 1094 ierr = PetscFree(l->name);CHKERRQ(ierr); 1095 ierr = PetscFree(l);CHKERRQ(ierr); 1096 } 1097 PetscFunctionReturn(0); 1098 } 1099 1100 /* put the solution callback into a functional callback */ 1101 static PetscErrorCode SolutionFunctional(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *modctx) 1102 { 1103 Model mod; 1104 PetscErrorCode ierr; 1105 PetscFunctionBegin; 1106 mod = (Model) modctx; 1107 ierr = (*mod->solution)(mod, time, x, u, mod->solutionctx);CHKERRQ(ierr); 1108 PetscFunctionReturn(0); 1109 } 1110 1111 PetscErrorCode SetInitialCondition(DM dm, Vec X, User user) 1112 { 1113 PetscErrorCode (*func[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 1114 void *ctx[1]; 1115 Model mod = user->model; 1116 PetscErrorCode ierr; 1117 1118 PetscFunctionBeginUser; 1119 func[0] = SolutionFunctional; 1120 ctx[0] = (void *) mod; 1121 ierr = DMProjectFunction(dm,0.0,func,ctx,INSERT_ALL_VALUES,X);CHKERRQ(ierr); 1122 PetscFunctionReturn(0); 1123 } 1124 1125 static PetscErrorCode OutputVTK(DM dm, const char *filename, PetscViewer *viewer) 1126 { 1127 PetscErrorCode ierr; 1128 1129 PetscFunctionBeginUser; 1130 ierr = PetscViewerCreate(PetscObjectComm((PetscObject)dm), viewer);CHKERRQ(ierr); 1131 ierr = PetscViewerSetType(*viewer, PETSCVIEWERVTK);CHKERRQ(ierr); 1132 ierr = PetscViewerFileSetName(*viewer, filename);CHKERRQ(ierr); 1133 PetscFunctionReturn(0); 1134 } 1135 1136 static PetscErrorCode MonitorVTK(TS ts,PetscInt stepnum,PetscReal time,Vec X,void *ctx) 1137 { 1138 User user = (User)ctx; 1139 DM dm, plex; 1140 PetscViewer viewer; 1141 char filename[PETSC_MAX_PATH_LEN],*ftable = NULL; 1142 PetscReal xnorm; 1143 PetscErrorCode ierr; 1144 1145 PetscFunctionBeginUser; 1146 ierr = PetscObjectSetName((PetscObject) X, "u");CHKERRQ(ierr); 1147 ierr = VecGetDM(X,&dm);CHKERRQ(ierr); 1148 ierr = VecNorm(X,NORM_INFINITY,&xnorm);CHKERRQ(ierr); 1149 1150 if (stepnum >= 0) { 1151 stepnum += user->monitorStepOffset; 1152 } 1153 if (stepnum >= 0) { /* No summary for final time */ 1154 Model mod = user->model; 1155 Vec cellgeom; 1156 PetscInt c,cStart,cEnd,fcount,i; 1157 size_t ftableused,ftablealloc; 1158 const PetscScalar *cgeom,*x; 1159 DM dmCell; 1160 DMLabel vtkLabel; 1161 PetscReal *fmin,*fmax,*fintegral,*ftmp; 1162 1163 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 1164 ierr = DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL);CHKERRQ(ierr); 1165 fcount = mod->maxComputed+1; 1166 ierr = PetscMalloc4(fcount,&fmin,fcount,&fmax,fcount,&fintegral,fcount,&ftmp);CHKERRQ(ierr); 1167 for (i=0; i<fcount; i++) { 1168 fmin[i] = PETSC_MAX_REAL; 1169 fmax[i] = PETSC_MIN_REAL; 1170 fintegral[i] = 0; 1171 } 1172 ierr = VecGetDM(cellgeom,&dmCell);CHKERRQ(ierr); 1173 ierr = DMPlexGetSimplexOrBoxCells(dmCell,0,&cStart,&cEnd);CHKERRQ(ierr); 1174 ierr = VecGetArrayRead(cellgeom,&cgeom);CHKERRQ(ierr); 1175 ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 1176 ierr = DMGetLabel(dm,"vtk",&vtkLabel);CHKERRQ(ierr); 1177 for (c = cStart; c < cEnd; ++c) { 1178 PetscFVCellGeom *cg; 1179 const PetscScalar *cx = NULL; 1180 PetscInt vtkVal = 0; 1181 1182 /* not that these two routines as currently implemented work for any dm with a 1183 * localSection/globalSection */ 1184 ierr = DMPlexPointLocalRead(dmCell,c,cgeom,&cg);CHKERRQ(ierr); 1185 ierr = DMPlexPointGlobalRead(dm,c,x,&cx);CHKERRQ(ierr); 1186 if (vtkLabel) {ierr = DMLabelGetValue(vtkLabel,c,&vtkVal);CHKERRQ(ierr);} 1187 if (!vtkVal || !cx) continue; /* ghost, or not a global cell */ 1188 for (i=0; i<mod->numCall; i++) { 1189 FunctionalLink flink = mod->functionalCall[i]; 1190 ierr = (*flink->func)(mod,time,cg->centroid,cx,ftmp,flink->ctx);CHKERRQ(ierr); 1191 } 1192 for (i=0; i<fcount; i++) { 1193 fmin[i] = PetscMin(fmin[i],ftmp[i]); 1194 fmax[i] = PetscMax(fmax[i],ftmp[i]); 1195 fintegral[i] += cg->volume * ftmp[i]; 1196 } 1197 } 1198 ierr = VecRestoreArrayRead(cellgeom,&cgeom);CHKERRQ(ierr); 1199 ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 1200 ierr = DMDestroy(&plex);CHKERRQ(ierr); 1201 ierr = MPI_Allreduce(MPI_IN_PLACE,fmin,fcount,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));CHKERRMPI(ierr); 1202 ierr = MPI_Allreduce(MPI_IN_PLACE,fmax,fcount,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));CHKERRMPI(ierr); 1203 ierr = MPI_Allreduce(MPI_IN_PLACE,fintegral,fcount,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));CHKERRMPI(ierr); 1204 1205 ftablealloc = fcount * 100; 1206 ftableused = 0; 1207 ierr = PetscMalloc1(ftablealloc,&ftable);CHKERRQ(ierr); 1208 for (i=0; i<mod->numMonitored; i++) { 1209 size_t countused; 1210 char buffer[256],*p; 1211 FunctionalLink flink = mod->functionalMonitored[i]; 1212 PetscInt id = flink->offset; 1213 if (i % 3) { 1214 ierr = PetscArraycpy(buffer," ",2);CHKERRQ(ierr); 1215 p = buffer + 2; 1216 } else if (i) { 1217 char newline[] = "\n"; 1218 ierr = PetscMemcpy(buffer,newline,sizeof(newline)-1);CHKERRQ(ierr); 1219 p = buffer + sizeof(newline) - 1; 1220 } else { 1221 p = buffer; 1222 } 1223 ierr = PetscSNPrintfCount(p,sizeof buffer-(p-buffer),"%12s [%10.7g,%10.7g] int %10.7g",&countused,flink->name,(double)fmin[id],(double)fmax[id],(double)fintegral[id]);CHKERRQ(ierr); 1224 countused--; 1225 countused += p - buffer; 1226 if (countused > ftablealloc-ftableused-1) { /* reallocate */ 1227 char *ftablenew; 1228 ftablealloc = 2*ftablealloc + countused; 1229 ierr = PetscMalloc(ftablealloc,&ftablenew);CHKERRQ(ierr); 1230 ierr = PetscArraycpy(ftablenew,ftable,ftableused);CHKERRQ(ierr); 1231 ierr = PetscFree(ftable);CHKERRQ(ierr); 1232 ftable = ftablenew; 1233 } 1234 ierr = PetscArraycpy(ftable+ftableused,buffer,countused);CHKERRQ(ierr); 1235 ftableused += countused; 1236 ftable[ftableused] = 0; 1237 } 1238 ierr = PetscFree4(fmin,fmax,fintegral,ftmp);CHKERRQ(ierr); 1239 1240 ierr = PetscPrintf(PetscObjectComm((PetscObject)ts),"% 3D time %8.4g |x| %8.4g %s\n",stepnum,(double)time,(double)xnorm,ftable ? ftable : "");CHKERRQ(ierr); 1241 ierr = PetscFree(ftable);CHKERRQ(ierr); 1242 } 1243 if (user->vtkInterval < 1) PetscFunctionReturn(0); 1244 if ((stepnum == -1) ^ (stepnum % user->vtkInterval == 0)) { 1245 if (stepnum == -1) { /* Final time is not multiple of normal time interval, write it anyway */ 1246 ierr = TSGetStepNumber(ts,&stepnum);CHKERRQ(ierr); 1247 } 1248 ierr = PetscSNPrintf(filename,sizeof filename,"%s-%03D.vtu",user->outputBasename,stepnum);CHKERRQ(ierr); 1249 ierr = OutputVTK(dm,filename,&viewer);CHKERRQ(ierr); 1250 ierr = VecView(X,viewer);CHKERRQ(ierr); 1251 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1252 } 1253 PetscFunctionReturn(0); 1254 } 1255 1256 static PetscErrorCode initializeTS(DM dm, User user, TS *ts) 1257 { 1258 PetscErrorCode ierr; 1259 1260 PetscFunctionBegin; 1261 ierr = TSCreate(PetscObjectComm((PetscObject)dm), ts);CHKERRQ(ierr); 1262 ierr = TSSetType(*ts, TSSSP);CHKERRQ(ierr); 1263 ierr = TSSetDM(*ts, dm);CHKERRQ(ierr); 1264 if (user->vtkmon) { 1265 ierr = TSMonitorSet(*ts,MonitorVTK,user,NULL);CHKERRQ(ierr); 1266 } 1267 ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, user);CHKERRQ(ierr); 1268 ierr = DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, user);CHKERRQ(ierr); 1269 ierr = TSSetMaxTime(*ts,2.0);CHKERRQ(ierr); 1270 ierr = TSSetExactFinalTime(*ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 1271 PetscFunctionReturn(0); 1272 } 1273 1274 static PetscErrorCode adaptToleranceFVM(PetscFV fvm, TS ts, Vec sol, VecTagger refineTag, VecTagger coarsenTag, User user, TS *tsNew, Vec *solNew) 1275 { 1276 DM dm, gradDM, plex, cellDM, adaptedDM = NULL; 1277 Vec cellGeom, faceGeom; 1278 PetscBool isForest, computeGradient; 1279 Vec grad, locGrad, locX, errVec; 1280 PetscInt cStart, cEnd, c, dim, nRefine, nCoarsen; 1281 PetscReal minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2], minInd, maxInd, time; 1282 PetscScalar *errArray; 1283 const PetscScalar *pointVals; 1284 const PetscScalar *pointGrads; 1285 const PetscScalar *pointGeom; 1286 DMLabel adaptLabel = NULL; 1287 IS refineIS, coarsenIS; 1288 PetscErrorCode ierr; 1289 1290 PetscFunctionBegin; 1291 ierr = TSGetTime(ts,&time);CHKERRQ(ierr); 1292 ierr = VecGetDM(sol, &dm);CHKERRQ(ierr); 1293 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 1294 ierr = PetscFVGetComputeGradients(fvm,&computeGradient);CHKERRQ(ierr); 1295 ierr = PetscFVSetComputeGradients(fvm,PETSC_TRUE);CHKERRQ(ierr); 1296 ierr = DMIsForest(dm, &isForest);CHKERRQ(ierr); 1297 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 1298 ierr = DMPlexGetDataFVM(plex, fvm, &cellGeom, &faceGeom, &gradDM);CHKERRQ(ierr); 1299 ierr = DMCreateLocalVector(plex,&locX);CHKERRQ(ierr); 1300 ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, faceGeom, cellGeom, NULL);CHKERRQ(ierr); 1301 ierr = DMGlobalToLocalBegin(plex, sol, INSERT_VALUES, locX);CHKERRQ(ierr); 1302 ierr = DMGlobalToLocalEnd (plex, sol, INSERT_VALUES, locX);CHKERRQ(ierr); 1303 ierr = DMCreateGlobalVector(gradDM, &grad);CHKERRQ(ierr); 1304 ierr = DMPlexReconstructGradientsFVM(plex, locX, grad);CHKERRQ(ierr); 1305 ierr = DMCreateLocalVector(gradDM, &locGrad);CHKERRQ(ierr); 1306 ierr = DMGlobalToLocalBegin(gradDM, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 1307 ierr = DMGlobalToLocalEnd(gradDM, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 1308 ierr = VecDestroy(&grad);CHKERRQ(ierr); 1309 ierr = DMPlexGetSimplexOrBoxCells(plex,0,&cStart,&cEnd);CHKERRQ(ierr); 1310 ierr = VecGetArrayRead(locGrad,&pointGrads);CHKERRQ(ierr); 1311 ierr = VecGetArrayRead(cellGeom,&pointGeom);CHKERRQ(ierr); 1312 ierr = VecGetArrayRead(locX,&pointVals);CHKERRQ(ierr); 1313 ierr = VecGetDM(cellGeom,&cellDM);CHKERRQ(ierr); 1314 ierr = DMLabelCreate(PETSC_COMM_SELF,"adapt",&adaptLabel);CHKERRQ(ierr); 1315 ierr = VecCreateMPI(PetscObjectComm((PetscObject)plex),cEnd-cStart,PETSC_DETERMINE,&errVec);CHKERRQ(ierr); 1316 ierr = VecSetUp(errVec);CHKERRQ(ierr); 1317 ierr = VecGetArray(errVec,&errArray);CHKERRQ(ierr); 1318 for (c = cStart; c < cEnd; c++) { 1319 PetscReal errInd = 0.; 1320 PetscScalar *pointGrad; 1321 PetscScalar *pointVal; 1322 PetscFVCellGeom *cg; 1323 1324 ierr = DMPlexPointLocalRead(gradDM,c,pointGrads,&pointGrad);CHKERRQ(ierr); 1325 ierr = DMPlexPointLocalRead(cellDM,c,pointGeom,&cg);CHKERRQ(ierr); 1326 ierr = DMPlexPointLocalRead(plex,c,pointVals,&pointVal);CHKERRQ(ierr); 1327 1328 ierr = (user->model->errorIndicator)(dim,cg->volume,user->model->physics->dof,pointVal,pointGrad,&errInd,user->model->errorCtx);CHKERRQ(ierr); 1329 errArray[c-cStart] = errInd; 1330 minMaxInd[0] = PetscMin(minMaxInd[0],errInd); 1331 minMaxInd[1] = PetscMax(minMaxInd[1],errInd); 1332 } 1333 ierr = VecRestoreArray(errVec,&errArray);CHKERRQ(ierr); 1334 ierr = VecRestoreArrayRead(locX,&pointVals);CHKERRQ(ierr); 1335 ierr = VecRestoreArrayRead(cellGeom,&pointGeom);CHKERRQ(ierr); 1336 ierr = VecRestoreArrayRead(locGrad,&pointGrads);CHKERRQ(ierr); 1337 ierr = VecDestroy(&locGrad);CHKERRQ(ierr); 1338 ierr = VecDestroy(&locX);CHKERRQ(ierr); 1339 ierr = DMDestroy(&plex);CHKERRQ(ierr); 1340 1341 ierr = VecTaggerComputeIS(refineTag,errVec,&refineIS);CHKERRQ(ierr); 1342 ierr = VecTaggerComputeIS(coarsenTag,errVec,&coarsenIS);CHKERRQ(ierr); 1343 ierr = ISGetSize(refineIS,&nRefine);CHKERRQ(ierr); 1344 ierr = ISGetSize(coarsenIS,&nCoarsen);CHKERRQ(ierr); 1345 if (nRefine) {ierr = DMLabelSetStratumIS(adaptLabel,DM_ADAPT_REFINE,refineIS);CHKERRQ(ierr);} 1346 if (nCoarsen) {ierr = DMLabelSetStratumIS(adaptLabel,DM_ADAPT_COARSEN,coarsenIS);CHKERRQ(ierr);} 1347 ierr = ISDestroy(&coarsenIS);CHKERRQ(ierr); 1348 ierr = ISDestroy(&refineIS);CHKERRQ(ierr); 1349 ierr = VecDestroy(&errVec);CHKERRQ(ierr); 1350 1351 ierr = PetscFVSetComputeGradients(fvm,computeGradient);CHKERRQ(ierr); 1352 minMaxInd[1] = -minMaxInd[1]; 1353 ierr = MPI_Allreduce(minMaxInd,minMaxIndGlobal,2,MPIU_REAL,MPI_MIN,PetscObjectComm((PetscObject)dm));CHKERRMPI(ierr); 1354 minInd = minMaxIndGlobal[0]; 1355 maxInd = -minMaxIndGlobal[1]; 1356 ierr = PetscInfo2(ts, "error indicator range (%E, %E)\n", minInd, maxInd);CHKERRQ(ierr); 1357 if (nRefine || nCoarsen) { /* at least one cell is over the refinement threshold */ 1358 ierr = DMAdaptLabel(dm,adaptLabel,&adaptedDM);CHKERRQ(ierr); 1359 } 1360 ierr = DMLabelDestroy(&adaptLabel);CHKERRQ(ierr); 1361 if (adaptedDM) { 1362 ierr = PetscInfo2(ts, "Adapted mesh, marking %D cells for refinement, and %D cells for coarsening\n", nRefine, nCoarsen);CHKERRQ(ierr); 1363 if (tsNew) {ierr = initializeTS(adaptedDM, user, tsNew);CHKERRQ(ierr);} 1364 if (solNew) { 1365 ierr = DMCreateGlobalVector(adaptedDM, solNew);CHKERRQ(ierr); 1366 ierr = PetscObjectSetName((PetscObject) *solNew, "solution");CHKERRQ(ierr); 1367 ierr = DMForestTransferVec(dm, sol, adaptedDM, *solNew, PETSC_TRUE, time);CHKERRQ(ierr); 1368 } 1369 if (isForest) {ierr = DMForestSetAdaptivityForest(adaptedDM,NULL);CHKERRQ(ierr);} /* clear internal references to the previous dm */ 1370 ierr = DMDestroy(&adaptedDM);CHKERRQ(ierr); 1371 } else { 1372 if (tsNew) *tsNew = NULL; 1373 if (solNew) *solNew = NULL; 1374 } 1375 PetscFunctionReturn(0); 1376 } 1377 1378 int main(int argc, char **argv) 1379 { 1380 MPI_Comm comm; 1381 PetscDS prob; 1382 PetscFV fvm; 1383 PetscLimiter limiter = NULL, noneLimiter = NULL; 1384 User user; 1385 Model mod; 1386 Physics phys; 1387 DM dm, plex; 1388 PetscReal ftime, cfl, dt, minRadius; 1389 PetscInt dim, nsteps; 1390 TS ts; 1391 TSConvergedReason reason; 1392 Vec X; 1393 PetscViewer viewer; 1394 PetscBool vtkCellGeom, useAMR; 1395 PetscInt adaptInterval; 1396 char physname[256] = "advect"; 1397 VecTagger refineTag = NULL, coarsenTag = NULL; 1398 PetscErrorCode ierr; 1399 1400 ierr = PetscInitialize(&argc, &argv, (char*) 0, help);if (ierr) return ierr; 1401 comm = PETSC_COMM_WORLD; 1402 1403 ierr = PetscNew(&user);CHKERRQ(ierr); 1404 ierr = PetscNew(&user->model);CHKERRQ(ierr); 1405 ierr = PetscNew(&user->model->physics);CHKERRQ(ierr); 1406 mod = user->model; 1407 phys = mod->physics; 1408 mod->comm = comm; 1409 useAMR = PETSC_FALSE; 1410 adaptInterval = 1; 1411 1412 /* Register physical models to be available on the command line */ 1413 ierr = PetscFunctionListAdd(&PhysicsList,"advect" ,PhysicsCreate_Advect);CHKERRQ(ierr); 1414 ierr = PetscFunctionListAdd(&PhysicsList,"sw" ,PhysicsCreate_SW);CHKERRQ(ierr); 1415 ierr = PetscFunctionListAdd(&PhysicsList,"euler" ,PhysicsCreate_Euler);CHKERRQ(ierr); 1416 1417 ierr = PetscOptionsBegin(comm,NULL,"Unstructured Finite Volume Mesh Options","");CHKERRQ(ierr); 1418 { 1419 cfl = 0.9 * 4; /* default SSPRKS2 with s=5 stages is stable for CFL number s-1 */ 1420 ierr = PetscOptionsReal("-ufv_cfl","CFL number per step","",cfl,&cfl,NULL);CHKERRQ(ierr); 1421 user->vtkInterval = 1; 1422 ierr = PetscOptionsInt("-ufv_vtk_interval","VTK output interval (0 to disable)","",user->vtkInterval,&user->vtkInterval,NULL);CHKERRQ(ierr); 1423 user->vtkmon = PETSC_TRUE; 1424 ierr = PetscOptionsBool("-ufv_vtk_monitor","Use VTKMonitor routine","",user->vtkmon,&user->vtkmon,NULL);CHKERRQ(ierr); 1425 vtkCellGeom = PETSC_FALSE; 1426 ierr = PetscStrcpy(user->outputBasename, "ex11");CHKERRQ(ierr); 1427 ierr = PetscOptionsString("-ufv_vtk_basename","VTK output basename","",user->outputBasename,user->outputBasename,sizeof(user->outputBasename),NULL);CHKERRQ(ierr); 1428 ierr = PetscOptionsBool("-ufv_vtk_cellgeom","Write cell geometry (for debugging)","",vtkCellGeom,&vtkCellGeom,NULL);CHKERRQ(ierr); 1429 ierr = PetscOptionsBool("-ufv_use_amr","use local adaptive mesh refinement","",useAMR,&useAMR,NULL);CHKERRQ(ierr); 1430 ierr = PetscOptionsInt("-ufv_adapt_interval","time steps between AMR","",adaptInterval,&adaptInterval,NULL);CHKERRQ(ierr); 1431 } 1432 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1433 1434 if (useAMR) { 1435 VecTaggerBox refineBox, coarsenBox; 1436 1437 refineBox.min = refineBox.max = PETSC_MAX_REAL; 1438 coarsenBox.min = coarsenBox.max = PETSC_MIN_REAL; 1439 1440 ierr = VecTaggerCreate(comm,&refineTag);CHKERRQ(ierr); 1441 ierr = PetscObjectSetOptionsPrefix((PetscObject)refineTag,"refine_");CHKERRQ(ierr); 1442 ierr = VecTaggerSetType(refineTag,VECTAGGERABSOLUTE);CHKERRQ(ierr); 1443 ierr = VecTaggerAbsoluteSetBox(refineTag,&refineBox);CHKERRQ(ierr); 1444 ierr = VecTaggerSetFromOptions(refineTag);CHKERRQ(ierr); 1445 ierr = VecTaggerSetUp(refineTag);CHKERRQ(ierr); 1446 ierr = PetscObjectViewFromOptions((PetscObject)refineTag,NULL,"-tag_view");CHKERRQ(ierr); 1447 1448 ierr = VecTaggerCreate(comm,&coarsenTag);CHKERRQ(ierr); 1449 ierr = PetscObjectSetOptionsPrefix((PetscObject)coarsenTag,"coarsen_");CHKERRQ(ierr); 1450 ierr = VecTaggerSetType(coarsenTag,VECTAGGERABSOLUTE);CHKERRQ(ierr); 1451 ierr = VecTaggerAbsoluteSetBox(coarsenTag,&coarsenBox);CHKERRQ(ierr); 1452 ierr = VecTaggerSetFromOptions(coarsenTag);CHKERRQ(ierr); 1453 ierr = VecTaggerSetUp(coarsenTag);CHKERRQ(ierr); 1454 ierr = PetscObjectViewFromOptions((PetscObject)coarsenTag,NULL,"-tag_view");CHKERRQ(ierr); 1455 } 1456 1457 ierr = PetscOptionsBegin(comm,NULL,"Unstructured Finite Volume Physics Options","");CHKERRQ(ierr); 1458 { 1459 PetscErrorCode (*physcreate)(Model,Physics,PetscOptionItems*); 1460 ierr = PetscOptionsFList("-physics","Physics module to solve","",PhysicsList,physname,physname,sizeof physname,NULL);CHKERRQ(ierr); 1461 ierr = PetscFunctionListFind(PhysicsList,physname,&physcreate);CHKERRQ(ierr); 1462 ierr = PetscMemzero(phys,sizeof(struct _n_Physics));CHKERRQ(ierr); 1463 ierr = (*physcreate)(mod,phys,PetscOptionsObject);CHKERRQ(ierr); 1464 /* Count number of fields and dofs */ 1465 for (phys->nfields=0,phys->dof=0; phys->field_desc[phys->nfields].name; phys->nfields++) phys->dof += phys->field_desc[phys->nfields].dof; 1466 if (phys->dof <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set dof",physname); 1467 ierr = ModelFunctionalSetFromOptions(mod,PetscOptionsObject);CHKERRQ(ierr); 1468 } 1469 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1470 1471 /* Create mesh */ 1472 { 1473 PetscInt i; 1474 1475 ierr = DMCreate(comm, &dm);CHKERRQ(ierr); 1476 ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr); 1477 ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 1478 for (i = 0; i < DIM; i++) { mod->bounds[2*i] = 0.; mod->bounds[2*i+1] = 1.;}; 1479 dim = DIM; 1480 { /* a null name means just do a hex box */ 1481 PetscInt cells[3] = {1, 1, 1}, n = 3; 1482 PetscBool flg2, skew = PETSC_FALSE; 1483 PetscInt nret2 = 2*DIM; 1484 ierr = PetscOptionsBegin(comm,NULL,"Rectangular mesh options","");CHKERRQ(ierr); 1485 ierr = PetscOptionsRealArray("-grid_bounds","bounds of the mesh in each direction (i.e., x_min,x_max,y_min,y_max","",mod->bounds,&nret2,&flg2);CHKERRQ(ierr); 1486 ierr = PetscOptionsBool("-grid_skew_60","Skew grid for 60 degree shock mesh","",skew,&skew,NULL);CHKERRQ(ierr); 1487 ierr = PetscOptionsIntArray("-dm_plex_box_faces", "Number of faces along each dimension", "", cells, &n, NULL);CHKERRQ(ierr); 1488 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1489 /* TODO Rewrite this with Mark, and remove grid_bounds at that time */ 1490 if (flg2) { 1491 PetscInt dimEmbed, i; 1492 PetscInt nCoords; 1493 PetscScalar *coords; 1494 Vec coordinates; 1495 1496 ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr); 1497 ierr = DMGetCoordinateDim(dm,&dimEmbed);CHKERRQ(ierr); 1498 ierr = VecGetLocalSize(coordinates,&nCoords);CHKERRQ(ierr); 1499 if (nCoords % dimEmbed) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Coordinate vector the wrong size"); 1500 ierr = VecGetArray(coordinates,&coords);CHKERRQ(ierr); 1501 for (i = 0; i < nCoords; i += dimEmbed) { 1502 PetscInt j; 1503 1504 PetscScalar *coord = &coords[i]; 1505 for (j = 0; j < dimEmbed; j++) { 1506 coord[j] = mod->bounds[2 * j] + coord[j] * (mod->bounds[2 * j + 1] - mod->bounds[2 * j]); 1507 if (dim==2 && cells[1]==1 && j==0 && skew) { 1508 if (cells[0] == 2 && i == 8) { 1509 coord[j] = .57735026918963; /* hack to get 60 deg skewed mesh */ 1510 } else if (cells[0] == 3) { 1511 if (i==2 || i==10) coord[j] = mod->bounds[1]/4.; 1512 else if (i==4) coord[j] = mod->bounds[1]/2.; 1513 else if (i==12) coord[j] = 1.57735026918963*mod->bounds[1]/2.; 1514 } 1515 } 1516 } 1517 } 1518 ierr = VecRestoreArray(coordinates,&coords);CHKERRQ(ierr); 1519 ierr = DMSetCoordinatesLocal(dm,coordinates);CHKERRQ(ierr); 1520 } 1521 } 1522 } 1523 ierr = DMViewFromOptions(dm, NULL, "-orig_dm_view");CHKERRQ(ierr); 1524 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1525 1526 /* set up BCs, functions, tags */ 1527 ierr = DMCreateLabel(dm, "Face Sets");CHKERRQ(ierr); 1528 mod->errorIndicator = ErrorIndicator_Simple; 1529 1530 { 1531 DM gdm; 1532 1533 ierr = DMPlexConstructGhostCells(dm, NULL, NULL, &gdm);CHKERRQ(ierr); 1534 ierr = DMDestroy(&dm);CHKERRQ(ierr); 1535 dm = gdm; 1536 ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr); 1537 } 1538 1539 ierr = PetscFVCreate(comm, &fvm);CHKERRQ(ierr); 1540 ierr = PetscFVSetFromOptions(fvm);CHKERRQ(ierr); 1541 ierr = PetscFVSetNumComponents(fvm, phys->dof);CHKERRQ(ierr); 1542 ierr = PetscFVSetSpatialDimension(fvm, dim);CHKERRQ(ierr); 1543 ierr = PetscObjectSetName((PetscObject) fvm,"");CHKERRQ(ierr); 1544 { 1545 PetscInt f, dof; 1546 for (f=0,dof=0; f < phys->nfields; f++) { 1547 PetscInt newDof = phys->field_desc[f].dof; 1548 1549 if (newDof == 1) { 1550 ierr = PetscFVSetComponentName(fvm,dof,phys->field_desc[f].name);CHKERRQ(ierr); 1551 } 1552 else { 1553 PetscInt j; 1554 1555 for (j = 0; j < newDof; j++) { 1556 char compName[256] = "Unknown"; 1557 1558 ierr = PetscSNPrintf(compName,sizeof(compName),"%s_%d",phys->field_desc[f].name,j);CHKERRQ(ierr); 1559 ierr = PetscFVSetComponentName(fvm,dof+j,compName);CHKERRQ(ierr); 1560 } 1561 } 1562 dof += newDof; 1563 } 1564 } 1565 /* FV is now structured with one field having all physics as components */ 1566 ierr = DMAddField(dm, NULL, (PetscObject) fvm);CHKERRQ(ierr); 1567 ierr = DMCreateDS(dm);CHKERRQ(ierr); 1568 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1569 ierr = PetscDSSetRiemannSolver(prob, 0, user->model->physics->riemann);CHKERRQ(ierr); 1570 ierr = PetscDSSetContext(prob, 0, user->model->physics);CHKERRQ(ierr); 1571 ierr = (*mod->setupbc)(dm, prob,phys);CHKERRQ(ierr); 1572 ierr = PetscDSSetFromOptions(prob);CHKERRQ(ierr); 1573 { 1574 char convType[256]; 1575 PetscBool flg; 1576 1577 ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr); 1578 ierr = PetscOptionsFList("-dm_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr); 1579 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1580 if (flg) { 1581 DM dmConv; 1582 1583 ierr = DMConvert(dm,convType,&dmConv);CHKERRQ(ierr); 1584 if (dmConv) { 1585 ierr = DMViewFromOptions(dmConv, NULL, "-dm_conv_view");CHKERRQ(ierr); 1586 ierr = DMDestroy(&dm);CHKERRQ(ierr); 1587 dm = dmConv; 1588 ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 1589 } 1590 } 1591 } 1592 1593 ierr = initializeTS(dm, user, &ts);CHKERRQ(ierr); 1594 1595 ierr = DMCreateGlobalVector(dm, &X);CHKERRQ(ierr); 1596 ierr = PetscObjectSetName((PetscObject) X, "solution");CHKERRQ(ierr); 1597 ierr = SetInitialCondition(dm, X, user);CHKERRQ(ierr); 1598 if (useAMR) { 1599 PetscInt adaptIter; 1600 1601 /* use no limiting when reconstructing gradients for adaptivity */ 1602 ierr = PetscFVGetLimiter(fvm, &limiter);CHKERRQ(ierr); 1603 ierr = PetscObjectReference((PetscObject) limiter);CHKERRQ(ierr); 1604 ierr = PetscLimiterCreate(PetscObjectComm((PetscObject) fvm), &noneLimiter);CHKERRQ(ierr); 1605 ierr = PetscLimiterSetType(noneLimiter, PETSCLIMITERNONE);CHKERRQ(ierr); 1606 1607 ierr = PetscFVSetLimiter(fvm, noneLimiter);CHKERRQ(ierr); 1608 for (adaptIter = 0; ; ++adaptIter) { 1609 PetscLogDouble bytes; 1610 TS tsNew = NULL; 1611 1612 ierr = PetscMemoryGetCurrentUsage(&bytes);CHKERRQ(ierr); 1613 ierr = PetscInfo2(ts, "refinement loop %D: memory used %g\n", adaptIter, bytes);CHKERRQ(ierr); 1614 ierr = DMViewFromOptions(dm, NULL, "-initial_dm_view");CHKERRQ(ierr); 1615 ierr = VecViewFromOptions(X, NULL, "-initial_vec_view");CHKERRQ(ierr); 1616 #if 0 1617 if (viewInitial) { 1618 PetscViewer viewer; 1619 char buf[256]; 1620 PetscBool isHDF5, isVTK; 1621 1622 ierr = PetscViewerCreate(comm,&viewer);CHKERRQ(ierr); 1623 ierr = PetscViewerSetType(viewer,PETSCVIEWERVTK);CHKERRQ(ierr); 1624 ierr = PetscViewerSetOptionsPrefix(viewer,"initial_");CHKERRQ(ierr); 1625 ierr = PetscViewerSetFromOptions(viewer);CHKERRQ(ierr); 1626 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&isHDF5);CHKERRQ(ierr); 1627 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isVTK);CHKERRQ(ierr); 1628 if (isHDF5) { 1629 ierr = PetscSNPrintf(buf, 256, "ex11-initial-%d.h5", adaptIter);CHKERRQ(ierr); 1630 } else if (isVTK) { 1631 ierr = PetscSNPrintf(buf, 256, "ex11-initial-%d.vtu", adaptIter);CHKERRQ(ierr); 1632 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTU);CHKERRQ(ierr); 1633 } 1634 ierr = PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);CHKERRQ(ierr); 1635 ierr = PetscViewerFileSetName(viewer,buf);CHKERRQ(ierr); 1636 if (isHDF5) { 1637 ierr = DMView(dm,viewer);CHKERRQ(ierr); 1638 ierr = PetscViewerFileSetMode(viewer,FILE_MODE_UPDATE);CHKERRQ(ierr); 1639 } 1640 ierr = VecView(X,viewer);CHKERRQ(ierr); 1641 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1642 } 1643 #endif 1644 1645 ierr = adaptToleranceFVM(fvm, ts, X, refineTag, coarsenTag, user, &tsNew, NULL);CHKERRQ(ierr); 1646 if (!tsNew) { 1647 break; 1648 } else { 1649 ierr = DMDestroy(&dm);CHKERRQ(ierr); 1650 ierr = VecDestroy(&X);CHKERRQ(ierr); 1651 ierr = TSDestroy(&ts);CHKERRQ(ierr); 1652 ts = tsNew; 1653 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1654 ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 1655 ierr = DMCreateGlobalVector(dm,&X);CHKERRQ(ierr); 1656 ierr = PetscObjectSetName((PetscObject) X, "solution");CHKERRQ(ierr); 1657 ierr = SetInitialCondition(dm, X, user);CHKERRQ(ierr); 1658 } 1659 } 1660 /* restore original limiter */ 1661 ierr = PetscFVSetLimiter(fvm, limiter);CHKERRQ(ierr); 1662 } 1663 1664 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 1665 if (vtkCellGeom) { 1666 DM dmCell; 1667 Vec cellgeom, partition; 1668 1669 ierr = DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL);CHKERRQ(ierr); 1670 ierr = OutputVTK(dm, "ex11-cellgeom.vtk", &viewer);CHKERRQ(ierr); 1671 ierr = VecView(cellgeom, viewer);CHKERRQ(ierr); 1672 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1673 ierr = CreatePartitionVec(dm, &dmCell, &partition);CHKERRQ(ierr); 1674 ierr = OutputVTK(dmCell, "ex11-partition.vtk", &viewer);CHKERRQ(ierr); 1675 ierr = VecView(partition, viewer);CHKERRQ(ierr); 1676 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1677 ierr = VecDestroy(&partition);CHKERRQ(ierr); 1678 ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 1679 } 1680 /* collect max maxspeed from all processes -- todo */ 1681 ierr = DMPlexGetGeometryFVM(plex, NULL, NULL, &minRadius);CHKERRQ(ierr); 1682 ierr = DMDestroy(&plex);CHKERRQ(ierr); 1683 ierr = MPI_Allreduce(&phys->maxspeed,&mod->maxspeed,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));CHKERRMPI(ierr); 1684 if (mod->maxspeed <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set maxspeed",physname); 1685 dt = cfl * minRadius / mod->maxspeed; 1686 ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 1687 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 1688 if (!useAMR) { 1689 ierr = TSSolve(ts,X);CHKERRQ(ierr); 1690 ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 1691 ierr = TSGetStepNumber(ts,&nsteps);CHKERRQ(ierr); 1692 } else { 1693 PetscReal finalTime; 1694 PetscInt adaptIter; 1695 TS tsNew = NULL; 1696 Vec solNew = NULL; 1697 1698 ierr = TSGetMaxTime(ts,&finalTime);CHKERRQ(ierr); 1699 ierr = TSSetMaxSteps(ts,adaptInterval);CHKERRQ(ierr); 1700 ierr = TSSolve(ts,X);CHKERRQ(ierr); 1701 ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 1702 ierr = TSGetStepNumber(ts,&nsteps);CHKERRQ(ierr); 1703 for (adaptIter = 0;ftime < finalTime;adaptIter++) { 1704 PetscLogDouble bytes; 1705 1706 ierr = PetscMemoryGetCurrentUsage(&bytes);CHKERRQ(ierr); 1707 ierr = PetscInfo2(ts, "AMR time step loop %D: memory used %g\n", adaptIter, bytes);CHKERRQ(ierr); 1708 ierr = PetscFVSetLimiter(fvm,noneLimiter);CHKERRQ(ierr); 1709 ierr = adaptToleranceFVM(fvm,ts,X,refineTag,coarsenTag,user,&tsNew,&solNew);CHKERRQ(ierr); 1710 ierr = PetscFVSetLimiter(fvm,limiter);CHKERRQ(ierr); 1711 if (tsNew) { 1712 ierr = PetscInfo(ts, "AMR used\n");CHKERRQ(ierr); 1713 ierr = DMDestroy(&dm);CHKERRQ(ierr); 1714 ierr = VecDestroy(&X);CHKERRQ(ierr); 1715 ierr = TSDestroy(&ts);CHKERRQ(ierr); 1716 ts = tsNew; 1717 X = solNew; 1718 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 1719 ierr = VecGetDM(X,&dm);CHKERRQ(ierr); 1720 ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 1721 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 1722 ierr = DMPlexGetGeometryFVM(dm, NULL, NULL, &minRadius);CHKERRQ(ierr); 1723 ierr = DMDestroy(&plex);CHKERRQ(ierr); 1724 ierr = MPI_Allreduce(&phys->maxspeed,&mod->maxspeed,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));CHKERRMPI(ierr); 1725 if (mod->maxspeed <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set maxspeed",physname); 1726 dt = cfl * minRadius / mod->maxspeed; 1727 ierr = TSSetStepNumber(ts,nsteps);CHKERRQ(ierr); 1728 ierr = TSSetTime(ts,ftime);CHKERRQ(ierr); 1729 ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 1730 } else { 1731 ierr = PetscInfo(ts, "AMR not used\n");CHKERRQ(ierr); 1732 } 1733 user->monitorStepOffset = nsteps; 1734 ierr = TSSetMaxSteps(ts,nsteps+adaptInterval);CHKERRQ(ierr); 1735 ierr = TSSolve(ts,X);CHKERRQ(ierr); 1736 ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 1737 ierr = TSGetStepNumber(ts,&nsteps);CHKERRQ(ierr); 1738 } 1739 } 1740 ierr = TSGetConvergedReason(ts,&reason);CHKERRQ(ierr); 1741 ierr = PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,nsteps);CHKERRQ(ierr); 1742 ierr = TSDestroy(&ts);CHKERRQ(ierr); 1743 1744 ierr = VecTaggerDestroy(&refineTag);CHKERRQ(ierr); 1745 ierr = VecTaggerDestroy(&coarsenTag);CHKERRQ(ierr); 1746 ierr = PetscFunctionListDestroy(&PhysicsList);CHKERRQ(ierr); 1747 ierr = PetscFunctionListDestroy(&PhysicsRiemannList_SW);CHKERRQ(ierr); 1748 ierr = FunctionalLinkDestroy(&user->model->functionalRegistry);CHKERRQ(ierr); 1749 ierr = PetscFree(user->model->functionalMonitored);CHKERRQ(ierr); 1750 ierr = PetscFree(user->model->functionalCall);CHKERRQ(ierr); 1751 ierr = PetscFree(user->model->physics->data);CHKERRQ(ierr); 1752 ierr = PetscFree(user->model->physics);CHKERRQ(ierr); 1753 ierr = PetscFree(user->model);CHKERRQ(ierr); 1754 ierr = PetscFree(user);CHKERRQ(ierr); 1755 ierr = VecDestroy(&X);CHKERRQ(ierr); 1756 ierr = PetscLimiterDestroy(&limiter);CHKERRQ(ierr); 1757 ierr = PetscLimiterDestroy(&noneLimiter);CHKERRQ(ierr); 1758 ierr = PetscFVDestroy(&fvm);CHKERRQ(ierr); 1759 ierr = DMDestroy(&dm);CHKERRQ(ierr); 1760 ierr = PetscFinalize(); 1761 return ierr; 1762 } 1763 1764 /* Godunov fluxs */ 1765 PetscScalar cvmgp_(PetscScalar *a, PetscScalar *b, PetscScalar *test) 1766 { 1767 /* System generated locals */ 1768 PetscScalar ret_val; 1769 1770 if (PetscRealPart(*test) > 0.) { 1771 goto L10; 1772 } 1773 ret_val = *b; 1774 return ret_val; 1775 L10: 1776 ret_val = *a; 1777 return ret_val; 1778 } /* cvmgp_ */ 1779 1780 PetscScalar cvmgm_(PetscScalar *a, PetscScalar *b, PetscScalar *test) 1781 { 1782 /* System generated locals */ 1783 PetscScalar ret_val; 1784 1785 if (PetscRealPart(*test) < 0.) { 1786 goto L10; 1787 } 1788 ret_val = *b; 1789 return ret_val; 1790 L10: 1791 ret_val = *a; 1792 return ret_val; 1793 } /* cvmgm_ */ 1794 1795 int riem1mdt( PetscScalar *gaml, PetscScalar *gamr, PetscScalar *rl, PetscScalar *pl, 1796 PetscScalar *uxl, PetscScalar *rr, PetscScalar *pr, 1797 PetscScalar *uxr, PetscScalar *rstarl, PetscScalar *rstarr, PetscScalar * 1798 pstar, PetscScalar *ustar) 1799 { 1800 /* Initialized data */ 1801 1802 static PetscScalar smallp = 1e-8; 1803 1804 /* System generated locals */ 1805 int i__1; 1806 PetscScalar d__1, d__2; 1807 1808 /* Local variables */ 1809 static int i0; 1810 static PetscScalar cl, cr, wl, zl, wr, zr, pst, durl, skpr1, skpr2; 1811 static int iwave; 1812 static PetscScalar gascl4, gascr4, cstarl, dpstar, cstarr; 1813 /* static PetscScalar csqrl, csqrr, gascl1, gascl2, gascl3, gascr1, gascr2, gascr3; */ 1814 static int iterno; 1815 static PetscScalar ustarl, ustarr, rarepr1, rarepr2; 1816 1817 /* gascl1 = *gaml - 1.; */ 1818 /* gascl2 = (*gaml + 1.) * .5; */ 1819 /* gascl3 = gascl2 / *gaml; */ 1820 gascl4 = 1. / (*gaml - 1.); 1821 1822 /* gascr1 = *gamr - 1.; */ 1823 /* gascr2 = (*gamr + 1.) * .5; */ 1824 /* gascr3 = gascr2 / *gamr; */ 1825 gascr4 = 1. / (*gamr - 1.); 1826 iterno = 10; 1827 /* find pstar: */ 1828 cl = PetscSqrtScalar(*gaml * *pl / *rl); 1829 cr = PetscSqrtScalar(*gamr * *pr / *rr); 1830 wl = *rl * cl; 1831 wr = *rr * cr; 1832 /* csqrl = wl * wl; */ 1833 /* csqrr = wr * wr; */ 1834 *pstar = (wl * *pr + wr * *pl) / (wl + wr); 1835 *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp)); 1836 pst = *pl / *pr; 1837 skpr1 = cr * (pst - 1.) * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst))); 1838 d__1 = (*gamr - 1.) / (*gamr * 2.); 1839 rarepr2 = gascr4 * 2. * cr * (1. - PetscPowScalar(pst, d__1)); 1840 pst = *pr / *pl; 1841 skpr2 = cl * (pst - 1.) * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst))); 1842 d__1 = (*gaml - 1.) / (*gaml * 2.); 1843 rarepr1 = gascl4 * 2. * cl * (1. - PetscPowScalar(pst, d__1)); 1844 durl = *uxr - *uxl; 1845 if (PetscRealPart(*pr) < PetscRealPart(*pl)) { 1846 if (PetscRealPart(durl) >= PetscRealPart(rarepr1)) { 1847 iwave = 100; 1848 } else if (PetscRealPart(durl) <= PetscRealPart(-skpr1)) { 1849 iwave = 300; 1850 } else { 1851 iwave = 400; 1852 } 1853 } else { 1854 if (PetscRealPart(durl) >= PetscRealPart(rarepr2)) { 1855 iwave = 100; 1856 } else if (PetscRealPart(durl) <= PetscRealPart(-skpr2)) { 1857 iwave = 300; 1858 } else { 1859 iwave = 200; 1860 } 1861 } 1862 if (iwave == 100) { 1863 /* 1-wave: rarefaction wave, 3-wave: rarefaction wave */ 1864 /* case (100) */ 1865 i__1 = iterno; 1866 for (i0 = 1; i0 <= i__1; ++i0) { 1867 d__1 = *pstar / *pl; 1868 d__2 = 1. / *gaml; 1869 *rstarl = *rl * PetscPowScalar(d__1, d__2); 1870 cstarl = PetscSqrtScalar(*gaml * *pstar / *rstarl); 1871 ustarl = *uxl - gascl4 * 2. * (cstarl - cl); 1872 zl = *rstarl * cstarl; 1873 d__1 = *pstar / *pr; 1874 d__2 = 1. / *gamr; 1875 *rstarr = *rr * PetscPowScalar(d__1, d__2); 1876 cstarr = PetscSqrtScalar(*gamr * *pstar / *rstarr); 1877 ustarr = *uxr + gascr4 * 2. * (cstarr - cr); 1878 zr = *rstarr * cstarr; 1879 dpstar = zl * zr * (ustarr - ustarl) / (zl + zr); 1880 *pstar -= dpstar; 1881 *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp)); 1882 if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) { 1883 #if 0 1884 break; 1885 #endif 1886 } 1887 } 1888 /* 1-wave: shock wave, 3-wave: rarefaction wave */ 1889 } else if (iwave == 200) { 1890 /* case (200) */ 1891 i__1 = iterno; 1892 for (i0 = 1; i0 <= i__1; ++i0) { 1893 pst = *pstar / *pl; 1894 ustarl = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst))); 1895 zl = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst); 1896 d__1 = *pstar / *pr; 1897 d__2 = 1. / *gamr; 1898 *rstarr = *rr * PetscPowScalar(d__1, d__2); 1899 cstarr = PetscSqrtScalar(*gamr * *pstar / *rstarr); 1900 zr = *rstarr * cstarr; 1901 ustarr = *uxr + gascr4 * 2. * (cstarr - cr); 1902 dpstar = zl * zr * (ustarr - ustarl) / (zl + zr); 1903 *pstar -= dpstar; 1904 *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp)); 1905 if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) { 1906 #if 0 1907 break; 1908 #endif 1909 } 1910 } 1911 /* 1-wave: shock wave, 3-wave: shock */ 1912 } else if (iwave == 300) { 1913 /* case (300) */ 1914 i__1 = iterno; 1915 for (i0 = 1; i0 <= i__1; ++i0) { 1916 pst = *pstar / *pl; 1917 ustarl = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst))); 1918 zl = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst); 1919 pst = *pstar / *pr; 1920 ustarr = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst))); 1921 zr = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst); 1922 dpstar = zl * zr * (ustarr - ustarl) / (zl + zr); 1923 *pstar -= dpstar; 1924 *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp)); 1925 if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) { 1926 #if 0 1927 break; 1928 #endif 1929 } 1930 } 1931 /* 1-wave: rarefaction wave, 3-wave: shock */ 1932 } else if (iwave == 400) { 1933 /* case (400) */ 1934 i__1 = iterno; 1935 for (i0 = 1; i0 <= i__1; ++i0) { 1936 d__1 = *pstar / *pl; 1937 d__2 = 1. / *gaml; 1938 *rstarl = *rl * PetscPowScalar(d__1, d__2); 1939 cstarl = PetscSqrtScalar(*gaml * *pstar / *rstarl); 1940 ustarl = *uxl - gascl4 * 2. * (cstarl - cl); 1941 zl = *rstarl * cstarl; 1942 pst = *pstar / *pr; 1943 ustarr = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst))); 1944 zr = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst); 1945 dpstar = zl * zr * (ustarr - ustarl) / (zl + zr); 1946 *pstar -= dpstar; 1947 *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp)); 1948 if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) { 1949 #if 0 1950 break; 1951 #endif 1952 } 1953 } 1954 } 1955 1956 *ustar = (zl * ustarr + zr * ustarl) / (zl + zr); 1957 if (PetscRealPart(*pstar) > PetscRealPart(*pl)) { 1958 pst = *pstar / *pl; 1959 *rstarl = ((*gaml + 1.) * pst + *gaml - 1.) / ((*gaml - 1.) * pst + * 1960 gaml + 1.) * *rl; 1961 } 1962 if (PetscRealPart(*pstar) > PetscRealPart(*pr)) { 1963 pst = *pstar / *pr; 1964 *rstarr = ((*gamr + 1.) * pst + *gamr - 1.) / ((*gamr - 1.) * pst + * 1965 gamr + 1.) * *rr; 1966 } 1967 return iwave; 1968 } 1969 1970 PetscScalar sign(PetscScalar x) 1971 { 1972 if (PetscRealPart(x) > 0) return 1.0; 1973 if (PetscRealPart(x) < 0) return -1.0; 1974 return 0.0; 1975 } 1976 /* Riemann Solver */ 1977 /* -------------------------------------------------------------------- */ 1978 int riemannsolver(PetscScalar *xcen, PetscScalar *xp, 1979 PetscScalar *dtt, PetscScalar *rl, PetscScalar *uxl, PetscScalar *pl, 1980 PetscScalar *utl, PetscScalar *ubl, PetscScalar *gaml, PetscScalar *rho1l, 1981 PetscScalar *rr, PetscScalar *uxr, PetscScalar *pr, PetscScalar *utr, 1982 PetscScalar *ubr, PetscScalar *gamr, PetscScalar *rho1r, PetscScalar *rx, 1983 PetscScalar *uxm, PetscScalar *px, PetscScalar *utx, PetscScalar *ubx, 1984 PetscScalar *gam, PetscScalar *rho1) 1985 { 1986 /* System generated locals */ 1987 PetscScalar d__1, d__2; 1988 1989 /* Local variables */ 1990 static PetscScalar s, c0, p0, r0, u0, w0, x0, x2, ri, cx, sgn0, wsp0, gasc1, gasc2, gasc3, gasc4; 1991 static PetscScalar cstar, pstar, rstar, ustar, xstar, wspst, ushock, streng, rstarl, rstarr, rstars; 1992 int iwave; 1993 1994 if (*rl == *rr && *pr == *pl && *uxl == *uxr && *gaml == *gamr) { 1995 *rx = *rl; 1996 *px = *pl; 1997 *uxm = *uxl; 1998 *gam = *gaml; 1999 x2 = *xcen + *uxm * *dtt; 2000 2001 if (PetscRealPart(*xp) >= PetscRealPart(x2)) { 2002 *utx = *utr; 2003 *ubx = *ubr; 2004 *rho1 = *rho1r; 2005 } else { 2006 *utx = *utl; 2007 *ubx = *ubl; 2008 *rho1 = *rho1l; 2009 } 2010 return 0; 2011 } 2012 iwave = riem1mdt(gaml, gamr, rl, pl, uxl, rr, pr, uxr, &rstarl, &rstarr, &pstar, &ustar); 2013 2014 x2 = *xcen + ustar * *dtt; 2015 d__1 = *xp - x2; 2016 sgn0 = sign(d__1); 2017 /* x is in 3-wave if sgn0 = 1 */ 2018 /* x is in 1-wave if sgn0 = -1 */ 2019 r0 = cvmgm_(rl, rr, &sgn0); 2020 p0 = cvmgm_(pl, pr, &sgn0); 2021 u0 = cvmgm_(uxl, uxr, &sgn0); 2022 *gam = cvmgm_(gaml, gamr, &sgn0); 2023 gasc1 = *gam - 1.; 2024 gasc2 = (*gam + 1.) * .5; 2025 gasc3 = gasc2 / *gam; 2026 gasc4 = 1. / (*gam - 1.); 2027 c0 = PetscSqrtScalar(*gam * p0 / r0); 2028 streng = pstar - p0; 2029 w0 = *gam * r0 * p0 * (gasc3 * streng / p0 + 1.); 2030 rstars = r0 / (1. - r0 * streng / w0); 2031 d__1 = p0 / pstar; 2032 d__2 = -1. / *gam; 2033 rstarr = r0 * PetscPowScalar(d__1, d__2); 2034 rstar = cvmgm_(&rstarr, &rstars, &streng); 2035 w0 = PetscSqrtScalar(w0); 2036 cstar = PetscSqrtScalar(*gam * pstar / rstar); 2037 wsp0 = u0 + sgn0 * c0; 2038 wspst = ustar + sgn0 * cstar; 2039 ushock = ustar + sgn0 * w0 / rstar; 2040 wspst = cvmgp_(&ushock, &wspst, &streng); 2041 wsp0 = cvmgp_(&ushock, &wsp0, &streng); 2042 x0 = *xcen + wsp0 * *dtt; 2043 xstar = *xcen + wspst * *dtt; 2044 /* using gas formula to evaluate rarefaction wave */ 2045 /* ri : reiman invariant */ 2046 ri = u0 - sgn0 * 2. * gasc4 * c0; 2047 cx = sgn0 * .5 * gasc1 / gasc2 * ((*xp - *xcen) / *dtt - ri); 2048 *uxm = ri + sgn0 * 2. * gasc4 * cx; 2049 s = p0 / PetscPowScalar(r0, *gam); 2050 d__1 = cx * cx / (*gam * s); 2051 *rx = PetscPowScalar(d__1, gasc4); 2052 *px = cx * cx * *rx / *gam; 2053 d__1 = sgn0 * (x0 - *xp); 2054 *rx = cvmgp_(rx, &r0, &d__1); 2055 d__1 = sgn0 * (x0 - *xp); 2056 *px = cvmgp_(px, &p0, &d__1); 2057 d__1 = sgn0 * (x0 - *xp); 2058 *uxm = cvmgp_(uxm, &u0, &d__1); 2059 d__1 = sgn0 * (xstar - *xp); 2060 *rx = cvmgm_(rx, &rstar, &d__1); 2061 d__1 = sgn0 * (xstar - *xp); 2062 *px = cvmgm_(px, &pstar, &d__1); 2063 d__1 = sgn0 * (xstar - *xp); 2064 *uxm = cvmgm_(uxm, &ustar, &d__1); 2065 if (PetscRealPart(*xp) >= PetscRealPart(x2)) { 2066 *utx = *utr; 2067 *ubx = *ubr; 2068 *rho1 = *rho1r; 2069 } else { 2070 *utx = *utl; 2071 *ubx = *ubl; 2072 *rho1 = *rho1l; 2073 } 2074 return iwave; 2075 } 2076 int godunovflux( const PetscScalar *ul, const PetscScalar *ur, 2077 PetscScalar *flux, const PetscReal *nn, const int *ndim, 2078 const PetscReal *gamma) 2079 { 2080 /* System generated locals */ 2081 int i__1,iwave; 2082 PetscScalar d__1, d__2, d__3; 2083 2084 /* Local variables */ 2085 static int k; 2086 static PetscScalar bn[3], fn, ft, tg[3], pl, rl, pm, pr, rr, xp, ubl, ubm, 2087 ubr, dtt, unm, tmp, utl, utm, uxl, utr, uxr, gaml, gamm, gamr, 2088 xcen, rhom, rho1l, rho1m, rho1r; 2089 2090 /* Function Body */ 2091 xcen = 0.; 2092 xp = 0.; 2093 i__1 = *ndim; 2094 for (k = 1; k <= i__1; ++k) { 2095 tg[k - 1] = 0.; 2096 bn[k - 1] = 0.; 2097 } 2098 dtt = 1.; 2099 if (*ndim == 3) { 2100 if (nn[0] == 0. && nn[1] == 0.) { 2101 tg[0] = 1.; 2102 } else { 2103 tg[0] = -nn[1]; 2104 tg[1] = nn[0]; 2105 } 2106 /* tmp=dsqrt(tg(1)**2+tg(2)**2) */ 2107 /* tg=tg/tmp */ 2108 bn[0] = -nn[2] * tg[1]; 2109 bn[1] = nn[2] * tg[0]; 2110 bn[2] = nn[0] * tg[1] - nn[1] * tg[0]; 2111 /* Computing 2nd power */ 2112 d__1 = bn[0]; 2113 /* Computing 2nd power */ 2114 d__2 = bn[1]; 2115 /* Computing 2nd power */ 2116 d__3 = bn[2]; 2117 tmp = PetscSqrtScalar(d__1 * d__1 + d__2 * d__2 + d__3 * d__3); 2118 i__1 = *ndim; 2119 for (k = 1; k <= i__1; ++k) { 2120 bn[k - 1] /= tmp; 2121 } 2122 } else if (*ndim == 2) { 2123 tg[0] = -nn[1]; 2124 tg[1] = nn[0]; 2125 /* tmp=dsqrt(tg(1)**2+tg(2)**2) */ 2126 /* tg=tg/tmp */ 2127 bn[0] = 0.; 2128 bn[1] = 0.; 2129 bn[2] = 1.; 2130 } 2131 rl = ul[0]; 2132 rr = ur[0]; 2133 uxl = 0.; 2134 uxr = 0.; 2135 utl = 0.; 2136 utr = 0.; 2137 ubl = 0.; 2138 ubr = 0.; 2139 i__1 = *ndim; 2140 for (k = 1; k <= i__1; ++k) { 2141 uxl += ul[k] * nn[k-1]; 2142 uxr += ur[k] * nn[k-1]; 2143 utl += ul[k] * tg[k - 1]; 2144 utr += ur[k] * tg[k - 1]; 2145 ubl += ul[k] * bn[k - 1]; 2146 ubr += ur[k] * bn[k - 1]; 2147 } 2148 uxl /= rl; 2149 uxr /= rr; 2150 utl /= rl; 2151 utr /= rr; 2152 ubl /= rl; 2153 ubr /= rr; 2154 2155 gaml = *gamma; 2156 gamr = *gamma; 2157 /* Computing 2nd power */ 2158 d__1 = uxl; 2159 /* Computing 2nd power */ 2160 d__2 = utl; 2161 /* Computing 2nd power */ 2162 d__3 = ubl; 2163 pl = (*gamma - 1.) * (ul[*ndim + 1] - rl * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3)); 2164 /* Computing 2nd power */ 2165 d__1 = uxr; 2166 /* Computing 2nd power */ 2167 d__2 = utr; 2168 /* Computing 2nd power */ 2169 d__3 = ubr; 2170 pr = (*gamma - 1.) * (ur[*ndim + 1] - rr * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3)); 2171 rho1l = rl; 2172 rho1r = rr; 2173 2174 iwave = riemannsolver(&xcen, &xp, &dtt, &rl, &uxl, &pl, &utl, &ubl, &gaml, & 2175 rho1l, &rr, &uxr, &pr, &utr, &ubr, &gamr, &rho1r, &rhom, &unm, & 2176 pm, &utm, &ubm, &gamm, &rho1m); 2177 2178 flux[0] = rhom * unm; 2179 fn = rhom * unm * unm + pm; 2180 ft = rhom * unm * utm; 2181 /* flux(2)=fn*nn(1)+ft*nn(2) */ 2182 /* flux(3)=fn*tg(1)+ft*tg(2) */ 2183 flux[1] = fn * nn[0] + ft * tg[0]; 2184 flux[2] = fn * nn[1] + ft * tg[1]; 2185 /* flux(2)=rhom*unm*(unm)+pm */ 2186 /* flux(3)=rhom*(unm)*utm */ 2187 if (*ndim == 3) { 2188 flux[3] = rhom * unm * ubm; 2189 } 2190 flux[*ndim + 1] = (rhom * .5 * (unm * unm + utm * utm + ubm * ubm) + gamm / (gamm - 1.) * pm) * unm; 2191 return iwave; 2192 } /* godunovflux_ */ 2193 2194 /* Subroutine to set up the initial conditions for the */ 2195 /* Shock Interface interaction or linear wave (Ravi Samtaney,Mark Adams). */ 2196 /* ----------------------------------------------------------------------- */ 2197 int projecteqstate(PetscReal wc[], const PetscReal ueq[], PetscReal lv[][3]) 2198 { 2199 int j,k; 2200 /* Wc=matmul(lv,Ueq) 3 vars */ 2201 for (k = 0; k < 3; ++k) { 2202 wc[k] = 0.; 2203 for (j = 0; j < 3; ++j) { 2204 wc[k] += lv[k][j]*ueq[j]; 2205 } 2206 } 2207 return 0; 2208 } 2209 /* ----------------------------------------------------------------------- */ 2210 int projecttoprim(PetscReal v[], const PetscReal wc[], PetscReal rv[][3]) 2211 { 2212 int k,j; 2213 /* V=matmul(rv,WC) 3 vars */ 2214 for (k = 0; k < 3; ++k) { 2215 v[k] = 0.; 2216 for (j = 0; j < 3; ++j) { 2217 v[k] += rv[k][j]*wc[j]; 2218 } 2219 } 2220 return 0; 2221 } 2222 /* ---------------------------------------------------------------------- */ 2223 int eigenvectors(PetscReal rv[][3], PetscReal lv[][3], const PetscReal ueq[], PetscReal gamma) 2224 { 2225 int j,k; 2226 PetscReal rho,csnd,p0; 2227 /* PetscScalar u; */ 2228 2229 for (k = 0; k < 3; ++k) for (j = 0; j < 3; ++j) { lv[k][j] = 0.; rv[k][j] = 0.; } 2230 rho = ueq[0]; 2231 /* u = ueq[1]; */ 2232 p0 = ueq[2]; 2233 csnd = PetscSqrtReal(gamma * p0 / rho); 2234 lv[0][1] = rho * .5; 2235 lv[0][2] = -.5 / csnd; 2236 lv[1][0] = csnd; 2237 lv[1][2] = -1. / csnd; 2238 lv[2][1] = rho * .5; 2239 lv[2][2] = .5 / csnd; 2240 rv[0][0] = -1. / csnd; 2241 rv[1][0] = 1. / rho; 2242 rv[2][0] = -csnd; 2243 rv[0][1] = 1. / csnd; 2244 rv[0][2] = 1. / csnd; 2245 rv[1][2] = 1. / rho; 2246 rv[2][2] = csnd; 2247 return 0; 2248 } 2249 2250 int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx) 2251 { 2252 PetscReal p0,u0,wcp[3],wc[3]; 2253 PetscReal lv[3][3]; 2254 PetscReal vp[3]; 2255 PetscReal rv[3][3]; 2256 PetscReal eps, ueq[3], rho0, twopi; 2257 2258 /* Function Body */ 2259 twopi = 2.*PETSC_PI; 2260 eps = 1e-4; /* perturbation */ 2261 rho0 = 1e3; /* density of water */ 2262 p0 = 101325.; /* init pressure of 1 atm (?) */ 2263 u0 = 0.; 2264 ueq[0] = rho0; 2265 ueq[1] = u0; 2266 ueq[2] = p0; 2267 /* Project initial state to characteristic variables */ 2268 eigenvectors(rv, lv, ueq, gamma); 2269 projecteqstate(wc, ueq, lv); 2270 wcp[0] = wc[0]; 2271 wcp[1] = wc[1]; 2272 wcp[2] = wc[2] + eps * PetscCosReal(coord[0] * 2. * twopi / Lx); 2273 projecttoprim(vp, wcp, rv); 2274 ux->r = vp[0]; /* density */ 2275 ux->ru[0] = vp[0] * vp[1]; /* x momentum */ 2276 ux->ru[1] = 0.; 2277 #if defined DIM > 2 2278 if (dim>2) ux->ru[2] = 0.; 2279 #endif 2280 /* E = rho * e + rho * v^2/2 = p/(gam-1) + rho*v^2/2 */ 2281 ux->E = vp[2]/(gamma - 1.) + 0.5*vp[0]*vp[1]*vp[1]; 2282 return 0; 2283 } 2284 2285 /*TEST 2286 2287 testset: 2288 args: -dm_plex_adj_cone -dm_plex_adj_closure 0 2289 2290 test: 2291 suffix: adv_2d_tri_0 2292 requires: triangle 2293 TODO: how did this ever get in main when there is no support for this 2294 args: -ufv_vtk_interval 0 -simplex -dm_refine 3 -dm_plex_faces 1,1 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3 2295 2296 test: 2297 suffix: adv_2d_tri_1 2298 requires: triangle 2299 TODO: how did this ever get in main when there is no support for this 2300 args: -ufv_vtk_interval 0 -simplex -dm_refine 5 -dm_plex_faces 1,1 -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1 2301 2302 test: 2303 suffix: tut_1 2304 requires: exodusii 2305 nsize: 1 2306 args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo 2307 2308 test: 2309 suffix: tut_2 2310 requires: exodusii 2311 nsize: 1 2312 args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw 2313 2314 test: 2315 suffix: tut_3 2316 requires: exodusii 2317 nsize: 4 2318 args: -dm_distribute -dm_distribute_overlap 1 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -monitor Error -advect_sol_type bump -petscfv_type leastsquares -petsclimiter_type sin 2319 2320 test: 2321 suffix: tut_4 2322 requires: exodusii 2323 nsize: 4 2324 args: -dm_distribute -dm_distribute_overlap 1 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -physics sw -monitor Height,Energy -petscfv_type leastsquares -petsclimiter_type minmod 2325 2326 testset: 2327 args: -dm_plex_adj_cone -dm_plex_adj_closure 0 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 2328 2329 # 2D Advection 0-10 2330 test: 2331 suffix: 0 2332 requires: exodusii 2333 args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo 2334 2335 test: 2336 suffix: 1 2337 requires: exodusii 2338 args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo 2339 2340 test: 2341 suffix: 2 2342 requires: exodusii 2343 nsize: 2 2344 args: -dm_distribute -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo 2345 2346 test: 2347 suffix: 3 2348 requires: exodusii 2349 nsize: 2 2350 args: -dm_distribute -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo 2351 2352 test: 2353 suffix: 4 2354 requires: exodusii 2355 nsize: 8 2356 args: -dm_distribute -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo 2357 2358 test: 2359 suffix: 5 2360 requires: exodusii 2361 args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw -ts_adapt_reject_safety 1 2362 2363 test: 2364 suffix: 7 2365 requires: exodusii 2366 args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1 2367 2368 test: 2369 suffix: 8 2370 requires: exodusii 2371 nsize: 2 2372 args: -dm_distribute -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1 2373 2374 test: 2375 suffix: 9 2376 requires: exodusii 2377 nsize: 8 2378 args: -dm_distribute -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1 2379 2380 test: 2381 suffix: 10 2382 requires: exodusii 2383 args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo 2384 2385 # 2D Shallow water 2386 test: 2387 suffix: sw_0 2388 requires: exodusii 2389 args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -bc_wall 100,101 -physics sw -ufv_cfl 5 -petscfv_type leastsquares -petsclimiter_type sin -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 -monitor height,energy 2390 2391 test: 2392 suffix: sw_hll 2393 args: -ufv_vtk_interval 0 -bc_wall 1,2,3,4 -physics sw -ufv_cfl 3 -petscfv_type leastsquares -petsclimiter_type sin -ts_max_steps 5 -ts_ssp_type rks2 -ts_ssp_nstages 10 -monitor height,energy -grid_bounds 0,5,0,5 -dm_plex_box_faces 25,25 -sw_riemann hll 2394 2395 # 2D Advection: p4est 2396 test: 2397 suffix: p4est_advec_2d 2398 requires: p4est 2399 args: -ufv_vtk_interval 0 -dm_type p4est -dm_forest_minimum_refinement 1 -dm_forest_initial_refinement 2 -dm_p4est_refine_pattern hash -dm_forest_maximum_refinement 5 2400 2401 # Advection in a box 2402 test: 2403 suffix: adv_2d_quad_0 2404 args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3 2405 2406 test: 2407 suffix: adv_2d_quad_1 2408 args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1 2409 timeoutfactor: 3 2410 2411 test: 2412 suffix: adv_2d_quad_p4est_0 2413 requires: p4est 2414 args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3 2415 2416 test: 2417 suffix: adv_2d_quad_p4est_1 2418 requires: p4est 2419 args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1 2420 timeoutfactor: 3 2421 2422 test: 2423 suffix: adv_2d_quad_p4est_adapt_0 2424 requires: p4est !__float128 #broken for quad precision 2425 args: -ufv_vtk_interval 0 -dm_refine 3 -dm_type p4est -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1 -ufv_use_amr -refine_vec_tagger_box 0.005,inf -coarsen_vec_tagger_box 0,1.e-5 -petscfv_type leastsquares -ts_max_time 0.01 2426 timeoutfactor: 3 2427 2428 test: 2429 suffix: adv_0 2430 requires: exodusii 2431 args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -bc_inflow 100,101,200 -bc_outflow 201 2432 2433 test: 2434 suffix: shock_0 2435 requires: p4est !single !complex 2436 args: -dm_plex_box_faces 2,1 -grid_bounds -1,1.,0.,1 -grid_skew_60 \ 2437 -dm_type p4est -dm_forest_partition_overlap 1 -dm_forest_maximum_refinement 6 -dm_forest_minimum_refinement 2 -dm_forest_initial_refinement 2 \ 2438 -ufv_use_amr -refine_vec_tagger_box 0.5,inf -coarsen_vec_tagger_box 0,1.e-2 -refine_tag_view -coarsen_tag_view \ 2439 -bc_wall 1,2,3,4 -physics euler -eu_type iv_shock -ufv_cfl 10 -eu_alpha 60. -eu_gamma 1.4 -eu_amach 2.02 -eu_rho2 3. \ 2440 -petscfv_type leastsquares -petsclimiter_type minmod -petscfv_compute_gradients 0 \ 2441 -ts_max_time 0.5 -ts_ssp_type rks2 -ts_ssp_nstages 10 \ 2442 -ufv_vtk_basename ${wPETSC_DIR}/ex11 -ufv_vtk_interval 0 -monitor density,energy 2443 timeoutfactor: 3 2444 2445 # Test GLVis visualization of PetscFV fields 2446 test: 2447 suffix: glvis_adv_2d_tet 2448 args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 \ 2449 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic 0 \ 2450 -ts_monitor_solution glvis: -ts_max_steps 0 2451 2452 test: 2453 suffix: glvis_adv_2d_quad 2454 args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 -bc_inflow 1,2,4 -bc_outflow 3 \ 2455 -dm_refine 5 -dm_plex_separate_marker \ 2456 -ts_monitor_solution glvis: -ts_max_steps 0 2457 2458 TEST*/ 2459