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