Branch data Line data Source code
1 : : // *****************************************************************************
2 : : /*!
3 : : \file src/PDE/MultiSpecies/DGMultiSpecies.hpp
4 : : \copyright 2012-2015 J. Bakosi,
5 : : 2016-2018 Los Alamos National Security, LLC.,
6 : : 2019-2021 Triad National Security, LLC.
7 : : All rights reserved. See the LICENSE file for details.
8 : : \brief Compressible multi-species flow using discontinuous Galerkin
9 : : finite elements
10 : : \details This file implements calls to the physics operators governing
11 : : compressible multi-species flow using discontinuous Galerkin discretization.
12 : : */
13 : : // *****************************************************************************
14 : : #ifndef DGMultiSpecies_h
15 : : #define DGMultiSpecies_h
16 : :
17 : : #include <cmath>
18 : : #include <algorithm>
19 : : #include <unordered_set>
20 : : #include <map>
21 : : #include <array>
22 : :
23 : : #include "Macro.hpp"
24 : : #include "Exception.hpp"
25 : : #include "Vector.hpp"
26 : : #include "ContainerUtil.hpp"
27 : : #include "UnsMesh.hpp"
28 : : #include "Inciter/InputDeck/InputDeck.hpp"
29 : : #include "Integrate/Basis.hpp"
30 : : #include "Integrate/Quadrature.hpp"
31 : : #include "Integrate/Initialize.hpp"
32 : : #include "Integrate/Mass.hpp"
33 : : #include "Integrate/Surface.hpp"
34 : : #include "Integrate/Boundary.hpp"
35 : : #include "Integrate/Volume.hpp"
36 : : #include "Integrate/Source.hpp"
37 : : #include "Integrate/SolidTerms.hpp"
38 : : #include "RiemannChoice.hpp"
39 : : #include "MultiSpecies/MultiSpeciesIndexing.hpp"
40 : : #include "Reconstruction.hpp"
41 : : #include "Limiter.hpp"
42 : : #include "Problem/FieldOutput.hpp"
43 : : #include "Problem/BoxInitialization.hpp"
44 : : #include "PrefIndicator.hpp"
45 : : #include "MultiSpecies/BCFunctions.hpp"
46 : : #include "MultiSpecies/MiscMultiSpeciesFns.hpp"
47 : : #include "EoS/GetMatProp.hpp"
48 : : #include "Mixture/Mixture.hpp"
49 : :
50 : : namespace inciter {
51 : :
52 : : extern ctr::InputDeck g_inputdeck;
53 : :
54 : : namespace dg {
55 : :
56 : : //! \brief MultiSpecies used polymorphically with tk::DGPDE
57 : : //! \details The template arguments specify policies and are used to configure
58 : : //! the behavior of the class. The policies are:
59 : : //! - Physics - physics configuration, see PDE/MultiSpecies/Physics.h
60 : : //! - Problem - problem configuration, see PDE/MultiSpecies/Problem.h
61 : : //! \note The default physics is Euler, which is set in
62 : : //! inciter::LuaParser::storeInputDeck()
63 : : template< class Physics, class Problem >
64 : : class MultiSpecies {
65 : :
66 : : private:
67 : : using eq = tag::multispecies;
68 : :
69 : : public:
70 : : //! Constructor
71 : 24 : explicit MultiSpecies() :
72 : : m_physics(),
73 : : m_ncomp( g_inputdeck.get< tag::ncomp >() ),
74 : : m_nprim(nprim()),
75 [ + - ]: 24 : m_riemann( multispeciesRiemannSolver( g_inputdeck.get< tag::flux >() ) )
76 : : {
77 : : // associate boundary condition configurations with state functions
78 [ + - ][ + - ]: 360 : brigand::for_each< ctr::bclist::Keys >( ConfigBC( m_bc,
[ + - ][ + + ]
[ + - ][ + + ]
[ + - ][ - - ]
[ - - ][ - - ]
[ - - ]
79 : : // BC State functions
80 : : { dirichlet
81 : : , symmetry
82 : : , invalidBC // Outlet BC not implemented
83 : : , farfield
84 : : , extrapolate
85 : : , noslipwall
86 : : , symmetry }, // Slip equivalent to symmetry without mesh motion
87 : : // BC Gradient functions
88 : : { noOpGrad
89 : : , symmetryGrad
90 : : , noOpGrad
91 : : , noOpGrad
92 : : , noOpGrad
93 : : , noOpGrad
94 : : , symmetryGrad }
95 : : ) );
96 : :
97 : : // EoS initialization
98 [ + - ]: 24 : initializeSpeciesEoS( m_mat_blk );
99 : 24 : }
100 : :
101 : : //! Find the number of primitive quantities required for this PDE system
102 : : //! \return The number of primitive quantities required to be stored for
103 : : //! this PDE system (zero for multi-species)
104 : : std::size_t nprim() const { return 0; }
105 : :
106 : : //! Find the number of materials set up for this PDE system
107 : : //! \return The number of materials set up for this PDE system
108 : : std::size_t nmat() const { return 1; }
109 : :
110 : : //! Assign number of DOFs per equation in the PDE system
111 : : //! \param[in,out] numEqDof Array storing number of Dofs for each PDE
112 : : //! equation
113 : : void numEquationDofs(std::vector< std::size_t >& numEqDof) const
114 : : {
115 : : // all equation-dofs initialized to ndofs
116 [ + + ]: 572 : for (std::size_t i=0; i<m_ncomp; ++i) {
117 : 484 : numEqDof.push_back(g_inputdeck.get< tag::ndof >());
118 : : }
119 : : }
120 : :
121 : : //! Determine elements that lie inside the user-defined IC box
122 : : //! \param[in] geoElem Element geometry array
123 : : //! \param[in] nielem Number of internal elements
124 : : //! \param[in,out] inbox List of nodes at which box user ICs are set for
125 : : //! each IC box
126 : : void IcBoxElems( const tk::Fields& geoElem,
127 : : std::size_t nielem,
128 : : std::vector< std::unordered_set< std::size_t > >& inbox ) const
129 : : {
130 : 88 : tk::BoxElems< eq >(geoElem, nielem, inbox);
131 : : }
132 : :
133 : : //! Find how many 'stiff equations' in this PDE system
134 : : //! \return number of stiff equations. Zero for now, but will need to change
135 : : //! as chemical non-equilibrium is added
136 : : std::size_t nstiffeq() const
137 : : {
138 : : return 0;
139 : : }
140 : :
141 : : //! Find how many 'non-stiff equations' in this PDE system
142 : : //! \return number of non-stiff equations
143 : : std::size_t nnonstiffeq() const
144 : : {
145 : : return m_ncomp-nstiffeq();
146 : : }
147 : :
148 : : //! Locate the stiff equations.
149 : : //! \param[out] stiffEqIdx list with pointers to stiff equations. Empty
150 : : //! for now but will have to index to chemical non-equilibrium when added
151 : : void setStiffEqIdx( std::vector< std::size_t >& stiffEqIdx ) const
152 : : {
153 : 0 : stiffEqIdx.resize(0);
154 : : }
155 : :
156 : : //! Locate the nonstiff equations.
157 : : //! \param[out] nonStiffEqIdx list with pointers to nonstiff equations
158 : : void setNonStiffEqIdx( std::vector< std::size_t >& nonStiffEqIdx ) const
159 : : {
160 : 0 : nonStiffEqIdx.resize(0);
161 : : }
162 : :
163 : : //! Initialize the compressible flow equations, prepare for time integration
164 : : //! \param[in] L Block diagonal mass matrix
165 : : //! \param[in] inpoel Element-node connectivity
166 : : //! \param[in] coord Array of nodal coordinates
167 : : //! \param[in] inbox List of elements at which box user ICs are set for
168 : : //! each IC box
169 : : //! \param[in] elemblkid Element ids associated with mesh block ids where
170 : : //! user ICs are set
171 : : //! \param[in,out] unk Array of unknowns
172 : : //! \param[in] t Physical time
173 : : //! \param[in] nielem Number of internal elements
174 [ + - ]: 88 : void initialize( const tk::Fields& L,
175 : : const std::vector< std::size_t >& inpoel,
176 : : const tk::UnsMesh::Coords& coord,
177 : : const std::vector< std::unordered_set< std::size_t > >& inbox,
178 : : const std::unordered_map< std::size_t, std::set< std::size_t > >&
179 : : elemblkid,
180 : : tk::Fields& unk,
181 : : tk::real t,
182 : : const std::size_t nielem ) const
183 : : {
184 [ + - ]: 88 : tk::initialize( m_ncomp, m_mat_blk, L, inpoel, coord,
185 : : Problem::initialize, unk, t, nielem );
186 : :
187 : 88 : const auto rdof = g_inputdeck.get< tag::rdof >();
188 : : const auto& ic = g_inputdeck.get< tag::ic >();
189 : : const auto& icbox = ic.get< tag::box >();
190 : : const auto& icmbk = ic.get< tag::meshblock >();
191 : :
192 : : // Set initial conditions inside user-defined IC boxes and mesh blocks
193 : 88 : std::vector< tk::real > s(m_ncomp, 0.0);
194 [ + + ]: 9184 : for (std::size_t e=0; e<nielem; ++e) {
195 : : // inside user-defined box
196 [ + - ]: 9096 : if (!icbox.empty()) {
197 : : std::size_t bcnt = 0;
198 [ + + ]: 18192 : for (const auto& b : icbox) { // for all boxes
199 [ + - ][ + + ]: 18192 : if (inbox.size() > bcnt && inbox[bcnt].find(e) != inbox[bcnt].end())
200 : : {
201 [ + - ][ - - ]: 4548 : std::vector< tk::real > box
202 : : { b.template get< tag::xmin >(), b.template get< tag::xmax >(),
203 : : b.template get< tag::ymin >(), b.template get< tag::ymax >(),
204 : : b.template get< tag::zmin >(), b.template get< tag::zmax >() };
205 : 4548 : auto V_ex = (box[1]-box[0]) * (box[3]-box[2]) * (box[5]-box[4]);
206 [ + + ]: 29562 : for (std::size_t c=0; c<m_ncomp; ++c) {
207 : 25014 : auto mark = c*rdof;
208 : 25014 : s[c] = unk(e,mark);
209 : : // set high-order DOFs to zero
210 [ + + ]: 100056 : for (std::size_t i=1; i<rdof; ++i)
211 : 75042 : unk(e,mark+i) = 0.0;
212 : : }
213 [ + - ]: 4548 : initializeBox<ctr::boxList>( m_mat_blk, V_ex, t, b, s );
214 : : // store box-initialization in solution vector
215 [ + + ]: 29562 : for (std::size_t c=0; c<m_ncomp; ++c) {
216 : 25014 : auto mark = c*rdof;
217 : 25014 : unk(e,mark) = s[c];
218 : : }
219 : : }
220 : 9096 : ++bcnt;
221 : : }
222 : : }
223 : :
224 : : // inside user-specified mesh blocks
225 [ - + ]: 9096 : if (!icmbk.empty()) {
226 [ - - ]: 0 : for (const auto& b : icmbk) { // for all blocks
227 : 0 : auto blid = b.get< tag::blockid >();
228 : 0 : auto V_ex = b.get< tag::volume >();
229 [ - - ]: 0 : if (elemblkid.find(blid) != elemblkid.end()) {
230 : : const auto& elset = tk::cref_find(elemblkid, blid);
231 [ - - ]: 0 : if (elset.find(e) != elset.end()) {
232 [ - - ]: 0 : initializeBox<ctr::meshblockList>( m_mat_blk, V_ex, t, b, s );
233 : : // store initialization in solution vector
234 [ - - ]: 0 : for (std::size_t c=0; c<m_ncomp; ++c) {
235 : 0 : auto mark = c*rdof;
236 : 0 : unk(e,mark) = s[c];
237 : : }
238 : : }
239 : : }
240 : : }
241 : : }
242 : : }
243 : 88 : }
244 : :
245 : : //! Compute density constraint for a given material. No-op
246 : : // //! \param[in] nelem Number of elements
247 : : // //! \param[in] unk Array of unknowns
248 : : //! \param[out] densityConstr Density Constraint: rho/(rho0*det(g))
249 : : void computeDensityConstr( std::size_t /*nelem*/,
250 : : tk::Fields& /*unk*/,
251 : : std::vector< tk::real >& densityConstr) const
252 : : {
253 : 176 : densityConstr.resize(0);
254 : : }
255 : :
256 : : //! Compute the left hand side block-diagonal mass matrix
257 : : //! \param[in] geoElem Element geometry array
258 : : //! \param[in,out] l Block diagonal mass matrix
259 : : void lhs( const tk::Fields& geoElem, tk::Fields& l ) const {
260 : 88 : const auto ndof = g_inputdeck.get< tag::ndof >();
261 : 88 : tk::mass( m_ncomp, ndof, geoElem, l );
262 : : }
263 : :
264 : : //! Update the interface cells to first order dofs. No-op.
265 : : // //! \param[in] unk Array of unknowns
266 : : // //! \param[in] nielem Number of internal elements
267 : : // //! \param[in,out] ndofel Array of dofs
268 : : // //! \param[in,out] interface Vector of interface marker
269 : : void updateInterfaceCells( tk::Fields& /*unk*/,
270 : : std::size_t /*nielem*/,
271 : : std::vector< std::size_t >& /*ndofel*/,
272 : : std::vector< std::size_t >& /*interface*/ ) const {}
273 : :
274 : : //! Update the primitives for this PDE system. No-op.
275 : : // //! \param[in] unk Array of unknowns
276 : : // //! \param[in] L The left hand side block-diagonal mass matrix
277 : : // //! \param[in] geoElem Element geometry array
278 : : // //! \param[in,out] prim Array of primitives
279 : : // //! \param[in] nielem Number of internal elements
280 : : // //! \param[in] ndofel Array of dofs
281 : : void updatePrimitives( const tk::Fields& /*unk*/,
282 : : const tk::Fields& /*L*/,
283 : : const tk::Fields& /*geoElem*/,
284 : : tk::Fields& /*prim*/,
285 : : std::size_t /*nielem*/,
286 : : std::vector< std::size_t >& /*ndofel*/ ) const {}
287 : :
288 : : //! Clean up the state of trace materials for this PDE system. No-op.
289 : : // //! \param[in] t Physical time
290 : : // //! \param[in] geoElem Element geometry array
291 : : // //! \param[in,out] unk Array of unknowns
292 : : // //! \param[in,out] prim Array of primitives
293 : : // //! \param[in] nielem Number of internal elements
294 : : void cleanTraceMaterial( tk::real /*t*/,
295 : : const tk::Fields& /*geoElem*/,
296 : : tk::Fields& /*unk*/,
297 : : tk::Fields& /*prim*/,
298 : : std::size_t /*nielem*/ ) const {}
299 : :
300 : : //! Reconstruct second-order solution from first-order
301 : : //! \param[in] geoElem Element geometry array
302 : : //! \param[in] fd Face connectivity and boundary conditions object
303 : : //! \param[in] esup Elements-surrounding-nodes connectivity
304 : : //! \param[in] inpoel Element-node connectivity
305 : : //! \param[in] coord Array of nodal coordinates
306 : : //! \param[in,out] U Solution vector at recent time step
307 : : // //! \param[in,out] P Vector of primitives at recent time step
308 : : //! \param[in] pref Indicator for p-adaptive algorithm
309 : : //! \param[in] ndofel Vector of local number of degrees of freedome
310 : 6600 : void reconstruct( tk::real,
311 : : const tk::Fields&,
312 : : const tk::Fields& geoElem,
313 : : const inciter::FaceData& fd,
314 : : const std::map< std::size_t, std::vector< std::size_t > >&
315 : : esup,
316 : : const std::vector< std::size_t >& inpoel,
317 : : const tk::UnsMesh::Coords& coord,
318 : : tk::Fields& U,
319 : : tk::Fields& /*P*/,
320 : : const bool pref,
321 : : const std::vector< std::size_t >& ndofel ) const
322 : : {
323 : 6600 : const auto rdof = g_inputdeck.get< tag::rdof >();
324 : 6600 : const auto ndof = g_inputdeck.get< tag::ndof >();
325 : :
326 : : bool is_p0p1(false);
327 [ + - ]: 6600 : if (rdof == 4 && ndof == 1)
328 : : is_p0p1 = true;
329 : :
330 : 6600 : const auto nelem = fd.Esuel().size()/4;
331 : :
332 : : Assert( U.nprop() == rdof*m_ncomp, "Number of components in solution "
333 : : "vector must equal "+ std::to_string(rdof*m_ncomp) );
334 : :
335 : : //----- reconstruction of conserved quantities -----
336 : : //--------------------------------------------------
337 : :
338 [ + + ]: 688800 : for (std::size_t e=0; e<nelem; ++e)
339 : : {
340 : : std::vector< std::size_t > vars;
341 : : // check if element is marked as p0p1
342 [ - + ][ - - ]: 682200 : if ( (pref && ndofel[e] == 1) || is_p0p1 ) {
[ + - ]
343 : : // 1. specify how many variables need to be reconstructed
344 [ + + ][ + - ]: 4434300 : for (std::size_t c=0; c<m_ncomp; ++c) vars.push_back(c);
345 : :
346 : : // 2. solve 3x3 least-squares system
347 : : // Reconstruct second-order dofs in Taylor space using nodal-stencils
348 [ + - ]: 682200 : tk::recoLeastSqExtStencil( rdof, e, esup, inpoel, geoElem, U, vars );
349 : :
350 : : // 3. transform reconstructed derivatives to Dubiner dofs
351 [ + - ]: 682200 : tk::transform_P0P1( rdof, e, inpoel, coord, U, vars );
352 : : }
353 : : }
354 : 6600 : }
355 : :
356 : : //! Limit second-order solution, and primitive quantities separately
357 : : // //! \param[in] pref Indicator for p-adaptive algorithm
358 : : //! \param[in] geoFace Face geometry array
359 : : //! \param[in] geoElem Element geometry array
360 : : //! \param[in] fd Face connectivity and boundary conditions object
361 : : //! \param[in] esup Elements-surrounding-nodes connectivity
362 : : //! \param[in] inpoel Element-node connectivity
363 : : //! \param[in] coord Array of nodal coordinates
364 : : //! \param[in] ndofel Vector of local number of degrees of freedome
365 : : // //! \param[in] gid Local->global node id map
366 : : // //! \param[in] bid Local chare-boundary node ids (value) associated to
367 : : // //! global node ids (key)
368 : : // //! \param[in] uNodalExtrm Chare-boundary nodal extrema for conservative
369 : : // //! variables
370 : : // //! \param[in] pNodalExtrm Chare-boundary nodal extrema for primitive
371 : : // //! variables
372 : : // //! \param[in] mtInv Inverse of Taylor mass matrix
373 : : //! \param[in,out] U Solution vector at recent time step
374 : : // //! \param[in,out] P Vector of primitives at recent time step
375 : : //! \param[in,out] shockmarker Vector of shock-marker values
376 : 6600 : void limit( [[maybe_unused]] tk::real,
377 : : const bool /*pref*/,
378 : : const tk::Fields& geoFace,
379 : : const tk::Fields& geoElem,
380 : : const inciter::FaceData& fd,
381 : : const std::map< std::size_t, std::vector< std::size_t > >& esup,
382 : : const std::vector< std::size_t >& inpoel,
383 : : const tk::UnsMesh::Coords& coord,
384 : : const std::vector< std::size_t >& ndofel,
385 : : const std::vector< std::size_t >& /*gid*/,
386 : : const std::unordered_map< std::size_t, std::size_t >& /*bid*/,
387 : : const std::vector< std::vector<tk::real> >& /*uNodalExtrm*/,
388 : : const std::vector< std::vector<tk::real> >& /*pNodalExtrm*/,
389 : : const std::vector< std::vector<tk::real> >& /*mtInv*/,
390 : : tk::Fields& U,
391 : : tk::Fields& /*P*/,
392 : : std::vector< std::size_t >& shockmarker ) const
393 : : {
394 : 6600 : const auto limiter = g_inputdeck.get< tag::limiter >();
395 : 6600 : auto nspec = g_inputdeck.get< tag::multispecies, tag::nspec >();
396 : 6600 : const auto rdof = g_inputdeck.get< tag::rdof >();
397 : : const auto& solidx = g_inputdeck.get< tag::matidxmap, tag::solidx >();
398 : :
399 : : // limit vectors of conserved and primitive quantities
400 [ + - ]: 6600 : if (limiter == ctr::LimiterType::VERTEXBASEDP1 && rdof == 4)
401 : : {
402 : 6600 : VertexBasedMultiSpecies_P1( esup, inpoel, ndofel, fd.Esuel().size()/4,
403 [ + - ]: 6600 : m_mat_blk, fd, geoFace, geoElem, coord, flux, solidx, U, nspec,
404 : : shockmarker );
405 : : }
406 [ - - ]: 0 : else if (limiter == ctr::LimiterType::VERTEXBASEDP1 && rdof == 10)
407 : : {
408 : 0 : VertexBasedMultiSpecies_P2( esup, inpoel, ndofel, fd.Esuel().size()/4,
409 [ - - ]: 0 : m_mat_blk, fd, geoFace, geoElem, coord, flux, solidx, U, nspec,
410 : : shockmarker );
411 : : }
412 [ - - ]: 0 : else if (limiter != ctr::LimiterType::NOLIMITER)
413 : : {
414 [ - - ][ - - ]: 0 : Throw("Limiter type not configured for multispecies.");
[ - - ][ - - ]
[ - - ][ - - ]
415 : : }
416 : 6600 : }
417 : :
418 : : //! Apply CPL to the conservative variable solution for this PDE system
419 : : // //! \param[in] prim Array of primitive variables
420 : : // //! \param[in] geoElem Element geometry array
421 : : // //! \param[in] inpoel Element-node connectivity
422 : : // //! \param[in] coord Array of nodal coordinates
423 : : // //! \param[in,out] unk Array of conservative variables
424 : : // //! \param[in] nielem Number of internal elements
425 : : //! \details This function applies CPL to obtain consistent dofs for
426 : : //! conservative quantities based on the limited primitive quantities.
427 : : //! No-op for now, but might need in the future, see appendix of paper.
428 : : //! See Pandare et al. (2023). On the Design of Stable,
429 : : //! Consistent, and Conservative High-Order Methods for Multi-Material
430 : : //! Hydrodynamics. J Comp Phys, 112313.
431 : : void CPL( const tk::Fields& /*prim*/,
432 : : const tk::Fields& /*geoElem*/,
433 : : const std::vector< std::size_t >& /*inpoel*/,
434 : : const tk::UnsMesh::Coords& /*coord*/,
435 : : tk::Fields& /*unk*/,
436 : : std::size_t /*nielem*/ ) const {}
437 : :
438 : : //! Return cell-average deformation gradient tensor. No-op.
439 : : std::array< std::vector< tk::real >, 9 > cellAvgDeformGrad(
440 : : const tk::Fields&,
441 : : std::size_t ) const
442 : : { return {}; }
443 : :
444 : : //! Reset the high order solution for p-adaptive scheme
445 : : //! \param[in] fd Face connectivity and boundary conditions object
446 : : //! \param[in,out] unk Solution vector at recent time step
447 : : //! \param[in,out] prim Primitive vector at recent time step
448 : : //! \param[in] ndofel Vector of local number of degrees of freedome
449 : : //! \details This function reset the high order coefficient for p-adaptive
450 : : //! solution polynomials. Unlike compflow class, the high order of fv
451 : : //! solution will not be reset since p0p1 is the base scheme for
452 : : //! multi-species p-adaptive DG method.
453 : 0 : void resetAdapSol( const inciter::FaceData& fd,
454 : : tk::Fields& unk,
455 : : tk::Fields& prim,
456 : : const std::vector< std::size_t >& ndofel ) const
457 : : {
458 : 0 : const auto rdof = g_inputdeck.get< tag::rdof >();
459 : 0 : const auto ncomp = unk.nprop() / rdof;
460 : 0 : const auto nprim = prim.nprop() / rdof;
461 : :
462 [ - - ]: 0 : for(std::size_t e = 0; e < fd.Esuel().size()/4; e++)
463 : : {
464 [ - - ]: 0 : if(ndofel[e] < 10)
465 : : {
466 [ - - ]: 0 : for (std::size_t c=0; c<ncomp; ++c)
467 : : {
468 : 0 : auto mark = c*rdof;
469 : 0 : unk(e, mark+4) = 0.0;
470 : 0 : unk(e, mark+5) = 0.0;
471 : 0 : unk(e, mark+6) = 0.0;
472 : 0 : unk(e, mark+7) = 0.0;
473 : 0 : unk(e, mark+8) = 0.0;
474 : 0 : unk(e, mark+9) = 0.0;
475 : : }
476 [ - - ]: 0 : for (std::size_t c=0; c<nprim; ++c)
477 : : {
478 : 0 : auto mark = c*rdof;
479 : 0 : prim(e, mark+4) = 0.0;
480 : 0 : prim(e, mark+5) = 0.0;
481 : 0 : prim(e, mark+6) = 0.0;
482 : 0 : prim(e, mark+7) = 0.0;
483 : 0 : prim(e, mark+8) = 0.0;
484 : 0 : prim(e, mark+9) = 0.0;
485 : : }
486 : : }
487 : : }
488 : 0 : }
489 : :
490 : : //! Compute right hand side
491 : : //! \param[in] t Physical time
492 : : //! \param[in] pref Indicator for p-adaptive algorithm
493 : : //! \param[in] geoFace Face geometry array
494 : : //! \param[in] geoElem Element geometry array
495 : : //! \param[in] fd Face connectivity and boundary conditions object
496 : : //! \param[in] inpoel Element-node connectivity
497 : : //! \param[in] coord Array of nodal coordinates
498 : : //! \param[in] U Solution vector at recent time step
499 : : //! \param[in] P Primitive vector at recent time step
500 : : //! \param[in] ndofel Vector of local number of degrees of freedom
501 : : //! \param[in] dt Delta time
502 : : //! \param[in,out] R Right-hand side vector computed
503 : 6600 : void rhs( tk::real t,
504 : : const bool pref,
505 : : const tk::Fields& geoFace,
506 : : const tk::Fields& geoElem,
507 : : const inciter::FaceData& fd,
508 : : const std::vector< std::size_t >& inpoel,
509 : : const std::vector< std::unordered_set< std::size_t > >&,
510 : : const tk::UnsMesh::Coords& coord,
511 : : const tk::Fields& U,
512 : : const tk::Fields& P,
513 : : const std::vector< std::size_t >& ndofel,
514 : : const tk::real dt,
515 : : tk::Fields& R ) const
516 : : {
517 : 6600 : const auto ndof = g_inputdeck.get< tag::ndof >();
518 : 6600 : const auto rdof = g_inputdeck.get< tag::rdof >();
519 : : const auto& solidx = g_inputdeck.get< tag::matidxmap, tag::solidx >();
520 : :
521 : 6600 : const auto nelem = fd.Esuel().size()/4;
522 : :
523 : : Assert( U.nunk() == P.nunk(), "Number of unknowns in solution "
524 : : "vector and primitive vector at recent time step incorrect" );
525 : : Assert( U.nunk() == R.nunk(), "Number of unknowns in solution "
526 : : "vector and right-hand side at recent time step incorrect" );
527 : : Assert( U.nprop() == rdof*m_ncomp, "Number of components in solution "
528 : : "vector must equal "+ std::to_string(rdof*m_ncomp) );
529 : : Assert( P.nprop() == rdof*m_nprim, "Number of components in primitive "
530 : : "vector must equal "+ std::to_string(rdof*m_nprim) );
531 : : Assert( R.nprop() == ndof*m_ncomp, "Number of components in right-hand "
532 : : "side vector must equal "+ std::to_string(ndof*m_ncomp) );
533 : : Assert( fd.Inpofa().size()/3 == fd.Esuf().size()/2,
534 : : "Mismatch in inpofa size" );
535 : :
536 : : // set rhs to zero
537 : : R.fill(0.0);
538 : :
539 : : // empty vector for non-conservative terms. This vector is unused for
540 : : // multi-species flow since, there are no non-conservative terms
541 : : // in the system of PDEs.
542 : 6600 : std::vector< std::vector< tk::real > > riemannDeriv;
543 : :
544 : 6600 : std::vector< std::vector< tk::real > > vriem;
545 : 6600 : std::vector< std::vector< tk::real > > riemannLoc;
546 : :
547 : : // configure a no-op lambda for prescribed velocity
548 : : auto velfn = []( ncomp_t, tk::real, tk::real, tk::real, tk::real ){
549 : : return tk::VelFn::result_type(); };
550 : :
551 : : // compute internal surface flux integrals
552 [ + - ]: 6600 : tk::surfInt( pref, 1, m_mat_blk, t, ndof, rdof, inpoel, solidx,
553 [ + - ]: 6600 : coord, fd, geoFace, geoElem, m_riemann, velfn, U, P, ndofel,
554 : : dt, R, riemannDeriv );
555 : :
556 : : // compute optional source term
557 [ + - ]: 6600 : tk::srcInt( m_mat_blk, t, ndof, fd.Esuel().size()/4, inpoel,
558 : : coord, geoElem, Problem::src, ndofel, R );
559 : :
560 [ - + ]: 6600 : if(ndof > 1)
561 : : // compute volume integrals
562 [ - - ][ - - ]: 0 : tk::volInt( 1, t, m_mat_blk, ndof, rdof, nelem, inpoel, coord, geoElem,
[ - - ]
563 : : flux, velfn, U, P, ndofel, R );
564 : :
565 : : // compute boundary surface flux integrals
566 [ + + ]: 52800 : for (const auto& b : m_bc)
567 [ + - ]: 92400 : tk::bndSurfInt( pref, 1, m_mat_blk, ndof, rdof, std::get<0>(b), fd,
568 : : geoFace, geoElem, inpoel, coord, t, m_riemann, velfn,
569 : : std::get<1>(b), U, P, ndofel, R, riemannDeriv );
570 : :
571 : : // compute external (energy) sources
572 : : //m_physics.physSrc(nspec, t, geoElem, {}, R, {});
573 : 6600 : }
574 : :
575 : : //! Evaluate the adaptive indicator and mark the ndof for each element
576 : : //! \param[in] nunk Number of unknowns
577 : : //! \param[in] coord Array of nodal coordinates
578 : : //! \param[in] inpoel Element-node connectivity
579 : : //! \param[in] fd Face connectivity and boundary conditions object
580 : : //! \param[in] unk Array of unknowns
581 : : // //! \param[in] prim Array of primitive quantities
582 : : //! \param[in] indicator p-refinement indicator type
583 : : //! \param[in] ndof Number of degrees of freedom in the solution
584 : : //! \param[in] ndofmax Max number of degrees of freedom for p-refinement
585 : : //! \param[in] tolref Tolerance for p-refinement
586 : : //! \param[in,out] ndofel Vector of local number of degrees of freedome
587 [ - - ]: 0 : void eval_ndof( std::size_t nunk,
588 : : [[maybe_unused]] const tk::UnsMesh::Coords& coord,
589 : : [[maybe_unused]] const std::vector< std::size_t >& inpoel,
590 : : const inciter::FaceData& fd,
591 : : const tk::Fields& unk,
592 : : const tk::Fields& /*prim*/,
593 : : inciter::ctr::PrefIndicatorType indicator,
594 : : std::size_t ndof,
595 : : std::size_t ndofmax,
596 : : tk::real tolref,
597 : : std::vector< std::size_t >& ndofel ) const
598 : : {
599 : : const auto& esuel = fd.Esuel();
600 : :
601 [ - - ]: 0 : if(indicator == inciter::ctr::PrefIndicatorType::SPECTRAL_DECAY)
602 : 0 : spectral_decay(1, nunk, esuel, unk, ndof, ndofmax, tolref, ndofel);
603 : : else
604 [ - - ][ - - ]: 0 : Throw( "No such adaptive indicator type" );
[ - - ][ - - ]
[ - - ][ - - ]
605 : 0 : }
606 : :
607 : : //! Compute the minimum time step size
608 : : //! \param[in] fd Face connectivity and boundary conditions object
609 : : //! \param[in] geoFace Face geometry array
610 : : //! \param[in] geoElem Element geometry array
611 : : // //! \param[in] ndofel Vector of local number of degrees of freedom
612 : : //! \param[in] U Solution vector at recent time step
613 : : //! \param[in] P Vector of primitive quantities at recent time step
614 : : //! \param[in] nielem Number of internal elements
615 : : //! \return Minimum time step size
616 : : //! \details The allowable dt is calculated by looking at the maximum
617 : : //! wave-speed in elements surrounding each face, times the area of that
618 : : //! face. Once the maximum of this quantity over the mesh is determined,
619 : : //! the volume of each cell is divided by this quantity. A minimum of this
620 : : //! ratio is found over the entire mesh, which gives the allowable dt.
621 : : tk::real dt( const std::array< std::vector< tk::real >, 3 >&,
622 : : const std::vector< std::size_t >&,
623 : : const inciter::FaceData& fd,
624 : : const tk::Fields& geoFace,
625 : : const tk::Fields& geoElem,
626 : : const std::vector< std::size_t >& /*ndofel*/,
627 : : const tk::Fields& U,
628 : : const tk::Fields& P,
629 : : const std::size_t nielem ) const
630 : : {
631 : : const auto ndof = g_inputdeck.get< tag::ndof >();
632 : : auto nspec = g_inputdeck.get< tag::multispecies, tag::nspec >();
633 : :
634 : : auto mindt = timeStepSizeMultiSpecies( m_mat_blk, fd.Esuf(), geoFace,
635 : : geoElem, nielem, nspec, U, P);
636 : :
637 : : //if (viscous)
638 : : // mindt = std::min(mindt, timeStepSizeViscousFV(geoElem, nielem, nspec, U));
639 : : //mindt = std::min(mindt, m_physics.dtRestriction(geoElem, nielem, {}));
640 : :
641 : : tk::real dgp = 0.0;
642 : : if (ndof == 4)
643 : : {
644 : : dgp = 1.0;
645 : : }
646 : : else if (ndof == 10)
647 : : {
648 : : dgp = 2.0;
649 : : }
650 : :
651 : : // Scale smallest dt with CFL coefficient and the CFL is scaled by (2*p+1)
652 : : // where p is the order of the DG polynomial by linear stability theory.
653 : : mindt /= (2.0*dgp + 1.0);
654 : : return mindt;
655 : : }
656 : :
657 : : //! Compute stiff terms for a single element. No-op until chem sources added
658 : : // //! \param[in] e Element number
659 : : // //! \param[in] geoElem Element geometry array
660 : : // //! \param[in] inpoel Element-node connectivity
661 : : // //! \param[in] coord Array of nodal coordinates
662 : : // //! \param[in] U Solution vector at recent time step
663 : : // //! \param[in] P Primitive vector at recent time step
664 : : // //! \param[in] ndofel Vector of local number of degrees of freedom
665 : : // //! \param[in,out] R Right-hand side vector computed
666 : : void stiff_rhs( std::size_t /*e*/,
667 : : const tk::Fields& /*geoElem*/,
668 : : const std::vector< std::size_t >& /*inpoel*/,
669 : : const tk::UnsMesh::Coords& /*coord*/,
670 : : const tk::Fields& /*U*/,
671 : : const tk::Fields& /*P*/,
672 : : const std::vector< std::size_t >& /*ndofel*/,
673 : : tk::Fields& /*R*/ ) const {}
674 : :
675 : : //! Extract the velocity field at cell nodes. Currently unused.
676 : : // //! \param[in] U Solution vector at recent time step
677 : : // //! \param[in] N Element node indices
678 : : //! \return Array of the four values of the velocity field
679 : : std::array< std::array< tk::real, 4 >, 3 >
680 : : velocity( const tk::Fields& /*U*/,
681 : : const std::array< std::vector< tk::real >, 3 >&,
682 : : const std::array< std::size_t, 4 >& /*N*/ ) const
683 : : {
684 : : std::array< std::array< tk::real, 4 >, 3 > v;
685 : : return v;
686 : : }
687 : :
688 : : //! Return a map that associates user-specified strings to functions
689 : : //! \return Map that associates user-specified strings to functions that
690 : : //! compute relevant quantities to be output to file
691 : : std::map< std::string, tk::GetVarFn > OutVarFn() const
692 : 352 : { return MultiSpeciesOutVarFn(); }
693 : :
694 : : //! Return analytic field names to be output to file
695 : : //! \return Vector of strings labelling analytic fields output in file
696 : : std::vector< std::string > analyticFieldNames() const {
697 : 0 : auto nspec = g_inputdeck.get< eq, tag::nspec >();
698 : :
699 : 0 : return MultiSpeciesFieldNames(nspec);
700 : : }
701 : :
702 : : //! Return time history field names to be output to file
703 : : //! \return Vector of strings labelling time history fields output in file
704 : : std::vector< std::string > histNames() const {
705 : 0 : return MultiSpeciesHistNames();
706 : : }
707 : :
708 : : //! Return surface field output going to file
709 : : std::vector< std::vector< tk::real > >
710 : : surfOutput( const std::map< int, std::vector< std::size_t > >&,
711 : : tk::Fields& ) const
712 : : {
713 : : std::vector< std::vector< tk::real > > s; // punt for now
714 : : return s;
715 : : }
716 : :
717 : : //! Return time history field output evaluated at time history points
718 : : //! \param[in] h History point data
719 : : //! \param[in] inpoel Element-node connectivity
720 : : //! \param[in] coord Array of nodal coordinates
721 : : //! \param[in] U Array of unknowns
722 : : // //! \param[in] P Array of primitive quantities
723 : : //! \return Vector of time history output of bulk flow quantities (density,
724 : : //! velocity, total energy, and pressure) evaluated at time history points
725 : : std::vector< std::vector< tk::real > >
726 : 0 : histOutput( const std::vector< HistData >& h,
727 : : const std::vector< std::size_t >& inpoel,
728 : : const tk::UnsMesh::Coords& coord,
729 : : const tk::Fields& U,
730 : : const tk::Fields& /*P*/ ) const
731 : : {
732 : 0 : const auto rdof = g_inputdeck.get< tag::rdof >();
733 : 0 : auto nspec = g_inputdeck.get< tag::multispecies, tag::nspec >();
734 : :
735 : : const auto& x = coord[0];
736 : : const auto& y = coord[1];
737 : : const auto& z = coord[2];
738 : :
739 : 0 : std::vector< std::vector< tk::real > > Up(h.size());
740 : :
741 : : std::size_t j = 0;
742 [ - - ]: 0 : for (const auto& p : h) {
743 : 0 : auto e = p.get< tag::elem >();
744 : 0 : auto chp = p.get< tag::coord >();
745 : :
746 : : // Evaluate inverse Jacobian
747 : 0 : std::array< std::array< tk::real, 3>, 4 > cp{{
748 : : {{ x[inpoel[4*e ]], y[inpoel[4*e ]], z[inpoel[4*e ]] }},
749 : : {{ x[inpoel[4*e+1]], y[inpoel[4*e+1]], z[inpoel[4*e+1]] }},
750 : : {{ x[inpoel[4*e+2]], y[inpoel[4*e+2]], z[inpoel[4*e+2]] }},
751 : : {{ x[inpoel[4*e+3]], y[inpoel[4*e+3]], z[inpoel[4*e+3]] }} }};
752 : 0 : auto J = tk::inverseJacobian( cp[0], cp[1], cp[2], cp[3] );
753 : :
754 : : // evaluate solution at history-point
755 : 0 : std::array< tk::real, 3 > dc{{chp[0]-cp[0][0], chp[1]-cp[0][1],
756 [ - - ]: 0 : chp[2]-cp[0][2]}};
757 [ - - ]: 0 : auto B = tk::eval_basis(rdof, tk::dot(J[0],dc), tk::dot(J[1],dc),
758 : : tk::dot(J[2],dc));
759 [ - - ]: 0 : auto uhp = eval_state(m_ncomp, rdof, rdof, e, U, B);
760 : :
761 : : // Mixture calculations, initialized
762 [ - - ]: 0 : Mixture mix(nspec, uhp, m_mat_blk);
763 : :
764 : : // store solution in history output vector
765 [ - - ][ - - ]: 0 : Up[j].resize(6+nspec, 0.0);
766 [ - - ]: 0 : Up[j][0] = mix.get_mix_density();
767 [ - - ]: 0 : Up[j][1] = uhp[multispecies::momentumIdx(nspec,0)]/Up[j][0];
768 : 0 : Up[j][2] = uhp[multispecies::momentumIdx(nspec,1)]/Up[j][0];
769 [ - - ]: 0 : Up[j][3] = uhp[multispecies::momentumIdx(nspec,2)]/Up[j][0];
770 : 0 : Up[j][4] = uhp[multispecies::energyIdx(nspec,0)];
771 [ - - ]: 0 : Up[j][5] = mix.pressure( Up[j][0], Up[j][1], Up[j][2], Up[j][3],
772 : : Up[j][4], m_mat_blk);
773 [ - - ]: 0 : for (std::size_t k=0; k<nspec; ++k) {
774 : 0 : Up[j][6+k] = uhp[multispecies::densityIdx(nspec,k)]/Up[j][0];
775 : : }
776 [ - - ]: 0 : ++j;
777 : : }
778 : :
779 : 0 : return Up;
780 : : }
781 : :
782 : : //! Return names of integral variables to be output to diagnostics file
783 : : //! \return Vector of strings labelling integral variables output
784 : : std::vector< std::string > names() const
785 : : {
786 : 6 : auto nspec = g_inputdeck.get< tag::multispecies, tag::nspec >();
787 : 6 : return MultiSpeciesDiagNames(nspec);
788 : : }
789 : :
790 : : //! Return analytic solution (if defined by Problem) at xi, yi, zi, t
791 : : //! \param[in] xi X-coordinate at which to evaluate the analytic solution
792 : : //! \param[in] yi Y-coordinate at which to evaluate the analytic solution
793 : : //! \param[in] zi Z-coordinate at which to evaluate the analytic solution
794 : : //! \param[in] t Physical time at which to evaluate the analytic solution
795 : : //! \return Vector of analytic solution at given location and time
796 : : std::vector< tk::real >
797 : : analyticSolution( tk::real xi, tk::real yi, tk::real zi, tk::real t ) const
798 : 0 : { return Problem::analyticSolution( m_ncomp, m_mat_blk, xi, yi, zi, t ); }
799 : :
800 : : //! Return analytic solution for conserved variables
801 : : //! \param[in] xi X-coordinate at which to evaluate the analytic solution
802 : : //! \param[in] yi Y-coordinate at which to evaluate the analytic solution
803 : : //! \param[in] zi Z-coordinate at which to evaluate the analytic solution
804 : : //! \param[in] t Physical time at which to evaluate the analytic solution
805 : : //! \return Vector of analytic solution at given location and time
806 : : std::vector< tk::real >
807 : : solution( tk::real xi, tk::real yi, tk::real zi, tk::real t ) const
808 : 227400 : { return Problem::initialize( m_ncomp, m_mat_blk, xi, yi, zi, t ); }
809 : :
810 : : //! Return cell-averaged specific total energy for an element
811 : : //! \param[in] e Element id for which total energy is required
812 : : //! \param[in] unk Vector of conserved quantities
813 : : //! \return Cell-averaged specific total energy for given element
814 : : tk::real sp_totalenergy(std::size_t e, const tk::Fields& unk) const
815 : : {
816 : 227400 : const auto rdof = g_inputdeck.get< tag::rdof >();
817 : 227400 : auto nspec = g_inputdeck.get< tag::multispecies, tag::nspec >();
818 : :
819 : 227400 : return unk(e, multispecies::energyDofIdx(nspec,0,rdof,0));
820 : : }
821 : :
822 : : private:
823 : : //! Physics policy
824 : : const Physics m_physics;
825 : : //! Number of components in this PDE system
826 : : const ncomp_t m_ncomp;
827 : : //! Number of primitive quantities stored in this PDE system
828 : : const ncomp_t m_nprim;
829 : : //! Riemann solver
830 : : tk::RiemannFluxFn m_riemann;
831 : : //! BC configuration
832 : : BCStateFn m_bc;
833 : : //! EOS material block
834 : : std::vector< EOS > m_mat_blk;
835 : :
836 : : //! Evaluate conservative part of physical flux function for this PDE system
837 : : //! \param[in] ncomp Number of scalar components in this PDE system
838 : : //! \param[in] ugp Numerical solution at the Gauss point at which to
839 : : //! evaluate the flux
840 : : //! \return Flux vectors for all components in this PDE system
841 : : //! \note The function signature must follow tk::FluxFn
842 : : static tk::FluxFn::result_type
843 : 0 : flux( [[maybe_unused]] ncomp_t ncomp,
844 : : const std::vector< EOS >& mat_blk,
845 : : const std::vector< tk::real >& ugp,
846 : : const std::vector< std::array< tk::real, 3 > >& )
847 : : {
848 : : Assert( ugp.size() == ncomp, "Size mismatch" );
849 : :
850 : 0 : auto nspec = g_inputdeck.get< tag::multispecies, tag::nspec >();
851 : :
852 : 0 : std::vector< std::array< tk::real, 3 > > fl( ugp.size() );
853 : :
854 [ - - ]: 0 : Mixture mix(nspec, ugp, mat_blk);
855 : 0 : auto rhob = mix.get_mix_density();
856 : :
857 [ - - ]: 0 : std::array< tk::real, 3 > u{{
858 [ - - ]: 0 : ugp[multispecies::momentumIdx(nspec,0)] / rhob,
859 : 0 : ugp[multispecies::momentumIdx(nspec,1)] / rhob,
860 [ - - ]: 0 : ugp[multispecies::momentumIdx(nspec,2)] / rhob }};
861 : 0 : auto rhoE0 = ugp[multispecies::energyIdx(nspec,0)];
862 [ - - ]: 0 : auto p = mix.pressure(rhob, u[0], u[1], u[2], rhoE0, mat_blk);
863 : :
864 : : // density flux
865 [ - - ]: 0 : for (std::size_t k=0; k<nspec; ++k) {
866 : : auto idx = multispecies::densityIdx(nspec, k);
867 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j) {
868 : 0 : fl[idx][j] = ugp[idx] * u[j];
869 : : }
870 : : }
871 : :
872 : : // momentum flux
873 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i) {
874 : : auto idx = multispecies::momentumIdx(nspec,i);
875 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j) {
876 [ - - ]: 0 : fl[idx][j] = ugp[idx] * u[j];
877 [ - - ]: 0 : if (i == j) fl[idx][j] += p;
878 : : }
879 : : }
880 : :
881 : : // energy flux
882 : : auto idx = multispecies::energyIdx(nspec,0);
883 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j) {
884 : 0 : fl[idx][j] = u[j] * (ugp[idx] + p);
885 : : }
886 : :
887 : 0 : return fl;
888 : : }
889 : :
890 : : //! \brief Boundary state function providing the left and right state of a
891 : : //! face at Dirichlet boundaries
892 : : //! \param[in] ncomp Number of scalar components in this PDE system
893 : : //! \param[in] mat_blk EOS material block
894 : : //! \param[in] ul Left (domain-internal) state
895 : : //! \param[in] x X-coordinate at which to compute the states
896 : : //! \param[in] y Y-coordinate at which to compute the states
897 : : //! \param[in] z Z-coordinate at which to compute the states
898 : : //! \param[in] t Physical time
899 : : //! \return Left and right states for all scalar components in this PDE
900 : : //! system
901 : : //! \note The function signature must follow tk::StateFn.
902 : : static tk::StateFn::result_type
903 : 2250 : dirichlet( ncomp_t ncomp,
904 : : const std::vector< EOS >& mat_blk,
905 : : const std::vector< tk::real >& ul, tk::real x, tk::real y,
906 : : tk::real z, tk::real t, const std::array< tk::real, 3 >& )
907 : : {
908 : 2250 : return {{ ul, Problem::initialize( ncomp, mat_blk, x, y, z, t ) }};
909 : : }
910 : :
911 : : // Other boundary condition types that do not depend on "Problem" should be
912 : : // added in BCFunctions.hpp
913 : : };
914 : :
915 : : } // dg::
916 : :
917 : : } // inciter::
918 : :
919 : : #endif // DGMultiSpecies_h
|