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