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