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