Actual source code: ex11.c
petsc-3.5.4 2015-05-23
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 <petscts.h>
38: #include <petscfv.h>
39: #include <petscdmplex.h>
40: #include <petscsf.h>
41: #include <petscblaslapack.h>
43: #define DIM 2 /* Geometric dimension */
44: #define ALEN(a) (sizeof(a)/sizeof((a)[0]))
46: static PetscFunctionList PhysicsList;
48: /* Represents continuum physical equations. */
49: typedef struct _n_Physics *Physics;
51: /* Physical model includes boundary conditions, initial conditions, and functionals of interest. It is
52: * discretization-independent, but its members depend on the scenario being solved. */
53: typedef struct _n_Model *Model;
55: /* 'User' implements a discretization of a continuous model. */
56: typedef struct _n_User *User;
58: typedef PetscErrorCode (*RiemannFunction)(const PetscReal *,const PetscReal *,const PetscScalar *,const PetscScalar *,PetscScalar *,void*) ;
59: typedef PetscErrorCode (*SolutionFunction)(Model,PetscReal ,const PetscReal *,PetscScalar *,void*) ;
60: typedef PetscErrorCode (*FunctionalFunction)(Model,PetscReal ,const PetscReal *,const PetscScalar *,PetscReal *,void*) ;
61: typedef PetscErrorCode (*SetupFields)(Physics,PetscSection ) ;
62: static PetscErrorCode ModelSolutionSetDefault(Model,SolutionFunction,void*) ;
63: static PetscErrorCode ModelFunctionalRegister(Model,const char*,PetscInt *,FunctionalFunction,void*) ;
64: static PetscErrorCode OutputVTK(DM,const char*,PetscViewer *) ;
66: struct FieldDescription {
67: const char *name;
68: PetscInt dof;
69: };
71: typedef struct _n_FunctionalLink *FunctionalLink;
72: struct _n_FunctionalLink {
73: char *name;
74: FunctionalFunction func;
75: void *ctx;
76: PetscInt offset;
77: FunctionalLink next;
78: };
80: struct _n_Physics {
81: RiemannFunction riemann;
82: PetscInt dof; /* number of degrees of freedom per cell */
83: PetscReal maxspeed; /* kludge to pick initial time step, need to add monitoring and step control */
84: void *data;
85: PetscInt nfields;
86: const struct FieldDescription *field_desc;
87: };
89: struct _n_Model {
90: MPI_Comm comm; /* Does not do collective communicaton, but some error conditions can be collective */
91: Physics physics;
92: FunctionalLink functionalRegistry;
93: PetscInt maxComputed;
94: PetscInt numMonitored;
95: FunctionalLink *functionalMonitored;
96: PetscInt numCall;
97: FunctionalLink *functionalCall;
98: SolutionFunction solution;
99: void *solutionctx;
100: PetscReal maxspeed; /* estimate of global maximum speed (for CFL calculation) */
101: };
103: struct _n_User {
104: PetscInt numSplitFaces;
105: PetscInt vtkInterval; /* For monitor */
106: Model model;
107: };
109: PETSC_STATIC_INLINE PetscScalar DotDIM(const PetscScalar *x,const PetscScalar *y)
110: {
111: PetscInt i;
112: PetscScalar prod=0.0;
114: for (i=0; i<DIM; i++) prod += x[i]*y[i];
115: return prod;
116: }
117: PETSC_STATIC_INLINE PetscReal NormDIM(const PetscScalar *x) { return PetscSqrtReal(PetscAbsScalar(DotDIM(x,x))); }
118: PETSC_STATIC_INLINE void axDIM(const PetscScalar a,PetscScalar *x)
119: {
120: PetscInt i;
121: for (i=0; i<DIM; i++) x[i] *= a;
122: }
123: PETSC_STATIC_INLINE void waxDIM(const PetscScalar a,const PetscScalar *x, PetscScalar *w)
124: {
125: PetscInt i;
126: for (i=0; i<DIM; i++) w[i] = x[i]*a;
127: }
128: PETSC_STATIC_INLINE void NormalSplitDIM(const PetscReal *n,const PetscScalar *x,PetscScalar *xn,PetscScalar *xt)
129: { /* Split x into normal and tangential components */
130: PetscInt i;
131: PetscScalar c;
132: c = DotDIM(x,n)/DotDIM(n,n);
133: for (i=0; i<DIM; i++) {
134: xn[i] = c*n[i];
135: xt[i] = x[i]-xn[i];
136: }
137: }
139: PETSC_STATIC_INLINE PetscScalar Dot2(const PetscScalar *x,const PetscScalar *y) { return x[0]*y[0] + x[1]*y[1];}
140: PETSC_STATIC_INLINE PetscReal Norm2(const PetscScalar *x) { return PetscSqrtReal(PetscAbsScalar(Dot2(x,x)));}
141: PETSC_STATIC_INLINE void Normalize2(PetscScalar *x) { PetscReal a = 1./Norm2(x); x[0] *= a; x[1] *= a; }
142: 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]; }
143: PETSC_STATIC_INLINE void Scale2(PetscScalar a,const PetscScalar *x,PetscScalar *y) { y[0] = a*x[0]; y[1] = a*x[1]; }
145: 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];}
146: 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;}
147: PETSC_STATIC_INLINE PetscReal NormD(PetscInt dim, const PetscScalar *x) {return PetscSqrtReal(PetscAbsScalar(DotD(dim,x,x)));}
149: PETSC_STATIC_INLINE void NormalSplit(const PetscReal *n,const PetscScalar *x,PetscScalar *xn,PetscScalar *xt)
150: { /* Split x into normal and tangential components */
151: Scale2(Dot2(x,n)/Dot2(n,n),n,xn);
152: Waxpy2(-1,xn,x,xt);
153: }
155: /******************* Advect ********************/
156: typedef enum {ADVECT_SOL_TILTED,ADVECT_SOL_BUMP} AdvectSolType;
157: static const char *const AdvectSolTypes[] = {"TILTED" ,"BUMP" ,"AdvectSolType" ,"ADVECT_SOL_" ,0};
158: typedef enum {ADVECT_SOL_BUMP_CONE,ADVECT_SOL_BUMP_COS} AdvectSolBumpType;
159: static const char *const AdvectSolBumpTypes[] = {"CONE" ,"COS" ,"AdvectSolBumpType" ,"ADVECT_SOL_BUMP_" ,0};
161: typedef struct {
162: PetscReal wind[DIM];
163: } Physics_Advect_Tilted;
164: typedef struct {
165: PetscReal center[DIM];
166: PetscReal radius;
167: AdvectSolBumpType type;
168: } Physics_Advect_Bump;
170: typedef struct {
171: PetscReal inflowState;
172: AdvectSolType soltype;
173: union {
174: Physics_Advect_Tilted tilted;
175: Physics_Advect_Bump bump;
176: } sol;
177: struct {
178: PetscInt Error;
179: } functional;
180: } Physics_Advect;
182: static const struct FieldDescription PhysicsFields_Advect[] = {{"U" ,1},{NULL,0}};
186: static PetscErrorCode PhysicsBoundary_Advect_Inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
187: {
188: Physics phys = (Physics)ctx;
189: Physics_Advect *advect = (Physics_Advect*)phys->data;
192: xG[0] = advect->inflowState;
193: return (0);
194: }
198: static PetscErrorCode PhysicsBoundary_Advect_Outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
199: {
201: xG[0] = xI[0];
202: return (0);
203: }
207: static PetscErrorCode PhysicsRiemann_Advect(const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscScalar *flux, Physics phys)
208: {
209: Physics_Advect *advect = (Physics_Advect*)phys->data;
210: PetscReal wind[DIM],wn;
213: switch (advect->soltype) {
214: case ADVECT_SOL_TILTED: {
215: Physics_Advect_Tilted *tilted = &advect->sol.tilted;
216: wind[0] = tilted->wind[0];
217: wind[1] = tilted->wind[1];
218: } break ;
219: case ADVECT_SOL_BUMP:
220: wind[0] = -qp[1];
221: wind[1] = qp[0];
222: break ;
223: default: SETERRQ1 (PETSC_COMM_SELF ,PETSC_ERR_SUP,"No support for solution type %s" ,AdvectSolBumpTypes[advect->soltype]);
224: }
225: wn = Dot2(wind, n);
226: flux[0] = (wn > 0 ? xL[0] : xR[0]) * wn;
227: return (0);
228: }
232: static PetscErrorCode PhysicsSolution_Advect(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx)
233: {
234: Physics phys = (Physics)ctx;
235: Physics_Advect *advect = (Physics_Advect*)phys->data;
238: switch (advect->soltype) {
239: case ADVECT_SOL_TILTED: {
240: PetscReal x0[DIM];
241: Physics_Advect_Tilted *tilted = &advect->sol.tilted;
242: Waxpy2(-time,tilted->wind,x,x0);
243: if (x0[1] > 0) u[0] = 1.*x[0] + 3.*x[1];
244: else u[0] = advect->inflowState;
245: } break ;
246: case ADVECT_SOL_BUMP: {
247: Physics_Advect_Bump *bump = &advect->sol.bump;
248: PetscReal x0[DIM],v[DIM],r,cost,sint;
249: cost = PetscCosReal(time);
250: sint = PetscSinReal(time);
251: x0[0] = cost*x[0] + sint*x[1];
252: x0[1] = -sint*x[0] + cost*x[1];
253: Waxpy2(-1,bump->center,x0,v);
254: r = Norm2(v);
255: switch (bump->type) {
256: case ADVECT_SOL_BUMP_CONE:
257: u[0] = PetscMax(1 - r/bump->radius,0);
258: break ;
259: case ADVECT_SOL_BUMP_COS:
260: u[0] = 0.5 + 0.5*PetscCosReal(PetscMin(r/bump->radius,1)*PETSC_PI);
261: break ;
262: }
263: } break ;
264: default: SETERRQ (PETSC_COMM_SELF ,PETSC_ERR_SUP,"Unknown solution type" );
265: }
266: return (0);
267: }
271: static PetscErrorCode PhysicsFunctional_Advect(Model mod,PetscReal time,const PetscScalar *x,const PetscScalar *y,PetscReal *f,void *ctx)
272: {
273: Physics phys = (Physics)ctx;
274: Physics_Advect *advect = (Physics_Advect*)phys->data;
275: PetscScalar yexact[1];
279: PhysicsSolution_Advect(mod,time,x,yexact,phys);
280: f[advect->functional.Error] = PetscAbsScalar(y[0]-yexact[0]);
281: return (0);
282: }
286: static PetscErrorCode PhysicsCreate_Advect(DM dm, Model mod,Physics phys)
287: {
288: Physics_Advect *advect;
292: phys->field_desc = PhysicsFields_Advect;
293: phys->riemann = (RiemannFunction) PhysicsRiemann_Advect;
294: PetscNew (&advect);
295: phys->data = advect;
296: PetscOptionsHead ("Advect options" );
297: {
298: PetscInt two = 2,dof = 1;
299: advect->soltype = ADVECT_SOL_TILTED;
300: PetscOptionsEnum ("-advect_sol_type" ,"solution type" ,"" ,AdvectSolTypes,(PetscEnum )advect->soltype,(PetscEnum *)&advect->soltype,NULL);
301: switch (advect->soltype) {
302: case ADVECT_SOL_TILTED: {
303: Physics_Advect_Tilted *tilted = &advect->sol.tilted;
304: two = 2;
305: tilted->wind[0] = 0.0;
306: tilted->wind[1] = 1.0;
307: PetscOptionsRealArray ("-advect_tilted_wind" ,"background wind vx,vy" ,"" ,tilted->wind,&two,NULL);
308: advect->inflowState = -2.0;
309: PetscOptionsRealArray ("-advect_tilted_inflow" ,"Inflow state" ,"" ,&advect->inflowState,&dof,NULL);
310: phys->maxspeed = Norm2(tilted->wind);
311: } break ;
312: case ADVECT_SOL_BUMP: {
313: Physics_Advect_Bump *bump = &advect->sol.bump;
314: two = 2;
315: bump->center[0] = 2.;
316: bump->center[1] = 0.;
317: PetscOptionsRealArray ("-advect_bump_center" ,"location of center of bump x,y" ,"" ,bump->center,&two,NULL);
318: bump->radius = 0.9;
319: PetscOptionsReal ("-advect_bump_radius" ,"radius of bump" ,"" ,bump->radius,&bump->radius,NULL);
320: bump->type = ADVECT_SOL_BUMP_CONE;
321: PetscOptionsEnum ("-advect_bump_type" ,"type of bump" ,"" ,AdvectSolBumpTypes,(PetscEnum )bump->type,(PetscEnum *)&bump->type,NULL);
322: phys->maxspeed = 3.; /* radius of mesh, kludge */
323: } break ;
324: }
325: }
326: PetscOptionsTail ();
327: {
328: const PetscInt inflowids[] = {100,200,300},outflowids[] = {101};
329: /* Register "canned" boundary conditions and defaults for where to apply. */
330: DMPlexAddBoundary(dm, PETSC_TRUE , "inflow" , "Face Sets" , 0, (void (*)()) PhysicsBoundary_Advect_Inflow, ALEN(inflowids), inflowids, phys);
331: DMPlexAddBoundary(dm, PETSC_TRUE , "outflow" , "Face Sets" , 0, (void (*)()) PhysicsBoundary_Advect_Outflow, ALEN(outflowids), outflowids, phys);
332: /* Initial/transient solution with default boundary conditions */
333: ModelSolutionSetDefault(mod,PhysicsSolution_Advect,phys);
334: /* Register "canned" functionals */
335: ModelFunctionalRegister(mod,"Error" ,&advect->functional.Error,PhysicsFunctional_Advect,phys);
336: }
337: return (0);
338: }
340: /******************* Shallow Water ********************/
341: typedef struct {
342: PetscReal gravity;
343: PetscReal boundaryHeight;
344: struct {
345: PetscInt Height;
346: PetscInt Speed;
347: PetscInt Energy;
348: } functional;
349: } Physics_SW;
350: typedef struct {
351: PetscScalar vals[0];
352: PetscScalar h;
353: PetscScalar uh[DIM];
354: } SWNode;
356: static const struct FieldDescription PhysicsFields_SW[] = {{"Height" ,1},{"Momentum" ,DIM},{NULL,0}};
360: /*
361: * h_t + div(uh) = 0
362: * (uh)_t + div (u\otimes uh + g h^2 / 2 I) = 0
363: *
364: * */
365: static PetscErrorCode SWFlux(Physics phys,const PetscReal *n,const SWNode *x,SWNode *f)
366: {
367: Physics_SW *sw = (Physics_SW*)phys->data;
368: PetscScalar uhn,u[DIM];
369: PetscInt i;
372: Scale2(1./x->h,x->uh,u);
373: uhn = Dot2(x->uh,n);
374: f->h = uhn;
375: for (i=0; i<DIM; i++) f->uh[i] = u[i] * uhn + sw->gravity * PetscSqr(x->h) * n[i];
376: return (0);
377: }
381: static PetscErrorCode PhysicsBoundary_SW_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
382: {
384: xG[0] = xI[0];
385: xG[1] = -xI[1];
386: xG[2] = -xI[2];
387: return (0);
388: }
392: static PetscErrorCode PhysicsRiemann_SW(const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscScalar *flux, Physics phys)
393: {
394: Physics_SW *sw = (Physics_SW*)phys->data;
395: PetscReal cL,cR,speed,nn[DIM];
396: const SWNode *uL = (const SWNode*)xL,*uR = (const SWNode*)xR;
397: SWNode fL,fR;
398: PetscInt i;
401: if (uL->h < 0 || uR->h < 0) SETERRQ (PETSC_COMM_SELF ,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative" );
402: nn[0] = n[0];
403: nn[1] = n[1];
404: Normalize2(nn);
405: SWFlux(phys,nn,uL,&fL);
406: SWFlux(phys,nn,uR,&fR);
407: cL = PetscSqrtReal(sw->gravity*PetscRealPart(uL->h));
408: cR = PetscSqrtReal(sw->gravity*PetscRealPart(uR->h)); /* gravity wave speed */
409: speed = PetscMax(PetscAbsScalar(Dot2(uL->uh,nn)/uL->h) + cL,PetscAbsScalar(Dot2(uR->uh,nn)/uR->h) + cR);
410: 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);
411: return (0);
412: }
416: static PetscErrorCode PhysicsSolution_SW(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx)
417: {
418: PetscReal dx[2],r,sigma;
421: if (time != 0.0) SETERRQ1 (mod->comm,PETSC_ERR_SUP,"No solution known for time %g" ,(double)time);
422: dx[0] = x[0] - 1.5;
423: dx[1] = x[1] - 1.0;
424: r = Norm2(dx);
425: sigma = 0.5;
426: u[0] = 1 + 2*PetscExpScalar(-PetscSqr(r)/(2*PetscSqr(sigma)));
427: u[1] = 0.0;
428: u[2] = 0.0;
429: return (0);
430: }
434: static PetscErrorCode PhysicsFunctional_SW(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
435: {
436: Physics phys = (Physics)ctx;
437: Physics_SW *sw = (Physics_SW*)phys->data;
438: const SWNode *x = (const SWNode*)xx;
439: PetscScalar u[2];
440: PetscReal h;
443: h = PetscRealPart(x->h);
444: Scale2(1./x->h,x->uh,u);
445: f[sw->functional.Height] = h;
446: f[sw->functional.Speed] = Norm2(u) + PetscSqrtReal(sw->gravity*h);
447: f[sw->functional.Energy] = 0.5*(Dot2(x->uh,u) + sw->gravity*PetscSqr(h));
448: return (0);
449: }
453: static PetscErrorCode PhysicsCreate_SW(DM dm, Model mod,Physics phys)
454: {
455: Physics_SW *sw;
459: phys->field_desc = PhysicsFields_SW;
460: phys->riemann = (RiemannFunction) PhysicsRiemann_SW;
461: PetscNew (&sw);
462: phys->data = sw;
463: PetscOptionsHead ("SW options" );
464: {
465: sw->gravity = 1.0;
466: PetscOptionsReal ("-sw_gravity" ,"Gravitational constant" ,"" ,sw->gravity,&sw->gravity,NULL);
467: }
468: PetscOptionsTail ();
469: phys->maxspeed = PetscSqrtReal(2.0*sw->gravity); /* Mach 1 for depth of 2 */
471: {
472: const PetscInt wallids[] = {100,101,200,300};
473: DMPlexAddBoundary(dm, PETSC_TRUE , "wall" , "Face Sets" , 0, (void (*)()) PhysicsBoundary_SW_Wall, ALEN(wallids), wallids, phys);
474: ModelSolutionSetDefault(mod,PhysicsSolution_SW,phys);
475: ModelFunctionalRegister(mod,"Height" ,&sw->functional.Height,PhysicsFunctional_SW,phys);
476: ModelFunctionalRegister(mod,"Speed" ,&sw->functional.Speed,PhysicsFunctional_SW,phys);
477: ModelFunctionalRegister(mod,"Energy" ,&sw->functional.Energy,PhysicsFunctional_SW,phys);
478: }
479: return (0);
480: }
482: /******************* Euler ********************/
483: typedef struct {
484: PetscScalar vals[0];
485: PetscScalar r;
486: PetscScalar ru[DIM];
487: PetscScalar e;
488: } EulerNode;
489: typedef PetscErrorCode (*EquationOfState)(const PetscReal *, const EulerNode*, PetscScalar *) ;
490: typedef struct {
491: PetscInt npars;
492: PetscReal pars[DIM];
493: EquationOfState pressure;
494: EquationOfState sound;
495: struct {
496: PetscInt Density;
497: PetscInt Momentum;
498: PetscInt Energy;
499: PetscInt Pressure;
500: PetscInt Speed;
501: } monitor;
502: } Physics_Euler;
504: static const struct FieldDescription PhysicsFields_Euler[] = {{"Density" ,1},{"Momentum" ,DIM},{"Energy" ,1},{NULL,0}};
508: static PetscErrorCode Pressure_PG(const PetscReal *pars,const EulerNode *x,PetscScalar *p)
509: {
510: PetscScalar ru2;
513: ru2 = DotDIM(x->ru,x->ru);
514: ru2 /= x->r;
515: /* kinematic dof = params[0] */
516: (*p)=2.0*(x->e-0.5*ru2)/pars[0];
517: return (0);
518: }
522: static PetscErrorCode SpeedOfSound_PG(const PetscReal *pars,const EulerNode *x,PetscScalar *c)
523: {
524: PetscScalar p;
527: /* TODO remove direct usage of Pressure_PG */
528: Pressure_PG(pars,x,&p);
529: /* TODO check the sign of p */
530: /* pars[1] = heat capacity ratio */
531: (*c)=PetscSqrtScalar(pars[1]*p/x->r);
532: return (0);
533: }
537: /*
538: * x = (rho,rho*(u_1),...,rho*e)^T
539: * x_t+div(f_1(x))+...+div(f_DIM(x)) = 0
540: *
541: * f_i(x) = u_i*x+(0,0,...,p,...,p*u_i)^T
542: *
543: * */
544: static PetscErrorCode EulerFlux(Physics phys,const PetscReal *n,const EulerNode *x,EulerNode *f)
545: {
546: Physics_Euler *eu = (Physics_Euler*)phys->data;
547: PetscScalar u,nu,p;
548: PetscInt i;
551: u = DotDIM(x->ru,x->ru);
552: u /= (x->r * x->r);
553: nu = DotDIM(x->ru,n);
554: /* TODO check the sign of p */
555: eu->pressure(eu->pars,x,&p);
556: f->r = nu * x->r;
557: for (i=0; i<DIM; i++) f->ru[i] = nu * x->ru[i] + n[i]*p;
558: f->e = nu*(x->e+p);
559: return (0);
560: }
562: /* PetscReal * => EulerNode* conversion */
565: static PetscErrorCode PhysicsBoundary_Euler_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
566: {
567: PetscInt i;
568: PetscScalar xn[DIM],xt[DIM];
571: xG[0] = xI[0];
572: NormalSplitDIM(n,xI+1,xn,xt);
573: for (i=0; i<DIM; i++) xG[i+1] = -xn[i]+xt[i];
574: xG[DIM+1] = xI[DIM+1];
575: return (0);
576: }
578: /* PetscReal * => EulerNode* conversion */
581: static PetscErrorCode PhysicsRiemann_Euler_Rusanov(const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscScalar *flux, Physics phys)
582: {
583: Physics_Euler *eu = (Physics_Euler*)phys->data;
584: PetscScalar cL,cR,speed;
585: const EulerNode *uL = (const EulerNode*)xL,*uR = (const EulerNode*)xR;
586: EulerNode fL,fR;
587: PetscInt i;
590: if (uL->r < 0 || uR->r < 0) SETERRQ (PETSC_COMM_SELF ,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed density is negative" );
591: EulerFlux(phys,n,uL,&fL);
592: EulerFlux(phys,n,uR,&fR);
593: eu->sound(eu->pars,uL,&cL);
594: eu->sound(eu->pars,uR,&cR);
595: speed = PetscMax(cL,cR)+PetscMax(PetscAbsScalar(DotDIM(uL->ru,n)/NormDIM(n)),PetscAbsScalar(DotDIM(uR->ru,n)/NormDIM(n)));
596: for (i=0; i<2+DIM; i++) flux[i] = 0.5*(fL.vals[i]+fR.vals[i])+0.5*speed*(xL[i]-xR[i]);
597: return (0);
598: }
602: static PetscErrorCode PhysicsSolution_Euler(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx)
603: {
604: PetscInt i;
607: if (time != 0.0) SETERRQ1 (mod->comm,PETSC_ERR_SUP,"No solution known for time %g" ,(double)time);
608: u[0] = 1.0;
609: u[DIM+1] = 1.0+PetscAbsReal(x[0]);
610: for (i=1; i<DIM+1; i++) u[i] = 0.0;
611: return (0);
612: }
616: static PetscErrorCode PhysicsFunctional_Euler(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
617: {
618: Physics phys = (Physics)ctx;
619: Physics_Euler *eu = (Physics_Euler*)phys->data;
620: const EulerNode *x = (const EulerNode*)xx;
621: PetscScalar p;
624: f[eu->monitor.Density] = x->r;
625: f[eu->monitor.Momentum] = NormDIM(x->ru);
626: f[eu->monitor.Energy] = x->e;
627: f[eu->monitor.Speed] = NormDIM(x->ru)/x->r;
628: eu->pressure(eu->pars, x, &p);
629: f[eu->monitor.Pressure] = p;
630: return (0);
631: }
635: static PetscErrorCode PhysicsCreate_Euler(DM dm, Model mod,Physics phys)
636: {
637: Physics_Euler *eu;
638: PetscErrorCode ierr;
641: phys->field_desc = PhysicsFields_Euler;
642: phys->riemann = (RiemannFunction) PhysicsRiemann_Euler_Rusanov;
643: PetscNew (&eu);
644: phys->data = eu;
645: PetscOptionsHead ("Euler options" );
646: {
647: eu->pars[0] = 3.0;
648: eu->pars[1] = 1.67;
649: PetscOptionsReal ("-eu_f" ,"Degrees of freedom" ,"" ,eu->pars[0],&eu->pars[0],NULL);
650: PetscOptionsReal ("-eu_gamma" ,"Heat capacity ratio" ,"" ,eu->pars[1],&eu->pars[1],NULL);
651: }
652: PetscOptionsTail ();
653: eu->pressure = Pressure_PG;
654: eu->sound = SpeedOfSound_PG;
655: phys->maxspeed = 1.0;
656: {
657: const PetscInt wallids[] = {100,101,200,300};
658: DMPlexAddBoundary(dm, PETSC_TRUE , "wall" , "Face Sets" , 0, (void (*)()) PhysicsBoundary_Euler_Wall, ALEN(wallids), wallids, phys);
659: ModelSolutionSetDefault(mod,PhysicsSolution_Euler,phys);
660: ModelFunctionalRegister(mod,"Speed" ,&eu->monitor.Speed,PhysicsFunctional_Euler,phys);
661: ModelFunctionalRegister(mod,"Energy" ,&eu->monitor.Energy,PhysicsFunctional_Euler,phys);
662: ModelFunctionalRegister(mod,"Density" ,&eu->monitor.Density,PhysicsFunctional_Euler,phys);
663: ModelFunctionalRegister(mod,"Momentum" ,&eu->monitor.Momentum,PhysicsFunctional_Euler,phys);
664: ModelFunctionalRegister(mod,"Pressure" ,&eu->monitor.Pressure,PhysicsFunctional_Euler,phys);
665: }
666: return (0);
667: }
671: PetscErrorCode ConstructCellBoundary(DM dm, User user)
672: {
673: const char *name = "Cell Sets" ;
674: const char *bdname = "split faces" ;
675: IS regionIS, innerIS;
676: const PetscInt *regions, *cells;
677: PetscInt numRegions, innerRegion, numCells, c;
678: PetscInt cStart, cEnd, cEndInterior, fStart, fEnd;
679: PetscBool hasLabel;
683: DMPlexGetHeightStratum (dm, 0, &cStart, &cEnd);
684: DMPlexGetHeightStratum (dm, 1, &fStart, &fEnd);
685: DMPlexGetHybridBounds (dm, &cEndInterior, NULL, NULL, NULL);
687: DMPlexHasLabel (dm, name, &hasLabel);
688: if (!hasLabel) return (0);
689: DMPlexGetLabelSize (dm, name, &numRegions);
690: if (numRegions != 2) return (0);
691: /* Get the inner id */
692: DMPlexGetLabelIdIS (dm, name, ®ionIS);
693: ISGetIndices (regionIS, ®ions);
694: innerRegion = regions[0];
695: ISRestoreIndices (regionIS, ®ions);
696: ISDestroy (®ionIS);
697: /* Find the faces between cells in different regions, could call DMPlexCreateNeighborCSR() */
698: DMPlexGetStratumIS (dm, name, innerRegion, &innerIS);
699: ISGetLocalSize (innerIS, &numCells);
700: ISGetIndices (innerIS, &cells);
701: DMPlexCreateLabel (dm, bdname);
702: for (c = 0; c < numCells; ++c) {
703: const PetscInt cell = cells[c];
704: const PetscInt *faces;
705: PetscInt numFaces, f;
707: if ((cell < cStart) || (cell >= cEnd)) SETERRQ1 (PETSC_COMM_SELF , PETSC_ERR_LIB, "Got invalid point %d which is not a cell" , cell);
708: DMPlexGetConeSize (dm, cell, &numFaces);
709: DMPlexGetCone (dm, cell, &faces);
710: for (f = 0; f < numFaces; ++f) {
711: const PetscInt face = faces[f];
712: const PetscInt *neighbors;
713: PetscInt nC, regionA, regionB;
715: if ((face < fStart) || (face >= fEnd)) SETERRQ1 (PETSC_COMM_SELF , PETSC_ERR_LIB, "Got invalid point %d which is not a face" , face);
716: DMPlexGetSupportSize (dm, face, &nC);
717: if (nC != 2) continue ;
718: DMPlexGetSupport (dm, face, &neighbors);
719: if ((neighbors[0] >= cEndInterior) || (neighbors[1] >= cEndInterior)) continue ;
720: 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]);
721: 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]);
722: DMPlexGetLabelValue (dm, name, neighbors[0], ®ionA);
723: DMPlexGetLabelValue (dm, name, neighbors[1], ®ionB);
724: if (regionA < 0) SETERRQ2 (PetscObjectComm ((PetscObject )dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value" , name, neighbors[0]);
725: if (regionB < 0) SETERRQ2 (PetscObjectComm ((PetscObject )dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value" , name, neighbors[1]);
726: if (regionA != regionB) {
727: DMPlexSetLabelValue (dm, bdname, faces[f], 1);
728: }
729: }
730: }
731: ISRestoreIndices (innerIS, &cells);
732: ISDestroy (&innerIS);
733: {
734: DMLabel label;
736: PetscViewerASCIISynchronizedAllow (PETSC_VIEWER_STDOUT_WORLD , PETSC_TRUE );
737: DMPlexGetLabel (dm, bdname, &label);
738: DMLabelView(label, PETSC_VIEWER_STDOUT_WORLD );
739: }
740: return (0);
741: }
745: /* Right now, I have just added duplicate faces, which see both cells. We can
746: - Add duplicate vertices and decouple the face cones
747: - Disconnect faces from cells across the rotation gap
748: */
749: PetscErrorCode SplitFaces(DM *dmSplit, const char labelName[], User user)
750: {
751: DM dm = *dmSplit, sdm;
752: PetscSF sfPoint, gsfPoint;
753: PetscSection coordSection, newCoordSection;
754: Vec coordinates;
755: IS idIS;
756: const PetscInt *ids;
757: PetscInt *newpoints;
758: PetscInt dim, depth, maxConeSize, maxSupportSize, numLabels, numGhostCells;
759: PetscInt numFS, fs, pStart, pEnd, p, cEnd, cEndInterior, vStart, vEnd, v, fStart, fEnd, newf, d, l;
760: PetscBool hasLabel;
764: DMPlexHasLabel (dm, labelName, &hasLabel);
765: if (!hasLabel) return (0);
766: DMCreate (PetscObjectComm ((PetscObject )dm), &sdm);
767: DMSetType (sdm, DMPLEX );
768: DMPlexGetDimension (dm, &dim);
769: DMPlexSetDimension (sdm, dim);
771: DMPlexGetLabelIdIS (dm, labelName, &idIS);
772: ISGetLocalSize (idIS, &numFS);
773: ISGetIndices (idIS, &ids);
775: user->numSplitFaces = 0;
776: for (fs = 0; fs < numFS; ++fs) {
777: PetscInt numBdFaces;
779: DMPlexGetStratumSize (dm, labelName, ids[fs], &numBdFaces);
780: user->numSplitFaces += numBdFaces;
781: }
782: DMPlexGetChart (dm, &pStart, &pEnd);
783: pEnd += user->numSplitFaces;
784: DMPlexSetChart (sdm, pStart, pEnd);
785: DMPlexGetHybridBounds (dm, &cEndInterior, NULL, NULL, NULL);
786: DMPlexSetHybridBounds (sdm, cEndInterior, PETSC_DETERMINE , PETSC_DETERMINE , PETSC_DETERMINE );
787: DMPlexGetHeightStratum (dm, 0, NULL, &cEnd);
788: numGhostCells = cEnd - cEndInterior;
789: /* Set cone and support sizes */
790: DMPlexGetDepth (dm, &depth);
791: for (d = 0; d <= depth; ++d) {
792: DMPlexGetDepthStratum (dm, d, &pStart, &pEnd);
793: for (p = pStart; p < pEnd; ++p) {
794: PetscInt newp = p;
795: PetscInt size;
797: DMPlexGetConeSize (dm, p, &size);
798: DMPlexSetConeSize (sdm, newp, size);
799: DMPlexGetSupportSize (dm, p, &size);
800: DMPlexSetSupportSize (sdm, newp, size);
801: }
802: }
803: DMPlexGetHeightStratum (dm, 1, &fStart, &fEnd);
804: for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
805: IS faceIS;
806: const PetscInt *faces;
807: PetscInt numFaces, f;
809: DMPlexGetStratumIS (dm, labelName, ids[fs], &faceIS);
810: ISGetLocalSize (faceIS, &numFaces);
811: ISGetIndices (faceIS, &faces);
812: for (f = 0; f < numFaces; ++f, ++newf) {
813: PetscInt size;
815: /* Right now I think that both faces should see both cells */
816: DMPlexGetConeSize (dm, faces[f], &size);
817: DMPlexSetConeSize (sdm, newf, size);
818: DMPlexGetSupportSize (dm, faces[f], &size);
819: DMPlexSetSupportSize (sdm, newf, size);
820: }
821: ISRestoreIndices (faceIS, &faces);
822: ISDestroy (&faceIS);
823: }
824: DMSetUp (sdm);
825: /* Set cones and supports */
826: DMPlexGetMaxSizes (dm, &maxConeSize, &maxSupportSize);
827: PetscMalloc (PetscMax(maxConeSize, maxSupportSize) * sizeof (PetscInt ), &newpoints);
828: DMPlexGetChart (dm, &pStart, &pEnd);
829: for (p = pStart; p < pEnd; ++p) {
830: const PetscInt *points, *orientations;
831: PetscInt size, i, newp = p;
833: DMPlexGetConeSize (dm, p, &size);
834: DMPlexGetCone (dm, p, &points);
835: DMPlexGetConeOrientation (dm, p, &orientations);
836: for (i = 0; i < size; ++i) newpoints[i] = points[i];
837: DMPlexSetCone (sdm, newp, newpoints);
838: DMPlexSetConeOrientation (sdm, newp, orientations);
839: DMPlexGetSupportSize (dm, p, &size);
840: DMPlexGetSupport (dm, p, &points);
841: for (i = 0; i < size; ++i) newpoints[i] = points[i];
842: DMPlexSetSupport (sdm, newp, newpoints);
843: }
844: PetscFree (newpoints);
845: for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
846: IS faceIS;
847: const PetscInt *faces;
848: PetscInt numFaces, f;
850: DMPlexGetStratumIS (dm, labelName, ids[fs], &faceIS);
851: ISGetLocalSize (faceIS, &numFaces);
852: ISGetIndices (faceIS, &faces);
853: for (f = 0; f < numFaces; ++f, ++newf) {
854: const PetscInt *points;
856: DMPlexGetCone (dm, faces[f], &points);
857: DMPlexSetCone (sdm, newf, points);
858: DMPlexGetSupport (dm, faces[f], &points);
859: DMPlexSetSupport (sdm, newf, points);
860: }
861: ISRestoreIndices (faceIS, &faces);
862: ISDestroy (&faceIS);
863: }
864: ISRestoreIndices (idIS, &ids);
865: ISDestroy (&idIS);
866: DMPlexStratify (sdm);
867: /* Convert coordinates */
868: DMPlexGetDepthStratum (dm, 0, &vStart, &vEnd);
869: DMGetCoordinateSection (dm, &coordSection);
870: PetscSectionCreate (PetscObjectComm ((PetscObject )dm), &newCoordSection);
871: PetscSectionSetNumFields (newCoordSection, 1);
872: PetscSectionSetFieldComponents (newCoordSection, 0, dim);
873: PetscSectionSetChart (newCoordSection, vStart, vEnd);
874: for (v = vStart; v < vEnd; ++v) {
875: PetscSectionSetDof (newCoordSection, v, dim);
876: PetscSectionSetFieldDof (newCoordSection, v, 0, dim);
877: }
878: PetscSectionSetUp (newCoordSection);
879: DMSetCoordinateSection (sdm, newCoordSection);
880: PetscSectionDestroy (&newCoordSection); /* relinquish our reference */
881: DMGetCoordinatesLocal (dm, &coordinates);
882: DMSetCoordinatesLocal (sdm, coordinates);
883: /* Convert labels */
884: DMPlexGetNumLabels (dm, &numLabels);
885: for (l = 0; l < numLabels; ++l) {
886: const char *lname;
887: PetscBool isDepth;
889: DMPlexGetLabelName (dm, l, &lname);
890: PetscStrcmp (lname, "depth" , &isDepth);
891: if (isDepth) continue ;
892: DMPlexCreateLabel (sdm, lname);
893: DMPlexGetLabelIdIS (dm, lname, &idIS);
894: ISGetLocalSize (idIS, &numFS);
895: ISGetIndices (idIS, &ids);
896: for (fs = 0; fs < numFS; ++fs) {
897: IS pointIS;
898: const PetscInt *points;
899: PetscInt numPoints;
901: DMPlexGetStratumIS (dm, lname, ids[fs], &pointIS);
902: ISGetLocalSize (pointIS, &numPoints);
903: ISGetIndices (pointIS, &points);
904: for (p = 0; p < numPoints; ++p) {
905: PetscInt newpoint = points[p];
907: DMPlexSetLabelValue (sdm, lname, newpoint, ids[fs]);
908: }
909: ISRestoreIndices (pointIS, &points);
910: ISDestroy (&pointIS);
911: }
912: ISRestoreIndices (idIS, &ids);
913: ISDestroy (&idIS);
914: }
915: /* Convert pointSF */
916: const PetscSFNode *remotePoints;
917: PetscSFNode *gremotePoints;
918: const PetscInt *localPoints;
919: PetscInt *glocalPoints,*newLocation,*newRemoteLocation;
920: PetscInt numRoots, numLeaves;
921: PetscMPIInt numProcs;
923: MPI_Comm_size (PetscObjectComm ((PetscObject )dm), &numProcs);
924: DMGetPointSF (dm, &sfPoint);
925: DMGetPointSF (sdm, &gsfPoint);
926: DMPlexGetChart (dm,&pStart,&pEnd);
927: PetscSFGetGraph (sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);
928: if (numRoots >= 0) {
929: PetscMalloc2 (numRoots,&newLocation,pEnd-pStart,&newRemoteLocation);
930: for (l=0; l<numRoots; l++) newLocation[l] = l; /* + (l >= cEnd ? numGhostCells : 0); */
931: PetscSFBcastBegin (sfPoint, MPIU_INT, newLocation, newRemoteLocation);
932: PetscSFBcastEnd (sfPoint, MPIU_INT, newLocation, newRemoteLocation);
933: PetscMalloc1 (numLeaves, &glocalPoints);
934: PetscMalloc1 (numLeaves, &gremotePoints);
935: for (l = 0; l < numLeaves; ++l) {
936: glocalPoints[l] = localPoints[l]; /* localPoints[l] >= cEnd ? localPoints[l] + numGhostCells : localPoints[l]; */
937: gremotePoints[l].rank = remotePoints[l].rank;
938: gremotePoints[l].index = newRemoteLocation[localPoints[l]];
939: }
940: PetscFree2 (newLocation,newRemoteLocation);
941: PetscSFSetGraph (gsfPoint, numRoots+numGhostCells, numLeaves, glocalPoints, PETSC_OWN_POINTER, gremotePoints, PETSC_OWN_POINTER);
942: }
943: DMDestroy (dmSplit);
944: *dmSplit = sdm;
945: return (0);
946: }
950: PetscErrorCode CreatePartitionVec(DM dm, DM *dmCell, Vec *partition)
951: {
952: PetscSF sfPoint;
953: PetscSection coordSection;
954: Vec coordinates;
955: PetscSection sectionCell;
956: PetscScalar *part;
957: PetscInt cStart, cEnd, c;
958: PetscMPIInt rank;
962: DMGetCoordinateSection (dm, &coordSection);
963: DMGetCoordinatesLocal (dm, &coordinates);
964: DMClone (dm, dmCell);
965: DMGetPointSF (dm, &sfPoint);
966: DMSetPointSF (*dmCell, sfPoint);
967: DMSetCoordinateSection (*dmCell, coordSection);
968: DMSetCoordinatesLocal (*dmCell, coordinates);
969: MPI_Comm_rank (PetscObjectComm ((PetscObject )dm), &rank);
970: PetscSectionCreate (PetscObjectComm ((PetscObject )dm), §ionCell);
971: DMPlexGetHeightStratum (*dmCell, 0, &cStart, &cEnd);
972: PetscSectionSetChart (sectionCell, cStart, cEnd);
973: for (c = cStart; c < cEnd; ++c) {
974: PetscSectionSetDof (sectionCell, c, 1);
975: }
976: PetscSectionSetUp (sectionCell);
977: DMSetDefaultSection (*dmCell, sectionCell);
978: PetscSectionDestroy (§ionCell);
979: DMCreateLocalVector (*dmCell, partition);
980: PetscObjectSetName ((PetscObject )*partition, "partition" );
981: VecGetArray (*partition, &part);
982: for (c = cStart; c < cEnd; ++c) {
983: PetscScalar *p;
985: DMPlexPointLocalRef (*dmCell, c, part, &p);
986: p[0] = rank;
987: }
988: VecRestoreArray (*partition, &part);
989: return (0);
990: }
994: PetscErrorCode CreateMassMatrix(DM dm, Vec *massMatrix, User user)
995: {
996: DM dmMass, dmFace, dmCell, dmCoord;
997: PetscSection coordSection;
998: Vec coordinates, facegeom, cellgeom;
999: PetscSection sectionMass;
1000: PetscScalar *m;
1001: const PetscScalar *fgeom, *cgeom, *coords;
1002: PetscInt vStart, vEnd, v;
1003: PetscErrorCode ierr;
1006: DMGetCoordinateSection (dm, &coordSection);
1007: DMGetCoordinatesLocal (dm, &coordinates);
1008: DMClone (dm, &dmMass);
1009: DMSetCoordinateSection (dmMass, coordSection);
1010: DMSetCoordinatesLocal (dmMass, coordinates);
1011: PetscSectionCreate (PetscObjectComm ((PetscObject )dm), §ionMass);
1012: DMPlexGetDepthStratum (dm, 0, &vStart, &vEnd);
1013: PetscSectionSetChart (sectionMass, vStart, vEnd);
1014: for (v = vStart; v < vEnd; ++v) {
1015: PetscInt numFaces;
1017: DMPlexGetSupportSize (dmMass, v, &numFaces);
1018: PetscSectionSetDof (sectionMass, v, numFaces*numFaces);
1019: }
1020: PetscSectionSetUp (sectionMass);
1021: DMSetDefaultSection (dmMass, sectionMass);
1022: PetscSectionDestroy (§ionMass);
1023: DMGetLocalVector (dmMass, massMatrix);
1024: VecGetArray (*massMatrix, &m);
1025: DMPlexTSGetGeometry (dm, &facegeom, &cellgeom, NULL);
1026: VecGetDM (facegeom, &dmFace);
1027: VecGetArrayRead (facegeom, &fgeom);
1028: VecGetDM (cellgeom, &dmCell);
1029: VecGetArrayRead (cellgeom, &cgeom);
1030: DMGetCoordinateDM (dm, &dmCoord);
1031: VecGetArrayRead (coordinates, &coords);
1032: for (v = vStart; v < vEnd; ++v) {
1033: const PetscInt *faces;
1034: const FaceGeom *fgA, *fgB, *cg;
1035: const PetscScalar *vertex;
1036: PetscInt numFaces, sides[2], f, g;
1038: DMPlexPointLocalRead (dmCoord, v, coords, &vertex);
1039: DMPlexGetSupportSize (dmMass, v, &numFaces);
1040: DMPlexGetSupport (dmMass, v, &faces);
1041: for (f = 0; f < numFaces; ++f) {
1042: sides[0] = faces[f];
1043: DMPlexPointLocalRead (dmFace, faces[f], fgeom, &fgA);
1044: for (g = 0; g < numFaces; ++g) {
1045: const PetscInt *cells = NULL;;
1046: PetscReal area = 0.0;
1047: PetscInt numCells;
1049: sides[1] = faces[g];
1050: DMPlexPointLocalRead (dmFace, faces[g], fgeom, &fgB);
1051: DMPlexGetJoin (dmMass, 2, sides, &numCells, &cells);
1052: if (numCells != 1) SETERRQ (PETSC_COMM_SELF , PETSC_ERR_LIB, "Invalid join for faces" );
1053: DMPlexPointLocalRead (dmCell, cells[0], cgeom, &cg);
1054: area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgA->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgA->centroid[0] - cg->centroid[0]));
1055: area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgB->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgB->centroid[0] - cg->centroid[0]));
1056: m[f*numFaces+g] = Dot2(fgA->normal, fgB->normal)*area*0.5;
1057: DMPlexRestoreJoin (dmMass, 2, sides, &numCells, &cells);
1058: }
1059: }
1060: }
1061: VecRestoreArrayRead (facegeom, &fgeom);
1062: VecRestoreArrayRead (cellgeom, &cgeom);
1063: VecRestoreArrayRead (coordinates, &coords);
1064: VecRestoreArray (*massMatrix, &m);
1065: DMDestroy (&dmMass);
1066: return (0);
1067: }
1071: PetscErrorCode SetUpLocalSpace(DM dm, User user)
1072: {
1073: PetscSection stateSection;
1074: Physics phys;
1075: PetscInt dof = user->model->physics->dof, *cind, d, stateSize, cStart, cEnd, cEndInterior, c, i;
1079: DMPlexGetHeightStratum (dm, 0, &cStart, &cEnd);
1080: DMPlexGetHybridBounds (dm, &cEndInterior, NULL, NULL, NULL);
1081: PetscSectionCreate (PetscObjectComm ((PetscObject )dm), &stateSection);
1082: phys = user->model->physics;
1083: PetscSectionSetNumFields (stateSection,phys->nfields);
1084: for (i=0; i<phys->nfields; i++) {
1085: PetscSectionSetFieldName (stateSection,i,phys->field_desc[i].name);
1086: PetscSectionSetFieldComponents (stateSection,i,phys->field_desc[i].dof);
1087: }
1088: PetscSectionSetChart (stateSection, cStart, cEnd);
1089: for (c = cStart; c < cEnd; ++c) {
1090: for (i=0; i<phys->nfields; i++) {
1091: PetscSectionSetFieldDof (stateSection,c,i,phys->field_desc[i].dof);
1092: }
1093: PetscSectionSetDof (stateSection, c, dof);
1094: }
1095: for (c = cEndInterior; c < cEnd; ++c) {
1096: PetscSectionSetConstraintDof (stateSection, c, dof);
1097: }
1098: PetscSectionSetUp (stateSection);
1099: PetscMalloc1 (dof, &cind);
1100: for (d = 0; d < dof; ++d) cind[d] = d;
1101: #if 0
1102: for (c = cStart; c < cEnd; ++c) {
1103: PetscInt val;
1105: DMPlexGetLabelValue (dm, "vtk" , c, &val);
1106: if (val < 0) {PetscSectionSetConstraintIndices (stateSection, c, cind);}
1107: }
1108: #endif
1109: for (c = cEndInterior; c < cEnd; ++c) {
1110: PetscSectionSetConstraintIndices (stateSection, c, cind);
1111: }
1112: PetscFree (cind);
1113: PetscSectionGetStorageSize (stateSection, &stateSize);
1114: DMSetDefaultSection (dm,stateSection);
1115: PetscSectionDestroy (&stateSection);
1116: return (0);
1117: }
1119: #if 0
1122: PetscErrorCode SetUpBoundaries(DM dm, User user)
1123: {
1124: Model mod = user->model;
1126: BoundaryLink b;
1129: PetscOptionsBegin (PetscObjectComm ((PetscObject )dm),NULL,"Boundary condition options" ,"" );
1130: for (b = mod->boundary; b; b=b->next) {
1131: char optname[512];
1132: PetscInt ids[512],len = 512;
1133: PetscBool flg;
1134: PetscSNPrintf (optname,sizeof optname,"-bc_%s" ,b->name);
1135: PetscMemzero (ids,sizeof (ids));
1136: PetscOptionsIntArray (optname,"List of boundary IDs" ,"" ,ids,&len,&flg);
1137: if (flg) {
1138: /* TODO: check all IDs to make sure they exist in the mesh */
1139: PetscFree (b->ids);
1140: b->numids = len;
1141: PetscMalloc1 (len,&b->ids);
1142: PetscMemcpy (b->ids,ids,len*sizeof (PetscInt ));
1143: }
1144: }
1145: PetscOptionsEnd ();
1146: return (0);
1147: }
1148: #endif
1152: /* Behavior will be different for multi-physics or when using non-default boundary conditions */
1153: static PetscErrorCode ModelSolutionSetDefault(Model mod,SolutionFunction func,void *ctx)
1154: {
1156: mod->solution = func;
1157: mod->solutionctx = ctx;
1158: return (0);
1159: }
1163: static PetscErrorCode ModelFunctionalRegister(Model mod,const char *name,PetscInt *offset,FunctionalFunction func,void *ctx)
1164: {
1166: FunctionalLink link,*ptr;
1167: PetscInt lastoffset = -1;
1170: for (ptr=&mod->functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset;
1171: PetscNew (&link);
1172: PetscStrallocpy (name,&link->name);
1173: link->offset = lastoffset + 1;
1174: link->func = func;
1175: link->ctx = ctx;
1176: link->next = NULL;
1177: *ptr = link;
1178: *offset = link->offset;
1179: return (0);
1180: }
1184: static PetscErrorCode ModelFunctionalSetFromOptions(Model mod)
1185: {
1187: PetscInt i,j;
1188: FunctionalLink link;
1189: char *names[256];
1192: mod->numMonitored = ALEN(names);
1193: PetscOptionsStringArray ("-monitor" ,"list of functionals to monitor" ,"" ,names,&mod->numMonitored,NULL);
1194: /* Create list of functionals that will be computed somehow */
1195: PetscMalloc1 (mod->numMonitored,&mod->functionalMonitored);
1196: /* Create index of calls that we will have to make to compute these functionals (over-allocation in general). */
1197: PetscMalloc1 (mod->numMonitored,&mod->functionalCall);
1198: mod->numCall = 0;
1199: for (i=0; i<mod->numMonitored; i++) {
1200: for (link=mod->functionalRegistry; link; link=link->next) {
1201: PetscBool match;
1202: PetscStrcasecmp (names[i],link->name,&match);
1203: if (match) break ;
1204: }
1205: if (!link) SETERRQ1 (mod->comm,PETSC_ERR_USER,"No known functional '%s'" ,names[i]);
1206: mod->functionalMonitored[i] = link;
1207: for (j=0; j<i; j++) {
1208: if (mod->functionalCall[j]->func == link->func && mod->functionalCall[j]->ctx == link->ctx) goto next_name;
1209: }
1210: mod->functionalCall[mod->numCall++] = link; /* Just points to the first link using the result. There may be more results. */
1211: next_name:
1212: PetscFree (names[i]);
1213: }
1215: /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */
1216: mod->maxComputed = -1;
1217: for (link=mod->functionalRegistry; link; link=link->next) {
1218: for (i=0; i<mod->numCall; i++) {
1219: FunctionalLink call = mod->functionalCall[i];
1220: if (link->func == call->func && link->ctx == call->ctx) {
1221: mod->maxComputed = PetscMax(mod->maxComputed,link->offset);
1222: }
1223: }
1224: }
1225: return (0);
1226: }
1230: static PetscErrorCode FunctionalLinkDestroy(FunctionalLink *link)
1231: {
1233: FunctionalLink l,next;
1236: if (!link) return (0);
1237: l = *link;
1238: *link = NULL;
1239: for (; l; l=next) {
1240: next = l->next;
1241: PetscFree (l->name);
1242: PetscFree (l);
1243: }
1244: return (0);
1245: }
1249: PetscErrorCode SetInitialCondition(DM dm, Vec X, User user)
1250: {
1251: DM dmCell;
1252: Model mod = user->model;
1253: Vec cellgeom;
1254: const PetscScalar *cgeom;
1255: PetscScalar *x;
1256: PetscInt cStart, cEnd, cEndInterior, c;
1257: PetscErrorCode ierr;
1260: DMPlexGetHybridBounds (dm, &cEndInterior, NULL, NULL, NULL);
1261: DMPlexTSGetGeometry (dm, NULL, &cellgeom, NULL);
1262: VecGetDM (cellgeom, &dmCell);
1263: DMPlexGetHeightStratum (dm, 0, &cStart, &cEnd);
1264: VecGetArrayRead (cellgeom, &cgeom);
1265: VecGetArray (X, &x);
1266: for (c = cStart; c < cEndInterior; ++c) {
1267: const CellGeom *cg;
1268: PetscScalar *xc;
1270: DMPlexPointLocalRead (dmCell,c,cgeom,&cg);
1271: DMPlexPointGlobalRef (dm,c,x,&xc);
1272: if (xc) {(*mod->solution)(mod,0.0,cg->centroid,xc,mod->solutionctx);}
1273: }
1274: VecRestoreArrayRead (cellgeom, &cgeom);
1275: VecRestoreArray (X, &x);
1276: return (0);
1277: }
1281: static PetscErrorCode OutputVTK(DM dm, const char *filename, PetscViewer *viewer)
1282: {
1286: PetscViewerCreate (PetscObjectComm ((PetscObject )dm), viewer);
1287: PetscViewerSetType (*viewer, PETSCVIEWERVTK);
1288: PetscViewerFileSetName (*viewer, filename);
1289: return (0);
1290: }
1294: static PetscErrorCode MonitorVTK(TS ts,PetscInt stepnum,PetscReal time,Vec X,void *ctx)
1295: {
1296: User user = (User)ctx;
1297: DM dm;
1298: Vec cellgeom;
1299: PetscViewer viewer;
1300: char filename[PETSC_MAX_PATH_LEN],*ftable = NULL;
1301: PetscReal xnorm;
1302: PetscInt cEndInterior;
1306: PetscObjectSetName ((PetscObject ) X, "solution" );
1307: VecGetDM (X,&dm);
1308: DMPlexTSGetGeometry (dm, NULL, &cellgeom, NULL);
1309: DMPlexGetHybridBounds (dm, &cEndInterior, NULL, NULL, NULL);
1310: VecNorm (X,NORM_INFINITY ,&xnorm);
1311: if (stepnum >= 0) { /* No summary for final time */
1312: Model mod = user->model;
1313: PetscInt c,cStart,cEnd,fcount,i;
1314: size_t ftableused,ftablealloc;
1315: const PetscScalar *cgeom,*x;
1316: DM dmCell;
1317: PetscReal *fmin,*fmax,*fintegral,*ftmp;
1318: fcount = mod->maxComputed+1;
1319: PetscMalloc4 (fcount,&fmin,fcount,&fmax,fcount,&fintegral,fcount,&ftmp);
1320: for (i=0; i<fcount; i++) {
1321: fmin[i] = PETSC_MAX_REAL;
1322: fmax[i] = PETSC_MIN_REAL;
1323: fintegral[i] = 0;
1324: }
1325: DMPlexGetHeightStratum (dm,0,&cStart,&cEnd);
1326: VecGetDM (cellgeom,&dmCell);
1327: VecGetArrayRead (cellgeom,&cgeom);
1328: VecGetArrayRead (X,&x);
1329: for (c = cStart; c < cEndInterior; ++c) {
1330: const CellGeom *cg;
1331: const PetscScalar *cx;
1332: DMPlexPointLocalRead (dmCell,c,cgeom,&cg);
1333: DMPlexPointGlobalRead (dm,c,x,&cx);
1334: if (!cx) continue ; /* not a global cell */
1335: for (i=0; i<mod->numCall; i++) {
1336: FunctionalLink flink = mod->functionalCall[i];
1337: (*flink->func)(mod,time,cg->centroid,cx,ftmp,flink->ctx);
1338: }
1339: for (i=0; i<fcount; i++) {
1340: fmin[i] = PetscMin(fmin[i],ftmp[i]);
1341: fmax[i] = PetscMax(fmax[i],ftmp[i]);
1342: fintegral[i] += cg->volume * ftmp[i];
1343: }
1344: }
1345: VecRestoreArrayRead (cellgeom,&cgeom);
1346: VecRestoreArrayRead (X,&x);
1347: MPI_Allreduce (MPI_IN_PLACE,fmin,fcount,MPIU_REAL,MPI_MIN,PetscObjectComm ((PetscObject )ts));
1348: MPI_Allreduce (MPI_IN_PLACE,fmax,fcount,MPIU_REAL,MPI_MAX,PetscObjectComm ((PetscObject )ts));
1349: MPI_Allreduce (MPI_IN_PLACE,fintegral,fcount,MPIU_REAL,MPI_SUM,PetscObjectComm ((PetscObject )ts));
1351: ftablealloc = fcount * 100;
1352: ftableused = 0;
1353: PetscMalloc1 (ftablealloc,&ftable);
1354: for (i=0; i<mod->numMonitored; i++) {
1355: size_t countused;
1356: char buffer[256],*p;
1357: FunctionalLink flink = mod->functionalMonitored[i];
1358: PetscInt id = flink->offset;
1359: if (i % 3) {
1360: PetscMemcpy (buffer," " ,2);
1361: p = buffer + 2;
1362: } else if (i) {
1363: char newline[] = "\n" ;
1364: PetscMemcpy (buffer,newline,sizeof newline-1);
1365: p = buffer + sizeof newline - 1;
1366: } else {
1367: p = buffer;
1368: }
1369: 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]);
1370: countused += p - buffer;
1371: if (countused > ftablealloc-ftableused-1) { /* reallocate */
1372: char *ftablenew;
1373: ftablealloc = 2*ftablealloc + countused;
1374: PetscMalloc (ftablealloc,&ftablenew);
1375: PetscMemcpy (ftablenew,ftable,ftableused);
1376: PetscFree (ftable);
1377: ftable = ftablenew;
1378: }
1379: PetscMemcpy (ftable+ftableused,buffer,countused);
1380: ftableused += countused;
1381: ftable[ftableused] = 0;
1382: }
1383: PetscFree4 (fmin,fmax,fintegral,ftmp);
1385: PetscPrintf (PetscObjectComm ((PetscObject )ts),"% 3D time %8.4g |x| %8.4g %s\n" ,stepnum,(double)time,(double)xnorm,ftable ? ftable : "" );
1386: PetscFree (ftable);
1387: }
1388: if (user->vtkInterval < 1) return (0);
1389: if ((stepnum == -1) ^ (stepnum % user->vtkInterval == 0)) {
1390: if (stepnum == -1) { /* Final time is not multiple of normal time interval, write it anyway */
1391: TSGetTimeStepNumber (ts,&stepnum);
1392: }
1393: PetscSNPrintf (filename,sizeof filename,"ex11-%03D.vtu" ,stepnum);
1394: OutputVTK(dm,filename,&viewer);
1395: VecView (X,viewer);
1396: PetscViewerDestroy (&viewer);
1397: }
1398: return (0);
1399: }
1403: int main(int argc, char **argv)
1404: {
1405: MPI_Comm comm;
1406: PetscFV fvm;
1407: User user;
1408: Model mod;
1409: Physics phys;
1410: DM dm;
1411: PetscReal ftime, cfl, dt, minRadius;
1412: PetscInt dim, nsteps;
1413: TS ts;
1414: TSConvergedReason reason;
1415: Vec X;
1416: PetscViewer viewer;
1417: PetscMPIInt rank;
1418: PetscBool vtkCellGeom, splitFaces;
1419: PetscInt overlap, f;
1420: char filename[PETSC_MAX_PATH_LEN] = "sevenside.exo" ;
1421: PetscErrorCode ierr;
1423: PetscInitialize (&argc, &argv, (char*) 0, help);
1424: comm = PETSC_COMM_WORLD ;
1425: MPI_Comm_rank (comm, &rank);
1427: PetscNew (&user);
1428: PetscNew (&user->model);
1429: PetscNew (&user->model->physics);
1430: mod = user->model;
1431: phys = mod->physics;
1432: mod->comm = comm;
1434: /* Register physical models to be available on the command line */
1435: PetscFunctionListAdd (&PhysicsList,"advect" ,PhysicsCreate_Advect);
1436: PetscFunctionListAdd (&PhysicsList,"sw" ,PhysicsCreate_SW);
1437: PetscFunctionListAdd (&PhysicsList,"euler" ,PhysicsCreate_Euler);
1440: PetscOptionsBegin (comm,NULL,"Unstructured Finite Volume Mesh Options" ,"" );
1441: {
1442: cfl = 0.9 * 4; /* default SSPRKS2 with s=5 stages is stable for CFL number s-1 */
1443: PetscOptionsReal ("-ufv_cfl" ,"CFL number per step" ,"" ,cfl,&cfl,NULL);
1444: PetscOptionsString ("-f" ,"Exodus.II filename to read" ,"" ,filename,filename,sizeof (filename),NULL);
1445: splitFaces = PETSC_FALSE ;
1446: PetscOptionsBool ("-ufv_split_faces" ,"Split faces between cell sets" ,"" ,splitFaces,&splitFaces,NULL);
1447: overlap = 1;
1448: PetscOptionsInt ("-ufv_mesh_overlap" ,"Number of cells to overlap partitions" ,"" ,overlap,&overlap,NULL);
1449: user->vtkInterval = 1;
1450: PetscOptionsInt ("-ufv_vtk_interval" ,"VTK output interval (0 to disable)" ,"" ,user->vtkInterval,&user->vtkInterval,NULL);
1451: vtkCellGeom = PETSC_FALSE ;
1452: PetscOptionsBool ("-ufv_vtk_cellgeom" ,"Write cell geometry (for debugging)" ,"" ,vtkCellGeom,&vtkCellGeom,NULL);
1453: }
1454: PetscOptionsEnd ();
1455: DMPlexCreateExodusFromFile(comm, filename, PETSC_TRUE , &dm);
1456: DMViewFromOptions(dm, NULL, "-dm_view" );
1457: DMPlexGetDimension (dm, &dim);
1459: PetscOptionsBegin (comm,NULL,"Unstructured Finite Volume Physics Options" ,"" );
1460: {
1461: PetscErrorCode (*physcreate)(DM,Model,Physics);
1462: char physname[256] = "advect" ;
1464: DMPlexCreateLabel (dm, "Face Sets" );
1465: PetscOptionsFList ("-physics" ,"Physics module to solve" ,"" ,PhysicsList,physname,physname,sizeof physname,NULL);
1466: PetscFunctionListFind (PhysicsList,physname,&physcreate);
1467: PetscMemzero (phys,sizeof (struct _n_Physics ));
1468: (*physcreate)(dm,mod,phys);
1469: mod->maxspeed = phys->maxspeed;
1470: /* Count number of fields and dofs */
1471: for (phys->nfields=0,phys->dof=0; phys->field_desc[phys->nfields].name; phys->nfields++) phys->dof += phys->field_desc[phys->nfields].dof;
1473: if (mod->maxspeed <= 0) SETERRQ1 (comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set maxspeed" ,physname);
1474: if (phys->dof <= 0) SETERRQ1 (comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set dof" ,physname);
1475: ModelFunctionalSetFromOptions(mod);
1476: }
1477: PetscOptionsEnd ();
1478: {
1479: DM dmDist;
1481: DMPlexSetAdjacencyUseCone (dm, PETSC_TRUE );
1482: DMPlexSetAdjacencyUseClosure (dm, PETSC_FALSE );
1483: DMPlexDistribute (dm, "chaco" , overlap, NULL, &dmDist);
1484: if (dmDist) {
1485: DMDestroy (&dm);
1486: dm = dmDist;
1487: }
1488: }
1489: DMSetFromOptions (dm);
1490: {
1491: DM gdm;
1493: DMPlexConstructGhostCells (dm, NULL, NULL, &gdm);
1494: DMDestroy (&dm);
1495: dm = gdm;
1496: DMViewFromOptions(dm, NULL, "-dm_view" );
1497: }
1498: if (splitFaces) {ConstructCellBoundary(dm, user);}
1499: SplitFaces(&dm, "split faces" , user);
1501: PetscFVCreate (comm, &fvm);
1502: PetscFVSetFromOptions (fvm);
1503: DMSetNumFields(dm, phys->nfields);
1504: for (f = 0; f < phys->nfields; ++f) {DMSetField (dm, f, (PetscObject ) fvm);}
1505: PetscFVSetNumComponents (fvm, phys->dof);
1506: PetscFVSetSpatialDimension (fvm, dim);
1508: /* Set up DM with section describing local vector and configure local vector. */
1509: SetUpLocalSpace(dm, user);
1511: TSCreate (comm, &ts);
1512: TSSetType (ts, TSSSP );
1513: TSSetDM (ts, dm);
1514: TSMonitorSet (ts,MonitorVTK,user,NULL);
1515: DMPlexTSSetRHSFunctionLocal (dm, user->model->physics->riemann, user->model->physics);
1517: DMCreateGlobalVector (dm, &X);
1518: PetscObjectSetName ((PetscObject ) X, "solution" );
1519: SetInitialCondition(dm, X, user);
1520: if (vtkCellGeom) {
1521: DM dmCell;
1522: Vec cellgeom, partition;
1524: DMPlexTSGetGeometry (dm, NULL, &cellgeom, NULL);
1525: OutputVTK(dm, "ex11-cellgeom.vtk" , &viewer);
1526: VecView (cellgeom, viewer);
1527: PetscViewerDestroy (&viewer);
1528: CreatePartitionVec(dm, &dmCell, &partition);
1529: OutputVTK(dmCell, "ex11-partition.vtk" , &viewer);
1530: VecView (partition, viewer);
1531: PetscViewerDestroy (&viewer);
1532: VecDestroy (&partition);
1533: DMDestroy (&dmCell);
1534: }
1536: DMPlexTSGetGeometry (dm, NULL, NULL, &minRadius);
1537: TSSetDuration (ts,1000,2.0);
1538: dt = cfl * minRadius / user->model->maxspeed;
1539: TSSetInitialTimeStep (ts,0.0,dt);
1540: TSSetFromOptions (ts);
1541: TSSolve (ts,X);
1542: TSGetSolveTime (ts,&ftime);
1543: TSGetTimeStepNumber (ts,&nsteps);
1544: TSGetConvergedReason (ts,&reason);
1545: PetscPrintf (PETSC_COMM_WORLD ,"%s at time %g after %D steps\n" ,TSConvergedReasons[reason],(double)ftime,nsteps);
1546: TSDestroy (&ts);
1548: PetscFunctionListDestroy (&PhysicsList);
1549: FunctionalLinkDestroy(&user->model->functionalRegistry);
1550: PetscFree (user->model->functionalMonitored);
1551: PetscFree (user->model->functionalCall);
1552: PetscFree (user->model->physics->data);
1553: PetscFree (user->model->physics);
1554: PetscFree (user->model);
1555: PetscFree (user);
1556: VecDestroy (&X);
1557: PetscFVDestroy (&fvm);
1558: DMDestroy (&dm);
1559: PetscFinalize ();
1560: return (0);
1561: }