Actual source code: ex11.c
petsc-3.6.1 2015-08-06
1: static char help[] = "Second Order TVD Finite Volume Example.\n" ;
We use a second order TVD finite volume method to evolve a system of PDEs. Our simple upwinded residual evaluation loops
over all mesh faces and uses a Riemann solver to produce the flux given the face geometry and cell values,
\begin{equation}
f_i = \mathrm{riemann}(\mathrm{phys}, p_\mathrm{centroid}, \hat n, x^L, x^R)
\end{equation}
and then update the cell values given the cell volume.
\begin{eqnarray}
f^L_i &-=& \frac{f_i}{vol^L} \\
f^R_i &+=& \frac{f_i}{vol^R}
\end{eqnarray}
As an example, we can consider the shallow water wave equation,
\begin{eqnarray}
h_t + \nabla\cdot \left( uh \right) &=& 0 \\
(uh)_t + \nabla\cdot \left( u\otimes uh + \frac{g h^2}{2} I \right) &=& 0
\end{eqnarray}
where $h$ is wave height, $u$ is wave velocity, and $g$ is the acceleration due to gravity.
A representative Riemann solver for the shallow water equations is given in the PhysicsRiemann_SW() function,
\begin{eqnarray}
f^{L,R}_h &=& uh^{L,R} \cdot \hat n \\
f^{L,R}_{uh} &=& \frac{f^{L,R}_h}{h^{L,R}} uh^{L,R} + g (h^{L,R})^2 \hat n \\
c^{L,R} &=& \sqrt{g h^{L,R}} \\
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) \\
f_i &=& \frac{A_\mathrm{face}}{2} \left( f^L_i + f^R_i + s \left( x^L_i - x^R_i \right) \right)
\end{eqnarray}
where $c$ is the local gravity wave speed and $f_i$ is a Rusanov flux.
The more sophisticated residual evaluation in RHSFunctionLocal_LS() uses a least-squares fit to a quadratic polynomial
over a neighborhood of the given element.
The mesh is read in from an ExodusII file, usually generated by Cubit.
37: #include <petscdmplex.h>
38: #include <petscds.h>
39: #include <petscts.h>
40: #include <petscsf.h> /* For SplitFaces() */
42: #define DIM 2 /* Geometric dimension */
43: #define ALEN(a) (sizeof(a)/sizeof((a)[0]))
45: static PetscFunctionList PhysicsList;
47: /* Represents continuum physical equations. */
48: typedef struct _n_Physics *Physics;
50: /* Physical model includes boundary conditions, initial conditions, and functionals of interest. It is
51: * discretization-independent, but its members depend on the scenario being solved. */
52: typedef struct _n_Model *Model;
54: /* 'User' implements a discretization of a continuous model. */
55: typedef struct _n_User *User;
56: typedef PetscErrorCode (*SolutionFunction)(Model,PetscReal ,const PetscReal *,PetscScalar *,void*) ;
57: typedef PetscErrorCode (*FunctionalFunction)(Model,PetscReal ,const PetscReal *,const PetscScalar *,PetscReal *,void*) ;
58: typedef PetscErrorCode (*SetupFields)(Physics,PetscSection ) ;
59: static PetscErrorCode ModelSolutionSetDefault(Model,SolutionFunction,void*) ;
60: static PetscErrorCode ModelFunctionalRegister(Model,const char*,PetscInt *,FunctionalFunction,void*) ;
61: static PetscErrorCode OutputVTK(DM ,const char*,PetscViewer *) ;
63: struct FieldDescription {
64: const char *name;
65: PetscInt dof;
66: };
68: typedef struct _n_FunctionalLink *FunctionalLink;
69: struct _n_FunctionalLink {
70: char *name;
71: FunctionalFunction func;
72: void *ctx;
73: PetscInt offset;
74: FunctionalLink next;
75: };
77: struct _n_Physics {
78: PetscRiemannFunc riemann;
79: PetscInt dof; /* number of degrees of freedom per cell */
80: PetscReal maxspeed; /* kludge to pick initial time step, need to add monitoring and step control */
81: void *data;
82: PetscInt nfields;
83: const struct FieldDescription *field_desc;
84: };
86: struct _n_Model {
87: MPI_Comm comm; /* Does not do collective communicaton, but some error conditions can be collective */
88: Physics physics;
89: FunctionalLink functionalRegistry;
90: PetscInt maxComputed;
91: PetscInt numMonitored;
92: FunctionalLink *functionalMonitored;
93: PetscInt numCall;
94: FunctionalLink *functionalCall;
95: SolutionFunction solution;
96: void *solutionctx;
97: PetscReal maxspeed; /* estimate of global maximum speed (for CFL calculation) */
98: };
100: struct _n_User {
101: PetscInt numSplitFaces;
102: PetscInt vtkInterval; /* For monitor */
103: Model model;
104: };
106: PETSC_STATIC_INLINE PetscScalar DotDIM(const PetscScalar *x,const PetscScalar *y)
107: {
108: PetscInt i;
109: PetscScalar prod=0.0;
111: for (i=0; i<DIM; i++) prod += x[i]*y[i];
112: return prod;
113: }
114: PETSC_STATIC_INLINE PetscReal NormDIM(const PetscScalar *x) { return PetscSqrtReal(PetscAbsScalar(DotDIM(x,x))); }
115: PETSC_STATIC_INLINE void axDIM(const PetscScalar a,PetscScalar *x)
116: {
117: PetscInt i;
118: for (i=0; i<DIM; i++) x[i] *= a;
119: }
120: PETSC_STATIC_INLINE void waxDIM(const PetscScalar a,const PetscScalar *x, PetscScalar *w)
121: {
122: PetscInt i;
123: for (i=0; i<DIM; i++) w[i] = x[i]*a;
124: }
125: PETSC_STATIC_INLINE void NormalSplitDIM(const PetscReal *n,const PetscScalar *x,PetscScalar *xn,PetscScalar *xt)
126: { /* Split x into normal and tangential components */
127: PetscInt i;
128: PetscScalar c;
129: c = DotDIM(x,n)/DotDIM(n,n);
130: for (i=0; i<DIM; i++) {
131: xn[i] = c*n[i];
132: xt[i] = x[i]-xn[i];
133: }
134: }
136: PETSC_STATIC_INLINE PetscScalar Dot2(const PetscScalar *x,const PetscScalar *y) { return x[0]*y[0] + x[1]*y[1];}
137: PETSC_STATIC_INLINE PetscReal Norm2(const PetscScalar *x) { return PetscSqrtReal(PetscAbsScalar(Dot2(x,x)));}
138: PETSC_STATIC_INLINE void Normalize2(PetscScalar *x) { PetscReal a = 1./Norm2(x); x[0] *= a; x[1] *= a; }
139: PETSC_STATIC_INLINE void Waxpy2(PetscScalar a,const PetscScalar *x,const PetscScalar *y,PetscScalar *w) { w[0] = a*x[0] + y[0]; w[1] = a*x[1] + y[1]; }
140: PETSC_STATIC_INLINE void Scale2(PetscScalar a,const PetscScalar *x,PetscScalar *y) { y[0] = a*x[0]; y[1] = a*x[1]; }
142: PETSC_STATIC_INLINE void WaxpyD(PetscInt dim, PetscScalar a, const PetscScalar *x, const PetscScalar *y, PetscScalar *w) {PetscInt d; for (d = 0; d < dim; ++d) w[d] = a*x[d] + y[d];}
143: PETSC_STATIC_INLINE PetscScalar DotD(PetscInt dim, const PetscScalar *x, const PetscScalar *y) {PetscScalar sum = 0.0; PetscInt d; for (d = 0; d < dim; ++d) sum += x[d]*y[d]; return sum;}
144: PETSC_STATIC_INLINE PetscReal NormD(PetscInt dim, const PetscScalar *x) {return PetscSqrtReal(PetscAbsScalar(DotD(dim,x,x)));}
146: PETSC_STATIC_INLINE void NormalSplit(const PetscReal *n,const PetscScalar *x,PetscScalar *xn,PetscScalar *xt)
147: { /* Split x into normal and tangential components */
148: Scale2(Dot2(x,n)/Dot2(n,n),n,xn);
149: Waxpy2(-1,xn,x,xt);
150: }
152: /******************* Advect ********************/
153: typedef enum {ADVECT_SOL_TILTED,ADVECT_SOL_BUMP} AdvectSolType;
154: static const char *const AdvectSolTypes[] = {"TILTED" ,"BUMP" ,"AdvectSolType" ,"ADVECT_SOL_" ,0};
155: typedef enum {ADVECT_SOL_BUMP_CONE,ADVECT_SOL_BUMP_COS} AdvectSolBumpType;
156: static const char *const AdvectSolBumpTypes[] = {"CONE" ,"COS" ,"AdvectSolBumpType" ,"ADVECT_SOL_BUMP_" ,0};
158: typedef struct {
159: PetscReal wind[DIM];
160: } Physics_Advect_Tilted;
161: typedef struct {
162: PetscReal center[DIM];
163: PetscReal radius;
164: AdvectSolBumpType type;
165: } Physics_Advect_Bump;
167: typedef struct {
168: PetscReal inflowState;
169: AdvectSolType soltype;
170: union {
171: Physics_Advect_Tilted tilted;
172: Physics_Advect_Bump bump;
173: } sol;
174: struct {
175: PetscInt Error;
176: } functional;
177: } Physics_Advect;
179: static const struct FieldDescription PhysicsFields_Advect[] = {{"U" ,1},{NULL,0}};
183: static PetscErrorCode PhysicsBoundary_Advect_Inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
184: {
185: Physics phys = (Physics)ctx;
186: Physics_Advect *advect = (Physics_Advect*)phys->data;
189: xG[0] = advect->inflowState;
190: return (0);
191: }
195: static PetscErrorCode PhysicsBoundary_Advect_Outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
196: {
198: xG[0] = xI[0];
199: return (0);
200: }
204: static void PhysicsRiemann_Advect(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscScalar *flux, Physics phys)
205: {
206: Physics_Advect *advect = (Physics_Advect*)phys->data;
207: PetscReal wind[DIM],wn;
209: switch (advect->soltype) {
210: case ADVECT_SOL_TILTED: {
211: Physics_Advect_Tilted *tilted = &advect->sol.tilted;
212: wind[0] = tilted->wind[0];
213: wind[1] = tilted->wind[1];
214: } break ;
215: case ADVECT_SOL_BUMP:
216: wind[0] = -qp[1];
217: wind[1] = qp[0];
218: break ;
219: /* default: SETERRQ1 (PETSC_COMM_SELF ,PETSC_ERR_SUP,"No support for solution type %s",AdvectSolBumpTypes[advect->soltype]); */
220: }
221: wn = Dot2(wind, n);
222: flux[0] = (wn > 0 ? xL[0] : xR[0]) * wn;
223: }
227: static PetscErrorCode PhysicsSolution_Advect(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx)
228: {
229: Physics phys = (Physics)ctx;
230: Physics_Advect *advect = (Physics_Advect*)phys->data;
233: switch (advect->soltype) {
234: case ADVECT_SOL_TILTED: {
235: PetscReal x0[DIM];
236: Physics_Advect_Tilted *tilted = &advect->sol.tilted;
237: Waxpy2(-time,tilted->wind,x,x0);
238: if (x0[1] > 0) u[0] = 1.*x[0] + 3.*x[1];
239: else u[0] = advect->inflowState;
240: } break ;
241: case ADVECT_SOL_BUMP: {
242: Physics_Advect_Bump *bump = &advect->sol.bump;
243: PetscReal x0[DIM],v[DIM],r,cost,sint;
244: cost = PetscCosReal(time);
245: sint = PetscSinReal(time);
246: x0[0] = cost*x[0] + sint*x[1];
247: x0[1] = -sint*x[0] + cost*x[1];
248: Waxpy2(-1,bump->center,x0,v);
249: r = Norm2(v);
250: switch (bump->type) {
251: case ADVECT_SOL_BUMP_CONE:
252: u[0] = PetscMax(1 - r/bump->radius,0);
253: break ;
254: case ADVECT_SOL_BUMP_COS:
255: u[0] = 0.5 + 0.5*PetscCosReal(PetscMin(r/bump->radius,1)*PETSC_PI);
256: break ;
257: }
258: } break ;
259: default: SETERRQ (PETSC_COMM_SELF ,PETSC_ERR_SUP,"Unknown solution type" );
260: }
261: return (0);
262: }
266: static PetscErrorCode PhysicsFunctional_Advect(Model mod,PetscReal time,const PetscScalar *x,const PetscScalar *y,PetscReal *f,void *ctx)
267: {
268: Physics phys = (Physics)ctx;
269: Physics_Advect *advect = (Physics_Advect*)phys->data;
270: PetscScalar yexact[1];
274: PhysicsSolution_Advect(mod,time,x,yexact,phys);
275: f[advect->functional.Error] = PetscAbsScalar(y[0]-yexact[0]);
276: return (0);
277: }
281: static PetscErrorCode PhysicsCreate_Advect(DM dm, Model mod,Physics phys,PetscOptions *PetscOptionsObject)
282: {
283: Physics_Advect *advect;
287: phys->field_desc = PhysicsFields_Advect;
288: phys->riemann = (PetscRiemannFunc)PhysicsRiemann_Advect;
289: PetscNew (&advect);
290: phys->data = advect;
291: PetscOptionsHead (PetscOptionsObject,"Advect options" );
292: {
293: PetscInt two = 2,dof = 1;
294: advect->soltype = ADVECT_SOL_TILTED;
295: PetscOptionsEnum ("-advect_sol_type" ,"solution type" ,"" ,AdvectSolTypes,(PetscEnum )advect->soltype,(PetscEnum *)&advect->soltype,NULL);
296: switch (advect->soltype) {
297: case ADVECT_SOL_TILTED: {
298: Physics_Advect_Tilted *tilted = &advect->sol.tilted;
299: two = 2;
300: tilted->wind[0] = 0.0;
301: tilted->wind[1] = 1.0;
302: PetscOptionsRealArray ("-advect_tilted_wind" ,"background wind vx,vy" ,"" ,tilted->wind,&two,NULL);
303: advect->inflowState = -2.0;
304: PetscOptionsRealArray ("-advect_tilted_inflow" ,"Inflow state" ,"" ,&advect->inflowState,&dof,NULL);
305: phys->maxspeed = Norm2(tilted->wind);
306: } break ;
307: case ADVECT_SOL_BUMP: {
308: Physics_Advect_Bump *bump = &advect->sol.bump;
309: two = 2;
310: bump->center[0] = 2.;
311: bump->center[1] = 0.;
312: PetscOptionsRealArray ("-advect_bump_center" ,"location of center of bump x,y" ,"" ,bump->center,&two,NULL);
313: bump->radius = 0.9;
314: PetscOptionsReal ("-advect_bump_radius" ,"radius of bump" ,"" ,bump->radius,&bump->radius,NULL);
315: bump->type = ADVECT_SOL_BUMP_CONE;
316: PetscOptionsEnum ("-advect_bump_type" ,"type of bump" ,"" ,AdvectSolBumpTypes,(PetscEnum )bump->type,(PetscEnum *)&bump->type,NULL);
317: phys->maxspeed = 3.; /* radius of mesh, kludge */
318: } break ;
319: }
320: }
321: PetscOptionsTail ();
322: {
323: const PetscInt inflowids[] = {100,200,300},outflowids[] = {101};
324: /* Register "canned" boundary conditions and defaults for where to apply. */
325: DMPlexAddBoundary (dm, PETSC_TRUE , "inflow" , "Face Sets" , 0, 0, NULL, (void (*)()) PhysicsBoundary_Advect_Inflow, ALEN(inflowids), inflowids, phys);
326: DMPlexAddBoundary (dm, PETSC_TRUE , "outflow" , "Face Sets" , 0, 0, NULL, (void (*)()) PhysicsBoundary_Advect_Outflow, ALEN(outflowids), outflowids, phys);
327: /* Initial/transient solution with default boundary conditions */
328: ModelSolutionSetDefault(mod,PhysicsSolution_Advect,phys);
329: /* Register "canned" functionals */
330: ModelFunctionalRegister(mod,"Error" ,&advect->functional.Error,PhysicsFunctional_Advect,phys);
331: }
332: return (0);
333: }
335: /******************* Shallow Water ********************/
336: typedef struct {
337: PetscReal gravity;
338: PetscReal boundaryHeight;
339: struct {
340: PetscInt Height;
341: PetscInt Speed;
342: PetscInt Energy;
343: } functional;
344: } Physics_SW;
345: typedef struct {
346: PetscScalar vals[0];
347: PetscScalar h;
348: PetscScalar uh[DIM];
349: } SWNode;
351: static const struct FieldDescription PhysicsFields_SW[] = {{"Height" ,1},{"Momentum" ,DIM},{NULL,0}};
355: /*
356: * h_t + div(uh) = 0
357: * (uh)_t + div (u\otimes uh + g h^2 / 2 I) = 0
358: *
359: * */
360: static PetscErrorCode SWFlux(Physics phys,const PetscReal *n,const SWNode *x,SWNode *f)
361: {
362: Physics_SW *sw = (Physics_SW*)phys->data;
363: PetscScalar uhn,u[DIM];
364: PetscInt i;
367: Scale2(1./x->h,x->uh,u);
368: uhn = Dot2(x->uh,n);
369: f->h = uhn;
370: for (i=0; i<DIM; i++) f->uh[i] = u[i] * uhn + sw->gravity * PetscSqr(x->h) * n[i];
371: return (0);
372: }
376: static PetscErrorCode PhysicsBoundary_SW_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
377: {
379: xG[0] = xI[0];
380: xG[1] = -xI[1];
381: xG[2] = -xI[2];
382: return (0);
383: }
387: static void PhysicsRiemann_SW(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscScalar *flux, Physics phys)
388: {
389: Physics_SW *sw = (Physics_SW*)phys->data;
390: PetscReal cL,cR,speed,nn[DIM];
391: const SWNode *uL = (const SWNode*)xL,*uR = (const SWNode*)xR;
392: SWNode fL,fR;
393: PetscInt i;
395: if (uL->h < 0 || uR->h < 0) {for (i=0; i<1+dim; i++) flux[i] = NAN; return ;} /* SETERRQ (PETSC_COMM_SELF ,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative"); */
396: nn[0] = n[0];
397: nn[1] = n[1];
398: Normalize2(nn);
399: SWFlux(phys,nn,uL,&fL);
400: SWFlux(phys,nn,uR,&fR);
401: cL = PetscSqrtReal(sw->gravity*PetscRealPart(uL->h));
402: cR = PetscSqrtReal(sw->gravity*PetscRealPart(uR->h)); /* gravity wave speed */
403: speed = PetscMax(PetscAbsScalar(Dot2(uL->uh,nn)/uL->h) + cL,PetscAbsScalar(Dot2(uR->uh,nn)/uR->h) + cR);
404: for (i=0; i<1+dim; i++) flux[i] = (0.5*(fL.vals[i] + fR.vals[i]) + 0.5*speed*(xL[i] - xR[i])) * Norm2(n);
405: }
409: static PetscErrorCode PhysicsSolution_SW(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx)
410: {
411: PetscReal dx[2],r,sigma;
414: if (time != 0.0) SETERRQ1 (mod->comm,PETSC_ERR_SUP,"No solution known for time %g" ,(double)time);
415: dx[0] = x[0] - 1.5;
416: dx[1] = x[1] - 1.0;
417: r = Norm2(dx);
418: sigma = 0.5;
419: u[0] = 1 + 2*PetscExpScalar(-PetscSqr(r)/(2*PetscSqr(sigma)));
420: u[1] = 0.0;
421: u[2] = 0.0;
422: return (0);
423: }
427: static PetscErrorCode PhysicsFunctional_SW(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
428: {
429: Physics phys = (Physics)ctx;
430: Physics_SW *sw = (Physics_SW*)phys->data;
431: const SWNode *x = (const SWNode*)xx;
432: PetscScalar u[2];
433: PetscReal h;
436: h = PetscRealPart(x->h);
437: Scale2(1./x->h,x->uh,u);
438: f[sw->functional.Height] = h;
439: f[sw->functional.Speed] = Norm2(u) + PetscSqrtReal(sw->gravity*h);
440: f[sw->functional.Energy] = 0.5*(Dot2(x->uh,u) + sw->gravity*PetscSqr(h));
441: return (0);
442: }
446: static PetscErrorCode PhysicsCreate_SW(DM dm, Model mod,Physics phys,PetscOptions *PetscOptionsObject)
447: {
448: Physics_SW *sw;
452: phys->field_desc = PhysicsFields_SW;
453: phys->riemann = (PetscRiemannFunc) PhysicsRiemann_SW;
454: PetscNew (&sw);
455: phys->data = sw;
456: PetscOptionsHead (PetscOptionsObject,"SW options" );
457: {
458: sw->gravity = 1.0;
459: PetscOptionsReal ("-sw_gravity" ,"Gravitational constant" ,"" ,sw->gravity,&sw->gravity,NULL);
460: }
461: PetscOptionsTail ();
462: phys->maxspeed = PetscSqrtReal(2.0*sw->gravity); /* Mach 1 for depth of 2 */
464: {
465: const PetscInt wallids[] = {100,101,200,300};
466: DMPlexAddBoundary (dm, PETSC_TRUE , "wall" , "Face Sets" , 0, 0, NULL, (void (*)()) PhysicsBoundary_SW_Wall, ALEN(wallids), wallids, phys);
467: ModelSolutionSetDefault(mod,PhysicsSolution_SW,phys);
468: ModelFunctionalRegister(mod,"Height" ,&sw->functional.Height,PhysicsFunctional_SW,phys);
469: ModelFunctionalRegister(mod,"Speed" ,&sw->functional.Speed,PhysicsFunctional_SW,phys);
470: ModelFunctionalRegister(mod,"Energy" ,&sw->functional.Energy,PhysicsFunctional_SW,phys);
471: }
472: return (0);
473: }
475: /******************* Euler ********************/
476: typedef struct {
477: PetscScalar vals[0];
478: PetscScalar r;
479: PetscScalar ru[DIM];
480: PetscScalar e;
481: } EulerNode;
482: typedef PetscErrorCode (*EquationOfState)(const PetscReal *, const EulerNode*, PetscScalar *) ;
483: typedef struct {
484: PetscInt npars;
485: PetscReal pars[DIM];
486: EquationOfState pressure;
487: EquationOfState sound;
488: struct {
489: PetscInt Density;
490: PetscInt Momentum;
491: PetscInt Energy;
492: PetscInt Pressure;
493: PetscInt Speed;
494: } monitor;
495: } Physics_Euler;
497: static const struct FieldDescription PhysicsFields_Euler[] = {{"Density" ,1},{"Momentum" ,DIM},{"Energy" ,1},{NULL,0}};
501: static PetscErrorCode Pressure_PG(const PetscReal *pars,const EulerNode *x,PetscScalar *p)
502: {
503: PetscScalar ru2;
506: ru2 = DotDIM(x->ru,x->ru);
507: ru2 /= x->r;
508: /* kinematic dof = params[0] */
509: (*p)=2.0*(x->e-0.5*ru2)/pars[0];
510: return (0);
511: }
515: static PetscErrorCode SpeedOfSound_PG(const PetscReal *pars,const EulerNode *x,PetscScalar *c)
516: {
517: PetscScalar p;
520: /* TODO remove direct usage of Pressure_PG */
521: Pressure_PG(pars,x,&p);
522: /* TODO check the sign of p */
523: /* pars[1] = heat capacity ratio */
524: (*c)=PetscSqrtScalar(pars[1]*p/x->r);
525: return (0);
526: }
530: /*
531: * x = (rho,rho*(u_1),...,rho*e)^T
532: * x_t+div(f_1(x))+...+div(f_DIM(x)) = 0
533: *
534: * f_i(x) = u_i*x+(0,0,...,p,...,p*u_i)^T
535: *
536: * */
537: static PetscErrorCode EulerFlux(Physics phys,const PetscReal *n,const EulerNode *x,EulerNode *f)
538: {
539: Physics_Euler *eu = (Physics_Euler*)phys->data;
540: PetscScalar u,nu,p;
541: PetscInt i;
544: u = DotDIM(x->ru,x->ru);
545: u /= (x->r * x->r);
546: nu = DotDIM(x->ru,n);
547: /* TODO check the sign of p */
548: eu->pressure(eu->pars,x,&p);
549: f->r = nu * x->r;
550: for (i=0; i<DIM; i++) f->ru[i] = nu * x->ru[i] + n[i]*p;
551: f->e = nu*(x->e+p);
552: return (0);
553: }
555: /* PetscReal * => EulerNode* conversion */
558: static PetscErrorCode PhysicsBoundary_Euler_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
559: {
560: PetscInt i;
561: PetscScalar xn[DIM],xt[DIM];
564: xG[0] = xI[0];
565: NormalSplitDIM(n,xI+1,xn,xt);
566: for (i=0; i<DIM; i++) xG[i+1] = -xn[i]+xt[i];
567: xG[DIM+1] = xI[DIM+1];
568: return (0);
569: }
571: /* PetscReal * => EulerNode* conversion */
574: static void PhysicsRiemann_Euler_Rusanov(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscScalar *flux, Physics phys)
575: {
576: Physics_Euler *eu = (Physics_Euler*)phys->data;
577: PetscScalar cL,cR,speed;
578: const EulerNode *uL = (const EulerNode*)xL,*uR = (const EulerNode*)xR;
579: EulerNode fL,fR;
580: PetscInt i;
582: if (uL->r < 0 || uR->r < 0) {for (i=0; i<2+dim; i++) flux[i] = NAN; return ;} /* SETERRQ (PETSC_COMM_SELF ,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed density is negative"); */
583: EulerFlux(phys,n,uL,&fL);
584: EulerFlux(phys,n,uR,&fR);
585: eu->sound(eu->pars,uL,&cL);
586: eu->sound(eu->pars,uR,&cR);
587: speed = PetscMax(cL,cR)+PetscMax(PetscAbsScalar(DotDIM(uL->ru,n)/NormDIM(n)),PetscAbsScalar(DotDIM(uR->ru,n)/NormDIM(n)));
588: for (i=0; i<2+dim; i++) flux[i] = 0.5*(fL.vals[i]+fR.vals[i])+0.5*speed*(xL[i]-xR[i]);
589: }
593: static PetscErrorCode PhysicsSolution_Euler(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx)
594: {
595: PetscInt i;
598: if (time != 0.0) SETERRQ1 (mod->comm,PETSC_ERR_SUP,"No solution known for time %g" ,(double)time);
599: u[0] = 1.0;
600: u[DIM+1] = 1.0+PetscAbsReal(x[0]);
601: for (i=1; i<DIM+1; i++) u[i] = 0.0;
602: return (0);
603: }
607: static PetscErrorCode PhysicsFunctional_Euler(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
608: {
609: Physics phys = (Physics)ctx;
610: Physics_Euler *eu = (Physics_Euler*)phys->data;
611: const EulerNode *x = (const EulerNode*)xx;
612: PetscScalar p;
615: f[eu->monitor.Density] = x->r;
616: f[eu->monitor.Momentum] = NormDIM(x->ru);
617: f[eu->monitor.Energy] = x->e;
618: f[eu->monitor.Speed] = NormDIM(x->ru)/x->r;
619: eu->pressure(eu->pars, x, &p);
620: f[eu->monitor.Pressure] = p;
621: return (0);
622: }
626: static PetscErrorCode PhysicsCreate_Euler(DM dm, Model mod,Physics phys,PetscOptions *PetscOptionsObject)
627: {
628: Physics_Euler *eu;
629: PetscErrorCode ierr;
632: phys->field_desc = PhysicsFields_Euler;
633: phys->riemann = (PetscRiemannFunc) PhysicsRiemann_Euler_Rusanov;
634: PetscNew (&eu);
635: phys->data = eu;
636: PetscOptionsHead (PetscOptionsObject,"Euler options" );
637: {
638: eu->pars[0] = 3.0;
639: eu->pars[1] = 1.67;
640: PetscOptionsReal ("-eu_f" ,"Degrees of freedom" ,"" ,eu->pars[0],&eu->pars[0],NULL);
641: PetscOptionsReal ("-eu_gamma" ,"Heat capacity ratio" ,"" ,eu->pars[1],&eu->pars[1],NULL);
642: }
643: PetscOptionsTail ();
644: eu->pressure = Pressure_PG;
645: eu->sound = SpeedOfSound_PG;
646: phys->maxspeed = 1.0;
647: {
648: const PetscInt wallids[] = {100,101,200,300};
649: DMPlexAddBoundary (dm, PETSC_TRUE , "wall" , "Face Sets" , 0, 0, NULL, (void (*)()) PhysicsBoundary_Euler_Wall, ALEN(wallids), wallids, phys);
650: ModelSolutionSetDefault(mod,PhysicsSolution_Euler,phys);
651: ModelFunctionalRegister(mod,"Speed" ,&eu->monitor.Speed,PhysicsFunctional_Euler,phys);
652: ModelFunctionalRegister(mod,"Energy" ,&eu->monitor.Energy,PhysicsFunctional_Euler,phys);
653: ModelFunctionalRegister(mod,"Density" ,&eu->monitor.Density,PhysicsFunctional_Euler,phys);
654: ModelFunctionalRegister(mod,"Momentum" ,&eu->monitor.Momentum,PhysicsFunctional_Euler,phys);
655: ModelFunctionalRegister(mod,"Pressure" ,&eu->monitor.Pressure,PhysicsFunctional_Euler,phys);
656: }
657: return (0);
658: }
662: PetscErrorCode ConstructCellBoundary(DM dm, User user)
663: {
664: const char *name = "Cell Sets" ;
665: const char *bdname = "split faces" ;
666: IS regionIS, innerIS;
667: const PetscInt *regions, *cells;
668: PetscInt numRegions, innerRegion, numCells, c;
669: PetscInt cStart, cEnd, cEndInterior, fStart, fEnd;
670: PetscBool hasLabel;
674: DMPlexGetHeightStratum (dm, 0, &cStart, &cEnd);
675: DMPlexGetHeightStratum (dm, 1, &fStart, &fEnd);
676: DMPlexGetHybridBounds (dm, &cEndInterior, NULL, NULL, NULL);
678: DMPlexHasLabel (dm, name, &hasLabel);
679: if (!hasLabel) return (0);
680: DMPlexGetLabelSize (dm, name, &numRegions);
681: if (numRegions != 2) return (0);
682: /* Get the inner id */
683: DMPlexGetLabelIdIS (dm, name, ®ionIS);
684: ISGetIndices (regionIS, ®ions);
685: innerRegion = regions[0];
686: ISRestoreIndices (regionIS, ®ions);
687: ISDestroy (®ionIS);
688: /* Find the faces between cells in different regions, could call DMPlexCreateNeighborCSR() */
689: DMPlexGetStratumIS (dm, name, innerRegion, &innerIS);
690: ISGetLocalSize (innerIS, &numCells);
691: ISGetIndices (innerIS, &cells);
692: DMPlexCreateLabel (dm, bdname);
693: for (c = 0; c < numCells; ++c) {
694: const PetscInt cell = cells[c];
695: const PetscInt *faces;
696: PetscInt numFaces, f;
698: if ((cell < cStart) || (cell >= cEnd)) SETERRQ1 (PETSC_COMM_SELF , PETSC_ERR_LIB, "Got invalid point %d which is not a cell" , cell);
699: DMPlexGetConeSize (dm, cell, &numFaces);
700: DMPlexGetCone (dm, cell, &faces);
701: for (f = 0; f < numFaces; ++f) {
702: const PetscInt face = faces[f];
703: const PetscInt *neighbors;
704: PetscInt nC, regionA, regionB;
706: if ((face < fStart) || (face >= fEnd)) SETERRQ1 (PETSC_COMM_SELF , PETSC_ERR_LIB, "Got invalid point %d which is not a face" , face);
707: DMPlexGetSupportSize (dm, face, &nC);
708: if (nC != 2) continue ;
709: DMPlexGetSupport (dm, face, &neighbors);
710: if ((neighbors[0] >= cEndInterior) || (neighbors[1] >= cEndInterior)) continue ;
711: 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]);
712: 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]);
713: DMPlexGetLabelValue (dm, name, neighbors[0], ®ionA);
714: DMPlexGetLabelValue (dm, name, neighbors[1], ®ionB);
715: if (regionA < 0) SETERRQ2 (PetscObjectComm ((PetscObject )dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value" , name, neighbors[0]);
716: if (regionB < 0) SETERRQ2 (PetscObjectComm ((PetscObject )dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value" , name, neighbors[1]);
717: if (regionA != regionB) {
718: DMPlexSetLabelValue (dm, bdname, faces[f], 1);
719: }
720: }
721: }
722: ISRestoreIndices (innerIS, &cells);
723: ISDestroy (&innerIS);
724: {
725: DMLabel label;
727: PetscViewerASCIISynchronizedAllow (PETSC_VIEWER_STDOUT_WORLD , PETSC_TRUE );
728: DMPlexGetLabel (dm, bdname, &label);
729: DMLabelView (label, PETSC_VIEWER_STDOUT_WORLD );
730: }
731: return (0);
732: }
736: /* Right now, I have just added duplicate faces, which see both cells. We can
737: - Add duplicate vertices and decouple the face cones
738: - Disconnect faces from cells across the rotation gap
739: */
740: PetscErrorCode SplitFaces(DM *dmSplit, const char labelName[], User user)
741: {
742: DM dm = *dmSplit, sdm;
743: PetscSF sfPoint, gsfPoint;
744: PetscSection coordSection, newCoordSection;
745: Vec coordinates;
746: IS idIS;
747: const PetscInt *ids;
748: PetscInt *newpoints;
749: PetscInt dim, depth, maxConeSize, maxSupportSize, numLabels, numGhostCells;
750: PetscInt numFS, fs, pStart, pEnd, p, cEnd, cEndInterior, vStart, vEnd, v, fStart, fEnd, newf, d, l;
751: PetscBool hasLabel;
755: DMPlexHasLabel (dm, labelName, &hasLabel);
756: if (!hasLabel) return (0);
757: DMCreate (PetscObjectComm ((PetscObject )dm), &sdm);
758: DMSetType (sdm, DMPLEX );
759: DMGetDimension (dm, &dim);
760: DMSetDimension (sdm, dim);
762: DMPlexGetLabelIdIS (dm, labelName, &idIS);
763: ISGetLocalSize (idIS, &numFS);
764: ISGetIndices (idIS, &ids);
766: user->numSplitFaces = 0;
767: for (fs = 0; fs < numFS; ++fs) {
768: PetscInt numBdFaces;
770: DMPlexGetStratumSize (dm, labelName, ids[fs], &numBdFaces);
771: user->numSplitFaces += numBdFaces;
772: }
773: DMPlexGetChart (dm, &pStart, &pEnd);
774: pEnd += user->numSplitFaces;
775: DMPlexSetChart (sdm, pStart, pEnd);
776: DMPlexGetHybridBounds (dm, &cEndInterior, NULL, NULL, NULL);
777: DMPlexSetHybridBounds (sdm, cEndInterior, PETSC_DETERMINE , PETSC_DETERMINE , PETSC_DETERMINE );
778: DMPlexGetHeightStratum (dm, 0, NULL, &cEnd);
779: numGhostCells = cEnd - cEndInterior;
780: /* Set cone and support sizes */
781: DMPlexGetDepth (dm, &depth);
782: for (d = 0; d <= depth; ++d) {
783: DMPlexGetDepthStratum (dm, d, &pStart, &pEnd);
784: for (p = pStart; p < pEnd; ++p) {
785: PetscInt newp = p;
786: PetscInt size;
788: DMPlexGetConeSize (dm, p, &size);
789: DMPlexSetConeSize (sdm, newp, size);
790: DMPlexGetSupportSize (dm, p, &size);
791: DMPlexSetSupportSize (sdm, newp, size);
792: }
793: }
794: DMPlexGetHeightStratum (dm, 1, &fStart, &fEnd);
795: for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
796: IS faceIS;
797: const PetscInt *faces;
798: PetscInt numFaces, f;
800: DMPlexGetStratumIS (dm, labelName, ids[fs], &faceIS);
801: ISGetLocalSize (faceIS, &numFaces);
802: ISGetIndices (faceIS, &faces);
803: for (f = 0; f < numFaces; ++f, ++newf) {
804: PetscInt size;
806: /* Right now I think that both faces should see both cells */
807: DMPlexGetConeSize (dm, faces[f], &size);
808: DMPlexSetConeSize (sdm, newf, size);
809: DMPlexGetSupportSize (dm, faces[f], &size);
810: DMPlexSetSupportSize (sdm, newf, size);
811: }
812: ISRestoreIndices (faceIS, &faces);
813: ISDestroy (&faceIS);
814: }
815: DMSetUp (sdm);
816: /* Set cones and supports */
817: DMPlexGetMaxSizes (dm, &maxConeSize, &maxSupportSize);
818: PetscMalloc (PetscMax(maxConeSize, maxSupportSize) * sizeof (PetscInt ), &newpoints);
819: DMPlexGetChart (dm, &pStart, &pEnd);
820: for (p = pStart; p < pEnd; ++p) {
821: const PetscInt *points, *orientations;
822: PetscInt size, i, newp = p;
824: DMPlexGetConeSize (dm, p, &size);
825: DMPlexGetCone (dm, p, &points);
826: DMPlexGetConeOrientation (dm, p, &orientations);
827: for (i = 0; i < size; ++i) newpoints[i] = points[i];
828: DMPlexSetCone (sdm, newp, newpoints);
829: DMPlexSetConeOrientation (sdm, newp, orientations);
830: DMPlexGetSupportSize (dm, p, &size);
831: DMPlexGetSupport (dm, p, &points);
832: for (i = 0; i < size; ++i) newpoints[i] = points[i];
833: DMPlexSetSupport (sdm, newp, newpoints);
834: }
835: PetscFree (newpoints);
836: for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
837: IS faceIS;
838: const PetscInt *faces;
839: PetscInt numFaces, f;
841: DMPlexGetStratumIS (dm, labelName, ids[fs], &faceIS);
842: ISGetLocalSize (faceIS, &numFaces);
843: ISGetIndices (faceIS, &faces);
844: for (f = 0; f < numFaces; ++f, ++newf) {
845: const PetscInt *points;
847: DMPlexGetCone (dm, faces[f], &points);
848: DMPlexSetCone (sdm, newf, points);
849: DMPlexGetSupport (dm, faces[f], &points);
850: DMPlexSetSupport (sdm, newf, points);
851: }
852: ISRestoreIndices (faceIS, &faces);
853: ISDestroy (&faceIS);
854: }
855: ISRestoreIndices (idIS, &ids);
856: ISDestroy (&idIS);
857: DMPlexStratify (sdm);
858: /* Convert coordinates */
859: DMPlexGetDepthStratum (dm, 0, &vStart, &vEnd);
860: DMGetCoordinateSection (dm, &coordSection);
861: PetscSectionCreate (PetscObjectComm ((PetscObject )dm), &newCoordSection);
862: PetscSectionSetNumFields (newCoordSection, 1);
863: PetscSectionSetFieldComponents (newCoordSection, 0, dim);
864: PetscSectionSetChart (newCoordSection, vStart, vEnd);
865: for (v = vStart; v < vEnd; ++v) {
866: PetscSectionSetDof (newCoordSection, v, dim);
867: PetscSectionSetFieldDof (newCoordSection, v, 0, dim);
868: }
869: PetscSectionSetUp (newCoordSection);
870: DMSetCoordinateSection (sdm, PETSC_DETERMINE , newCoordSection);
871: PetscSectionDestroy (&newCoordSection); /* relinquish our reference */
872: DMGetCoordinatesLocal (dm, &coordinates);
873: DMSetCoordinatesLocal (sdm, coordinates);
874: /* Convert labels */
875: DMPlexGetNumLabels (dm, &numLabels);
876: for (l = 0; l < numLabels; ++l) {
877: const char *lname;
878: PetscBool isDepth;
880: DMPlexGetLabelName (dm, l, &lname);
881: PetscStrcmp (lname, "depth" , &isDepth);
882: if (isDepth) continue ;
883: DMPlexCreateLabel (sdm, lname);
884: DMPlexGetLabelIdIS (dm, lname, &idIS);
885: ISGetLocalSize (idIS, &numFS);
886: ISGetIndices (idIS, &ids);
887: for (fs = 0; fs < numFS; ++fs) {
888: IS pointIS;
889: const PetscInt *points;
890: PetscInt numPoints;
892: DMPlexGetStratumIS (dm, lname, ids[fs], &pointIS);
893: ISGetLocalSize (pointIS, &numPoints);
894: ISGetIndices (pointIS, &points);
895: for (p = 0; p < numPoints; ++p) {
896: PetscInt newpoint = points[p];
898: DMPlexSetLabelValue (sdm, lname, newpoint, ids[fs]);
899: }
900: ISRestoreIndices (pointIS, &points);
901: ISDestroy (&pointIS);
902: }
903: ISRestoreIndices (idIS, &ids);
904: ISDestroy (&idIS);
905: }
906: /* Convert pointSF */
907: const PetscSFNode *remotePoints;
908: PetscSFNode *gremotePoints;
909: const PetscInt *localPoints;
910: PetscInt *glocalPoints,*newLocation,*newRemoteLocation;
911: PetscInt numRoots, numLeaves;
912: PetscMPIInt numProcs;
914: MPI_Comm_size (PetscObjectComm ((PetscObject )dm), &numProcs);
915: DMGetPointSF (dm, &sfPoint);
916: DMGetPointSF (sdm, &gsfPoint);
917: DMPlexGetChart (dm,&pStart,&pEnd);
918: PetscSFGetGraph (sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);
919: if (numRoots >= 0) {
920: PetscMalloc2 (numRoots,&newLocation,pEnd-pStart,&newRemoteLocation);
921: for (l=0; l<numRoots; l++) newLocation[l] = l; /* + (l >= cEnd ? numGhostCells : 0); */
922: PetscSFBcastBegin (sfPoint, MPIU_INT, newLocation, newRemoteLocation);
923: PetscSFBcastEnd (sfPoint, MPIU_INT, newLocation, newRemoteLocation);
924: PetscMalloc1 (numLeaves, &glocalPoints);
925: PetscMalloc1 (numLeaves, &gremotePoints);
926: for (l = 0; l < numLeaves; ++l) {
927: glocalPoints[l] = localPoints[l]; /* localPoints[l] >= cEnd ? localPoints[l] + numGhostCells : localPoints[l]; */
928: gremotePoints[l].rank = remotePoints[l].rank;
929: gremotePoints[l].index = newRemoteLocation[localPoints[l]];
930: }
931: PetscFree2 (newLocation,newRemoteLocation);
932: PetscSFSetGraph (gsfPoint, numRoots+numGhostCells, numLeaves, glocalPoints, PETSC_OWN_POINTER, gremotePoints, PETSC_OWN_POINTER);
933: }
934: DMDestroy (dmSplit);
935: *dmSplit = sdm;
936: return (0);
937: }
941: PetscErrorCode CreatePartitionVec(DM dm, DM *dmCell, Vec *partition)
942: {
943: PetscSF sfPoint;
944: PetscSection coordSection;
945: Vec coordinates;
946: PetscSection sectionCell;
947: PetscScalar *part;
948: PetscInt cStart, cEnd, c;
949: PetscMPIInt rank;
953: DMGetCoordinateSection (dm, &coordSection);
954: DMGetCoordinatesLocal (dm, &coordinates);
955: DMClone (dm, dmCell);
956: DMGetPointSF (dm, &sfPoint);
957: DMSetPointSF (*dmCell, sfPoint);
958: DMSetCoordinateSection (*dmCell, PETSC_DETERMINE , coordSection);
959: DMSetCoordinatesLocal (*dmCell, coordinates);
960: MPI_Comm_rank (PetscObjectComm ((PetscObject )dm), &rank);
961: PetscSectionCreate (PetscObjectComm ((PetscObject )dm), §ionCell);
962: DMPlexGetHeightStratum (*dmCell, 0, &cStart, &cEnd);
963: PetscSectionSetChart (sectionCell, cStart, cEnd);
964: for (c = cStart; c < cEnd; ++c) {
965: PetscSectionSetDof (sectionCell, c, 1);
966: }
967: PetscSectionSetUp (sectionCell);
968: DMSetDefaultSection (*dmCell, sectionCell);
969: PetscSectionDestroy (§ionCell);
970: DMCreateLocalVector (*dmCell, partition);
971: PetscObjectSetName ((PetscObject )*partition, "partition" );
972: VecGetArray (*partition, &part);
973: for (c = cStart; c < cEnd; ++c) {
974: PetscScalar *p;
976: DMPlexPointLocalRef (*dmCell, c, part, &p);
977: p[0] = rank;
978: }
979: VecRestoreArray (*partition, &part);
980: return (0);
981: }
985: PetscErrorCode CreateMassMatrix(DM dm, Vec *massMatrix, User user)
986: {
987: DM dmMass, dmFace, dmCell, dmCoord;
988: PetscSection coordSection;
989: Vec coordinates, facegeom, cellgeom;
990: PetscSection sectionMass;
991: PetscScalar *m;
992: const PetscScalar *fgeom, *cgeom, *coords;
993: PetscInt vStart, vEnd, v;
994: PetscErrorCode ierr;
997: DMGetCoordinateSection (dm, &coordSection);
998: DMGetCoordinatesLocal (dm, &coordinates);
999: DMClone (dm, &dmMass);
1000: DMSetCoordinateSection (dmMass, PETSC_DETERMINE , coordSection);
1001: DMSetCoordinatesLocal (dmMass, coordinates);
1002: PetscSectionCreate (PetscObjectComm ((PetscObject )dm), §ionMass);
1003: DMPlexGetDepthStratum (dm, 0, &vStart, &vEnd);
1004: PetscSectionSetChart (sectionMass, vStart, vEnd);
1005: for (v = vStart; v < vEnd; ++v) {
1006: PetscInt numFaces;
1008: DMPlexGetSupportSize (dmMass, v, &numFaces);
1009: PetscSectionSetDof (sectionMass, v, numFaces*numFaces);
1010: }
1011: PetscSectionSetUp (sectionMass);
1012: DMSetDefaultSection (dmMass, sectionMass);
1013: PetscSectionDestroy (§ionMass);
1014: DMGetLocalVector (dmMass, massMatrix);
1015: VecGetArray (*massMatrix, &m);
1016: DMPlexTSGetGeometryFVM (dm, &facegeom, &cellgeom, NULL);
1017: VecGetDM (facegeom, &dmFace);
1018: VecGetArrayRead (facegeom, &fgeom);
1019: VecGetDM (cellgeom, &dmCell);
1020: VecGetArrayRead (cellgeom, &cgeom);
1021: DMGetCoordinateDM (dm, &dmCoord);
1022: VecGetArrayRead (coordinates, &coords);
1023: for (v = vStart; v < vEnd; ++v) {
1024: const PetscInt *faces;
1025: const PetscFVFaceGeom *fgA, *fgB, *cg;
1026: const PetscScalar *vertex;
1027: PetscInt numFaces, sides[2], f, g;
1029: DMPlexPointLocalRead (dmCoord, v, coords, &vertex);
1030: DMPlexGetSupportSize (dmMass, v, &numFaces);
1031: DMPlexGetSupport (dmMass, v, &faces);
1032: for (f = 0; f < numFaces; ++f) {
1033: sides[0] = faces[f];
1034: DMPlexPointLocalRead (dmFace, faces[f], fgeom, &fgA);
1035: for (g = 0; g < numFaces; ++g) {
1036: const PetscInt *cells = NULL;;
1037: PetscReal area = 0.0;
1038: PetscInt numCells;
1040: sides[1] = faces[g];
1041: DMPlexPointLocalRead (dmFace, faces[g], fgeom, &fgB);
1042: DMPlexGetJoin (dmMass, 2, sides, &numCells, &cells);
1043: if (numCells != 1) SETERRQ (PETSC_COMM_SELF , PETSC_ERR_LIB, "Invalid join for faces" );
1044: DMPlexPointLocalRead (dmCell, cells[0], cgeom, &cg);
1045: area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgA->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgA->centroid[0] - cg->centroid[0]));
1046: area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgB->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgB->centroid[0] - cg->centroid[0]));
1047: m[f*numFaces+g] = Dot2(fgA->normal, fgB->normal)*area*0.5;
1048: DMPlexRestoreJoin (dmMass, 2, sides, &numCells, &cells);
1049: }
1050: }
1051: }
1052: VecRestoreArrayRead (facegeom, &fgeom);
1053: VecRestoreArrayRead (cellgeom, &cgeom);
1054: VecRestoreArrayRead (coordinates, &coords);
1055: VecRestoreArray (*massMatrix, &m);
1056: DMDestroy (&dmMass);
1057: return (0);
1058: }
1062: /* Behavior will be different for multi-physics or when using non-default boundary conditions */
1063: static PetscErrorCode ModelSolutionSetDefault(Model mod,SolutionFunction func,void *ctx)
1064: {
1066: mod->solution = func;
1067: mod->solutionctx = ctx;
1068: return (0);
1069: }
1073: static PetscErrorCode ModelFunctionalRegister(Model mod,const char *name,PetscInt *offset,FunctionalFunction func,void *ctx)
1074: {
1076: FunctionalLink link,*ptr;
1077: PetscInt lastoffset = -1;
1080: for (ptr=&mod->functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset;
1081: PetscNew (&link);
1082: PetscStrallocpy (name,&link->name);
1083: link->offset = lastoffset + 1;
1084: link->func = func;
1085: link->ctx = ctx;
1086: link->next = NULL;
1087: *ptr = link;
1088: *offset = link->offset;
1089: return (0);
1090: }
1094: static PetscErrorCode ModelFunctionalSetFromOptions(Model mod,PetscOptions *PetscOptionsObject)
1095: {
1097: PetscInt i,j;
1098: FunctionalLink link;
1099: char *names[256];
1102: mod->numMonitored = ALEN(names);
1103: PetscOptionsStringArray ("-monitor" ,"list of functionals to monitor" ,"" ,names,&mod->numMonitored,NULL);
1104: /* Create list of functionals that will be computed somehow */
1105: PetscMalloc1 (mod->numMonitored,&mod->functionalMonitored);
1106: /* Create index of calls that we will have to make to compute these functionals (over-allocation in general). */
1107: PetscMalloc1 (mod->numMonitored,&mod->functionalCall);
1108: mod->numCall = 0;
1109: for (i=0; i<mod->numMonitored; i++) {
1110: for (link=mod->functionalRegistry; link; link=link->next) {
1111: PetscBool match;
1112: PetscStrcasecmp (names[i],link->name,&match);
1113: if (match) break ;
1114: }
1115: if (!link) SETERRQ1 (mod->comm,PETSC_ERR_USER,"No known functional '%s'" ,names[i]);
1116: mod->functionalMonitored[i] = link;
1117: for (j=0; j<i; j++) {
1118: if (mod->functionalCall[j]->func == link->func && mod->functionalCall[j]->ctx == link->ctx) goto next_name;
1119: }
1120: mod->functionalCall[mod->numCall++] = link; /* Just points to the first link using the result. There may be more results. */
1121: next_name:
1122: PetscFree (names[i]);
1123: }
1125: /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */
1126: mod->maxComputed = -1;
1127: for (link=mod->functionalRegistry; link; link=link->next) {
1128: for (i=0; i<mod->numCall; i++) {
1129: FunctionalLink call = mod->functionalCall[i];
1130: if (link->func == call->func && link->ctx == call->ctx) {
1131: mod->maxComputed = PetscMax(mod->maxComputed,link->offset);
1132: }
1133: }
1134: }
1135: return (0);
1136: }
1140: static PetscErrorCode FunctionalLinkDestroy(FunctionalLink *link)
1141: {
1143: FunctionalLink l,next;
1146: if (!link) return (0);
1147: l = *link;
1148: *link = NULL;
1149: for (; l; l=next) {
1150: next = l->next;
1151: PetscFree (l->name);
1152: PetscFree (l);
1153: }
1154: return (0);
1155: }
1159: PetscErrorCode SetInitialCondition(DM dm, Vec X, User user)
1160: {
1161: DM dmCell;
1162: Model mod = user->model;
1163: Vec cellgeom;
1164: const PetscScalar *cgeom;
1165: PetscScalar *x;
1166: PetscInt cStart, cEnd, cEndInterior, c;
1167: PetscErrorCode ierr;
1170: DMPlexGetHybridBounds (dm, &cEndInterior, NULL, NULL, NULL);
1171: DMPlexTSGetGeometryFVM (dm, NULL, &cellgeom, NULL);
1172: VecGetDM (cellgeom, &dmCell);
1173: DMPlexGetHeightStratum (dm, 0, &cStart, &cEnd);
1174: VecGetArrayRead (cellgeom, &cgeom);
1175: VecGetArray (X, &x);
1176: for (c = cStart; c < cEndInterior; ++c) {
1177: const PetscFVCellGeom *cg;
1178: PetscScalar *xc;
1180: DMPlexPointLocalRead (dmCell,c,cgeom,&cg);
1181: DMPlexPointGlobalRef (dm,c,x,&xc);
1182: if (xc) {(*mod->solution)(mod,0.0,cg->centroid,xc,mod->solutionctx);}
1183: }
1184: VecRestoreArrayRead (cellgeom, &cgeom);
1185: VecRestoreArray (X, &x);
1186: return (0);
1187: }
1191: static PetscErrorCode OutputVTK(DM dm, const char *filename, PetscViewer *viewer)
1192: {
1196: PetscViewerCreate (PetscObjectComm ((PetscObject )dm), viewer);
1197: PetscViewerSetType (*viewer, PETSCVIEWERVTK);
1198: PetscViewerFileSetName (*viewer, filename);
1199: return (0);
1200: }
1204: static PetscErrorCode MonitorVTK(TS ts,PetscInt stepnum,PetscReal time,Vec X,void *ctx)
1205: {
1206: User user = (User)ctx;
1207: DM dm;
1208: Vec cellgeom;
1209: PetscViewer viewer;
1210: char filename[PETSC_MAX_PATH_LEN],*ftable = NULL;
1211: PetscReal xnorm;
1212: PetscInt cEndInterior;
1216: PetscObjectSetName ((PetscObject ) X, "solution" );
1217: VecGetDM (X,&dm);
1218: DMPlexTSGetGeometryFVM (dm, NULL, &cellgeom, NULL);
1219: DMPlexGetHybridBounds (dm, &cEndInterior, NULL, NULL, NULL);
1220: VecNorm (X,NORM_INFINITY ,&xnorm);
1221: if (stepnum >= 0) { /* No summary for final time */
1222: Model mod = user->model;
1223: PetscInt c,cStart,cEnd,fcount,i;
1224: size_t ftableused,ftablealloc;
1225: const PetscScalar *cgeom,*x;
1226: DM dmCell;
1227: PetscReal *fmin,*fmax,*fintegral,*ftmp;
1228: fcount = mod->maxComputed+1;
1229: PetscMalloc4 (fcount,&fmin,fcount,&fmax,fcount,&fintegral,fcount,&ftmp);
1230: for (i=0; i<fcount; i++) {
1231: fmin[i] = PETSC_MAX_REAL;
1232: fmax[i] = PETSC_MIN_REAL;
1233: fintegral[i] = 0;
1234: }
1235: DMPlexGetHeightStratum (dm,0,&cStart,&cEnd);
1236: VecGetDM (cellgeom,&dmCell);
1237: VecGetArrayRead (cellgeom,&cgeom);
1238: VecGetArrayRead (X,&x);
1239: for (c = cStart; c < cEndInterior; ++c) {
1240: const PetscFVCellGeom *cg;
1241: const PetscScalar *cx;
1242: DMPlexPointLocalRead (dmCell,c,cgeom,&cg);
1243: DMPlexPointGlobalRead (dm,c,x,&cx);
1244: if (!cx) continue ; /* not a global cell */
1245: for (i=0; i<mod->numCall; i++) {
1246: FunctionalLink flink = mod->functionalCall[i];
1247: (*flink->func)(mod,time,cg->centroid,cx,ftmp,flink->ctx);
1248: }
1249: for (i=0; i<fcount; i++) {
1250: fmin[i] = PetscMin(fmin[i],ftmp[i]);
1251: fmax[i] = PetscMax(fmax[i],ftmp[i]);
1252: fintegral[i] += cg->volume * ftmp[i];
1253: }
1254: }
1255: VecRestoreArrayRead (cellgeom,&cgeom);
1256: VecRestoreArrayRead (X,&x);
1257: MPI_Allreduce (MPI_IN_PLACE,fmin,fcount,MPIU_REAL,MPIU_MIN,PetscObjectComm ((PetscObject )ts));
1258: MPI_Allreduce (MPI_IN_PLACE,fmax,fcount,MPIU_REAL,MPIU_MAX,PetscObjectComm ((PetscObject )ts));
1259: MPI_Allreduce (MPI_IN_PLACE,fintegral,fcount,MPIU_REAL,MPIU_SUM,PetscObjectComm ((PetscObject )ts));
1261: ftablealloc = fcount * 100;
1262: ftableused = 0;
1263: PetscMalloc1 (ftablealloc,&ftable);
1264: for (i=0; i<mod->numMonitored; i++) {
1265: size_t countused;
1266: char buffer[256],*p;
1267: FunctionalLink flink = mod->functionalMonitored[i];
1268: PetscInt id = flink->offset;
1269: if (i % 3) {
1270: PetscMemcpy (buffer," " ,2);
1271: p = buffer + 2;
1272: } else if (i) {
1273: char newline[] = "\n" ;
1274: PetscMemcpy (buffer,newline,sizeof newline-1);
1275: p = buffer + sizeof newline - 1;
1276: } else {
1277: p = buffer;
1278: }
1279: 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]);
1280: countused += p - buffer;
1281: if (countused > ftablealloc-ftableused-1) { /* reallocate */
1282: char *ftablenew;
1283: ftablealloc = 2*ftablealloc + countused;
1284: PetscMalloc (ftablealloc,&ftablenew);
1285: PetscMemcpy (ftablenew,ftable,ftableused);
1286: PetscFree (ftable);
1287: ftable = ftablenew;
1288: }
1289: PetscMemcpy (ftable+ftableused,buffer,countused);
1290: ftableused += countused;
1291: ftable[ftableused] = 0;
1292: }
1293: PetscFree4 (fmin,fmax,fintegral,ftmp);
1295: PetscPrintf (PetscObjectComm ((PetscObject )ts),"% 3D time %8.4g |x| %8.4g %s\n" ,stepnum,(double)time,(double)xnorm,ftable ? ftable : "" );
1296: PetscFree (ftable);
1297: }
1298: if (user->vtkInterval < 1) return (0);
1299: if ((stepnum == -1) ^ (stepnum % user->vtkInterval == 0)) {
1300: if (stepnum == -1) { /* Final time is not multiple of normal time interval, write it anyway */
1301: TSGetTimeStepNumber (ts,&stepnum);
1302: }
1303: PetscSNPrintf (filename,sizeof filename,"ex11-%03D.vtu" ,stepnum);
1304: OutputVTK(dm,filename,&viewer);
1305: VecView (X,viewer);
1306: PetscViewerDestroy (&viewer);
1307: }
1308: return (0);
1309: }
1313: int main(int argc, char **argv)
1314: {
1315: MPI_Comm comm;
1316: PetscDS prob;
1317: PetscFV fvm;
1318: User user;
1319: Model mod;
1320: Physics phys;
1321: DM dm;
1322: PetscReal ftime, cfl, dt, minRadius;
1323: PetscInt dim, nsteps;
1324: TS ts;
1325: TSConvergedReason reason;
1326: Vec X;
1327: PetscViewer viewer;
1328: PetscBool vtkCellGeom, splitFaces;
1329: PetscInt overlap;
1330: char filename[PETSC_MAX_PATH_LEN] = "sevenside.exo" ;
1331: PetscErrorCode ierr;
1333: PetscInitialize (&argc, &argv, (char*) 0, help);
1334: comm = PETSC_COMM_WORLD ;
1336: PetscNew (&user);
1337: PetscNew (&user->model);
1338: PetscNew (&user->model->physics);
1339: mod = user->model;
1340: phys = mod->physics;
1341: mod->comm = comm;
1343: /* Register physical models to be available on the command line */
1344: PetscFunctionListAdd (&PhysicsList,"advect" ,PhysicsCreate_Advect);
1345: PetscFunctionListAdd (&PhysicsList,"sw" ,PhysicsCreate_SW);
1346: PetscFunctionListAdd (&PhysicsList,"euler" ,PhysicsCreate_Euler);
1349: PetscOptionsBegin (comm,NULL,"Unstructured Finite Volume Mesh Options" ,"" );
1350: {
1351: cfl = 0.9 * 4; /* default SSPRKS2 with s=5 stages is stable for CFL number s-1 */
1352: PetscOptionsReal ("-ufv_cfl" ,"CFL number per step" ,"" ,cfl,&cfl,NULL);
1353: PetscOptionsString ("-f" ,"Exodus.II filename to read" ,"" ,filename,filename,sizeof (filename),NULL);
1354: splitFaces = PETSC_FALSE ;
1355: PetscOptionsBool ("-ufv_split_faces" ,"Split faces between cell sets" ,"" ,splitFaces,&splitFaces,NULL);
1356: overlap = 1;
1357: PetscOptionsInt ("-ufv_mesh_overlap" ,"Number of cells to overlap partitions" ,"" ,overlap,&overlap,NULL);
1358: user->vtkInterval = 1;
1359: PetscOptionsInt ("-ufv_vtk_interval" ,"VTK output interval (0 to disable)" ,"" ,user->vtkInterval,&user->vtkInterval,NULL);
1360: vtkCellGeom = PETSC_FALSE ;
1361: PetscOptionsBool ("-ufv_vtk_cellgeom" ,"Write cell geometry (for debugging)" ,"" ,vtkCellGeom,&vtkCellGeom,NULL);
1362: }
1363: PetscOptionsEnd ();
1364: DMPlexCreateFromFile (comm, filename, PETSC_TRUE , &dm);
1365: DMViewFromOptions(dm, NULL, "-dm_view" );
1366: DMGetDimension (dm, &dim);
1368: PetscOptionsBegin (comm,NULL,"Unstructured Finite Volume Physics Options" ,"" );
1369: {
1370: PetscErrorCode (*physcreate)(DM ,Model,Physics,PetscOptions*);
1371: char physname[256] = "advect" ;
1373: DMPlexCreateLabel (dm, "Face Sets" );
1374: PetscOptionsFList ("-physics" ,"Physics module to solve" ,"" ,PhysicsList,physname,physname,sizeof physname,NULL);
1375: PetscFunctionListFind (PhysicsList,physname,&physcreate);
1376: PetscMemzero (phys,sizeof (struct _n_Physics ));
1377: (*physcreate)(dm,mod,phys,PetscOptionsObject);
1378: mod->maxspeed = phys->maxspeed;
1379: /* Count number of fields and dofs */
1380: for (phys->nfields=0,phys->dof=0; phys->field_desc[phys->nfields].name; phys->nfields++) phys->dof += phys->field_desc[phys->nfields].dof;
1382: if (mod->maxspeed <= 0) SETERRQ1 (comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set maxspeed" ,physname);
1383: if (phys->dof <= 0) SETERRQ1 (comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set dof" ,physname);
1384: ModelFunctionalSetFromOptions(mod,PetscOptionsObject);
1385: }
1386: PetscOptionsEnd ();
1387: {
1388: DM dmDist;
1390: DMPlexSetAdjacencyUseCone (dm, PETSC_TRUE );
1391: DMPlexSetAdjacencyUseClosure (dm, PETSC_FALSE );
1392: DMPlexDistribute (dm, overlap, NULL, &dmDist);
1393: if (dmDist) {
1394: DMDestroy (&dm);
1395: dm = dmDist;
1396: }
1397: }
1398: DMSetFromOptions (dm);
1399: {
1400: DM gdm;
1402: DMPlexConstructGhostCells (dm, NULL, NULL, &gdm);
1403: DMDestroy (&dm);
1404: dm = gdm;
1405: DMViewFromOptions(dm, NULL, "-dm_view" );
1406: }
1407: if (splitFaces) {ConstructCellBoundary(dm, user);}
1408: SplitFaces(&dm, "split faces" , user);
1410: PetscFVCreate (comm, &fvm);
1411: PetscFVSetFromOptions (fvm);
1412: PetscFVSetNumComponents (fvm, phys->dof);
1413: PetscFVSetSpatialDimension (fvm, dim);
1414: DMGetDS (dm, &prob);
1415: /* FV is now structured with one field having all physics as components */
1416: PetscDSAddDiscretization (prob, (PetscObject ) fvm);
1417: PetscDSSetRiemannSolver (prob, 0, user->model->physics->riemann);
1418: PetscDSSetContext(prob, 0, user->model->physics);
1420: TSCreate (comm, &ts);
1421: TSSetType (ts, TSSSP );
1422: TSSetDM (ts, dm);
1423: TSMonitorSet (ts,MonitorVTK,user,NULL);
1424: DMTSSetRHSFunctionLocal (dm, DMPlexTSComputeRHSFunctionFVM , user);
1426: DMCreateGlobalVector (dm, &X);
1427: PetscObjectSetName ((PetscObject ) X, "solution" );
1428: SetInitialCondition(dm, X, user);
1429: if (vtkCellGeom) {
1430: DM dmCell;
1431: Vec cellgeom, partition;
1433: DMPlexTSGetGeometryFVM (dm, NULL, &cellgeom, NULL);
1434: OutputVTK(dm, "ex11-cellgeom.vtk" , &viewer);
1435: VecView (cellgeom, viewer);
1436: PetscViewerDestroy (&viewer);
1437: CreatePartitionVec(dm, &dmCell, &partition);
1438: OutputVTK(dmCell, "ex11-partition.vtk" , &viewer);
1439: VecView (partition, viewer);
1440: PetscViewerDestroy (&viewer);
1441: VecDestroy (&partition);
1442: DMDestroy (&dmCell);
1443: }
1445: DMPlexTSGetGeometryFVM (dm, NULL, NULL, &minRadius);
1446: TSSetDuration (ts,1000,2.0);
1447: dt = cfl * minRadius / user->model->maxspeed;
1448: TSSetInitialTimeStep (ts,0.0,dt);
1449: TSSetFromOptions (ts);
1450: TSSolve (ts,X);
1451: TSGetSolveTime (ts,&ftime);
1452: TSGetTimeStepNumber (ts,&nsteps);
1453: TSGetConvergedReason (ts,&reason);
1454: PetscPrintf (PETSC_COMM_WORLD ,"%s at time %g after %D steps\n" ,TSConvergedReasons[reason],(double)ftime,nsteps);
1455: TSDestroy (&ts);
1457: PetscFunctionListDestroy (&PhysicsList);
1458: FunctionalLinkDestroy(&user->model->functionalRegistry);
1459: PetscFree (user->model->functionalMonitored);
1460: PetscFree (user->model->functionalCall);
1461: PetscFree (user->model->physics->data);
1462: PetscFree (user->model->physics);
1463: PetscFree (user->model);
1464: PetscFree (user);
1465: VecDestroy (&X);
1466: PetscFVDestroy (&fvm);
1467: DMDestroy (&dm);
1468: PetscFinalize ();
1469: return (0);
1470: }