Actual source code: ex11.c
petsc-3.6.0 2015-06-09
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;
57: typedef void (*RiemannFunction)(const PetscReal *,const PetscReal *,const PetscScalar *,const PetscScalar *,PetscScalar *,void*) ;
58: typedef PetscErrorCode (*SolutionFunction)(Model,PetscReal ,const PetscReal *,PetscScalar *,void*) ;
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 *) ;
65: struct FieldDescription {
66: const char *name;
67: PetscInt dof;
68: };
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: };
79: struct _n_Physics {
80: RiemannFunction 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: };
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: void *solutionctx;
99: PetscReal maxspeed; /* estimate of global maximum speed (for CFL calculation) */
100: };
102: struct _n_User {
103: PetscInt numSplitFaces;
104: PetscInt vtkInterval; /* For monitor */
105: Model model;
106: };
108: PETSC_STATIC_INLINE PetscScalar DotDIM(const PetscScalar *x,const PetscScalar *y)
109: {
110: PetscInt i;
111: PetscScalar prod=0.0;
113: for (i=0; i<DIM; i++) prod += x[i]*y[i];
114: return prod;
115: }
116: PETSC_STATIC_INLINE PetscReal NormDIM(const PetscScalar *x) { return PetscSqrtReal(PetscAbsScalar(DotDIM(x,x))); }
117: PETSC_STATIC_INLINE void axDIM(const PetscScalar a,PetscScalar *x)
118: {
119: PetscInt i;
120: for (i=0; i<DIM; i++) x[i] *= a;
121: }
122: PETSC_STATIC_INLINE void waxDIM(const PetscScalar a,const PetscScalar *x, PetscScalar *w)
123: {
124: PetscInt i;
125: for (i=0; i<DIM; i++) w[i] = x[i]*a;
126: }
127: PETSC_STATIC_INLINE void NormalSplitDIM(const PetscReal *n,const PetscScalar *x,PetscScalar *xn,PetscScalar *xt)
128: { /* Split x into normal and tangential components */
129: PetscInt i;
130: PetscScalar c;
131: c = DotDIM(x,n)/DotDIM(n,n);
132: for (i=0; i<DIM; i++) {
133: xn[i] = c*n[i];
134: xt[i] = x[i]-xn[i];
135: }
136: }
138: PETSC_STATIC_INLINE PetscScalar Dot2(const PetscScalar *x,const PetscScalar *y) { return x[0]*y[0] + x[1]*y[1];}
139: PETSC_STATIC_INLINE PetscReal Norm2(const PetscScalar *x) { return PetscSqrtReal(PetscAbsScalar(Dot2(x,x)));}
140: PETSC_STATIC_INLINE void Normalize2(PetscScalar *x) { PetscReal a = 1./Norm2(x); x[0] *= a; x[1] *= a; }
141: 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]; }
142: PETSC_STATIC_INLINE void Scale2(PetscScalar a,const PetscScalar *x,PetscScalar *y) { y[0] = a*x[0]; y[1] = a*x[1]; }
144: 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];}
145: 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;}
146: PETSC_STATIC_INLINE PetscReal NormD(PetscInt dim, const PetscScalar *x) {return PetscSqrtReal(PetscAbsScalar(DotD(dim,x,x)));}
148: PETSC_STATIC_INLINE void NormalSplit(const PetscReal *n,const PetscScalar *x,PetscScalar *xn,PetscScalar *xt)
149: { /* Split x into normal and tangential components */
150: Scale2(Dot2(x,n)/Dot2(n,n),n,xn);
151: Waxpy2(-1,xn,x,xt);
152: }
154: /******************* Advect ********************/
155: typedef enum {ADVECT_SOL_TILTED,ADVECT_SOL_BUMP} AdvectSolType;
156: static const char *const AdvectSolTypes[] = {"TILTED" ,"BUMP" ,"AdvectSolType" ,"ADVECT_SOL_" ,0};
157: typedef enum {ADVECT_SOL_BUMP_CONE,ADVECT_SOL_BUMP_COS} AdvectSolBumpType;
158: static const char *const AdvectSolBumpTypes[] = {"CONE" ,"COS" ,"AdvectSolBumpType" ,"ADVECT_SOL_BUMP_" ,0};
160: typedef struct {
161: PetscReal wind[DIM];
162: } Physics_Advect_Tilted;
163: typedef struct {
164: PetscReal center[DIM];
165: PetscReal radius;
166: AdvectSolBumpType type;
167: } Physics_Advect_Bump;
169: typedef struct {
170: PetscReal inflowState;
171: AdvectSolType soltype;
172: union {
173: Physics_Advect_Tilted tilted;
174: Physics_Advect_Bump bump;
175: } sol;
176: struct {
177: PetscInt Error;
178: } functional;
179: } Physics_Advect;
181: static const struct FieldDescription PhysicsFields_Advect[] = {{"U" ,1},{NULL,0}};
185: static PetscErrorCode PhysicsBoundary_Advect_Inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
186: {
187: Physics phys = (Physics)ctx;
188: Physics_Advect *advect = (Physics_Advect*)phys->data;
191: xG[0] = advect->inflowState;
192: return (0);
193: }
197: static PetscErrorCode PhysicsBoundary_Advect_Outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
198: {
200: xG[0] = xI[0];
201: return (0);
202: }
206: static void PhysicsRiemann_Advect(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscScalar *flux, Physics phys)
207: {
208: Physics_Advect *advect = (Physics_Advect*)phys->data;
209: PetscReal wind[DIM],wn;
211: switch (advect->soltype) {
212: case ADVECT_SOL_TILTED: {
213: Physics_Advect_Tilted *tilted = &advect->sol.tilted;
214: wind[0] = tilted->wind[0];
215: wind[1] = tilted->wind[1];
216: } break ;
217: case ADVECT_SOL_BUMP:
218: wind[0] = -qp[1];
219: wind[1] = qp[0];
220: break ;
221: /* default: SETERRQ1 (PETSC_COMM_SELF ,PETSC_ERR_SUP,"No support for solution type %s",AdvectSolBumpTypes[advect->soltype]); */
222: }
223: wn = Dot2(wind, n);
224: flux[0] = (wn > 0 ? xL[0] : xR[0]) * wn;
225: }
229: static PetscErrorCode PhysicsSolution_Advect(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx)
230: {
231: Physics phys = (Physics)ctx;
232: Physics_Advect *advect = (Physics_Advect*)phys->data;
235: switch (advect->soltype) {
236: case ADVECT_SOL_TILTED: {
237: PetscReal x0[DIM];
238: Physics_Advect_Tilted *tilted = &advect->sol.tilted;
239: Waxpy2(-time,tilted->wind,x,x0);
240: if (x0[1] > 0) u[0] = 1.*x[0] + 3.*x[1];
241: else u[0] = advect->inflowState;
242: } break ;
243: case ADVECT_SOL_BUMP: {
244: Physics_Advect_Bump *bump = &advect->sol.bump;
245: PetscReal x0[DIM],v[DIM],r,cost,sint;
246: cost = PetscCosReal(time);
247: sint = PetscSinReal(time);
248: x0[0] = cost*x[0] + sint*x[1];
249: x0[1] = -sint*x[0] + cost*x[1];
250: Waxpy2(-1,bump->center,x0,v);
251: r = Norm2(v);
252: switch (bump->type) {
253: case ADVECT_SOL_BUMP_CONE:
254: u[0] = PetscMax(1 - r/bump->radius,0);
255: break ;
256: case ADVECT_SOL_BUMP_COS:
257: u[0] = 0.5 + 0.5*PetscCosReal(PetscMin(r/bump->radius,1)*PETSC_PI);
258: break ;
259: }
260: } break ;
261: default: SETERRQ (PETSC_COMM_SELF ,PETSC_ERR_SUP,"Unknown solution type" );
262: }
263: return (0);
264: }
268: static PetscErrorCode PhysicsFunctional_Advect(Model mod,PetscReal time,const PetscScalar *x,const PetscScalar *y,PetscReal *f,void *ctx)
269: {
270: Physics phys = (Physics)ctx;
271: Physics_Advect *advect = (Physics_Advect*)phys->data;
272: PetscScalar yexact[1];
276: PhysicsSolution_Advect(mod,time,x,yexact,phys);
277: f[advect->functional.Error] = PetscAbsScalar(y[0]-yexact[0]);
278: return (0);
279: }
283: static PetscErrorCode PhysicsCreate_Advect(DM dm, Model mod,Physics phys,PetscOptions *PetscOptionsObject)
284: {
285: Physics_Advect *advect;
289: phys->field_desc = PhysicsFields_Advect;
290: phys->riemann = PhysicsRiemann_Advect;
291: PetscNew (&advect);
292: phys->data = advect;
293: PetscOptionsHead (PetscOptionsObject,"Advect options" );
294: {
295: PetscInt two = 2,dof = 1;
296: advect->soltype = ADVECT_SOL_TILTED;
297: PetscOptionsEnum ("-advect_sol_type" ,"solution type" ,"" ,AdvectSolTypes,(PetscEnum )advect->soltype,(PetscEnum *)&advect->soltype,NULL);
298: switch (advect->soltype) {
299: case ADVECT_SOL_TILTED: {
300: Physics_Advect_Tilted *tilted = &advect->sol.tilted;
301: two = 2;
302: tilted->wind[0] = 0.0;
303: tilted->wind[1] = 1.0;
304: PetscOptionsRealArray ("-advect_tilted_wind" ,"background wind vx,vy" ,"" ,tilted->wind,&two,NULL);
305: advect->inflowState = -2.0;
306: PetscOptionsRealArray ("-advect_tilted_inflow" ,"Inflow state" ,"" ,&advect->inflowState,&dof,NULL);
307: phys->maxspeed = Norm2(tilted->wind);
308: } break ;
309: case ADVECT_SOL_BUMP: {
310: Physics_Advect_Bump *bump = &advect->sol.bump;
311: two = 2;
312: bump->center[0] = 2.;
313: bump->center[1] = 0.;
314: PetscOptionsRealArray ("-advect_bump_center" ,"location of center of bump x,y" ,"" ,bump->center,&two,NULL);
315: bump->radius = 0.9;
316: PetscOptionsReal ("-advect_bump_radius" ,"radius of bump" ,"" ,bump->radius,&bump->radius,NULL);
317: bump->type = ADVECT_SOL_BUMP_CONE;
318: PetscOptionsEnum ("-advect_bump_type" ,"type of bump" ,"" ,AdvectSolBumpTypes,(PetscEnum )bump->type,(PetscEnum *)&bump->type,NULL);
319: phys->maxspeed = 3.; /* radius of mesh, kludge */
320: } break ;
321: }
322: }
323: PetscOptionsTail ();
324: {
325: const PetscInt inflowids[] = {100,200,300},outflowids[] = {101};
326: /* Register "canned" boundary conditions and defaults for where to apply. */
327: DMPlexAddBoundary (dm, PETSC_TRUE , "inflow" , "Face Sets" , 0, 0, NULL, (void (*)()) PhysicsBoundary_Advect_Inflow, ALEN(inflowids), inflowids, phys);
328: DMPlexAddBoundary (dm, PETSC_TRUE , "outflow" , "Face Sets" , 0, 0, NULL, (void (*)()) PhysicsBoundary_Advect_Outflow, ALEN(outflowids), outflowids, phys);
329: /* Initial/transient solution with default boundary conditions */
330: ModelSolutionSetDefault(mod,PhysicsSolution_Advect,phys);
331: /* Register "canned" functionals */
332: ModelFunctionalRegister(mod,"Error" ,&advect->functional.Error,PhysicsFunctional_Advect,phys);
333: }
334: return (0);
335: }
337: /******************* Shallow Water ********************/
338: typedef struct {
339: PetscReal gravity;
340: PetscReal boundaryHeight;
341: struct {
342: PetscInt Height;
343: PetscInt Speed;
344: PetscInt Energy;
345: } functional;
346: } Physics_SW;
347: typedef struct {
348: PetscScalar vals[0];
349: PetscScalar h;
350: PetscScalar uh[DIM];
351: } SWNode;
353: static const struct FieldDescription PhysicsFields_SW[] = {{"Height" ,1},{"Momentum" ,DIM},{NULL,0}};
357: /*
358: * h_t + div(uh) = 0
359: * (uh)_t + div (u\otimes uh + g h^2 / 2 I) = 0
360: *
361: * */
362: static PetscErrorCode SWFlux(Physics phys,const PetscReal *n,const SWNode *x,SWNode *f)
363: {
364: Physics_SW *sw = (Physics_SW*)phys->data;
365: PetscScalar uhn,u[DIM];
366: PetscInt i;
369: Scale2(1./x->h,x->uh,u);
370: uhn = Dot2(x->uh,n);
371: f->h = uhn;
372: for (i=0; i<DIM; i++) f->uh[i] = u[i] * uhn + sw->gravity * PetscSqr(x->h) * n[i];
373: return (0);
374: }
378: static PetscErrorCode PhysicsBoundary_SW_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
379: {
381: xG[0] = xI[0];
382: xG[1] = -xI[1];
383: xG[2] = -xI[2];
384: return (0);
385: }
389: static void PhysicsRiemann_SW(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscScalar *flux, Physics phys)
390: {
391: Physics_SW *sw = (Physics_SW*)phys->data;
392: PetscReal cL,cR,speed,nn[DIM];
393: const SWNode *uL = (const SWNode*)xL,*uR = (const SWNode*)xR;
394: SWNode fL,fR;
395: PetscInt i;
397: 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"); */
398: nn[0] = n[0];
399: nn[1] = n[1];
400: Normalize2(nn);
401: SWFlux(phys,nn,uL,&fL);
402: SWFlux(phys,nn,uR,&fR);
403: cL = PetscSqrtReal(sw->gravity*PetscRealPart(uL->h));
404: cR = PetscSqrtReal(sw->gravity*PetscRealPart(uR->h)); /* gravity wave speed */
405: speed = PetscMax(PetscAbsScalar(Dot2(uL->uh,nn)/uL->h) + cL,PetscAbsScalar(Dot2(uR->uh,nn)/uR->h) + cR);
406: 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);
407: }
411: static PetscErrorCode PhysicsSolution_SW(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx)
412: {
413: PetscReal dx[2],r,sigma;
416: if (time != 0.0) SETERRQ1 (mod->comm,PETSC_ERR_SUP,"No solution known for time %g" ,(double)time);
417: dx[0] = x[0] - 1.5;
418: dx[1] = x[1] - 1.0;
419: r = Norm2(dx);
420: sigma = 0.5;
421: u[0] = 1 + 2*PetscExpScalar(-PetscSqr(r)/(2*PetscSqr(sigma)));
422: u[1] = 0.0;
423: u[2] = 0.0;
424: return (0);
425: }
429: static PetscErrorCode PhysicsFunctional_SW(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
430: {
431: Physics phys = (Physics)ctx;
432: Physics_SW *sw = (Physics_SW*)phys->data;
433: const SWNode *x = (const SWNode*)xx;
434: PetscScalar u[2];
435: PetscReal h;
438: h = PetscRealPart(x->h);
439: Scale2(1./x->h,x->uh,u);
440: f[sw->functional.Height] = h;
441: f[sw->functional.Speed] = Norm2(u) + PetscSqrtReal(sw->gravity*h);
442: f[sw->functional.Energy] = 0.5*(Dot2(x->uh,u) + sw->gravity*PetscSqr(h));
443: return (0);
444: }
448: static PetscErrorCode PhysicsCreate_SW(DM dm, Model mod,Physics phys,PetscOptions *PetscOptionsObject)
449: {
450: Physics_SW *sw;
454: phys->field_desc = PhysicsFields_SW;
455: phys->riemann = (RiemannFunction) PhysicsRiemann_SW;
456: PetscNew (&sw);
457: phys->data = sw;
458: PetscOptionsHead (PetscOptionsObject,"SW options" );
459: {
460: sw->gravity = 1.0;
461: PetscOptionsReal ("-sw_gravity" ,"Gravitational constant" ,"" ,sw->gravity,&sw->gravity,NULL);
462: }
463: PetscOptionsTail ();
464: phys->maxspeed = PetscSqrtReal(2.0*sw->gravity); /* Mach 1 for depth of 2 */
466: {
467: const PetscInt wallids[] = {100,101,200,300};
468: DMPlexAddBoundary (dm, PETSC_TRUE , "wall" , "Face Sets" , 0, 0, NULL, (void (*)()) PhysicsBoundary_SW_Wall, ALEN(wallids), wallids, phys);
469: ModelSolutionSetDefault(mod,PhysicsSolution_SW,phys);
470: ModelFunctionalRegister(mod,"Height" ,&sw->functional.Height,PhysicsFunctional_SW,phys);
471: ModelFunctionalRegister(mod,"Speed" ,&sw->functional.Speed,PhysicsFunctional_SW,phys);
472: ModelFunctionalRegister(mod,"Energy" ,&sw->functional.Energy,PhysicsFunctional_SW,phys);
473: }
474: return (0);
475: }
477: /******************* Euler ********************/
478: typedef struct {
479: PetscScalar vals[0];
480: PetscScalar r;
481: PetscScalar ru[DIM];
482: PetscScalar e;
483: } EulerNode;
484: typedef PetscErrorCode (*EquationOfState)(const PetscReal *, const EulerNode*, PetscScalar *) ;
485: typedef struct {
486: PetscInt npars;
487: PetscReal pars[DIM];
488: EquationOfState pressure;
489: EquationOfState sound;
490: struct {
491: PetscInt Density;
492: PetscInt Momentum;
493: PetscInt Energy;
494: PetscInt Pressure;
495: PetscInt Speed;
496: } monitor;
497: } Physics_Euler;
499: static const struct FieldDescription PhysicsFields_Euler[] = {{"Density" ,1},{"Momentum" ,DIM},{"Energy" ,1},{NULL,0}};
503: static PetscErrorCode Pressure_PG(const PetscReal *pars,const EulerNode *x,PetscScalar *p)
504: {
505: PetscScalar ru2;
508: ru2 = DotDIM(x->ru,x->ru);
509: ru2 /= x->r;
510: /* kinematic dof = params[0] */
511: (*p)=2.0*(x->e-0.5*ru2)/pars[0];
512: return (0);
513: }
517: static PetscErrorCode SpeedOfSound_PG(const PetscReal *pars,const EulerNode *x,PetscScalar *c)
518: {
519: PetscScalar p;
522: /* TODO remove direct usage of Pressure_PG */
523: Pressure_PG(pars,x,&p);
524: /* TODO check the sign of p */
525: /* pars[1] = heat capacity ratio */
526: (*c)=PetscSqrtScalar(pars[1]*p/x->r);
527: return (0);
528: }
532: /*
533: * x = (rho,rho*(u_1),...,rho*e)^T
534: * x_t+div(f_1(x))+...+div(f_DIM(x)) = 0
535: *
536: * f_i(x) = u_i*x+(0,0,...,p,...,p*u_i)^T
537: *
538: * */
539: static PetscErrorCode EulerFlux(Physics phys,const PetscReal *n,const EulerNode *x,EulerNode *f)
540: {
541: Physics_Euler *eu = (Physics_Euler*)phys->data;
542: PetscScalar u,nu,p;
543: PetscInt i;
546: u = DotDIM(x->ru,x->ru);
547: u /= (x->r * x->r);
548: nu = DotDIM(x->ru,n);
549: /* TODO check the sign of p */
550: eu->pressure(eu->pars,x,&p);
551: f->r = nu * x->r;
552: for (i=0; i<DIM; i++) f->ru[i] = nu * x->ru[i] + n[i]*p;
553: f->e = nu*(x->e+p);
554: return (0);
555: }
557: /* PetscReal * => EulerNode* conversion */
560: static PetscErrorCode PhysicsBoundary_Euler_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
561: {
562: PetscInt i;
563: PetscScalar xn[DIM],xt[DIM];
566: xG[0] = xI[0];
567: NormalSplitDIM(n,xI+1,xn,xt);
568: for (i=0; i<DIM; i++) xG[i+1] = -xn[i]+xt[i];
569: xG[DIM+1] = xI[DIM+1];
570: return (0);
571: }
573: /* PetscReal * => EulerNode* conversion */
576: 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)
577: {
578: Physics_Euler *eu = (Physics_Euler*)phys->data;
579: PetscScalar cL,cR,speed;
580: const EulerNode *uL = (const EulerNode*)xL,*uR = (const EulerNode*)xR;
581: EulerNode fL,fR;
582: PetscInt i;
584: 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"); */
585: EulerFlux(phys,n,uL,&fL);
586: EulerFlux(phys,n,uR,&fR);
587: eu->sound(eu->pars,uL,&cL);
588: eu->sound(eu->pars,uR,&cR);
589: speed = PetscMax(cL,cR)+PetscMax(PetscAbsScalar(DotDIM(uL->ru,n)/NormDIM(n)),PetscAbsScalar(DotDIM(uR->ru,n)/NormDIM(n)));
590: for (i=0; i<2+dim; i++) flux[i] = 0.5*(fL.vals[i]+fR.vals[i])+0.5*speed*(xL[i]-xR[i]);
591: }
595: static PetscErrorCode PhysicsSolution_Euler(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx)
596: {
597: PetscInt i;
600: if (time != 0.0) SETERRQ1 (mod->comm,PETSC_ERR_SUP,"No solution known for time %g" ,(double)time);
601: u[0] = 1.0;
602: u[DIM+1] = 1.0+PetscAbsReal(x[0]);
603: for (i=1; i<DIM+1; i++) u[i] = 0.0;
604: return (0);
605: }
609: static PetscErrorCode PhysicsFunctional_Euler(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
610: {
611: Physics phys = (Physics)ctx;
612: Physics_Euler *eu = (Physics_Euler*)phys->data;
613: const EulerNode *x = (const EulerNode*)xx;
614: PetscScalar p;
617: f[eu->monitor.Density] = x->r;
618: f[eu->monitor.Momentum] = NormDIM(x->ru);
619: f[eu->monitor.Energy] = x->e;
620: f[eu->monitor.Speed] = NormDIM(x->ru)/x->r;
621: eu->pressure(eu->pars, x, &p);
622: f[eu->monitor.Pressure] = p;
623: return (0);
624: }
628: static PetscErrorCode PhysicsCreate_Euler(DM dm, Model mod,Physics phys,PetscOptions *PetscOptionsObject)
629: {
630: Physics_Euler *eu;
631: PetscErrorCode ierr;
634: phys->field_desc = PhysicsFields_Euler;
635: phys->riemann = (RiemannFunction) PhysicsRiemann_Euler_Rusanov;
636: PetscNew (&eu);
637: phys->data = eu;
638: PetscOptionsHead (PetscOptionsObject,"Euler options" );
639: {
640: eu->pars[0] = 3.0;
641: eu->pars[1] = 1.67;
642: PetscOptionsReal ("-eu_f" ,"Degrees of freedom" ,"" ,eu->pars[0],&eu->pars[0],NULL);
643: PetscOptionsReal ("-eu_gamma" ,"Heat capacity ratio" ,"" ,eu->pars[1],&eu->pars[1],NULL);
644: }
645: PetscOptionsTail ();
646: eu->pressure = Pressure_PG;
647: eu->sound = SpeedOfSound_PG;
648: phys->maxspeed = 1.0;
649: {
650: const PetscInt wallids[] = {100,101,200,300};
651: DMPlexAddBoundary (dm, PETSC_TRUE , "wall" , "Face Sets" , 0, 0, NULL, (void (*)()) PhysicsBoundary_Euler_Wall, ALEN(wallids), wallids, phys);
652: ModelSolutionSetDefault(mod,PhysicsSolution_Euler,phys);
653: ModelFunctionalRegister(mod,"Speed" ,&eu->monitor.Speed,PhysicsFunctional_Euler,phys);
654: ModelFunctionalRegister(mod,"Energy" ,&eu->monitor.Energy,PhysicsFunctional_Euler,phys);
655: ModelFunctionalRegister(mod,"Density" ,&eu->monitor.Density,PhysicsFunctional_Euler,phys);
656: ModelFunctionalRegister(mod,"Momentum" ,&eu->monitor.Momentum,PhysicsFunctional_Euler,phys);
657: ModelFunctionalRegister(mod,"Pressure" ,&eu->monitor.Pressure,PhysicsFunctional_Euler,phys);
658: }
659: return (0);
660: }
664: PetscErrorCode ConstructCellBoundary(DM dm, User user)
665: {
666: const char *name = "Cell Sets" ;
667: const char *bdname = "split faces" ;
668: IS regionIS, innerIS;
669: const PetscInt *regions, *cells;
670: PetscInt numRegions, innerRegion, numCells, c;
671: PetscInt cStart, cEnd, cEndInterior, fStart, fEnd;
672: PetscBool hasLabel;
676: DMPlexGetHeightStratum (dm, 0, &cStart, &cEnd);
677: DMPlexGetHeightStratum (dm, 1, &fStart, &fEnd);
678: DMPlexGetHybridBounds (dm, &cEndInterior, NULL, NULL, NULL);
680: DMPlexHasLabel (dm, name, &hasLabel);
681: if (!hasLabel) return (0);
682: DMPlexGetLabelSize (dm, name, &numRegions);
683: if (numRegions != 2) return (0);
684: /* Get the inner id */
685: DMPlexGetLabelIdIS (dm, name, ®ionIS);
686: ISGetIndices (regionIS, ®ions);
687: innerRegion = regions[0];
688: ISRestoreIndices (regionIS, ®ions);
689: ISDestroy (®ionIS);
690: /* Find the faces between cells in different regions, could call DMPlexCreateNeighborCSR() */
691: DMPlexGetStratumIS (dm, name, innerRegion, &innerIS);
692: ISGetLocalSize (innerIS, &numCells);
693: ISGetIndices (innerIS, &cells);
694: DMPlexCreateLabel (dm, bdname);
695: for (c = 0; c < numCells; ++c) {
696: const PetscInt cell = cells[c];
697: const PetscInt *faces;
698: PetscInt numFaces, f;
700: if ((cell < cStart) || (cell >= cEnd)) SETERRQ1 (PETSC_COMM_SELF , PETSC_ERR_LIB, "Got invalid point %d which is not a cell" , cell);
701: DMPlexGetConeSize (dm, cell, &numFaces);
702: DMPlexGetCone (dm, cell, &faces);
703: for (f = 0; f < numFaces; ++f) {
704: const PetscInt face = faces[f];
705: const PetscInt *neighbors;
706: PetscInt nC, regionA, regionB;
708: if ((face < fStart) || (face >= fEnd)) SETERRQ1 (PETSC_COMM_SELF , PETSC_ERR_LIB, "Got invalid point %d which is not a face" , face);
709: DMPlexGetSupportSize (dm, face, &nC);
710: if (nC != 2) continue ;
711: DMPlexGetSupport (dm, face, &neighbors);
712: if ((neighbors[0] >= cEndInterior) || (neighbors[1] >= cEndInterior)) continue ;
713: 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]);
714: 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]);
715: DMPlexGetLabelValue (dm, name, neighbors[0], ®ionA);
716: DMPlexGetLabelValue (dm, name, neighbors[1], ®ionB);
717: if (regionA < 0) SETERRQ2 (PetscObjectComm ((PetscObject )dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value" , name, neighbors[0]);
718: if (regionB < 0) SETERRQ2 (PetscObjectComm ((PetscObject )dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value" , name, neighbors[1]);
719: if (regionA != regionB) {
720: DMPlexSetLabelValue (dm, bdname, faces[f], 1);
721: }
722: }
723: }
724: ISRestoreIndices (innerIS, &cells);
725: ISDestroy (&innerIS);
726: {
727: DMLabel label;
729: PetscViewerASCIISynchronizedAllow (PETSC_VIEWER_STDOUT_WORLD , PETSC_TRUE );
730: DMPlexGetLabel (dm, bdname, &label);
731: DMLabelView (label, PETSC_VIEWER_STDOUT_WORLD );
732: }
733: return (0);
734: }
738: /* Right now, I have just added duplicate faces, which see both cells. We can
739: - Add duplicate vertices and decouple the face cones
740: - Disconnect faces from cells across the rotation gap
741: */
742: PetscErrorCode SplitFaces(DM *dmSplit, const char labelName[], User user)
743: {
744: DM dm = *dmSplit, sdm;
745: PetscSF sfPoint, gsfPoint;
746: PetscSection coordSection, newCoordSection;
747: Vec coordinates;
748: IS idIS;
749: const PetscInt *ids;
750: PetscInt *newpoints;
751: PetscInt dim, depth, maxConeSize, maxSupportSize, numLabels, numGhostCells;
752: PetscInt numFS, fs, pStart, pEnd, p, cEnd, cEndInterior, vStart, vEnd, v, fStart, fEnd, newf, d, l;
753: PetscBool hasLabel;
757: DMPlexHasLabel (dm, labelName, &hasLabel);
758: if (!hasLabel) return (0);
759: DMCreate (PetscObjectComm ((PetscObject )dm), &sdm);
760: DMSetType (sdm, DMPLEX );
761: DMGetDimension (dm, &dim);
762: DMSetDimension (sdm, dim);
764: DMPlexGetLabelIdIS (dm, labelName, &idIS);
765: ISGetLocalSize (idIS, &numFS);
766: ISGetIndices (idIS, &ids);
768: user->numSplitFaces = 0;
769: for (fs = 0; fs < numFS; ++fs) {
770: PetscInt numBdFaces;
772: DMPlexGetStratumSize (dm, labelName, ids[fs], &numBdFaces);
773: user->numSplitFaces += numBdFaces;
774: }
775: DMPlexGetChart (dm, &pStart, &pEnd);
776: pEnd += user->numSplitFaces;
777: DMPlexSetChart (sdm, pStart, pEnd);
778: DMPlexGetHybridBounds (dm, &cEndInterior, NULL, NULL, NULL);
779: DMPlexSetHybridBounds (sdm, cEndInterior, PETSC_DETERMINE , PETSC_DETERMINE , PETSC_DETERMINE );
780: DMPlexGetHeightStratum (dm, 0, NULL, &cEnd);
781: numGhostCells = cEnd - cEndInterior;
782: /* Set cone and support sizes */
783: DMPlexGetDepth (dm, &depth);
784: for (d = 0; d <= depth; ++d) {
785: DMPlexGetDepthStratum (dm, d, &pStart, &pEnd);
786: for (p = pStart; p < pEnd; ++p) {
787: PetscInt newp = p;
788: PetscInt size;
790: DMPlexGetConeSize (dm, p, &size);
791: DMPlexSetConeSize (sdm, newp, size);
792: DMPlexGetSupportSize (dm, p, &size);
793: DMPlexSetSupportSize (sdm, newp, size);
794: }
795: }
796: DMPlexGetHeightStratum (dm, 1, &fStart, &fEnd);
797: for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
798: IS faceIS;
799: const PetscInt *faces;
800: PetscInt numFaces, f;
802: DMPlexGetStratumIS (dm, labelName, ids[fs], &faceIS);
803: ISGetLocalSize (faceIS, &numFaces);
804: ISGetIndices (faceIS, &faces);
805: for (f = 0; f < numFaces; ++f, ++newf) {
806: PetscInt size;
808: /* Right now I think that both faces should see both cells */
809: DMPlexGetConeSize (dm, faces[f], &size);
810: DMPlexSetConeSize (sdm, newf, size);
811: DMPlexGetSupportSize (dm, faces[f], &size);
812: DMPlexSetSupportSize (sdm, newf, size);
813: }
814: ISRestoreIndices (faceIS, &faces);
815: ISDestroy (&faceIS);
816: }
817: DMSetUp (sdm);
818: /* Set cones and supports */
819: DMPlexGetMaxSizes (dm, &maxConeSize, &maxSupportSize);
820: PetscMalloc (PetscMax(maxConeSize, maxSupportSize) * sizeof (PetscInt ), &newpoints);
821: DMPlexGetChart (dm, &pStart, &pEnd);
822: for (p = pStart; p < pEnd; ++p) {
823: const PetscInt *points, *orientations;
824: PetscInt size, i, newp = p;
826: DMPlexGetConeSize (dm, p, &size);
827: DMPlexGetCone (dm, p, &points);
828: DMPlexGetConeOrientation (dm, p, &orientations);
829: for (i = 0; i < size; ++i) newpoints[i] = points[i];
830: DMPlexSetCone (sdm, newp, newpoints);
831: DMPlexSetConeOrientation (sdm, newp, orientations);
832: DMPlexGetSupportSize (dm, p, &size);
833: DMPlexGetSupport (dm, p, &points);
834: for (i = 0; i < size; ++i) newpoints[i] = points[i];
835: DMPlexSetSupport (sdm, newp, newpoints);
836: }
837: PetscFree (newpoints);
838: for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
839: IS faceIS;
840: const PetscInt *faces;
841: PetscInt numFaces, f;
843: DMPlexGetStratumIS (dm, labelName, ids[fs], &faceIS);
844: ISGetLocalSize (faceIS, &numFaces);
845: ISGetIndices (faceIS, &faces);
846: for (f = 0; f < numFaces; ++f, ++newf) {
847: const PetscInt *points;
849: DMPlexGetCone (dm, faces[f], &points);
850: DMPlexSetCone (sdm, newf, points);
851: DMPlexGetSupport (dm, faces[f], &points);
852: DMPlexSetSupport (sdm, newf, points);
853: }
854: ISRestoreIndices (faceIS, &faces);
855: ISDestroy (&faceIS);
856: }
857: ISRestoreIndices (idIS, &ids);
858: ISDestroy (&idIS);
859: DMPlexStratify (sdm);
860: /* Convert coordinates */
861: DMPlexGetDepthStratum (dm, 0, &vStart, &vEnd);
862: DMGetCoordinateSection (dm, &coordSection);
863: PetscSectionCreate (PetscObjectComm ((PetscObject )dm), &newCoordSection);
864: PetscSectionSetNumFields (newCoordSection, 1);
865: PetscSectionSetFieldComponents (newCoordSection, 0, dim);
866: PetscSectionSetChart (newCoordSection, vStart, vEnd);
867: for (v = vStart; v < vEnd; ++v) {
868: PetscSectionSetDof (newCoordSection, v, dim);
869: PetscSectionSetFieldDof (newCoordSection, v, 0, dim);
870: }
871: PetscSectionSetUp (newCoordSection);
872: DMSetCoordinateSection (sdm, PETSC_DETERMINE , newCoordSection);
873: PetscSectionDestroy (&newCoordSection); /* relinquish our reference */
874: DMGetCoordinatesLocal (dm, &coordinates);
875: DMSetCoordinatesLocal (sdm, coordinates);
876: /* Convert labels */
877: DMPlexGetNumLabels (dm, &numLabels);
878: for (l = 0; l < numLabels; ++l) {
879: const char *lname;
880: PetscBool isDepth;
882: DMPlexGetLabelName (dm, l, &lname);
883: PetscStrcmp (lname, "depth" , &isDepth);
884: if (isDepth) continue ;
885: DMPlexCreateLabel (sdm, lname);
886: DMPlexGetLabelIdIS (dm, lname, &idIS);
887: ISGetLocalSize (idIS, &numFS);
888: ISGetIndices (idIS, &ids);
889: for (fs = 0; fs < numFS; ++fs) {
890: IS pointIS;
891: const PetscInt *points;
892: PetscInt numPoints;
894: DMPlexGetStratumIS (dm, lname, ids[fs], &pointIS);
895: ISGetLocalSize (pointIS, &numPoints);
896: ISGetIndices (pointIS, &points);
897: for (p = 0; p < numPoints; ++p) {
898: PetscInt newpoint = points[p];
900: DMPlexSetLabelValue (sdm, lname, newpoint, ids[fs]);
901: }
902: ISRestoreIndices (pointIS, &points);
903: ISDestroy (&pointIS);
904: }
905: ISRestoreIndices (idIS, &ids);
906: ISDestroy (&idIS);
907: }
908: /* Convert pointSF */
909: const PetscSFNode *remotePoints;
910: PetscSFNode *gremotePoints;
911: const PetscInt *localPoints;
912: PetscInt *glocalPoints,*newLocation,*newRemoteLocation;
913: PetscInt numRoots, numLeaves;
914: PetscMPIInt numProcs;
916: MPI_Comm_size (PetscObjectComm ((PetscObject )dm), &numProcs);
917: DMGetPointSF (dm, &sfPoint);
918: DMGetPointSF (sdm, &gsfPoint);
919: DMPlexGetChart (dm,&pStart,&pEnd);
920: PetscSFGetGraph (sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);
921: if (numRoots >= 0) {
922: PetscMalloc2 (numRoots,&newLocation,pEnd-pStart,&newRemoteLocation);
923: for (l=0; l<numRoots; l++) newLocation[l] = l; /* + (l >= cEnd ? numGhostCells : 0); */
924: PetscSFBcastBegin (sfPoint, MPIU_INT, newLocation, newRemoteLocation);
925: PetscSFBcastEnd (sfPoint, MPIU_INT, newLocation, newRemoteLocation);
926: PetscMalloc1 (numLeaves, &glocalPoints);
927: PetscMalloc1 (numLeaves, &gremotePoints);
928: for (l = 0; l < numLeaves; ++l) {
929: glocalPoints[l] = localPoints[l]; /* localPoints[l] >= cEnd ? localPoints[l] + numGhostCells : localPoints[l]; */
930: gremotePoints[l].rank = remotePoints[l].rank;
931: gremotePoints[l].index = newRemoteLocation[localPoints[l]];
932: }
933: PetscFree2 (newLocation,newRemoteLocation);
934: PetscSFSetGraph (gsfPoint, numRoots+numGhostCells, numLeaves, glocalPoints, PETSC_OWN_POINTER, gremotePoints, PETSC_OWN_POINTER);
935: }
936: DMDestroy (dmSplit);
937: *dmSplit = sdm;
938: return (0);
939: }
943: PetscErrorCode CreatePartitionVec(DM dm, DM *dmCell, Vec *partition)
944: {
945: PetscSF sfPoint;
946: PetscSection coordSection;
947: Vec coordinates;
948: PetscSection sectionCell;
949: PetscScalar *part;
950: PetscInt cStart, cEnd, c;
951: PetscMPIInt rank;
955: DMGetCoordinateSection (dm, &coordSection);
956: DMGetCoordinatesLocal (dm, &coordinates);
957: DMClone (dm, dmCell);
958: DMGetPointSF (dm, &sfPoint);
959: DMSetPointSF (*dmCell, sfPoint);
960: DMSetCoordinateSection (*dmCell, PETSC_DETERMINE , coordSection);
961: DMSetCoordinatesLocal (*dmCell, coordinates);
962: MPI_Comm_rank (PetscObjectComm ((PetscObject )dm), &rank);
963: PetscSectionCreate (PetscObjectComm ((PetscObject )dm), §ionCell);
964: DMPlexGetHeightStratum (*dmCell, 0, &cStart, &cEnd);
965: PetscSectionSetChart (sectionCell, cStart, cEnd);
966: for (c = cStart; c < cEnd; ++c) {
967: PetscSectionSetDof (sectionCell, c, 1);
968: }
969: PetscSectionSetUp (sectionCell);
970: DMSetDefaultSection (*dmCell, sectionCell);
971: PetscSectionDestroy (§ionCell);
972: DMCreateLocalVector (*dmCell, partition);
973: PetscObjectSetName ((PetscObject )*partition, "partition" );
974: VecGetArray (*partition, &part);
975: for (c = cStart; c < cEnd; ++c) {
976: PetscScalar *p;
978: DMPlexPointLocalRef (*dmCell, c, part, &p);
979: p[0] = rank;
980: }
981: VecRestoreArray (*partition, &part);
982: return (0);
983: }
987: PetscErrorCode CreateMassMatrix(DM dm, Vec *massMatrix, User user)
988: {
989: DM dmMass, dmFace, dmCell, dmCoord;
990: PetscSection coordSection;
991: Vec coordinates, facegeom, cellgeom;
992: PetscSection sectionMass;
993: PetscScalar *m;
994: const PetscScalar *fgeom, *cgeom, *coords;
995: PetscInt vStart, vEnd, v;
996: PetscErrorCode ierr;
999: DMGetCoordinateSection (dm, &coordSection);
1000: DMGetCoordinatesLocal (dm, &coordinates);
1001: DMClone (dm, &dmMass);
1002: DMSetCoordinateSection (dmMass, PETSC_DETERMINE , coordSection);
1003: DMSetCoordinatesLocal (dmMass, coordinates);
1004: PetscSectionCreate (PetscObjectComm ((PetscObject )dm), §ionMass);
1005: DMPlexGetDepthStratum (dm, 0, &vStart, &vEnd);
1006: PetscSectionSetChart (sectionMass, vStart, vEnd);
1007: for (v = vStart; v < vEnd; ++v) {
1008: PetscInt numFaces;
1010: DMPlexGetSupportSize (dmMass, v, &numFaces);
1011: PetscSectionSetDof (sectionMass, v, numFaces*numFaces);
1012: }
1013: PetscSectionSetUp (sectionMass);
1014: DMSetDefaultSection (dmMass, sectionMass);
1015: PetscSectionDestroy (§ionMass);
1016: DMGetLocalVector (dmMass, massMatrix);
1017: VecGetArray (*massMatrix, &m);
1018: DMPlexTSGetGeometryFVM (dm, &facegeom, &cellgeom, NULL);
1019: VecGetDM (facegeom, &dmFace);
1020: VecGetArrayRead (facegeom, &fgeom);
1021: VecGetDM (cellgeom, &dmCell);
1022: VecGetArrayRead (cellgeom, &cgeom);
1023: DMGetCoordinateDM (dm, &dmCoord);
1024: VecGetArrayRead (coordinates, &coords);
1025: for (v = vStart; v < vEnd; ++v) {
1026: const PetscInt *faces;
1027: const PetscFVFaceGeom *fgA, *fgB, *cg;
1028: const PetscScalar *vertex;
1029: PetscInt numFaces, sides[2], f, g;
1031: DMPlexPointLocalRead (dmCoord, v, coords, &vertex);
1032: DMPlexGetSupportSize (dmMass, v, &numFaces);
1033: DMPlexGetSupport (dmMass, v, &faces);
1034: for (f = 0; f < numFaces; ++f) {
1035: sides[0] = faces[f];
1036: DMPlexPointLocalRead (dmFace, faces[f], fgeom, &fgA);
1037: for (g = 0; g < numFaces; ++g) {
1038: const PetscInt *cells = NULL;;
1039: PetscReal area = 0.0;
1040: PetscInt numCells;
1042: sides[1] = faces[g];
1043: DMPlexPointLocalRead (dmFace, faces[g], fgeom, &fgB);
1044: DMPlexGetJoin (dmMass, 2, sides, &numCells, &cells);
1045: if (numCells != 1) SETERRQ (PETSC_COMM_SELF , PETSC_ERR_LIB, "Invalid join for faces" );
1046: DMPlexPointLocalRead (dmCell, cells[0], cgeom, &cg);
1047: area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgA->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgA->centroid[0] - cg->centroid[0]));
1048: area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgB->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgB->centroid[0] - cg->centroid[0]));
1049: m[f*numFaces+g] = Dot2(fgA->normal, fgB->normal)*area*0.5;
1050: DMPlexRestoreJoin (dmMass, 2, sides, &numCells, &cells);
1051: }
1052: }
1053: }
1054: VecRestoreArrayRead (facegeom, &fgeom);
1055: VecRestoreArrayRead (cellgeom, &cgeom);
1056: VecRestoreArrayRead (coordinates, &coords);
1057: VecRestoreArray (*massMatrix, &m);
1058: DMDestroy (&dmMass);
1059: return (0);
1060: }
1064: /* Behavior will be different for multi-physics or when using non-default boundary conditions */
1065: static PetscErrorCode ModelSolutionSetDefault(Model mod,SolutionFunction func,void *ctx)
1066: {
1068: mod->solution = func;
1069: mod->solutionctx = ctx;
1070: return (0);
1071: }
1075: static PetscErrorCode ModelFunctionalRegister(Model mod,const char *name,PetscInt *offset,FunctionalFunction func,void *ctx)
1076: {
1078: FunctionalLink link,*ptr;
1079: PetscInt lastoffset = -1;
1082: for (ptr=&mod->functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset;
1083: PetscNew (&link);
1084: PetscStrallocpy (name,&link->name);
1085: link->offset = lastoffset + 1;
1086: link->func = func;
1087: link->ctx = ctx;
1088: link->next = NULL;
1089: *ptr = link;
1090: *offset = link->offset;
1091: return (0);
1092: }
1096: static PetscErrorCode ModelFunctionalSetFromOptions(Model mod,PetscOptions *PetscOptionsObject)
1097: {
1099: PetscInt i,j;
1100: FunctionalLink link;
1101: char *names[256];
1104: mod->numMonitored = ALEN(names);
1105: PetscOptionsStringArray ("-monitor" ,"list of functionals to monitor" ,"" ,names,&mod->numMonitored,NULL);
1106: /* Create list of functionals that will be computed somehow */
1107: PetscMalloc1 (mod->numMonitored,&mod->functionalMonitored);
1108: /* Create index of calls that we will have to make to compute these functionals (over-allocation in general). */
1109: PetscMalloc1 (mod->numMonitored,&mod->functionalCall);
1110: mod->numCall = 0;
1111: for (i=0; i<mod->numMonitored; i++) {
1112: for (link=mod->functionalRegistry; link; link=link->next) {
1113: PetscBool match;
1114: PetscStrcasecmp (names[i],link->name,&match);
1115: if (match) break ;
1116: }
1117: if (!link) SETERRQ1 (mod->comm,PETSC_ERR_USER,"No known functional '%s'" ,names[i]);
1118: mod->functionalMonitored[i] = link;
1119: for (j=0; j<i; j++) {
1120: if (mod->functionalCall[j]->func == link->func && mod->functionalCall[j]->ctx == link->ctx) goto next_name;
1121: }
1122: mod->functionalCall[mod->numCall++] = link; /* Just points to the first link using the result. There may be more results. */
1123: next_name:
1124: PetscFree (names[i]);
1125: }
1127: /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */
1128: mod->maxComputed = -1;
1129: for (link=mod->functionalRegistry; link; link=link->next) {
1130: for (i=0; i<mod->numCall; i++) {
1131: FunctionalLink call = mod->functionalCall[i];
1132: if (link->func == call->func && link->ctx == call->ctx) {
1133: mod->maxComputed = PetscMax(mod->maxComputed,link->offset);
1134: }
1135: }
1136: }
1137: return (0);
1138: }
1142: static PetscErrorCode FunctionalLinkDestroy(FunctionalLink *link)
1143: {
1145: FunctionalLink l,next;
1148: if (!link) return (0);
1149: l = *link;
1150: *link = NULL;
1151: for (; l; l=next) {
1152: next = l->next;
1153: PetscFree (l->name);
1154: PetscFree (l);
1155: }
1156: return (0);
1157: }
1161: PetscErrorCode SetInitialCondition(DM dm, Vec X, User user)
1162: {
1163: DM dmCell;
1164: Model mod = user->model;
1165: Vec cellgeom;
1166: const PetscScalar *cgeom;
1167: PetscScalar *x;
1168: PetscInt cStart, cEnd, cEndInterior, c;
1169: PetscErrorCode ierr;
1172: DMPlexGetHybridBounds (dm, &cEndInterior, NULL, NULL, NULL);
1173: DMPlexTSGetGeometryFVM (dm, NULL, &cellgeom, NULL);
1174: VecGetDM (cellgeom, &dmCell);
1175: DMPlexGetHeightStratum (dm, 0, &cStart, &cEnd);
1176: VecGetArrayRead (cellgeom, &cgeom);
1177: VecGetArray (X, &x);
1178: for (c = cStart; c < cEndInterior; ++c) {
1179: const PetscFVCellGeom *cg;
1180: PetscScalar *xc;
1182: DMPlexPointLocalRead (dmCell,c,cgeom,&cg);
1183: DMPlexPointGlobalRef (dm,c,x,&xc);
1184: if (xc) {(*mod->solution)(mod,0.0,cg->centroid,xc,mod->solutionctx);}
1185: }
1186: VecRestoreArrayRead (cellgeom, &cgeom);
1187: VecRestoreArray (X, &x);
1188: return (0);
1189: }
1193: static PetscErrorCode OutputVTK(DM dm, const char *filename, PetscViewer *viewer)
1194: {
1198: PetscViewerCreate (PetscObjectComm ((PetscObject )dm), viewer);
1199: PetscViewerSetType (*viewer, PETSCVIEWERVTK);
1200: PetscViewerFileSetName (*viewer, filename);
1201: return (0);
1202: }
1206: static PetscErrorCode MonitorVTK(TS ts,PetscInt stepnum,PetscReal time,Vec X,void *ctx)
1207: {
1208: User user = (User)ctx;
1209: DM dm;
1210: Vec cellgeom;
1211: PetscViewer viewer;
1212: char filename[PETSC_MAX_PATH_LEN],*ftable = NULL;
1213: PetscReal xnorm;
1214: PetscInt cEndInterior;
1218: PetscObjectSetName ((PetscObject ) X, "solution" );
1219: VecGetDM (X,&dm);
1220: DMPlexTSGetGeometryFVM (dm, NULL, &cellgeom, NULL);
1221: DMPlexGetHybridBounds (dm, &cEndInterior, NULL, NULL, NULL);
1222: VecNorm (X,NORM_INFINITY ,&xnorm);
1223: if (stepnum >= 0) { /* No summary for final time */
1224: Model mod = user->model;
1225: PetscInt c,cStart,cEnd,fcount,i;
1226: size_t ftableused,ftablealloc;
1227: const PetscScalar *cgeom,*x;
1228: DM dmCell;
1229: PetscReal *fmin,*fmax,*fintegral,*ftmp;
1230: fcount = mod->maxComputed+1;
1231: PetscMalloc4 (fcount,&fmin,fcount,&fmax,fcount,&fintegral,fcount,&ftmp);
1232: for (i=0; i<fcount; i++) {
1233: fmin[i] = PETSC_MAX_REAL;
1234: fmax[i] = PETSC_MIN_REAL;
1235: fintegral[i] = 0;
1236: }
1237: DMPlexGetHeightStratum (dm,0,&cStart,&cEnd);
1238: VecGetDM (cellgeom,&dmCell);
1239: VecGetArrayRead (cellgeom,&cgeom);
1240: VecGetArrayRead (X,&x);
1241: for (c = cStart; c < cEndInterior; ++c) {
1242: const PetscFVCellGeom *cg;
1243: const PetscScalar *cx;
1244: DMPlexPointLocalRead (dmCell,c,cgeom,&cg);
1245: DMPlexPointGlobalRead (dm,c,x,&cx);
1246: if (!cx) continue ; /* not a global cell */
1247: for (i=0; i<mod->numCall; i++) {
1248: FunctionalLink flink = mod->functionalCall[i];
1249: (*flink->func)(mod,time,cg->centroid,cx,ftmp,flink->ctx);
1250: }
1251: for (i=0; i<fcount; i++) {
1252: fmin[i] = PetscMin(fmin[i],ftmp[i]);
1253: fmax[i] = PetscMax(fmax[i],ftmp[i]);
1254: fintegral[i] += cg->volume * ftmp[i];
1255: }
1256: }
1257: VecRestoreArrayRead (cellgeom,&cgeom);
1258: VecRestoreArrayRead (X,&x);
1259: MPI_Allreduce (MPI_IN_PLACE,fmin,fcount,MPIU_REAL,MPIU_MIN,PetscObjectComm ((PetscObject )ts));
1260: MPI_Allreduce (MPI_IN_PLACE,fmax,fcount,MPIU_REAL,MPIU_MAX,PetscObjectComm ((PetscObject )ts));
1261: MPI_Allreduce (MPI_IN_PLACE,fintegral,fcount,MPIU_REAL,MPIU_SUM,PetscObjectComm ((PetscObject )ts));
1263: ftablealloc = fcount * 100;
1264: ftableused = 0;
1265: PetscMalloc1 (ftablealloc,&ftable);
1266: for (i=0; i<mod->numMonitored; i++) {
1267: size_t countused;
1268: char buffer[256],*p;
1269: FunctionalLink flink = mod->functionalMonitored[i];
1270: PetscInt id = flink->offset;
1271: if (i % 3) {
1272: PetscMemcpy (buffer," " ,2);
1273: p = buffer + 2;
1274: } else if (i) {
1275: char newline[] = "\n" ;
1276: PetscMemcpy (buffer,newline,sizeof newline-1);
1277: p = buffer + sizeof newline - 1;
1278: } else {
1279: p = buffer;
1280: }
1281: 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]);
1282: countused += p - buffer;
1283: if (countused > ftablealloc-ftableused-1) { /* reallocate */
1284: char *ftablenew;
1285: ftablealloc = 2*ftablealloc + countused;
1286: PetscMalloc (ftablealloc,&ftablenew);
1287: PetscMemcpy (ftablenew,ftable,ftableused);
1288: PetscFree (ftable);
1289: ftable = ftablenew;
1290: }
1291: PetscMemcpy (ftable+ftableused,buffer,countused);
1292: ftableused += countused;
1293: ftable[ftableused] = 0;
1294: }
1295: PetscFree4 (fmin,fmax,fintegral,ftmp);
1297: PetscPrintf (PetscObjectComm ((PetscObject )ts),"% 3D time %8.4g |x| %8.4g %s\n" ,stepnum,(double)time,(double)xnorm,ftable ? ftable : "" );
1298: PetscFree (ftable);
1299: }
1300: if (user->vtkInterval < 1) return (0);
1301: if ((stepnum == -1) ^ (stepnum % user->vtkInterval == 0)) {
1302: if (stepnum == -1) { /* Final time is not multiple of normal time interval, write it anyway */
1303: TSGetTimeStepNumber (ts,&stepnum);
1304: }
1305: PetscSNPrintf (filename,sizeof filename,"ex11-%03D.vtu" ,stepnum);
1306: OutputVTK(dm,filename,&viewer);
1307: VecView (X,viewer);
1308: PetscViewerDestroy (&viewer);
1309: }
1310: return (0);
1311: }
1315: int main(int argc, char **argv)
1316: {
1317: MPI_Comm comm;
1318: PetscDS prob;
1319: PetscFV fvm;
1320: User user;
1321: Model mod;
1322: Physics phys;
1323: DM dm;
1324: PetscReal ftime, cfl, dt, minRadius;
1325: PetscInt dim, nsteps;
1326: TS ts;
1327: TSConvergedReason reason;
1328: Vec X;
1329: PetscViewer viewer;
1330: PetscBool vtkCellGeom, splitFaces;
1331: PetscInt overlap;
1332: char filename[PETSC_MAX_PATH_LEN] = "sevenside.exo" ;
1333: PetscErrorCode ierr;
1335: PetscInitialize (&argc, &argv, (char*) 0, help);
1336: comm = PETSC_COMM_WORLD ;
1338: PetscNew (&user);
1339: PetscNew (&user->model);
1340: PetscNew (&user->model->physics);
1341: mod = user->model;
1342: phys = mod->physics;
1343: mod->comm = comm;
1345: /* Register physical models to be available on the command line */
1346: PetscFunctionListAdd (&PhysicsList,"advect" ,PhysicsCreate_Advect);
1347: PetscFunctionListAdd (&PhysicsList,"sw" ,PhysicsCreate_SW);
1348: PetscFunctionListAdd (&PhysicsList,"euler" ,PhysicsCreate_Euler);
1351: PetscOptionsBegin (comm,NULL,"Unstructured Finite Volume Mesh Options" ,"" );
1352: {
1353: cfl = 0.9 * 4; /* default SSPRKS2 with s=5 stages is stable for CFL number s-1 */
1354: PetscOptionsReal ("-ufv_cfl" ,"CFL number per step" ,"" ,cfl,&cfl,NULL);
1355: PetscOptionsString ("-f" ,"Exodus.II filename to read" ,"" ,filename,filename,sizeof (filename),NULL);
1356: splitFaces = PETSC_FALSE ;
1357: PetscOptionsBool ("-ufv_split_faces" ,"Split faces between cell sets" ,"" ,splitFaces,&splitFaces,NULL);
1358: overlap = 1;
1359: PetscOptionsInt ("-ufv_mesh_overlap" ,"Number of cells to overlap partitions" ,"" ,overlap,&overlap,NULL);
1360: user->vtkInterval = 1;
1361: PetscOptionsInt ("-ufv_vtk_interval" ,"VTK output interval (0 to disable)" ,"" ,user->vtkInterval,&user->vtkInterval,NULL);
1362: vtkCellGeom = PETSC_FALSE ;
1363: PetscOptionsBool ("-ufv_vtk_cellgeom" ,"Write cell geometry (for debugging)" ,"" ,vtkCellGeom,&vtkCellGeom,NULL);
1364: }
1365: PetscOptionsEnd ();
1366: DMPlexCreateFromFile (comm, filename, PETSC_TRUE , &dm);
1367: DMViewFromOptions(dm, NULL, "-dm_view" );
1368: DMGetDimension (dm, &dim);
1370: PetscOptionsBegin (comm,NULL,"Unstructured Finite Volume Physics Options" ,"" );
1371: {
1372: PetscErrorCode (*physcreate)(DM ,Model,Physics,PetscOptions*);
1373: char physname[256] = "advect" ;
1375: DMPlexCreateLabel (dm, "Face Sets" );
1376: PetscOptionsFList ("-physics" ,"Physics module to solve" ,"" ,PhysicsList,physname,physname,sizeof physname,NULL);
1377: PetscFunctionListFind (PhysicsList,physname,&physcreate);
1378: PetscMemzero (phys,sizeof (struct _n_Physics ));
1379: (*physcreate)(dm,mod,phys,PetscOptionsObject);
1380: mod->maxspeed = phys->maxspeed;
1381: /* Count number of fields and dofs */
1382: for (phys->nfields=0,phys->dof=0; phys->field_desc[phys->nfields].name; phys->nfields++) phys->dof += phys->field_desc[phys->nfields].dof;
1384: if (mod->maxspeed <= 0) SETERRQ1 (comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set maxspeed" ,physname);
1385: if (phys->dof <= 0) SETERRQ1 (comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set dof" ,physname);
1386: ModelFunctionalSetFromOptions(mod,PetscOptionsObject);
1387: }
1388: PetscOptionsEnd ();
1389: {
1390: DM dmDist;
1392: DMPlexSetAdjacencyUseCone (dm, PETSC_TRUE );
1393: DMPlexSetAdjacencyUseClosure (dm, PETSC_FALSE );
1394: DMPlexDistribute (dm, overlap, NULL, &dmDist);
1395: if (dmDist) {
1396: DMDestroy (&dm);
1397: dm = dmDist;
1398: }
1399: }
1400: DMSetFromOptions (dm);
1401: {
1402: DM gdm;
1404: DMPlexConstructGhostCells (dm, NULL, NULL, &gdm);
1405: DMDestroy (&dm);
1406: dm = gdm;
1407: DMViewFromOptions(dm, NULL, "-dm_view" );
1408: }
1409: if (splitFaces) {ConstructCellBoundary(dm, user);}
1410: SplitFaces(&dm, "split faces" , user);
1412: PetscFVCreate (comm, &fvm);
1413: PetscFVSetFromOptions (fvm);
1414: PetscFVSetNumComponents (fvm, phys->dof);
1415: PetscFVSetSpatialDimension (fvm, dim);
1416: DMGetDS (dm, &prob);
1417: /* FV is now structured with one field having all physics as components */
1418: PetscDSAddDiscretization (prob, (PetscObject ) fvm);
1419: PetscDSSetRiemannSolver (prob, 0, user->model->physics->riemann);
1420: PetscDSSetContext(prob, 0, user->model->physics);
1422: TSCreate (comm, &ts);
1423: TSSetType (ts, TSSSP );
1424: TSSetDM (ts, dm);
1425: TSMonitorSet (ts,MonitorVTK,user,NULL);
1426: DMTSSetRHSFunctionLocal (dm, DMPlexTSComputeRHSFunctionFVM , user);
1428: DMCreateGlobalVector (dm, &X);
1429: PetscObjectSetName ((PetscObject ) X, "solution" );
1430: SetInitialCondition(dm, X, user);
1431: if (vtkCellGeom) {
1432: DM dmCell;
1433: Vec cellgeom, partition;
1435: DMPlexTSGetGeometryFVM (dm, NULL, &cellgeom, NULL);
1436: OutputVTK(dm, "ex11-cellgeom.vtk" , &viewer);
1437: VecView (cellgeom, viewer);
1438: PetscViewerDestroy (&viewer);
1439: CreatePartitionVec(dm, &dmCell, &partition);
1440: OutputVTK(dmCell, "ex11-partition.vtk" , &viewer);
1441: VecView (partition, viewer);
1442: PetscViewerDestroy (&viewer);
1443: VecDestroy (&partition);
1444: DMDestroy (&dmCell);
1445: }
1447: DMPlexTSGetGeometryFVM (dm, NULL, NULL, &minRadius);
1448: TSSetDuration (ts,1000,2.0);
1449: dt = cfl * minRadius / user->model->maxspeed;
1450: TSSetInitialTimeStep (ts,0.0,dt);
1451: TSSetFromOptions (ts);
1452: TSSolve (ts,X);
1453: TSGetSolveTime (ts,&ftime);
1454: TSGetTimeStepNumber (ts,&nsteps);
1455: TSGetConvergedReason (ts,&reason);
1456: PetscPrintf (PETSC_COMM_WORLD ,"%s at time %g after %D steps\n" ,TSConvergedReasons[reason],(double)ftime,nsteps);
1457: TSDestroy (&ts);
1459: PetscFunctionListDestroy (&PhysicsList);
1460: FunctionalLinkDestroy(&user->model->functionalRegistry);
1461: PetscFree (user->model->functionalMonitored);
1462: PetscFree (user->model->functionalCall);
1463: PetscFree (user->model->physics->data);
1464: PetscFree (user->model->physics);
1465: PetscFree (user->model);
1466: PetscFree (user);
1467: VecDestroy (&X);
1468: PetscFVDestroy (&fvm);
1469: DMDestroy (&dm);
1470: PetscFinalize ();
1471: return (0);
1472: }