Branch data Line data Source code
1 : : // *****************************************************************************
2 : : /*!
3 : : \file src/PDE/MultiMat/FVMultiMat.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-material flow using finite volumes
9 : : \details This file implements calls to the physics operators governing
10 : : compressible multi-material flow (with velocity equilibrium) using finite
11 : : volume discretizations.
12 : : */
13 : : // *****************************************************************************
14 : : #ifndef FVMultiMat_h
15 : : #define FVMultiMat_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/MultiMatTerms.hpp"
37 : : #include "Integrate/Source.hpp"
38 : : #include "RiemannChoice.hpp"
39 : : #include "MultiMat/MultiMatIndexing.hpp"
40 : : #include "Reconstruction.hpp"
41 : : #include "Limiter.hpp"
42 : : #include "Problem/FieldOutput.hpp"
43 : : #include "Problem/BoxInitialization.hpp"
44 : : #include "MultiMat/BCFunctions.hpp"
45 : : #include "MultiMat/MiscMultiMatFns.hpp"
46 : :
47 : : namespace inciter {
48 : :
49 : : extern ctr::InputDeck g_inputdeck;
50 : :
51 : : namespace fv {
52 : :
53 : : //! \brief MultiMat used polymorphically with tk::FVPDE
54 : : //! \details The template arguments specify policies and are used to configure
55 : : //! the behavior of the class. The policies are:
56 : : //! - Physics - physics configuration, see PDE/MultiMat/Physics.h
57 : : //! - Problem - problem configuration, see PDE/MultiMat/Problem.h
58 : : //! \note The default physics is Euler, set in inciter::deck::check_multimat()
59 : : template< class Physics, class Problem >
60 : : class MultiMat {
61 : :
62 : : private:
63 : : using eq = tag::multimat;
64 : :
65 : : public:
66 : : //! Constructor
67 : 81 : explicit MultiMat() :
68 : : m_physics(),
69 : 81 : m_ncomp( g_inputdeck.get< tag::ncomp >() ),
70 : : m_riemann( multimatRiemannSolver(
71 : 81 : g_inputdeck.get< tag::flux >() ) )
72 : : {
73 : : // associate boundary condition configurations with state functions
74 [ + - ][ + - ]: 1215 : brigand::for_each< ctr::bclist::Keys >( ConfigBC( m_bc,
[ + - ][ + - ]
[ + - ][ + - ]
[ + - ][ + - ]
[ + - ][ + - ]
[ + - ][ + - ]
[ + - ][ + - ]
[ + - ][ + - ]
[ + - ][ + + ]
[ + + ][ - - ]
[ - - ]
75 : : // BC State functions
76 : : { dirichlet
77 : : , symmetry
78 : : , invalidBC // Inlet BC not implemented
79 : : , invalidBC // Outlet BC not implemented
80 : : , farfield
81 : : , extrapolate
82 : : , noslipwall },
83 : : // BC Gradient functions
84 : : { noOpGrad
85 : : , symmetryGrad
86 : : , noOpGrad
87 : : , noOpGrad
88 : : , noOpGrad
89 : : , noOpGrad
90 : : , noOpGrad }
91 : : ) );
92 : :
93 : : // EoS initialization
94 [ + - ]: 81 : initializeMaterialEoS( m_mat_blk );
95 : 81 : }
96 : :
97 : : //! Find the number of primitive quantities required for this PDE system
98 : : //! \return The number of primitive quantities required to be stored for
99 : : //! this PDE system
100 : 16776 : std::size_t nprim() const
101 : : {
102 : 16776 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
103 : : // multimat needs individual material pressures and velocities currently
104 : 16776 : return (nmat+3);
105 : : }
106 : :
107 : : //! Find the number of materials set up for this PDE system
108 : : //! \return The number of materials set up for this PDE system
109 : 0 : std::size_t nmat() const
110 : : {
111 : 0 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
112 : 0 : return nmat;
113 : : }
114 : :
115 : : //! Determine elements that lie inside the user-defined IC box
116 : : //! \param[in] geoElem Element geometry array
117 : : //! \param[in] nielem Number of internal elements
118 : : //! \param[in,out] inbox List of nodes at which box user ICs are set for
119 : : //! each IC box
120 : 298 : void IcBoxElems( const tk::Fields& geoElem,
121 : : std::size_t nielem,
122 : : std::vector< std::unordered_set< std::size_t > >& inbox ) const
123 : : {
124 : 298 : tk::BoxElems< eq >(geoElem, nielem, inbox);
125 : 298 : }
126 : :
127 : : //! Initalize the compressible flow equations, prepare for time integration
128 : : //! \param[in] L Block diagonal mass matrix
129 : : //! \param[in] inpoel Element-node connectivity
130 : : //! \param[in] coord Array of nodal coordinates
131 : : //! \param[in] inbox List of elements at which box user ICs are set for
132 : : //! each IC box
133 : : //! \param[in] elemblkid Element ids associated with mesh block ids where
134 : : //! user ICs are set
135 : : //! \param[in,out] unk Array of unknowns
136 : : //! \param[in] t Physical time
137 : : //! \param[in] nielem Number of internal elements
138 : 298 : void initialize( const tk::Fields& L,
139 : : const std::vector< std::size_t >& inpoel,
140 : : const tk::UnsMesh::Coords& coord,
141 : : const std::vector< std::unordered_set< std::size_t > >& inbox,
142 : : const std::unordered_map< std::size_t, std::set< std::size_t > >&
143 : : elemblkid,
144 : : tk::Fields& unk,
145 : : tk::real t,
146 : : const std::size_t nielem ) const
147 : : {
148 [ + - ][ + - ]: 298 : tk::initialize( m_ncomp, m_mat_blk, L, inpoel, coord,
149 : : Problem::initialize, unk, t, nielem );
150 : :
151 : 298 : const auto rdof = g_inputdeck.get< tag::rdof >();
152 : 298 : const auto& ic = g_inputdeck.get< tag::ic >();
153 : 298 : const auto& icbox = ic.get< tag::box >();
154 : 298 : const auto& icmbk = ic.get< tag::meshblock >();
155 : :
156 : 298 : const auto& bgpre = ic.get< tag::pressure >();
157 : 298 : const auto& bgtemp = ic.get< tag::temperature >();
158 : :
159 : : // Set initial conditions inside user-defined IC boxes and mesh blocks
160 [ + - ]: 596 : std::vector< tk::real > s(m_ncomp, 0.0);
161 [ + + ]: 17394 : for (std::size_t e=0; e<nielem; ++e) {
162 : : // inside user-defined box
163 [ - + ]: 17096 : if (!icbox.empty()) {
164 : 0 : std::size_t bcnt = 0;
165 [ - - ]: 0 : for (const auto& b : icbox) { // for all boxes
166 [ - - ][ - - ]: 0 : if (inbox.size() > bcnt && inbox[bcnt].find(e) != inbox[bcnt].end())
[ - - ][ - - ]
167 : : {
168 [ - - ]: 0 : std::vector< tk::real > box
169 : 0 : { b.template get< tag::xmin >(), b.template get< tag::xmax >(),
170 : 0 : b.template get< tag::ymin >(), b.template get< tag::ymax >(),
171 : 0 : b.template get< tag::zmin >(), b.template get< tag::zmax >() };
172 : 0 : auto V_ex = (box[1]-box[0]) * (box[3]-box[2]) * (box[5]-box[4]);
173 [ - - ]: 0 : for (std::size_t c=0; c<m_ncomp; ++c) {
174 : 0 : auto mark = c*rdof;
175 [ - - ]: 0 : s[c] = unk(e,mark);
176 : : // set high-order DOFs to zero
177 [ - - ]: 0 : for (std::size_t i=1; i<rdof; ++i)
178 [ - - ]: 0 : unk(e,mark+i) = 0.0;
179 : : }
180 [ - - ]: 0 : initializeBox<ctr::boxList>( m_mat_blk, V_ex, t, b, bgpre,
181 : : bgtemp, s );
182 : : // store box-initialization in solution vector
183 [ - - ]: 0 : for (std::size_t c=0; c<m_ncomp; ++c) {
184 : 0 : auto mark = c*rdof;
185 [ - - ]: 0 : unk(e,mark) = s[c];
186 : : }
187 : : }
188 : 0 : ++bcnt;
189 : : }
190 : : }
191 : :
192 : : // inside user-specified mesh blocks
193 [ + + ]: 18020 : for (const auto& b : icmbk) { // for all blocks
194 : 924 : auto blid = b.get< tag::blockid >();
195 : 924 : auto V_ex = b.get< tag::volume >();
196 [ + - ][ + - ]: 924 : if (elemblkid.find(blid) != elemblkid.end()) {
197 [ + - ]: 924 : const auto& elset = tk::cref_find(elemblkid, blid);
198 [ + - ][ + + ]: 924 : if (elset.find(e) != elset.end()) {
199 [ + - ]: 104 : initializeBox<ctr::meshblockList>( m_mat_blk, V_ex, t, b,
200 : : bgpre, bgtemp, s );
201 : : // store initialization in solution vector
202 [ + + ]: 1040 : for (std::size_t c=0; c<m_ncomp; ++c) {
203 : 936 : auto mark = c*rdof;
204 [ + - ]: 936 : unk(e,mark) = s[c];
205 : : }
206 : : }
207 : : }
208 : : }
209 : : }
210 : 298 : }
211 : :
212 : : //! Compute the left hand side block-diagonal mass matrix
213 : : //! \param[in] geoElem Element geometry array
214 : : //! \param[in,out] l Block diagonal mass matrix
215 : 298 : void lhs( const tk::Fields& geoElem, tk::Fields& l ) const {
216 : 298 : const auto nelem = geoElem.nunk();
217 [ + + ]: 55093 : for (std::size_t e=0; e<nelem; ++e)
218 [ + + ]: 694308 : for (ncomp_t c=0; c<m_ncomp; ++c)
219 : 639513 : l(e, c) = geoElem(e,0);
220 : 298 : }
221 : :
222 : : //! Update the primitives for this PDE system
223 : : //! \param[in] unk Array of unknowns
224 : : //! \param[in,out] prim Array of primitives
225 : : //! \param[in] nielem Number of internal elements
226 : : //! \details This function computes and stores the dofs for primitive
227 : : //! quantities, which are required for obtaining reconstructed states used
228 : : //! in the Riemann solver. See /PDE/Riemann/AUSM.hpp, where the
229 : : //! normal velocity for advection is calculated from independently
230 : : //! reconstructed velocities.
231 : 1966 : void updatePrimitives( const tk::Fields& unk,
232 : : tk::Fields& prim,
233 : : std::size_t nielem ) const
234 : : {
235 : 1966 : const auto rdof = g_inputdeck.get< tag::rdof >();
236 : 1966 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
237 : :
238 [ - + ][ - - ]: 1966 : Assert( unk.nunk() == prim.nunk(), "Number of unknowns in solution "
[ - - ][ - - ]
239 : : "vector and primitive vector at recent time step incorrect" );
240 [ - + ][ - - ]: 1966 : Assert( unk.nprop() == rdof*m_ncomp, "Number of components in solution "
[ - - ][ - - ]
[ - - ]
241 : : "vector must equal "+ std::to_string(rdof*m_ncomp) );
242 [ - + ][ - - ]: 1966 : Assert( prim.nprop() == rdof*nprim(), "Number of components in vector of "
[ - - ][ - - ]
[ - - ]
243 : : "primitive quantities must equal "+ std::to_string(rdof*nprim()) );
244 : :
245 [ + + ]: 315622 : for (std::size_t e=0; e<nielem; ++e)
246 : : {
247 : : // cell-average bulk density
248 : 313656 : tk::real rhob(0.0);
249 [ + + ]: 1006668 : for (std::size_t k=0; k<nmat; ++k)
250 : : {
251 [ + - ]: 693012 : rhob += unk(e, densityDofIdx(nmat, k, rdof, 0));
252 : : }
253 : :
254 : : // cell-average velocity
255 : : std::array< tk::real, 3 >
256 [ + - ]: 313656 : vel{{ unk(e, momentumDofIdx(nmat, 0, rdof, 0))/rhob,
257 [ + - ]: 313656 : unk(e, momentumDofIdx(nmat, 1, rdof, 0))/rhob,
258 [ + - ]: 313656 : unk(e, momentumDofIdx(nmat, 2, rdof, 0))/rhob }};
259 : :
260 [ + + ]: 1254624 : for (std::size_t idir=0; idir<3; ++idir)
261 : : {
262 [ + - ]: 940968 : prim(e, velocityDofIdx(nmat, idir, rdof, 0)) = vel[idir];
263 [ + + ]: 3763872 : for (std::size_t idof=1; idof<rdof; ++idof)
264 [ + - ]: 2822904 : prim(e, velocityDofIdx(nmat, idir, rdof, idof)) = 0.0;
265 : : }
266 : :
267 : : // cell-average material pressure
268 [ + + ]: 1006668 : for (std::size_t k=0; k<nmat; ++k)
269 : : {
270 [ + - ]: 693012 : tk::real arhomat = unk(e, densityDofIdx(nmat, k, rdof, 0));
271 [ + - ]: 693012 : tk::real arhoemat = unk(e, energyDofIdx(nmat, k, rdof, 0));
272 [ + - ]: 693012 : tk::real alphamat = unk(e, volfracDofIdx(nmat, k, rdof, 0));
273 [ + - ][ + - ]: 693012 : auto gmat = getDeformGrad(nmat, k, unk.extract(e));
274 [ + - ]: 693012 : prim(e, pressureDofIdx(nmat, k, rdof, 0)) =
275 [ + - ]: 693012 : m_mat_blk[k].compute< EOS::pressure >( arhomat, vel[0], vel[1],
276 : : vel[2], arhoemat, alphamat, k, gmat );
277 [ + - ]: 693012 : prim(e, pressureDofIdx(nmat, k, rdof, 0)) =
278 : 693012 : constrain_pressure( m_mat_blk,
279 [ + - ][ + - ]: 693012 : prim(e, pressureDofIdx(nmat, k, rdof, 0)), arhomat, alphamat, k);
280 [ + + ]: 2772048 : for (std::size_t idof=1; idof<rdof; ++idof)
281 [ + - ]: 2079036 : prim(e, pressureDofIdx(nmat, k, rdof, idof)) = 0.0;
282 : : }
283 : : }
284 : 1966 : }
285 : :
286 : : //! Clean up the state of trace materials for this PDE system
287 : : //! \param[in] t Physical time
288 : : //! \param[in] geoElem Element geometry array
289 : : //! \param[in,out] unk Array of unknowns
290 : : //! \param[in,out] prim Array of primitives
291 : : //! \param[in] nielem Number of internal elements
292 : : //! \details This function cleans up the state of materials present in trace
293 : : //! quantities in each cell. Specifically, the state of materials with
294 : : //! very low volume-fractions in a cell is replaced by the state of the
295 : : //! material which is present in the largest quantity in that cell. This
296 : : //! becomes necessary when shocks pass through cells which contain a very
297 : : //! small amount of material. The state of that tiny material might
298 : : //! become unphysical and cause solution to diverge; thus requiring such
299 : : //! a "reset".
300 : 1668 : void cleanTraceMaterial( tk::real t,
301 : : const tk::Fields& geoElem,
302 : : tk::Fields& unk,
303 : : tk::Fields& prim,
304 : : std::size_t nielem ) const
305 : : {
306 : 1668 : [[maybe_unused]] const auto rdof = g_inputdeck.get< tag::rdof >();
307 : 1668 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
308 : :
309 [ - + ][ - - ]: 1668 : Assert( unk.nunk() == prim.nunk(), "Number of unknowns in solution "
[ - - ][ - - ]
310 : : "vector and primitive vector at recent time step incorrect" );
311 [ - + ][ - - ]: 1668 : Assert( unk.nprop() == rdof*m_ncomp, "Number of components in solution "
[ - - ][ - - ]
[ - - ]
312 : : "vector must equal "+ std::to_string(rdof*m_ncomp) );
313 [ - + ][ - - ]: 1668 : Assert( prim.nprop() == rdof*nprim(), "Number of components in vector of "
[ - - ][ - - ]
[ - - ]
314 : : "primitive quantities must equal "+ std::to_string(rdof*nprim()) );
315 [ - + ][ - - ]: 1668 : Assert( (g_inputdeck.get< tag::ndof >()) <= 4, "High-order "
[ - - ][ - - ]
316 : : "discretizations not set up for multimat cleanTraceMaterial()" );
317 : :
318 : 1668 : auto neg_density = cleanTraceMultiMat(t, nielem, m_mat_blk, geoElem, nmat,
319 : : unk, prim);
320 : :
321 [ - + ][ - - ]: 1668 : if (neg_density) Throw("Negative partial density.");
[ - - ][ - - ]
322 : 1668 : }
323 : :
324 : : //! Reconstruct second-order solution from first-order
325 : : //! \param[in] geoElem Element geometry array
326 : : //! \param[in] fd Face connectivity and boundary conditions object
327 : : //! \param[in] esup Elements-surrounding-nodes connectivity
328 : : //! \param[in] inpoel Element-node connectivity
329 : : //! \param[in] coord Array of nodal coordinates
330 : : //! \param[in,out] U Solution vector at recent time step
331 : : //! \param[in,out] P Vector of primitives at recent time step
332 : 1668 : void reconstruct( const tk::Fields& geoElem,
333 : : const inciter::FaceData& fd,
334 : : const std::map< std::size_t, std::vector< std::size_t > >&
335 : : esup,
336 : : const std::vector< std::size_t >& inpoel,
337 : : const tk::UnsMesh::Coords& coord,
338 : : tk::Fields& U,
339 : : tk::Fields& P ) const
340 : : {
341 : 1668 : const auto rdof = g_inputdeck.get< tag::rdof >();
342 : 1668 : const auto nelem = fd.Esuel().size()/4;
343 : 1668 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
344 : :
345 [ - + ][ - - ]: 1668 : Assert( U.nprop() == rdof*m_ncomp, "Number of components in solution "
[ - - ][ - - ]
[ - - ]
346 : : "vector must equal "+ std::to_string(rdof*m_ncomp) );
347 : :
348 : : //----- reconstruction of conserved quantities -----
349 : : //--------------------------------------------------
350 : : // specify how many variables need to be reconstructed
351 : 3336 : std::vector< std::size_t > vars;
352 [ + + ]: 6172 : for (std::size_t k=0; k<nmat; ++k) {
353 [ + - ]: 4504 : vars.push_back(volfracIdx(nmat,k));
354 [ + - ]: 4504 : vars.push_back(densityIdx(nmat,k));
355 : : }
356 : :
357 [ + + ]: 298228 : for (std::size_t e=0; e<nelem; ++e)
358 : : {
359 : : // 1. solve 3x3 least-squares system
360 : : // Reconstruct second-order dofs of volume-fractions in Taylor space
361 : : // using nodal-stencils, for a good interface-normal estimate
362 [ + - ]: 296560 : tk::recoLeastSqExtStencil( rdof, e, esup, inpoel, geoElem, U, vars );
363 : :
364 : : // 2. transform reconstructed derivatives to Dubiner dofs
365 [ + - ]: 296560 : tk::transform_P0P1(rdof, e, inpoel, coord, U, vars);
366 : : }
367 : :
368 : : //----- reconstruction of primitive quantities -----
369 : : //--------------------------------------------------
370 : : // For multimat, conserved and primitive quantities are reconstructed
371 : : // separately.
372 : 1668 : vars.clear();
373 [ + + ][ + - ]: 11176 : for (std::size_t c=0; c<nprim(); ++c) vars.push_back(c);
374 [ + + ]: 298228 : for (std::size_t e=0; e<nelem; ++e)
375 : : {
376 : : // 1.
377 : : // Reconstruct second-order dofs of volume-fractions in Taylor space
378 : : // using nodal-stencils, for a good interface-normal estimate
379 [ + - ]: 296560 : tk::recoLeastSqExtStencil( rdof, e, esup, inpoel, geoElem, P, vars );
380 : :
381 : : // 2.
382 [ + - ]: 296560 : tk::transform_P0P1(rdof, e, inpoel, coord, P, vars );
383 : : }
384 : 1668 : }
385 : :
386 : : //! Limit second-order solution, and primitive quantities separately
387 : : //! \param[in] geoFace Face geometry array
388 : : //! \param[in] fd Face connectivity and boundary conditions object
389 : : //! \param[in] esup Elements-surrounding-nodes connectivity
390 : : //! \param[in] inpoel Element-node connectivity
391 : : //! \param[in] coord Array of nodal coordinates
392 : : //! \param[in] srcFlag Whether the energy source was added
393 : : //! \param[in,out] U Solution vector at recent time step
394 : : //! \param[in,out] P Vector of primitives at recent time step
395 : 1668 : void limit( const tk::Fields& geoFace,
396 : : const inciter::FaceData& fd,
397 : : const std::map< std::size_t, std::vector< std::size_t > >& esup,
398 : : const std::vector< std::size_t >& inpoel,
399 : : const tk::UnsMesh::Coords& coord,
400 : : const std::vector< int >& srcFlag,
401 : : tk::Fields& U,
402 : : tk::Fields& P ) const
403 : : {
404 [ - + ][ - - ]: 1668 : Assert( U.nunk() == P.nunk(), "Number of unknowns in solution "
[ - - ][ - - ]
405 : : "vector and primitive vector at recent time step incorrect" );
406 : :
407 : 1668 : const auto limiter = g_inputdeck.get< tag::limiter >();
408 : 1668 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
409 : : const auto& solidx = g_inputdeck.get<
410 : 1668 : tag::matidxmap, tag::solidx >();
411 : :
412 : : // limit vectors of conserved and primitive quantities
413 [ + - ]: 1668 : if (limiter == ctr::LimiterType::VERTEXBASEDP1)
414 : : {
415 : 1668 : VertexBasedMultiMat_FV( esup, inpoel, fd.Esuel().size()/4,
416 : : coord, srcFlag, solidx, U, P, nmat );
417 : 1668 : PositivityPreservingMultiMat_FV( inpoel, fd.Esuel().size()/4, nmat,
418 : 1668 : m_mat_blk, coord, geoFace, U, P );
419 : : }
420 [ - - ]: 0 : else if (limiter != ctr::LimiterType::NOLIMITER)
421 : : {
422 [ - - ][ - - ]: 0 : Throw("Limiter type not configured for multimat.");
[ - - ]
423 : : }
424 : 1668 : }
425 : :
426 : : //! Apply CPL to the conservative variable solution for this PDE system
427 : : //! \param[in] prim Array of primitive variables
428 : : //! \param[in] geoElem Element geometry array
429 : : //! \param[in] inpoel Element-node connectivity
430 : : //! \param[in] coord Array of nodal coordinates
431 : : //! \param[in,out] unk Array of conservative variables
432 : : //! \param[in] nielem Number of internal elements
433 : : //! \details This function applies CPL to obtain consistent dofs for
434 : : //! conservative quantities based on the limited primitive quantities.
435 : : //! See Pandare et al. (2023). On the Design of Stable,
436 : : //! Consistent, and Conservative High-Order Methods for Multi-Material
437 : : //! Hydrodynamics. J Comp Phys, 112313.
438 : : void CPL( const tk::Fields& prim,
439 : : const tk::Fields& geoElem,
440 : : const std::vector< std::size_t >& inpoel,
441 : : const tk::UnsMesh::Coords& coord,
442 : : tk::Fields& unk,
443 : : std::size_t nielem ) const
444 : : {
445 : : [[maybe_unused]] const auto rdof = g_inputdeck.get< tag::rdof >();
446 : : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
447 : :
448 : : Assert( unk.nunk() == prim.nunk(), "Number of unknowns in solution "
449 : : "vector and primitive vector at recent time step incorrect" );
450 : : Assert( unk.nprop() == rdof*m_ncomp, "Number of components in solution "
451 : : "vector must equal "+ std::to_string(rdof*m_ncomp) );
452 : : Assert( prim.nprop() == rdof*nprim(), "Number of components in vector of "
453 : : "primitive quantities must equal "+ std::to_string(rdof*nprim()) );
454 : :
455 : : correctLimConservMultiMat(nielem, m_mat_blk, nmat, inpoel,
456 : : coord, geoElem, prim, unk);
457 : : }
458 : :
459 : : //! Compute right hand side
460 : : //! \param[in] t Physical time
461 : : //! \param[in] geoFace Face geometry array
462 : : //! \param[in] geoElem Element geometry array
463 : : //! \param[in] fd Face connectivity and boundary conditions object
464 : : //! \param[in] inpoel Element-node connectivity
465 : : //! \param[in] coord Array of nodal coordinates
466 : : //! \param[in] elemblkid Element ids associated with mesh block ids where
467 : : //! user ICs are set
468 : : //! \param[in] U Solution vector at recent time step
469 : : //! \param[in] P Primitive vector at recent time step
470 : : //! \param[in,out] R Right-hand side vector computed
471 : : //! \param[in,out] srcFlag Whether the energy source was added
472 : 1668 : void rhs( tk::real t,
473 : : const tk::Fields& geoFace,
474 : : const tk::Fields& geoElem,
475 : : const inciter::FaceData& fd,
476 : : const std::vector< std::size_t >& inpoel,
477 : : const tk::UnsMesh::Coords& coord,
478 : : const std::unordered_map< std::size_t, std::set< std::size_t > >&
479 : : elemblkid,
480 : : const tk::Fields& U,
481 : : const tk::Fields& P,
482 : : tk::Fields& R,
483 : : std::vector< int >& srcFlag ) const
484 : : {
485 : 1668 : const auto rdof = g_inputdeck.get< tag::rdof >();
486 : 1668 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
487 : 1668 : const auto intsharp =
488 : 1668 : g_inputdeck.get< tag::multimat, tag::intsharp >();
489 : 1668 : auto viscous = g_inputdeck.get< tag::multimat, tag::viscous >();
490 : :
491 : 1668 : const auto nelem = fd.Esuel().size()/4;
492 : :
493 [ - + ][ - - ]: 1668 : Assert( U.nunk() == P.nunk(), "Number of unknowns in solution "
[ - - ][ - - ]
494 : : "vector and primitive vector at recent time step incorrect" );
495 [ - + ][ - - ]: 1668 : Assert( U.nunk() == R.nunk(), "Number of unknowns in solution "
[ - - ][ - - ]
496 : : "vector and right-hand side at recent time step incorrect" );
497 [ - + ][ - - ]: 1668 : Assert( U.nprop() == rdof*m_ncomp, "Number of components in solution "
[ - - ][ - - ]
[ - - ]
498 : : "vector must equal "+ std::to_string(rdof*m_ncomp) );
499 [ - + ][ - - ]: 1668 : Assert( P.nprop() == rdof*nprim(), "Number of components in primitive "
[ - - ][ - - ]
[ - - ]
500 : : "vector must equal "+ std::to_string(rdof*nprim()) );
501 [ - + ][ - - ]: 1668 : Assert( fd.Inpofa().size()/3 == fd.Esuf().size()/2,
[ - - ][ - - ]
502 : : "Mismatch in inpofa size" );
503 : :
504 : : // set rhs to zero
505 [ + - ]: 1668 : R.fill(0.0);
506 : :
507 : : // configure a no-op lambda for prescribed velocity
508 : 722616 : auto velfn = []( ncomp_t, tk::real, tk::real, tk::real, tk::real ){
509 : 722616 : return tk::VelFn::result_type(); };
510 : :
511 : : // compute internal surface flux (including non-conservative) integrals
512 [ + - ]: 3336 : tk::surfIntFV( nmat, m_mat_blk, t, rdof, inpoel,
513 [ + - ]: 1668 : coord, fd, geoFace, geoElem, m_riemann, velfn, U, P,
514 : : srcFlag, R, intsharp );
515 : : // compute internal surface viscous flux integrals
516 [ + - ]: 3336 : tk::Fields T( U.nunk(), rdof*nmat );
517 [ - + ]: 1668 : if (viscous) {
518 [ - - ]: 0 : computeTemperaturesFV( m_mat_blk, nmat, inpoel, coord, geoElem,
519 : : U, P, srcFlag, T );
520 : :
521 [ - - ]: 0 : tk::surfIntViscousFV( nmat, m_mat_blk, rdof, inpoel,
522 : : coord, fd, geoFace, geoElem, U, P, T,
523 : : srcFlag, R, intsharp );
524 : : }
525 : :
526 : : // compute boundary surface flux (including non-conservative) integrals
527 [ + + ]: 13344 : for (const auto& b : m_bc) {
528 [ + - ][ + - ]: 23352 : tk::bndSurfIntFV( nmat, m_mat_blk, rdof, std::get<0>(b),
529 : 11676 : fd, geoFace, geoElem, inpoel, coord, t, m_riemann,
530 : 11676 : velfn, std::get<1>(b), U, P, srcFlag, R, intsharp );
531 [ - + ]: 11676 : if (viscous)
532 [ - - ]: 0 : tk::bndSurfIntViscousFV( nmat, m_mat_blk, rdof, std::get<0>(b),
533 : : fd, geoFace, geoElem, inpoel, coord, t,
534 : 0 : std::get<1>(b), std::get<2>(b), U, P, T,
535 : : srcFlag, R, intsharp );
536 : : }
537 : :
538 : : // compute optional source term
539 [ + - ][ + - ]: 1668 : tk::srcIntFV( m_mat_blk, t, fd.Esuel().size()/4,
540 : : geoElem, Problem::src, R, nmat );
541 : :
542 : : // compute finite pressure relaxation terms
543 [ + + ]: 1668 : if (g_inputdeck.get< tag::multimat, tag::prelax >())
544 : : {
545 : 1618 : const auto ct = g_inputdeck.get< tag::multimat,
546 : 1618 : tag::prelax_timescale >();
547 [ + - ]: 1618 : tk::pressureRelaxationIntFV( nmat, m_mat_blk, rdof,
548 : : nelem, inpoel, coord, geoElem, U, P, ct,
549 : : R );
550 : : }
551 : :
552 : : // compute external (energy) sources
553 [ + - ]: 1668 : m_physics.physSrc(nmat, t, geoElem, elemblkid, R, srcFlag);
554 : 1668 : }
555 : :
556 : : //! Compute the minimum time step size
557 : : // //! \param[in] fd Face connectivity and boundary conditions object
558 : : // //! \param[in] geoFace Face geometry array
559 : : //! \param[in] geoElem Element geometry array
560 : : //! \param[in] U Solution vector at recent time step
561 : : //! \param[in] P Vector of primitive quantities at recent time step
562 : : //! \param[in] nielem Number of internal elements
563 : : //! \param[in] srcFlag Whether the energy source was added
564 : : //! \param[in,out] local_dte Time step size for each element (for local
565 : : //! time stepping)
566 : : //! \return Minimum time step size
567 : : //! \details The allowable dt is calculated by looking at the maximum
568 : : //! wave-speed in elements surrounding each face, times the area of that
569 : : //! face. Once the maximum of this quantity over the mesh is determined,
570 : : //! the volume of each cell is divided by this quantity. A minimum of this
571 : : //! ratio is found over the entire mesh, which gives the allowable dt.
572 : 250 : tk::real dt( const inciter::FaceData& /*fd*/,
573 : : const tk::Fields& /*geoFace*/,
574 : : const tk::Fields& geoElem,
575 : : const tk::Fields& U,
576 : : const tk::Fields& P,
577 : : const std::size_t nielem,
578 : : const std::vector< int >& srcFlag,
579 : : std::vector< tk::real >& local_dte ) const
580 : : {
581 : 250 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
582 : 250 : auto viscous = g_inputdeck.get< tag::multimat, tag::viscous >();
583 : :
584 : : // obtain dt restrictions from all physics
585 [ + - ]: 250 : auto dt_e = timeStepSizeMultiMatFV(m_mat_blk, geoElem, nielem, nmat, U,
586 : : P, local_dte);
587 [ - + ]: 250 : if (viscous)
588 [ - - ]: 0 : dt_e = std::min(dt_e, timeStepSizeViscousFV(geoElem, nielem, nmat, U));
589 [ + - ]: 250 : auto dt_p = m_physics.dtRestriction(geoElem, nielem, srcFlag);
590 : :
591 : 250 : return std::min(dt_e, dt_p);
592 : : }
593 : :
594 : : //! Extract the velocity field at cell nodes. Currently unused.
595 : : //! \param[in] U Solution vector at recent time step
596 : : //! \param[in] N Element node indices
597 : : //! \return Array of the four values of the velocity field
598 : : std::array< std::array< tk::real, 4 >, 3 >
599 : : velocity( const tk::Fields& U,
600 : : const std::array< std::vector< tk::real >, 3 >&,
601 : : const std::array< std::size_t, 4 >& N ) const
602 : : {
603 : : const auto rdof = g_inputdeck.get< tag::rdof >();
604 : : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
605 : :
606 : : std::array< std::array< tk::real, 4 >, 3 > v;
607 : : v[0] = U.extract( momentumDofIdx(nmat, 0, rdof, 0), N );
608 : : v[1] = U.extract( momentumDofIdx(nmat, 1, rdof, 0), N );
609 : : v[2] = U.extract( momentumDofIdx(nmat, 2, rdof, 0), N );
610 : :
611 : : std::vector< std::array< tk::real, 4 > > ar;
612 : : ar.resize(nmat);
613 : : for (std::size_t k=0; k<nmat; ++k)
614 : : ar[k] = U.extract( densityDofIdx(nmat, k, rdof, 0), N );
615 : :
616 : : std::array< tk::real, 4 > r{{ 0.0, 0.0, 0.0, 0.0 }};
617 : : for (std::size_t i=0; i<r.size(); ++i) {
618 : : for (std::size_t k=0; k<nmat; ++k)
619 : : r[i] += ar[k][i];
620 : : }
621 : :
622 : : std::transform( r.begin(), r.end(), v[0].begin(), v[0].begin(),
623 : : []( tk::real s, tk::real& d ){ return d /= s; } );
624 : : std::transform( r.begin(), r.end(), v[1].begin(), v[1].begin(),
625 : : []( tk::real s, tk::real& d ){ return d /= s; } );
626 : : std::transform( r.begin(), r.end(), v[2].begin(), v[2].begin(),
627 : : []( tk::real s, tk::real& d ){ return d /= s; } );
628 : : return v;
629 : : }
630 : :
631 : : //! Return a map that associates user-specified strings to functions
632 : : //! \return Map that associates user-specified strings to functions that
633 : : //! compute relevant quantities to be output to file
634 : 56 : std::map< std::string, tk::GetVarFn > OutVarFn() const
635 : 56 : { return MultiMatOutVarFn(); }
636 : :
637 : : //! Return analytic field names to be output to file
638 : : //! \return Vector of strings labelling analytic fields output in file
639 : 0 : std::vector< std::string > analyticFieldNames() const {
640 : 0 : auto nmat = g_inputdeck.get< eq, tag::nmat >();
641 : :
642 : 0 : return MultiMatFieldNames(nmat);
643 : : }
644 : :
645 : : //! Return surface field names to be output to file
646 : : //! \return Vector of strings labelling surface fields output in file
647 : 28 : std::vector< std::string > surfNames() const
648 : 28 : { return MultiMatSurfNames(); }
649 : :
650 : : //! Return time history field names to be output to file
651 : : //! \return Vector of strings labelling time history fields output in file
652 : 0 : std::vector< std::string > histNames() const {
653 : 0 : return MultiMatHistNames();
654 : : }
655 : :
656 : : //! Return surface field output going to file
657 : : std::vector< std::vector< tk::real > >
658 : 28 : surfOutput( const inciter::FaceData& fd,
659 : : const tk::Fields& U,
660 : : const tk::Fields& P ) const
661 : : {
662 : 28 : const auto rdof = g_inputdeck.get< tag::rdof >();
663 : 28 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
664 : :
665 : 28 : return MultiMatSurfOutput( nmat, rdof, fd, U, P );
666 : : }
667 : :
668 : : //! Return time history field output evaluated at time history points
669 : : //! \param[in] h History point data
670 : : //! \param[in] inpoel Element-node connectivity
671 : : //! \param[in] coord Array of nodal coordinates
672 : : //! \param[in] U Array of unknowns
673 : : //! \param[in] P Array of primitive quantities
674 : : //! \return Vector of time history output of bulk flow quantities (density,
675 : : //! velocity, total energy, pressure, and volume fraction) evaluated at
676 : : //! time history points
677 : : std::vector< std::vector< tk::real > >
678 : 0 : histOutput( const std::vector< HistData >& h,
679 : : const std::vector< std::size_t >& inpoel,
680 : : const tk::UnsMesh::Coords& coord,
681 : : const tk::Fields& U,
682 : : const tk::Fields& P ) const
683 : : {
684 : 0 : const auto rdof = g_inputdeck.get< tag::rdof >();
685 : 0 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
686 : :
687 : 0 : const auto& x = coord[0];
688 : 0 : const auto& y = coord[1];
689 : 0 : const auto& z = coord[2];
690 : :
691 [ - - ]: 0 : std::vector< std::vector< tk::real > > Up(h.size());
692 : :
693 : 0 : std::size_t j = 0;
694 [ - - ]: 0 : for (const auto& p : h) {
695 : 0 : auto e = p.get< tag::elem >();
696 : 0 : auto chp = p.get< tag::coord >();
697 : :
698 : : // Evaluate inverse Jacobian
699 : 0 : std::array< std::array< tk::real, 3>, 4 > cp{{
700 : : {{ x[inpoel[4*e ]], y[inpoel[4*e ]], z[inpoel[4*e ]] }},
701 : 0 : {{ x[inpoel[4*e+1]], y[inpoel[4*e+1]], z[inpoel[4*e+1]] }},
702 : 0 : {{ x[inpoel[4*e+2]], y[inpoel[4*e+2]], z[inpoel[4*e+2]] }},
703 : 0 : {{ x[inpoel[4*e+3]], y[inpoel[4*e+3]], z[inpoel[4*e+3]] }} }};
704 : 0 : auto J = tk::inverseJacobian( cp[0], cp[1], cp[2], cp[3] );
705 : :
706 : : // evaluate solution at history-point
707 : 0 : std::array< tk::real, 3 > dc{{chp[0]-cp[0][0], chp[1]-cp[0][1],
708 : 0 : chp[2]-cp[0][2]}};
709 [ - - ]: 0 : auto B = tk::eval_basis(rdof, tk::dot(J[0],dc), tk::dot(J[1],dc),
710 : 0 : tk::dot(J[2],dc));
711 [ - - ]: 0 : auto uhp = eval_state(m_ncomp, rdof, rdof, e, U, B);
712 [ - - ]: 0 : auto php = eval_state(nprim(), rdof, rdof, e, P, B);
713 : :
714 : : // store solution in history output vector
715 [ - - ]: 0 : Up[j].resize(6+nmat, 0.0);
716 [ - - ]: 0 : for (std::size_t k=0; k<nmat; ++k) {
717 : 0 : Up[j][0] += uhp[densityIdx(nmat,k)];
718 : 0 : Up[j][4] += uhp[energyIdx(nmat,k)];
719 : 0 : Up[j][5] += php[pressureIdx(nmat,k)];
720 : 0 : Up[j][6+k] = uhp[volfracIdx(nmat,k)];
721 : : }
722 : 0 : Up[j][1] = php[velocityIdx(nmat,0)];
723 : 0 : Up[j][2] = php[velocityIdx(nmat,1)];
724 : 0 : Up[j][3] = php[velocityIdx(nmat,2)];
725 : 0 : ++j;
726 : : }
727 : :
728 : 0 : return Up;
729 : : }
730 : :
731 : : //! Return names of integral variables to be output to diagnostics file
732 : : //! \return Vector of strings labelling integral variables output
733 : 21 : std::vector< std::string > names() const
734 : : {
735 : 21 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
736 : 21 : return MultiMatDiagNames(nmat);
737 : : }
738 : :
739 : : //! Return analytic solution (if defined by Problem) at xi, yi, zi, t
740 : : //! \param[in] xi X-coordinate at which to evaluate the analytic solution
741 : : //! \param[in] yi Y-coordinate at which to evaluate the analytic solution
742 : : //! \param[in] zi Z-coordinate at which to evaluate the analytic solution
743 : : //! \param[in] t Physical time at which to evaluate the analytic solution
744 : : //! \return Vector of analytic solution at given location and time
745 : : std::vector< tk::real >
746 : 0 : analyticSolution( tk::real xi, tk::real yi, tk::real zi, tk::real t ) const
747 : 0 : { return Problem::analyticSolution( m_ncomp, m_mat_blk, xi, yi, zi, t ); }
748 : :
749 : : //! Return analytic solution for conserved variables
750 : : //! \param[in] xi X-coordinate at which to evaluate the analytic solution
751 : : //! \param[in] yi Y-coordinate at which to evaluate the analytic solution
752 : : //! \param[in] zi Z-coordinate at which to evaluate the analytic solution
753 : : //! \param[in] t Physical time at which to evaluate the analytic solution
754 : : //! \return Vector of analytic solution at given location and time
755 : : std::vector< tk::real >
756 : 0 : solution( tk::real xi, tk::real yi, tk::real zi, tk::real t ) const
757 : 0 : { return Problem::initialize( m_ncomp, m_mat_blk, xi, yi, zi, t ); }
758 : :
759 : : //! Return cell-averaged specific total energy for an element
760 : : //! \param[in] e Element id for which total energy is required
761 : : //! \param[in] unk Vector of conserved quantities
762 : : //! \return Cell-averaged specific total energy for given element
763 : 0 : tk::real sp_totalenergy(std::size_t e, const tk::Fields& unk) const
764 : : {
765 : 0 : const auto rdof = g_inputdeck.get< tag::rdof >();
766 : 0 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
767 : :
768 : 0 : tk::real sp_te(0.0);
769 : : // sum each material total energy
770 [ - - ]: 0 : for (std::size_t k=0; k<nmat; ++k) {
771 : 0 : sp_te += unk(e, energyDofIdx(nmat,k,rdof,0));
772 : : }
773 : 0 : return sp_te;
774 : : }
775 : :
776 : : //! Compute relevant sound speed for output
777 : : //! \param[in] nielem Number of internal elements
778 : : //! \param[in] U Solution vector at recent time step
779 : : //! \param[in] P Primitive vector at recent time step
780 : : //! \param[in,out] ss Sound speed vector
781 : 28 : void soundspeed(
782 : : std::size_t nielem,
783 : : const tk::Fields& U,
784 : : const tk::Fields& P,
785 : : std::vector< tk::real >& ss) const
786 : : {
787 [ - + ][ - - ]: 28 : Assert( ss.size() == nielem, "Size of sound speed vector incorrect " );
[ - - ][ - - ]
788 : :
789 : 28 : const auto ndof = g_inputdeck.get< tag::ndof >();
790 : 28 : const auto rdof = g_inputdeck.get< tag::rdof >();
791 : 28 : const auto use_mass_avg =
792 : 28 : g_inputdeck.get< tag::multimat, tag::dt_sos_massavg >();
793 : 28 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
794 : 28 : std::size_t ncomp = U.nprop()/rdof;
795 : 28 : std::size_t nprim = P.nprop()/rdof;
796 : :
797 [ + - ][ + - ]: 56 : std::vector< tk::real > ugp(ncomp, 0.0), pgp(nprim, 0.0);
798 : :
799 [ + + ]: 11636 : for (std::size_t e=0; e<nielem; ++e) {
800 : : // basis function at centroid
801 [ + - ]: 23216 : std::vector< tk::real > B(rdof, 0.0);
802 : 11608 : B[0] = 1.0;
803 : :
804 : : // get conserved quantities
805 [ + - ]: 11608 : ugp = eval_state(ncomp, rdof, ndof, e, U, B);
806 : : // get primitive quantities
807 [ + - ]: 11608 : pgp = eval_state(nprim, rdof, ndof, e, P, B);
808 : :
809 : : // acoustic speed (this should be consistent with time-step calculation)
810 : 11608 : ss[e] = 0.0;
811 : 11608 : tk::real mixtureDensity = 0.0;
812 [ + + ]: 34824 : for (std::size_t k=0; k<nmat; ++k)
813 : : {
814 [ - + ]: 23216 : if (use_mass_avg > 0)
815 : : {
816 : : // mass averaging SoS
817 : 0 : ss[e] += ugp[densityIdx(nmat,k)]*
818 [ - - ]: 0 : m_mat_blk[k].compute< EOS::soundspeed >(
819 : : ugp[densityIdx(nmat, k)], pgp[pressureIdx(nmat, k)],
820 : : ugp[volfracIdx(nmat, k)], k );
821 : :
822 : 0 : mixtureDensity += ugp[densityIdx(nmat,k)];
823 : : }
824 : : else
825 : : {
826 [ + + ]: 23216 : if (ugp[volfracIdx(nmat, k)] > 1.0e-04)
827 : : {
828 [ + - ]: 27458 : ss[e] = std::max( ss[e], m_mat_blk[k].compute< EOS::soundspeed >(
829 : : ugp[densityIdx(nmat, k)], pgp[pressureIdx(nmat, k)],
830 : 13729 : ugp[volfracIdx(nmat, k)], k ) );
831 : : }
832 : : }
833 : : }
834 [ - + ]: 11608 : if (use_mass_avg > 0) ss[e] /= mixtureDensity;
835 : : }
836 : 28 : }
837 : :
838 : : private:
839 : : //! Physics policy
840 : : const Physics m_physics;
841 : : //! Number of components in this PDE system
842 : : const ncomp_t m_ncomp;
843 : : //! Riemann solver
844 : : tk::RiemannFluxFn m_riemann;
845 : : //! BC configuration
846 : : BCStateFn m_bc;
847 : : //! EOS material block
848 : : std::vector< EOS > m_mat_blk;
849 : :
850 : : //! Evaluate conservative part of physical flux function for this PDE system
851 : : //! \param[in] ncomp Number of scalar components in this PDE system
852 : : //! \param[in] ugp Numerical solution at the Gauss point at which to
853 : : //! evaluate the flux
854 : : //! \return Flux vectors for all components in this PDE system
855 : : //! \note The function signature must follow tk::FluxFn
856 : : static tk::FluxFn::result_type
857 : : flux( ncomp_t ncomp,
858 : : const std::vector< EOS >& mat_blk,
859 : : const std::vector< tk::real >& ugp,
860 : : const std::vector< std::array< tk::real, 3 > >& )
861 : : {
862 : : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
863 : :
864 : : return tk::fluxTerms(ncomp, nmat, mat_blk, ugp);
865 : : }
866 : :
867 : : //! \brief Boundary state function providing the left and right state of a
868 : : //! face at Dirichlet boundaries
869 : : //! \param[in] ncomp Number of scalar components in this PDE system
870 : : //! \param[in] ul Left (domain-internal) state
871 : : //! \param[in] x X-coordinate at which to compute the states
872 : : //! \param[in] y Y-coordinate at which to compute the states
873 : : //! \param[in] z Z-coordinate at which to compute the states
874 : : //! \param[in] t Physical time
875 : : //! \return Left and right states for all scalar components in this PDE
876 : : //! system
877 : : //! \note The function signature must follow tk::StateFn. For multimat, the
878 : : //! left or right state is the vector of conserved quantities, followed by
879 : : //! the vector of primitive quantities appended to it.
880 : : static tk::StateFn::result_type
881 : 0 : dirichlet( ncomp_t ncomp,
882 : : const std::vector< EOS >& mat_blk,
883 : : const std::vector< tk::real >& ul, tk::real x, tk::real y,
884 : : tk::real z, tk::real t, const std::array< tk::real, 3 >& )
885 : : {
886 : 0 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
887 : :
888 [ - - ]: 0 : auto ur = Problem::initialize( ncomp, mat_blk, x, y, z, t );
889 [ - - ][ - - ]: 0 : Assert( ur.size() == ncomp, "Incorrect size for boundary state vector" );
[ - - ][ - - ]
890 : :
891 [ - - ]: 0 : ur.resize(ul.size());
892 : :
893 : 0 : tk::real rho(0.0);
894 [ - - ]: 0 : for (std::size_t k=0; k<nmat; ++k)
895 : 0 : rho += ur[densityIdx(nmat, k)];
896 : :
897 : : // get primitives in boundary state
898 : :
899 : : // velocity
900 : 0 : ur[ncomp+velocityIdx(nmat, 0)] = ur[momentumIdx(nmat, 0)] / rho;
901 : 0 : ur[ncomp+velocityIdx(nmat, 1)] = ur[momentumIdx(nmat, 1)] / rho;
902 : 0 : ur[ncomp+velocityIdx(nmat, 2)] = ur[momentumIdx(nmat, 2)] / rho;
903 : :
904 : : // material pressures
905 [ - - ]: 0 : for (std::size_t k=0; k<nmat; ++k)
906 : : {
907 [ - - ]: 0 : auto gk = getDeformGrad(nmat, k, ur);
908 [ - - ]: 0 : ur[ncomp+pressureIdx(nmat, k)] = mat_blk[k].compute< EOS::pressure >(
909 : 0 : ur[densityIdx(nmat, k)], ur[ncomp+velocityIdx(nmat, 0)],
910 : 0 : ur[ncomp+velocityIdx(nmat, 1)], ur[ncomp+velocityIdx(nmat, 2)],
911 : 0 : ur[energyIdx(nmat, k)], ur[volfracIdx(nmat, k)], k, gk );
912 : : }
913 : :
914 [ - - ][ - - ]: 0 : Assert( ur.size() == ncomp+nmat+3, "Incorrect size for appended "
[ - - ][ - - ]
915 : : "boundary state vector" );
916 : :
917 [ - - ]: 0 : return {{ std::move(ul), std::move(ur) }};
918 : : }
919 : :
920 : : // Other boundary condition types that do not depend on "Problem" should be
921 : : // added in BCFunctions.hpp
922 : : };
923 : :
924 : : } // fv::
925 : :
926 : : } // inciter::
927 : :
928 : : #endif // FVMultiMat_h
|