Branch data Line data Source code
1 : : // *****************************************************************************
2 : : /*!
3 : : \file src/PDE/MultiMat/DGMultiMat.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 discontinuous Galerkin
9 : : finite elements
10 : : \details This file implements calls to the physics operators governing
11 : : compressible multi-material flow (with velocity equilibrium) using
12 : : discontinuous Galerkin discretizations.
13 : : */
14 : : // *****************************************************************************
15 : : #ifndef DGMultiMat_h
16 : : #define DGMultiMat_h
17 : :
18 : : #include <cmath>
19 : : #include <algorithm>
20 : : #include <unordered_set>
21 : : #include <map>
22 : : #include <array>
23 : :
24 : : #include "Macro.hpp"
25 : : #include "Exception.hpp"
26 : : #include "Vector.hpp"
27 : : #include "ContainerUtil.hpp"
28 : : #include "UnsMesh.hpp"
29 : : #include "Inciter/InputDeck/InputDeck.hpp"
30 : : #include "Integrate/Basis.hpp"
31 : : #include "Integrate/Quadrature.hpp"
32 : : #include "Integrate/Initialize.hpp"
33 : : #include "Integrate/Mass.hpp"
34 : : #include "Integrate/Surface.hpp"
35 : : #include "Integrate/Boundary.hpp"
36 : : #include "Integrate/Volume.hpp"
37 : : #include "Integrate/MultiMatTerms.hpp"
38 : : #include "Integrate/Source.hpp"
39 : : #include "Integrate/SolidTerms.hpp"
40 : : #include "RiemannChoice.hpp"
41 : : #include "MultiMat/MultiMatIndexing.hpp"
42 : : #include "Reconstruction.hpp"
43 : : #include "Limiter.hpp"
44 : : #include "Problem/FieldOutput.hpp"
45 : : #include "Problem/BoxInitialization.hpp"
46 : : #include "PrefIndicator.hpp"
47 : : #include "MultiMat/BCFunctions.hpp"
48 : : #include "MultiMat/MiscMultiMatFns.hpp"
49 : : #include "EoS/GetMatProp.hpp"
50 : :
51 : : namespace inciter {
52 : :
53 : : extern ctr::InputDeck g_inputdeck;
54 : :
55 : : namespace dg {
56 : :
57 : : //! \brief MultiMat used polymorphically with tk::DGPDE
58 : : //! \details The template arguments specify policies and are used to configure
59 : : //! the behavior of the class. The policies are:
60 : : //! - Physics - physics configuration, see PDE/MultiMat/Physics.h
61 : : //! - Problem - problem configuration, see PDE/MultiMat/Problem.h
62 : : //! \note The default physics is Euler, set in inciter::deck::check_multimat()
63 : : template< class Physics, class Problem >
64 : : class MultiMat {
65 : :
66 : : private:
67 : : using eq = tag::multimat;
68 : :
69 : : public:
70 : : //! Constructor
71 : 118 : explicit MultiMat() :
72 : : m_ncomp( g_inputdeck.get< tag::ncomp >() ),
73 : : m_nprim(nprim()),
74 : : m_riemann( multimatRiemannSolver(
75 [ + - ]: 118 : g_inputdeck.get< tag::flux >() ) )
76 : : {
77 : : // associate boundary condition configurations with state functions
78 [ + - ][ + - ]: 1770 : 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 : : // Inlet BC has a different structure than above BCs, so it must be
98 : : // handled differently than with ConfigBC
99 [ + - ][ + - ]: 354 : ConfigInletBC(m_bc, inlet, zeroGrad);
[ + - ][ - - ]
100 : :
101 : : // Back pressure BC has a different structure than above BCs, so it must
102 : : // be handled differently than with ConfigBC
103 [ + - ][ + - ]: 236 : ConfigBackPressureBC(m_bc, back_pressure, noOpGrad);
[ - - ]
104 : :
105 : : // EoS initialization
106 [ + - ]: 118 : initializeMaterialEoS( m_mat_blk );
107 : 118 : }
108 : :
109 : : //! Find the number of primitive quantities required for this PDE system
110 : : //! \return The number of primitive quantities required to be stored for
111 : : //! this PDE system
112 : : std::size_t nprim() const
113 : : {
114 : 183 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
115 : : const auto& solidx = inciter::g_inputdeck.get<
116 : : tag::matidxmap, tag::solidx >();
117 : :
118 : : // individual material pressures and three velocity components
119 : 183 : std::size_t np(nmat+3);
120 : :
121 [ - - ][ + + ]: 633 : for (std::size_t k=0; k<nmat; ++k) {
[ + + ][ + + ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ + + ][ + + ]
[ + + ][ - - ]
122 [ - - ][ - + ]: 450 : if (solidx[k] > 0) {
[ - + ][ - + ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - + ][ - + ]
[ - + ][ - - ]
123 : : // individual material Cauchy stress tensor components
124 : 0 : np += 6;
125 : : }
126 : : }
127 : :
128 : : return np;
129 : : }
130 : :
131 : : //! Find the number of materials set up for this PDE system
132 : : //! \return The number of materials set up for this PDE system
133 : : std::size_t nmat() const
134 : : {
135 : 0 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
136 : : return nmat;
137 : : }
138 : :
139 : : //! Assign number of DOFs per equation in the PDE system
140 : : //! \param[in,out] numEqDof Array storing number of Dofs for each PDE
141 : : //! equation
142 : : void numEquationDofs(std::vector< std::size_t >& numEqDof) const
143 : : {
144 : : // all equation-dofs initialized to ndofs first
145 : : for (std::size_t i=0; i<m_ncomp; ++i) {
146 : : numEqDof.push_back(g_inputdeck.get< tag::ndof >());
147 : : }
148 : :
149 : : // volume fractions are P0Pm (ndof = 1) for multi-material simulations
150 : : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
151 : : if(nmat > 1)
152 : : for (std::size_t k=0; k<nmat; ++k)
153 : : numEqDof[volfracIdx(nmat, k)] = 1;
154 : : }
155 : :
156 : : //! Determine elements that lie inside the user-defined IC box
157 : : //! \param[in] geoElem Element geometry array
158 : : //! \param[in] nielem Number of internal elements
159 : : //! \param[in,out] inbox List of nodes at which box user ICs are set for
160 : : //! each IC box
161 : : void IcBoxElems( const tk::Fields& geoElem,
162 : : std::size_t nielem,
163 : : std::vector< std::unordered_set< std::size_t > >& inbox ) const
164 : : {
165 : 65 : tk::BoxElems< eq >(geoElem, nielem, inbox);
166 : : }
167 : :
168 : : //! Find how many 'stiff equations', which are the inverse
169 : : //! deformation equations because of plasticity
170 : : //! \return number of stiff equations
171 : : std::size_t nstiffeq() const
172 : : {
173 : : const auto& solidx = g_inputdeck.get< tag::matidxmap, tag::solidx >();
174 : 390 : std::size_t nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
175 : 390 : return 9*numSolids(nmat, solidx);
176 : : }
177 : :
178 : : //! Find how many 'non-stiff equations', which are the inverse
179 : : //! deformation equations because of plasticity
180 : : //! \return number of stiff equations
181 : : std::size_t nnonstiffeq() const
182 : : {
183 : 130 : return m_ncomp-nstiffeq();
184 : : }
185 : :
186 : : //! Locate the stiff equations.
187 : : //! \param[out] stiffEqIdx list with pointers to stiff equations
188 : 0 : void setStiffEqIdx( std::vector< std::size_t >& stiffEqIdx ) const
189 : : {
190 : 0 : stiffEqIdx.resize(nstiffeq(), 0);
191 : : const auto& solidx = g_inputdeck.get< tag::matidxmap, tag::solidx >();
192 : 0 : std::size_t nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
193 : : std::size_t icnt = 0;
194 [ - - ]: 0 : for (std::size_t k=0; k<nmat; ++k)
195 [ - - ]: 0 : if (solidx[k] > 0)
196 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
197 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
198 : : {
199 : 0 : stiffEqIdx[icnt] =
200 : 0 : inciter::deformIdx(nmat, solidx[k], i, j);
201 : 0 : icnt++;
202 : : }
203 : 0 : }
204 : :
205 : : //! Locate the nonstiff equations.
206 : : //! \param[out] nonStiffEqIdx list with pointers to nonstiff equations
207 : 0 : void setNonStiffEqIdx( std::vector< std::size_t >& nonStiffEqIdx ) const
208 : : {
209 : 0 : nonStiffEqIdx.resize(nnonstiffeq(), 0);
210 [ - - ]: 0 : for (std::size_t icomp=0; icomp<nnonstiffeq(); icomp++)
211 : 0 : nonStiffEqIdx[icomp] = icomp;
212 : 0 : }
213 : :
214 : : //! Initialize the compressible flow equations, prepare for time integration
215 : : //! \param[in] L Block diagonal mass matrix
216 : : //! \param[in] inpoel Element-node connectivity
217 : : //! \param[in] coord Array of nodal coordinates
218 : : //! \param[in] inbox List of elements at which box user ICs are set for
219 : : //! each IC box
220 : : //! \param[in] elemblkid Element ids associated with mesh block ids where
221 : : //! user ICs are set
222 : : //! \param[in,out] unk Array of unknowns
223 : : //! \param[in] t Physical time
224 : : //! \param[in] nielem Number of internal elements
225 [ + - ]: 65 : void initialize( const tk::Fields& L,
226 : : const std::vector< std::size_t >& inpoel,
227 : : const tk::UnsMesh::Coords& coord,
228 : : const std::vector< std::unordered_set< std::size_t > >& inbox,
229 : : const std::unordered_map< std::size_t, std::set< std::size_t > >&
230 : : elemblkid,
231 : : tk::Fields& unk,
232 : : tk::real t,
233 : : const std::size_t nielem ) const
234 : : {
235 [ + - ]: 65 : tk::initialize( m_ncomp, m_mat_blk, L, inpoel, coord,
236 : : Problem::initialize, unk, t, nielem );
237 : :
238 : 65 : const auto rdof = g_inputdeck.get< tag::rdof >();
239 : : const auto& ic = g_inputdeck.get< tag::ic >();
240 : : const auto& icbox = ic.get< tag::box >();
241 : : const auto& icmbk = ic.get< tag::meshblock >();
242 : :
243 : : const auto& bgpre = ic.get< tag::pressure >();
244 : : const auto& bgtemp = ic.get< tag::temperature >();
245 : :
246 : : // Set initial conditions inside user-defined IC boxes and mesh blocks
247 : 65 : std::vector< tk::real > s(m_ncomp, 0.0);
248 [ + + ]: 22511 : for (std::size_t e=0; e<nielem; ++e) {
249 : : // inside user-defined box
250 [ - + ]: 22446 : if (!icbox.empty()) {
251 : : std::size_t bcnt = 0;
252 [ - - ]: 0 : for (const auto& b : icbox) { // for all boxes
253 [ - - ][ - - ]: 0 : if (inbox.size() > bcnt && inbox[bcnt].find(e) != inbox[bcnt].end())
254 : : {
255 [ - - ][ - - ]: 0 : std::vector< tk::real > box
256 : : { b.template get< tag::xmin >(), b.template get< tag::xmax >(),
257 : : b.template get< tag::ymin >(), b.template get< tag::ymax >(),
258 : : b.template get< tag::zmin >(), b.template get< tag::zmax >() };
259 : 0 : auto V_ex = (box[1]-box[0]) * (box[3]-box[2]) * (box[5]-box[4]);
260 [ - - ]: 0 : for (std::size_t c=0; c<m_ncomp; ++c) {
261 : 0 : auto mark = c*rdof;
262 : 0 : s[c] = unk(e,mark);
263 : : // set high-order DOFs to zero
264 [ - - ]: 0 : for (std::size_t i=1; i<rdof; ++i)
265 : 0 : unk(e,mark+i) = 0.0;
266 : : }
267 [ - - ]: 0 : initializeBox<ctr::boxList>( m_mat_blk, V_ex, t, b, bgpre,
268 : : bgtemp, s );
269 : : // store box-initialization in solution vector
270 [ - - ]: 0 : for (std::size_t c=0; c<m_ncomp; ++c) {
271 : 0 : auto mark = c*rdof;
272 : 0 : unk(e,mark) = s[c];
273 : : }
274 : : }
275 : 0 : ++bcnt;
276 : : }
277 : : }
278 : :
279 : : // inside user-specified mesh blocks
280 [ - + ]: 22446 : if (!icmbk.empty()) {
281 [ - - ]: 0 : for (const auto& b : icmbk) { // for all blocks
282 : 0 : auto blid = b.get< tag::blockid >();
283 : 0 : auto V_ex = b.get< tag::volume >();
284 [ - - ]: 0 : if (elemblkid.find(blid) != elemblkid.end()) {
285 : : const auto& elset = tk::cref_find(elemblkid, blid);
286 [ - - ]: 0 : if (elset.find(e) != elset.end()) {
287 [ - - ]: 0 : initializeBox<ctr::meshblockList>( m_mat_blk, V_ex, t, b,
288 : : bgpre, bgtemp, s );
289 : : // store initialization in solution vector
290 [ - - ]: 0 : for (std::size_t c=0; c<m_ncomp; ++c) {
291 : 0 : auto mark = c*rdof;
292 : 0 : unk(e,mark) = s[c];
293 : : }
294 : : }
295 : : }
296 : : }
297 : : }
298 : : }
299 : 65 : }
300 : :
301 : : //! Compute density constraint for a given material
302 : : //! \param[in] nelem Number of elements
303 : : //! \param[in] unk Array of unknowns
304 : : //! \param[out] densityConstr Density Constraint: rho/(rho0*det(g))
305 : 155 : void computeDensityConstr( std::size_t nelem,
306 : : tk::Fields& unk,
307 : : std::vector< tk::real >& densityConstr) const
308 : : {
309 : : const auto& solidx = g_inputdeck.get< tag::matidxmap, tag::solidx >();
310 : 155 : std::size_t rdof = g_inputdeck.get< tag::rdof >();
311 : 155 : std::size_t nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
312 [ + + ]: 64461 : for (std::size_t e=0; e<nelem; ++e)
313 : 64306 : densityConstr[e] = 0.0;
314 [ + + ]: 480 : for (std::size_t imat=0; imat<nmat; ++imat)
315 [ - + ]: 325 : if (solidx[imat] > 0)
316 : : {
317 [ - - ]: 0 : for (std::size_t e=0; e<nelem; ++e)
318 : : {
319 : : // Retrieve unknowns
320 : 0 : tk::real arho = unk(e, densityDofIdx(nmat, imat, rdof, 0));
321 : : std::array< std::array< tk::real, 3 >, 3 > g;
322 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
323 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
324 : 0 : g[i][j] = unk(e, deformDofIdx(nmat, solidx[imat], i, j, rdof, 0));
325 : : // Compute determinant of g
326 : 0 : tk::real detg = tk::determinant(g);
327 : : // Compute constraint measure
328 : 0 : densityConstr[e] += arho/(m_mat_blk[imat].compute< EOS::rho0 >()*detg);
329 : : }
330 : : }
331 : : else
332 : : {
333 [ + + ]: 150795 : for (std::size_t e=0; e<nelem; ++e)
334 : : {
335 : : // Retrieve alpha and add it to the constraint measure
336 : 150470 : tk::real alpha = unk(e, volfracDofIdx(nmat, imat, rdof, 0));
337 : 150470 : densityConstr[e] += alpha;
338 : : }
339 : : }
340 : 155 : }
341 : :
342 : : //! Compute the left hand side block-diagonal mass matrix
343 : : //! \param[in] geoElem Element geometry array
344 : : //! \param[in,out] l Block diagonal mass matrix
345 : : void lhs( const tk::Fields& geoElem, tk::Fields& l ) const {
346 : 65 : const auto ndof = g_inputdeck.get< tag::ndof >();
347 : : // Unlike Compflow and Transport, there is a weak reconstruction about
348 : : // conservative variable after limiting function which will require the
349 : : // size of left hand side vector to be rdof
350 : 65 : tk::mass( m_ncomp, ndof, geoElem, l );
351 : : }
352 : :
353 : : //! Update the interface cells to first order dofs
354 : : //! \param[in] unk Array of unknowns
355 : : //! \param[in] nielem Number of internal elements
356 : : //! \param[in,out] ndofel Array of dofs
357 : : //! \param[in,out] interface Vector of interface marker
358 : : //! \details This function resets the high-order terms in interface cells.
359 : 5205 : void updateInterfaceCells( tk::Fields& unk,
360 : : std::size_t nielem,
361 : : std::vector< std::size_t >& ndofel,
362 : : std::vector< std::size_t >& interface ) const
363 : : {
364 : 5205 : auto intsharp =
365 : : g_inputdeck.get< tag::multimat, tag::intsharp >();
366 : : // If this cell is not material interface, return this function
367 [ + + ]: 5205 : if(not intsharp) return;
368 : :
369 : 375 : auto rdof = g_inputdeck.get< tag::rdof >();
370 : 375 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
371 : : const auto& solidx = g_inputdeck.get<
372 : : tag::matidxmap, tag::solidx >();
373 : :
374 [ + + ]: 227775 : for (std::size_t e=0; e<nielem; ++e) {
375 : 227400 : std::vector< std::size_t > matInt(nmat, 0);
376 [ + - ][ - - ]: 227400 : std::vector< tk::real > alAvg(nmat, 0.0);
377 [ + + ]: 682200 : for (std::size_t k=0; k<nmat; ++k)
378 : 454800 : alAvg[k] = unk(e, volfracDofIdx(nmat,k,rdof,0));
379 [ + - ]: 227400 : auto intInd = interfaceIndicator(nmat, alAvg, matInt);
380 : :
381 : : // interface cells cannot be high-order
382 [ + + ]: 227400 : if (intInd) {
383 : 10834 : interface[e] = 1;
384 [ + + ]: 32502 : for (std::size_t k=0; k<nmat; ++k) {
385 [ + - ]: 21668 : if (matInt[k]) {
386 [ + + ]: 86672 : for (std::size_t i=1; i<rdof; ++i) {
387 : 65004 : unk(e, densityDofIdx(nmat,k,rdof,i)) = 0.0;
388 : 65004 : unk(e, energyDofIdx(nmat,k,rdof,i)) = 0.0;
389 : : }
390 [ - + ]: 21668 : if (solidx[k] > 0) {
391 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
392 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
393 [ - - ]: 0 : for (std::size_t idof=1; idof<rdof; ++idof) {
394 : 0 : unk(e, deformDofIdx(nmat,solidx[k],i,j,rdof,idof)) = 0.0;
395 : : }
396 : : }
397 : : }
398 : : }
399 [ + + ]: 43336 : for (std::size_t idir=0; idir<3; ++idir) {
400 [ + + ]: 130008 : for (std::size_t i=1; i<rdof; ++i) {
401 : 97506 : unk(e, momentumDofIdx(nmat,idir,rdof,i)) = 0.0;
402 : : }
403 : : }
404 : : } else {
405 : : // If the cell is marked as interface cell in the previous timestep
406 : : // and does not marked as interface for the current timestep, DGP2
407 : : // will be applied for the current timestep in p-adaptive process
408 : : // Please note this block is added since the spectral decay indicator
409 : : // does not applied to P0 cells.
410 [ + + ]: 216566 : if (interface[e] == 1) {
411 [ + - ][ - + ]: 6 : if(ndofel[e] < 10 && rdof == 10) {
412 : 0 : ndofel[e] = 10;
413 [ - - ]: 0 : for (std::size_t k=0; k<nmat; ++k) {
414 [ - - ]: 0 : for (std::size_t i=1; i<rdof; ++i) {
415 : 0 : unk(e, densityDofIdx(nmat,k,rdof,i)) = 0.0;
416 : 0 : unk(e, energyDofIdx(nmat,k,rdof,i)) = 0.0;
417 : : }
418 : : }
419 [ - - ]: 0 : for (std::size_t idir=0; idir<3; ++idir) {
420 [ - - ]: 0 : for (std::size_t i=1; i<rdof; ++i) {
421 : 0 : unk(e, momentumDofIdx(nmat,idir,rdof,i)) = 0.0;
422 : : }
423 : : }
424 : : }
425 : : }
426 : 216566 : interface[e] = 0;
427 : : }
428 : : }
429 : : }
430 : :
431 : : //! Update the primitives for this PDE system
432 : : //! \param[in] unk Array of unknowns
433 : : //! \param[in] L The left hand side block-diagonal mass matrix
434 : : //! \param[in] geoElem Element geometry array
435 : : //! \param[in,out] prim Array of primitives
436 : : //! \param[in] nielem Number of internal elements
437 : : //! \param[in] ndofel Array of dofs
438 : : //! \details This function computes and stores the dofs for primitive
439 : : //! quantities, which are required for obtaining reconstructed states used
440 : : //! in the Riemann solver. See /PDE/Riemann/AUSM.hpp, where the
441 : : //! normal velocity for advection is calculated from independently
442 : : //! reconstructed velocities.
443 : 5270 : void updatePrimitives( const tk::Fields& unk,
444 : : const tk::Fields& L,
445 : : const tk::Fields& geoElem,
446 : : tk::Fields& prim,
447 : : std::size_t nielem,
448 : : std::vector< std::size_t >& ndofel ) const
449 : : {
450 : 5270 : const auto rdof = g_inputdeck.get< tag::rdof >();
451 : 5270 : const auto ndof = g_inputdeck.get< tag::ndof >();
452 : 5270 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
453 : : const auto& solidx = g_inputdeck.get< tag::matidxmap, tag::solidx >();
454 : :
455 : : Assert( unk.nunk() == prim.nunk(), "Number of unknowns in solution "
456 : : "vector and primitive vector at recent time step incorrect" );
457 : : Assert( unk.nprop() == rdof*m_ncomp, "Number of components in solution "
458 : : "vector must equal "+ std::to_string(rdof*m_ncomp) );
459 : : Assert( prim.nprop() == rdof*m_nprim, "Number of components in vector of "
460 : : "primitive quantities must equal "+ std::to_string(rdof*m_nprim) );
461 : :
462 [ + + ]: 2189396 : for (std::size_t e=0; e<nielem; ++e)
463 : : {
464 : 2184126 : std::vector< tk::real > R(m_nprim*ndof, 0.0);
465 : :
466 [ + - ]: 2184126 : auto ng = tk::NGvol(ndof);
467 : :
468 : : // arrays for quadrature points
469 : : std::array< std::vector< tk::real >, 3 > coordgp;
470 : : std::vector< tk::real > wgp;
471 : :
472 [ + - ]: 2184126 : coordgp[0].resize( ng );
473 [ + - ]: 2184126 : coordgp[1].resize( ng );
474 [ + - ]: 2184126 : coordgp[2].resize( ng );
475 [ + - ]: 2184126 : wgp.resize( ng );
476 : :
477 [ + - ]: 2184126 : tk::GaussQuadratureTet( ng, coordgp, wgp );
478 : :
479 : : // Local degree of freedom
480 : 2184126 : auto dof_el = ndofel[e];
481 : :
482 : : // Loop over quadrature points in element e
483 [ + + ]: 6399692 : for (std::size_t igp=0; igp<ng; ++igp)
484 : : {
485 : : // Compute the basis function
486 [ + - ]: 4215566 : auto B =
487 [ + - ]: 4215566 : tk::eval_basis( dof_el, coordgp[0][igp], coordgp[1][igp], coordgp[2][igp] );
488 : :
489 [ + - ]: 4215566 : auto w = wgp[igp] * geoElem(e, 0);
490 : :
491 [ + - ]: 4215566 : auto state = tk::eval_state( m_ncomp, rdof, dof_el, e, unk, B );
492 : :
493 : : // bulk density at quadrature point
494 : : tk::real rhob(0.0);
495 [ + + ]: 13746884 : for (std::size_t k=0; k<nmat; ++k)
496 : 9531318 : rhob += state[densityIdx(nmat, k)];
497 : :
498 : : // velocity vector at quadrature point
499 : : std::array< tk::real, 3 >
500 [ + - ]: 4215566 : vel{ state[momentumIdx(nmat, 0)]/rhob,
501 : 4215566 : state[momentumIdx(nmat, 1)]/rhob,
502 : 4215566 : state[momentumIdx(nmat, 2)]/rhob };
503 : :
504 [ + - ][ - - ]: 4215566 : std::vector< tk::real > pri(m_nprim, 0.0);
505 : :
506 : : // Evaluate material pressure at quadrature point
507 [ + + ]: 13746884 : for(std::size_t imat = 0; imat < nmat; imat++)
508 : : {
509 [ + - ]: 9531318 : auto alphamat = state[volfracIdx(nmat, imat)];
510 [ + - ]: 9531318 : auto arhomat = state[densityIdx(nmat, imat)];
511 : 9531318 : auto arhoemat = state[energyIdx(nmat, imat)];
512 [ + - ]: 9531318 : auto gmat = getDeformGrad(nmat, imat, state);
513 [ + - ]: 9531318 : pri[pressureIdx(nmat,imat)] = m_mat_blk[imat].compute<
514 [ + - ]: 9531318 : EOS::pressure >( arhomat, vel[0], vel[1], vel[2], arhoemat,
515 : : alphamat, imat, gmat );
516 : :
517 [ + - ][ - + ]: 9531318 : pri[pressureIdx(nmat,imat)] = constrain_pressure( m_mat_blk,
518 : : pri[pressureIdx(nmat,imat)], arhomat, alphamat, imat);
519 : :
520 [ - + ]: 9531318 : if (solidx[imat] > 0) {
521 [ - - ]: 0 : auto asigmat = m_mat_blk[imat].computeTensor< EOS::CauchyStress >(
522 : : arhomat, vel[0], vel[1], vel[2], arhoemat,
523 : : alphamat, imat, gmat );
524 : :
525 : 0 : pri[stressIdx(nmat,solidx[imat],0)] = asigmat[0][0];
526 : 0 : pri[stressIdx(nmat,solidx[imat],1)] = asigmat[1][1];
527 : 0 : pri[stressIdx(nmat,solidx[imat],2)] = asigmat[2][2];
528 : 0 : pri[stressIdx(nmat,solidx[imat],3)] = asigmat[0][1];
529 : 0 : pri[stressIdx(nmat,solidx[imat],4)] = asigmat[0][2];
530 : 0 : pri[stressIdx(nmat,solidx[imat],5)] = asigmat[1][2];
531 : : }
532 : : }
533 : :
534 : : // Evaluate bulk velocity at quadrature point
535 [ + + ]: 16862264 : for (std::size_t idir=0; idir<3; ++idir) {
536 : 12646698 : pri[velocityIdx(nmat,idir)] = vel[idir];
537 : : }
538 : :
539 [ + + ]: 26393582 : for(std::size_t k = 0; k < m_nprim; k++)
540 : : {
541 : 22178016 : auto mark = k * ndof;
542 [ + + ]: 82445532 : for(std::size_t idof = 0; idof < dof_el; idof++)
543 : 60267516 : R[mark+idof] += w * pri[k] * B[idof];
544 : : }
545 : : }
546 : :
547 : : // Update the DG solution of primitive variables
548 [ + + ]: 14204942 : for(std::size_t k = 0; k < m_nprim; k++)
549 : : {
550 : 12020816 : auto mark = k * ndof;
551 : 12020816 : auto rmark = k * rdof;
552 [ + + ]: 31659532 : for(std::size_t idof = 0; idof < dof_el; idof++)
553 : : {
554 [ + + ]: 19638716 : prim(e, rmark+idof) = R[mark+idof] / L(e, mark+idof);
555 [ + + ]: 19638716 : if(fabs(prim(e, rmark+idof)) < 1e-16)
556 : 5168405 : prim(e, rmark+idof) = 0;
557 : : }
558 : : }
559 : : }
560 : 5270 : }
561 : :
562 : : //! Clean up the state of trace materials for this PDE system
563 : : //! \param[in] t Physical time
564 : : //! \param[in] geoElem Element geometry array
565 : : //! \param[in,out] unk Array of unknowns
566 : : //! \param[in,out] prim Array of primitives
567 : : //! \param[in] nielem Number of internal elements
568 : : //! \details This function cleans up the state of materials present in trace
569 : : //! quantities in each cell. Specifically, the state of materials with
570 : : //! very low volume-fractions in a cell is replaced by the state of the
571 : : //! material which is present in the largest quantity in that cell. This
572 : : //! becomes necessary when shocks pass through cells which contain a very
573 : : //! small amount of material. The state of that tiny material might
574 : : //! become unphysical and cause solution to diverge; thus requiring such
575 : : //! a "reset".
576 : 5205 : void cleanTraceMaterial( tk::real t,
577 : : const tk::Fields& geoElem,
578 : : tk::Fields& unk,
579 : : tk::Fields& prim,
580 : : std::size_t nielem ) const
581 : : {
582 : : [[maybe_unused]] const auto rdof = g_inputdeck.get< tag::rdof >();
583 : 5205 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
584 : :
585 : : Assert( unk.nunk() == prim.nunk(), "Number of unknowns in solution "
586 : : "vector and primitive vector at recent time step incorrect" );
587 : : Assert( unk.nprop() == rdof*m_ncomp, "Number of components in solution "
588 : : "vector must equal "+ std::to_string(rdof*m_ncomp) );
589 : : Assert( prim.nprop() == rdof*m_nprim, "Number of components in vector of "
590 : : "primitive quantities must equal "+ std::to_string(rdof*m_nprim) );
591 : :
592 : 5205 : auto neg_density = cleanTraceMultiMat(t, nielem, m_mat_blk, geoElem, nmat,
593 : : unk, prim);
594 : :
595 [ - + ][ - - ]: 5205 : if (neg_density) Throw("Negative partial density.");
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ]
596 : 5205 : }
597 : :
598 : : //! Reconstruct second-order solution from first-order
599 : : //! \param[in] geoElem Element geometry array
600 : : //! \param[in] fd Face connectivity and boundary conditions object
601 : : //! \param[in] esup Elements-surrounding-nodes connectivity
602 : : //! \param[in] inpoel Element-node connectivity
603 : : //! \param[in] coord Array of nodal coordinates
604 : : //! \param[in,out] U Solution vector at recent time step
605 : : //! \param[in,out] P Vector of primitives at recent time step
606 : : //! \param[in] pref Indicator for p-adaptive algorithm
607 : : //! \param[in] ndofel Vector of local number of degrees of freedome
608 : 4080 : void reconstruct( tk::real,
609 : : const tk::Fields&,
610 : : const tk::Fields& geoElem,
611 : : const inciter::FaceData& fd,
612 : : const std::map< std::size_t, std::vector< std::size_t > >&
613 : : esup,
614 : : const std::vector< std::size_t >& inpoel,
615 : : const tk::UnsMesh::Coords& coord,
616 : : tk::Fields& U,
617 : : tk::Fields& P,
618 : : const bool pref,
619 : : const std::vector< std::size_t >& ndofel ) const
620 : : {
621 : 4080 : const auto rdof = g_inputdeck.get< tag::rdof >();
622 : 4080 : const auto ndof = g_inputdeck.get< tag::ndof >();
623 : :
624 : : bool is_p0p1(false);
625 [ + + ]: 4080 : if (rdof == 4 && ndof == 1)
626 : : is_p0p1 = true;
627 : :
628 : 4080 : const auto nelem = fd.Esuel().size()/4;
629 : 4080 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
630 : :
631 : : Assert( U.nprop() == rdof*m_ncomp, "Number of components in solution "
632 : : "vector must equal "+ std::to_string(rdof*m_ncomp) );
633 : :
634 : : //----- reconstruction of conserved quantities -----
635 : : //--------------------------------------------------
636 : :
637 [ + + ]: 845460 : for (std::size_t e=0; e<nelem; ++e)
638 : : {
639 : : // 1. specify how many variables need to be reconstructed
640 : : std::vector< std::size_t > vars;
641 : : // for p-adaptive DG
642 [ - + ]: 841380 : if (pref) {
643 : : // If DG is applied, reconstruct only volume fractions
644 [ - - ]: 0 : if(ndofel[e] > 1) {
645 [ - - ][ - - ]: 0 : for (std::size_t k=0; k<nmat; ++k) vars.push_back(volfracIdx(nmat, k));
646 : : }
647 : : else // If P0P1 is applied for this element
648 [ - - ][ - - ]: 0 : for (std::size_t c=0; c<m_ncomp; ++c) vars.push_back(c);
649 : : }
650 : : else {
651 : : // for P0P1, reconstruct all variables
652 [ + + ]: 841380 : if (is_p0p1)
653 [ + + ][ + - ]: 3411000 : for (std::size_t c=0; c<m_ncomp; ++c) vars.push_back(c);
654 : : // for high-order DG, reconstruct only volume fractions
655 [ + - ]: 500280 : else if (ndof > 1)
656 [ + + ][ + - ]: 1500840 : for (std::size_t k=0; k<nmat; ++k) vars.push_back(volfracIdx(nmat, k));
[ - - ]
657 : : }
658 : :
659 : : // 2. solve 3x3 least-squares system
660 : : // Reconstruct second-order dofs of volume-fractions in Taylor space
661 : : // using nodal-stencils, for a good interface-normal estimate
662 [ + - ]: 841380 : tk::recoLeastSqExtStencil( rdof, e, esup, inpoel, geoElem, U, vars );
663 : :
664 : : // 3. transform reconstructed derivatives to Dubiner dofs
665 [ + - ]: 841380 : tk::transform_P0P1( rdof, e, inpoel, coord, U, vars );
666 : : }
667 : :
668 : : //----- reconstruction of primitive quantities -----
669 : : //--------------------------------------------------
670 : : // For multimat, conserved and primitive quantities are reconstructed
671 : : // separately.
672 : :
673 [ + + ]: 845460 : for (std::size_t e=0; e<nelem; ++e)
674 : : {
675 : : // There are two conditions that requires the reconstruction of the
676 : : // primitive variables:
677 : : // 1. p-adaptive is triggered and P0P1 scheme is applied to specific
678 : : // elements
679 : : // 2. p-adaptive is not triggered and P0P1 scheme is applied to the
680 : : // whole computation domain
681 [ - + ][ - - ]: 841380 : if ((pref && ndofel[e] == 1) || (!pref && is_p0p1)) {
[ + + ]
682 : : std::vector< std::size_t > vars;
683 [ + + ][ + - ]: 2046600 : for (std::size_t c=0; c<m_nprim; ++c) vars.push_back(c);
684 : :
685 : : // 1.
686 : : // Reconstruct second-order dofs of volume-fractions in Taylor space
687 : : // using nodal-stencils, for a good interface-normal estimate
688 [ + - ]: 341100 : tk::recoLeastSqExtStencil( rdof, e, esup, inpoel, geoElem, P, vars );
689 : :
690 : : // 2.
691 [ + - ]: 341100 : tk::transform_P0P1(rdof, e, inpoel, coord, P, vars);
692 : : }
693 : : }
694 : :
695 : 4080 : }
696 : :
697 : : //! Limit second-order solution, and primitive quantities separately
698 : : //! \param[in] t Physical time
699 : : //! \param[in] pref Indicator for p-adaptive algorithm
700 : : //! \param[in] geoFace Face geometry array
701 : : //! \param[in] geoElem Element geometry array
702 : : //! \param[in] fd Face connectivity and boundary conditions object
703 : : //! \param[in] esup Elements-surrounding-nodes connectivity
704 : : //! \param[in] inpoel Element-node connectivity
705 : : //! \param[in] coord Array of nodal coordinates
706 : : //! \param[in] ndofel Vector of local number of degrees of freedome
707 : : //! \param[in] gid Local->global node id map
708 : : //! \param[in] bid Local chare-boundary node ids (value) associated to
709 : : //! global node ids (key)
710 : : //! \param[in] uNodalExtrm Chare-boundary nodal extrema for conservative
711 : : //! variables
712 : : //! \param[in] pNodalExtrm Chare-boundary nodal extrema for primitive
713 : : //! variables
714 : : //! \param[in] mtInv Inverse of Taylor mass matrix
715 : : //! \param[in,out] U Solution vector at recent time step
716 : : //! \param[in,out] P Vector of primitives at recent time step
717 : : //! \param[in,out] shockmarker Vector of shock-marker values
718 : 4080 : void limit( [[maybe_unused]] tk::real t,
719 : : const bool pref,
720 : : const tk::Fields& geoFace,
721 : : const tk::Fields& geoElem,
722 : : const inciter::FaceData& fd,
723 : : const std::map< std::size_t, std::vector< std::size_t > >& esup,
724 : : const std::vector< std::size_t >& inpoel,
725 : : const tk::UnsMesh::Coords& coord,
726 : : const std::vector< std::size_t >& ndofel,
727 : : const std::vector< std::size_t >& gid,
728 : : const std::unordered_map< std::size_t, std::size_t >& bid,
729 : : const std::vector< std::vector<tk::real> >& uNodalExtrm,
730 : : const std::vector< std::vector<tk::real> >& pNodalExtrm,
731 : : const std::vector< std::vector<tk::real> >& mtInv,
732 : : tk::Fields& U,
733 : : tk::Fields& P,
734 : : std::vector< std::size_t >& shockmarker ) const
735 : : {
736 : : Assert( U.nunk() == P.nunk(), "Number of unknowns in solution "
737 : : "vector and primitive vector at recent time step incorrect" );
738 : :
739 : 4080 : const auto limiter = g_inputdeck.get< tag::limiter >();
740 : 4080 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
741 : 4080 : const auto rdof = g_inputdeck.get< tag::rdof >();
742 : : const auto& solidx = g_inputdeck.get<
743 : : tag::matidxmap, tag::solidx >();
744 : :
745 : : // limit vectors of conserved and primitive quantities
746 [ - + ]: 4080 : if (limiter == ctr::LimiterType::SUPERBEEP1)
747 : : {
748 : 0 : SuperbeeMultiMat_P1( fd.Esuel(), inpoel, ndofel,
749 : : coord, solidx, U, P, nmat );
750 : : }
751 [ + - ]: 4080 : else if (limiter == ctr::LimiterType::VERTEXBASEDP1 && rdof == 4)
752 : : {
753 : 4080 : VertexBasedMultiMat_P1( esup, inpoel, ndofel, fd.Esuel().size()/4,
754 [ + - ]: 4080 : m_mat_blk, fd, geoFace, geoElem, coord, flux, solidx, U, P,
755 : : nmat, shockmarker );
756 : : }
757 [ - - ]: 0 : else if (limiter == ctr::LimiterType::VERTEXBASEDP1 && rdof == 10)
758 : : {
759 : 0 : VertexBasedMultiMat_P2( pref, esup, inpoel, ndofel, fd.Esuel().size()/4,
760 [ - - ]: 0 : m_mat_blk, fd, geoFace, geoElem, coord, gid, bid,
761 : : uNodalExtrm, pNodalExtrm, mtInv, flux, solidx, U, P, nmat,
762 : : shockmarker );
763 : : }
764 [ - - ]: 0 : else if (limiter != ctr::LimiterType::NOLIMITER)
765 : : {
766 [ - - ][ - - ]: 0 : Throw("Limiter type not configured for multimat.");
[ - - ][ - - ]
[ - - ][ - - ]
767 : : }
768 : 4080 : }
769 : :
770 : : //! Apply CPL to the conservative variable solution for this PDE system
771 : : //! \param[in] prim Array of primitive variables
772 : : //! \param[in] geoElem Element geometry array
773 : : //! \param[in] inpoel Element-node connectivity
774 : : //! \param[in] coord Array of nodal coordinates
775 : : //! \param[in,out] unk Array of conservative variables
776 : : //! \param[in] nielem Number of internal elements
777 : : //! \details This function applies CPL to obtain consistent dofs for
778 : : //! conservative quantities based on the limited primitive quantities.
779 : : //! See Pandare et al. (2023). On the Design of Stable,
780 : : //! Consistent, and Conservative High-Order Methods for Multi-Material
781 : : //! Hydrodynamics. J Comp Phys, 112313.
782 : : void CPL( const tk::Fields& prim,
783 : : const tk::Fields& geoElem,
784 : : const std::vector< std::size_t >& inpoel,
785 : : const tk::UnsMesh::Coords& coord,
786 : : tk::Fields& unk,
787 : : std::size_t nielem ) const
788 : : {
789 : : [[maybe_unused]] const auto rdof = g_inputdeck.get< tag::rdof >();
790 : 30 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
791 : :
792 : : Assert( unk.nunk() == prim.nunk(), "Number of unknowns in solution "
793 : : "vector and primitive vector at recent time step incorrect" );
794 : : Assert( unk.nprop() == rdof*m_ncomp, "Number of components in solution "
795 : : "vector must equal "+ std::to_string(rdof*m_ncomp) );
796 : : Assert( prim.nprop() == rdof*m_nprim, "Number of components in vector of "
797 : : "primitive quantities must equal "+ std::to_string(rdof*m_nprim) );
798 : :
799 : 30 : correctLimConservMultiMat(nielem, m_mat_blk, nmat, inpoel,
800 : : coord, geoElem, prim, unk);
801 : : }
802 : :
803 : : //! Return cell-average deformation gradient tensor
804 : : //! \param[in] unk Solution vector at recent time step
805 : : //! \param[in] nielem Number of internal elements
806 : : //! \details This function returns the bulk cell-average inverse
807 : : //! deformation gradient tensor
808 : 0 : std::array< std::vector< tk::real >, 9 > cellAvgDeformGrad(
809 : : const tk::Fields& unk,
810 : : std::size_t nielem ) const
811 : : {
812 : 0 : const auto rdof = g_inputdeck.get< tag::rdof >();
813 : 0 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
814 : : const auto& solidx = g_inputdeck.get< tag::matidxmap, tag::solidx >();
815 : :
816 : : std::array< std::vector< tk::real >, 9 > gb;
817 [ - - ][ - - ]: 0 : if (inciter::haveSolid(nmat, solidx)) {
818 [ - - ]: 0 : for (auto& gij : gb)
819 [ - - ]: 0 : gij.resize(nielem, 0.0);
820 [ - - ]: 0 : for (std::size_t e=0; e<nielem; ++e) {
821 [ - - ]: 0 : for (std::size_t k=0; k<nmat; ++k) {
822 [ - - ]: 0 : if (solidx[k] > 0) {
823 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
824 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
825 : 0 : gb[3*i+j][e] += unk(e, volfracDofIdx(nmat,k,rdof,0)) *
826 : : unk(e,deformDofIdx(nmat,solidx[k],i,j,rdof,0));
827 : : }
828 : : }
829 : : }
830 : : }
831 : :
832 : 0 : return gb;
833 : : }
834 : :
835 : :
836 : : //! Reset the high order solution for p-adaptive scheme
837 : : //! \param[in] fd Face connectivity and boundary conditions object
838 : : //! \param[in,out] unk Solution vector at recent time step
839 : : //! \param[in,out] prim Primitive vector at recent time step
840 : : //! \param[in] ndofel Vector of local number of degrees of freedome
841 : : //! \details This function reset the high order coefficient for p-adaptive
842 : : //! solution polynomials. Unlike compflow class, the high order of fv
843 : : //! solution will not be reset since p0p1 is the base scheme for
844 : : //! multi-material p-adaptive DG method.
845 : 0 : void resetAdapSol( const inciter::FaceData& fd,
846 : : tk::Fields& unk,
847 : : tk::Fields& prim,
848 : : const std::vector< std::size_t >& ndofel ) const
849 : : {
850 : 0 : const auto rdof = g_inputdeck.get< tag::rdof >();
851 : 0 : const auto ncomp = unk.nprop() / rdof;
852 : 0 : const auto nprim = prim.nprop() / rdof;
853 : :
854 [ - - ]: 0 : for(std::size_t e = 0; e < fd.Esuel().size()/4; e++)
855 : : {
856 [ - - ]: 0 : if(ndofel[e] < 10)
857 : : {
858 [ - - ]: 0 : for (std::size_t c=0; c<ncomp; ++c)
859 : : {
860 : 0 : auto mark = c*rdof;
861 : 0 : unk(e, mark+4) = 0.0;
862 : 0 : unk(e, mark+5) = 0.0;
863 : 0 : unk(e, mark+6) = 0.0;
864 : 0 : unk(e, mark+7) = 0.0;
865 : 0 : unk(e, mark+8) = 0.0;
866 : 0 : unk(e, mark+9) = 0.0;
867 : : }
868 [ - - ]: 0 : for (std::size_t c=0; c<nprim; ++c)
869 : : {
870 : 0 : auto mark = c*rdof;
871 : 0 : prim(e, mark+4) = 0.0;
872 : 0 : prim(e, mark+5) = 0.0;
873 : 0 : prim(e, mark+6) = 0.0;
874 : 0 : prim(e, mark+7) = 0.0;
875 : 0 : prim(e, mark+8) = 0.0;
876 : 0 : prim(e, mark+9) = 0.0;
877 : : }
878 : : }
879 : : }
880 : 0 : }
881 : :
882 : : //! Compute right hand side
883 : : //! \param[in] t Physical time
884 : : //! \param[in] pref Indicator for p-adaptive algorithm
885 : : //! \param[in] geoFace Face geometry array
886 : : //! \param[in] geoElem Element geometry array
887 : : //! \param[in] fd Face connectivity and boundary conditions object
888 : : //! \param[in] inpoel Element-node connectivity
889 : : //! \param[in] coord Array of nodal coordinates
890 : : //! \param[in] U Solution vector at recent time step
891 : : //! \param[in] P Primitive vector at recent time step
892 : : //! \param[in] ndofel Vector of local number of degrees of freedom
893 : : //! \param[in] dt Delta time
894 : : //! \param[in,out] R Right-hand side vector computed
895 : 5205 : void rhs( tk::real t,
896 : : const bool pref,
897 : : const tk::Fields& geoFace,
898 : : const tk::Fields& geoElem,
899 : : const inciter::FaceData& fd,
900 : : const std::vector< std::size_t >& inpoel,
901 : : const std::vector< std::unordered_set< std::size_t > >&,
902 : : const tk::UnsMesh::Coords& coord,
903 : : const tk::Fields& U,
904 : : const tk::Fields& P,
905 : : const std::vector< std::size_t >& ndofel,
906 : : const tk::real dt,
907 : : tk::Fields& R ) const
908 : : {
909 : 5205 : const auto ndof = g_inputdeck.get< tag::ndof >();
910 : 5205 : const auto rdof = g_inputdeck.get< tag::rdof >();
911 : 5205 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
912 : 5205 : const auto intsharp =
913 : : g_inputdeck.get< tag::multimat, tag::intsharp >();
914 : : const auto& solidx = inciter::g_inputdeck.get<
915 : : tag::matidxmap, tag::solidx >();
916 : 5205 : auto nsld = numSolids(nmat, solidx);
917 : :
918 : 5205 : const auto nelem = fd.Esuel().size()/4;
919 : :
920 : : Assert( U.nunk() == P.nunk(), "Number of unknowns in solution "
921 : : "vector and primitive vector at recent time step incorrect" );
922 : : Assert( U.nunk() == R.nunk(), "Number of unknowns in solution "
923 : : "vector and right-hand side at recent time step incorrect" );
924 : : Assert( U.nprop() == rdof*m_ncomp, "Number of components in solution "
925 : : "vector must equal "+ std::to_string(rdof*m_ncomp) );
926 : : Assert( P.nprop() == rdof*m_nprim, "Number of components in primitive "
927 : : "vector must equal "+ std::to_string(rdof*m_nprim) );
928 : : Assert( R.nprop() == ndof*m_ncomp, "Number of components in right-hand "
929 : : "side vector must equal "+ std::to_string(ndof*m_ncomp) );
930 : : Assert( fd.Inpofa().size()/3 == fd.Esuf().size()/2,
931 : : "Mismatch in inpofa size" );
932 : :
933 : : // set rhs to zero
934 : : R.fill(0.0);
935 : :
936 : : // Allocate space for Riemann derivatives used in non-conservative terms.
937 : : // The following Riemann derivatives are stored, in order:
938 : : // 1) 3*nmat terms: derivatives of partial pressure of each material,
939 : : // for the energy equations.
940 : : // 2) ndof terms: derivatives of Riemann velocity times the basis
941 : : // function, for the volume fraction equations.
942 : : // 3) nmat*3*3*9 terms: 3 derivatives of u_l*g_ij for each material, for
943 : : // the deformation gradient equations.
944 : : // 4) 3*nsld terms: 3 derivatives of \alpha \sigma_ij for each solid
945 : : // material, for the energy equations.
946 : : std::vector< std::vector< tk::real > >
947 [ + - ][ + - ]: 15615 : riemannDeriv(3*nmat+ndof+3*nsld, std::vector<tk::real>(U.nunk(),0.0));
[ + - ]
948 : :
949 : : // configure a no-op lambda for prescribed velocity
950 : : auto velfn = []( ncomp_t, tk::real, tk::real, tk::real, tk::real ){
951 : : return tk::VelFn::result_type(); };
952 : :
953 : : // compute internal surface flux integrals
954 [ + - ]: 5205 : tk::surfInt( pref, nmat, m_mat_blk, t, ndof, rdof, inpoel, solidx,
955 [ + - ]: 5205 : coord, fd, geoFace, geoElem, m_riemann, velfn, U, P, ndofel,
956 : : dt, R, riemannDeriv, intsharp );
957 : :
958 : : // compute optional source term
959 [ + - ]: 5205 : tk::srcInt( m_mat_blk, t, ndof, fd.Esuel().size()/4, inpoel,
960 : : coord, geoElem, Problem::src, ndofel, R, nmat );
961 : :
962 [ + + ]: 5205 : if(ndof > 1)
963 : : // compute volume integrals
964 [ + - ][ + - ]: 2340 : tk::volInt( nmat, t, m_mat_blk, ndof, rdof, nelem,
[ - - ]
965 : : inpoel, coord, geoElem, flux, velfn, U, P, ndofel, R,
966 : : intsharp );
967 : :
968 : : // compute boundary surface flux integrals
969 [ + + ]: 46845 : for (const auto& b : m_bc)
970 [ + - ]: 83280 : tk::bndSurfInt( pref, nmat, m_mat_blk, ndof, rdof,
971 : : std::get<0>(b), fd, geoFace, geoElem, inpoel, coord, t,
972 : : m_riemann, velfn, std::get<1>(b), U, P, ndofel, R,
973 : : riemannDeriv, intsharp );
974 : :
975 : : Assert( riemannDeriv.size() == 3*nmat+ndof+3*nsld, "Size of "
976 : : "Riemann derivative vector incorrect" );
977 : :
978 : : // get derivatives from riemannDeriv
979 [ + + ]: 46230 : for (std::size_t k=0; k<riemannDeriv.size(); ++k)
980 : : {
981 : : Assert( riemannDeriv[k].size() == U.nunk(), "Riemann derivative vector "
982 : : "for non-conservative terms has incorrect size" );
983 [ + + ]: 28165725 : for (std::size_t e=0; e<U.nunk(); ++e)
984 : 28124700 : riemannDeriv[k][e] /= geoElem(e, 0);
985 : : }
986 : :
987 : : // compute volume integrals of non-conservative terms
988 [ + - ]: 5205 : tk::nonConservativeInt( pref, nmat, m_mat_blk, ndof, rdof, nelem,
989 : : inpoel, coord, geoElem, U, P, riemannDeriv,
990 : : ndofel, R, intsharp );
991 : :
992 : : // Compute integrals for inverse deformation correction in solid materials
993 [ + - ][ - + ]: 5205 : if (inciter::haveSolid(nmat, solidx) &&
[ - - ]
994 : : g_inputdeck.get< tag::multimat, tag::rho0constraint >())
995 [ - - ]: 0 : tk::solidTermsVolInt( nmat, m_mat_blk, ndof, rdof, nelem,
996 : : inpoel, coord, geoElem, U, P, ndofel,
997 : : dt, R);
998 : :
999 : : // compute finite pressure relaxation terms
1000 [ + + ]: 5205 : if (g_inputdeck.get< tag::multimat, tag::prelax >())
1001 : : {
1002 : 375 : const auto ct = g_inputdeck.get< tag::multimat,
1003 : : tag::prelax_timescale >();
1004 [ + - ]: 375 : tk::pressureRelaxationInt( pref, nmat, m_mat_blk, ndof,
1005 : : rdof, nelem, inpoel, coord, geoElem, U, P,
1006 : : ndofel, ct, R, intsharp );
1007 : : }
1008 : 5205 : }
1009 : :
1010 : : //! Evaluate the adaptive indicator and mark the ndof for each element
1011 : : //! \param[in] nunk Number of unknowns
1012 : : //! \param[in] coord Array of nodal coordinates
1013 : : //! \param[in] inpoel Element-node connectivity
1014 : : //! \param[in] fd Face connectivity and boundary conditions object
1015 : : //! \param[in] unk Array of unknowns
1016 : : //! \param[in] prim Array of primitive quantities
1017 : : //! \param[in] indicator p-refinement indicator type
1018 : : //! \param[in] ndof Number of degrees of freedom in the solution
1019 : : //! \param[in] ndofmax Max number of degrees of freedom for p-refinement
1020 : : //! \param[in] tolref Tolerance for p-refinement
1021 : : //! \param[in,out] ndofel Vector of local number of degrees of freedome
1022 [ - - ]: 0 : void eval_ndof( std::size_t nunk,
1023 : : [[maybe_unused]] const tk::UnsMesh::Coords& coord,
1024 : : [[maybe_unused]] const std::vector< std::size_t >& inpoel,
1025 : : const inciter::FaceData& fd,
1026 : : const tk::Fields& unk,
1027 : : [[maybe_unused]] const tk::Fields& prim,
1028 : : inciter::ctr::PrefIndicatorType indicator,
1029 : : std::size_t ndof,
1030 : : std::size_t ndofmax,
1031 : : tk::real tolref,
1032 : : std::vector< std::size_t >& ndofel ) const
1033 : : {
1034 : : const auto& esuel = fd.Esuel();
1035 : 0 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1036 : :
1037 [ - - ]: 0 : if(indicator == inciter::ctr::PrefIndicatorType::SPECTRAL_DECAY)
1038 : 0 : spectral_decay(nmat, nunk, esuel, unk, ndof, ndofmax, tolref, ndofel);
1039 : : else
1040 [ - - ][ - - ]: 0 : Throw( "No such adaptive indicator type" );
[ - - ][ - - ]
[ - - ][ - - ]
1041 : 0 : }
1042 : :
1043 : : //! Compute the minimum time step size
1044 : : //! \param[in] fd Face connectivity and boundary conditions object
1045 : : //! \param[in] geoFace Face geometry array
1046 : : //! \param[in] geoElem Element geometry array
1047 : : // //! \param[in] ndofel Vector of local number of degrees of freedom
1048 : : //! \param[in] U Solution vector at recent time step
1049 : : //! \param[in] P Vector of primitive quantities at recent time step
1050 : : //! \param[in] nielem Number of internal elements
1051 : : //! \return Minimum time step size
1052 : : //! \details The allowable dt is calculated by looking at the maximum
1053 : : //! wave-speed in elements surrounding each face, times the area of that
1054 : : //! face. Once the maximum of this quantity over the mesh is determined,
1055 : : //! the volume of each cell is divided by this quantity. A minimum of this
1056 : : //! ratio is found over the entire mesh, which gives the allowable dt.
1057 : : tk::real dt( const std::array< std::vector< tk::real >, 3 >&,
1058 : : const std::vector< std::size_t >&,
1059 : : const inciter::FaceData& fd,
1060 : : const tk::Fields& geoFace,
1061 : : const tk::Fields& geoElem,
1062 : : const std::vector< std::size_t >& /*ndofel*/,
1063 : : const tk::Fields& U,
1064 : : const tk::Fields& P,
1065 : : const std::size_t nielem ) const
1066 : : {
1067 : : const auto ndof = g_inputdeck.get< tag::ndof >();
1068 : : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1069 : :
1070 : : auto mindt = timeStepSizeMultiMat( m_mat_blk, fd.Esuf(), geoFace, geoElem,
1071 : : nielem, nmat, U, P);
1072 : :
1073 : : tk::real dgp = 0.0;
1074 : : if (ndof == 4)
1075 : : {
1076 : : dgp = 1.0;
1077 : : }
1078 : : else if (ndof == 10)
1079 : : {
1080 : : dgp = 2.0;
1081 : : }
1082 : :
1083 : : // Scale smallest dt with CFL coefficient and the CFL is scaled by (2*p+1)
1084 : : // where p is the order of the DG polynomial by linear stability theory.
1085 : : mindt /= (2.0*dgp + 1.0);
1086 : : return mindt;
1087 : : }
1088 : :
1089 : : //! Compute stiff terms for a single element
1090 : : //! \param[in] e Element number
1091 : : //! \param[in] geoElem Element geometry array
1092 : : //! \param[in] inpoel Element-node connectivity
1093 : : //! \param[in] coord Array of nodal coordinates
1094 : : //! \param[in] U Solution vector at recent time step
1095 : : //! \param[in] P Primitive vector at recent time step
1096 : : //! \param[in] ndofel Vector of local number of degrees of freedom
1097 : : //! \param[in,out] R Right-hand side vector computed
1098 : 0 : void stiff_rhs( std::size_t e,
1099 : : const tk::Fields& geoElem,
1100 : : const std::vector< std::size_t >& inpoel,
1101 : : const tk::UnsMesh::Coords& coord,
1102 : : const tk::Fields& U,
1103 : : const tk::Fields& P,
1104 : : const std::vector< std::size_t >& ndofel,
1105 : : tk::Fields& R ) const
1106 : : {
1107 : 0 : const auto ndof = g_inputdeck.get< tag::ndof >();
1108 : 0 : const auto rdof = g_inputdeck.get< tag::rdof >();
1109 : 0 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1110 : 0 : const auto intsharp =
1111 : : g_inputdeck.get< tag::multimat, tag::intsharp >();
1112 : : const auto& solidx = inciter::g_inputdeck.get<
1113 : : tag::matidxmap, tag::solidx >();
1114 : :
1115 : : Assert( U.nunk() == P.nunk(), "Number of unknowns in solution "
1116 : : "vector and primitive vector at recent time step incorrect" );
1117 : : Assert( U.nprop() == rdof*m_ncomp, "Number of components in solution "
1118 : : "vector must equal "+ std::to_string(rdof*m_ncomp) );
1119 : : Assert( P.nprop() == rdof*m_nprim, "Number of components in primitive "
1120 : : "vector must equal "+ std::to_string(rdof*m_nprim) );
1121 : : Assert( R.nprop() == ndof*nstiffeq(), "Number of components in "
1122 : : "right-hand side must equal "+ std::to_string(ndof*nstiffeq()) );
1123 : :
1124 : : // set rhs to zero for element e
1125 [ - - ]: 0 : for (std::size_t i=0; i<ndof*nstiffeq(); ++i)
1126 : 0 : R(e, i) = 0.0;
1127 : :
1128 : : const auto& cx = coord[0];
1129 : : const auto& cy = coord[1];
1130 : : const auto& cz = coord[2];
1131 : :
1132 : 0 : auto ncomp = U.nprop()/rdof;
1133 : 0 : auto nprim = P.nprop()/rdof;
1134 : :
1135 : 0 : auto ng = tk::NGvol(ndofel[e]);
1136 : :
1137 : : // arrays for quadrature points
1138 : : std::array< std::vector< tk::real >, 3 > coordgp;
1139 : : std::vector< tk::real > wgp;
1140 : :
1141 [ - - ]: 0 : coordgp[0].resize( ng );
1142 [ - - ]: 0 : coordgp[1].resize( ng );
1143 [ - - ]: 0 : coordgp[2].resize( ng );
1144 [ - - ]: 0 : wgp.resize( ng );
1145 : :
1146 [ - - ]: 0 : tk::GaussQuadratureTet( ng, coordgp, wgp );
1147 : :
1148 : : // Extract the element coordinates
1149 : 0 : std::array< std::array< tk::real, 3>, 4 > coordel {{
1150 : : {{ cx[ inpoel[4*e ] ], cy[ inpoel[4*e ] ], cz[ inpoel[4*e ] ] }},
1151 : : {{ cx[ inpoel[4*e+1] ], cy[ inpoel[4*e+1] ], cz[ inpoel[4*e+1] ] }},
1152 : : {{ cx[ inpoel[4*e+2] ], cy[ inpoel[4*e+2] ], cz[ inpoel[4*e+2] ] }},
1153 : : {{ cx[ inpoel[4*e+3] ], cy[ inpoel[4*e+3] ], cz[ inpoel[4*e+3] ] }}
1154 : : }};
1155 : :
1156 : : // Gaussian quadrature
1157 [ - - ]: 0 : for (std::size_t igp=0; igp<ng; ++igp)
1158 : : {
1159 : : // Compute the coordinates of quadrature point at physical domain
1160 [ - - ]: 0 : auto gp = tk::eval_gp( igp, coordel, coordgp );
1161 : :
1162 : : // Compute the basis function
1163 [ - - ]: 0 : auto B = tk::eval_basis( ndofel[e], coordgp[0][igp], coordgp[1][igp],
1164 : : coordgp[2][igp] );
1165 : :
1166 [ - - ]: 0 : auto state = tk::evalPolynomialSol(m_mat_blk, intsharp, ncomp, nprim,
1167 : : rdof, nmat, e, ndofel[e], inpoel, coord, geoElem, gp, B, U, P);
1168 : :
1169 : : // compute source
1170 : : // Loop through materials
1171 : : std::size_t ksld = 0;
1172 [ - - ]: 0 : for (std::size_t k=0; k<nmat; ++k)
1173 : : {
1174 [ - - ]: 0 : if (solidx[k] > 0)
1175 : : {
1176 : 0 : tk::real alpha = state[inciter::volfracIdx(nmat, k)];
1177 : : std::array< std::array< tk::real, 3 >, 3 > g;
1178 : : // Compute the source terms
1179 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1180 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1181 : 0 : g[i][j] = state[inciter::deformIdx(nmat,solidx[k],i,j)];
1182 : :
1183 : : // Compute Lp
1184 : : // Reference: Ortega, A. L., Lombardini, M., Pullin, D. I., &
1185 : : // Meiron, D. I. (2014). Numerical simulation of elastic–plastic
1186 : : // solid mechanics using an Eulerian stretch tensor approach and
1187 : : // HLLD Riemann solver. Journal of Computational Physics, 257,
1188 : : // 414-441
1189 : : std::array< std::array< tk::real, 3 >, 3 > Lp;
1190 : :
1191 : : // 1. Compute dev(sigma)
1192 [ - - ]: 0 : auto sigma_dev = m_mat_blk[k].computeTensor< EOS::CauchyStress >(
1193 [ - - ]: 0 : state[inciter::densityIdx(nmat, k)], 0.0, 0.0, 0.0, 0.0, alpha, k,
1194 : : g );
1195 : 0 : tk::real apr = state[ncomp+inciter::pressureIdx(nmat, k)];
1196 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i) sigma_dev[i][i] -= apr;
1197 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1198 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1199 : 0 : sigma_dev[i][j] /= alpha;
1200 : 0 : tk::real sigma_trace =
1201 : 0 : sigma_dev[0][0]+sigma_dev[1][1]+sigma_dev[2][2];
1202 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1203 : 0 : sigma_dev[i][i] -= sigma_trace/3.0;
1204 : :
1205 : : // 2. Compute inv(g)
1206 : : double ginv[9];
1207 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1208 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1209 : 0 : ginv[3*i+j] = g[i][j];
1210 : : lapack_int ipiv[3];
1211 : : #ifndef NDEBUG
1212 : : lapack_int ierr =
1213 : : #endif
1214 [ - - ]: 0 : LAPACKE_dgetrf(LAPACK_ROW_MAJOR, 3, 3, ginv, 3, ipiv);
1215 : : Assert(ierr==0, "Lapack error in LU factorization of g");
1216 : : #ifndef NDEBUG
1217 : : lapack_int jerr =
1218 : : #endif
1219 [ - - ]: 0 : LAPACKE_dgetri(LAPACK_ROW_MAJOR, 3, ginv, 3, ipiv);
1220 : : Assert(jerr==0, "Lapack error in inverting g");
1221 : :
1222 : : // 3. Compute dev(sigma)*inv(g)
1223 : : std::array< std::array< tk::real, 3 >, 3 > aux_mat;
1224 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1225 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1226 : : {
1227 : : tk::real sum = 0.0;
1228 [ - - ]: 0 : for (std::size_t l=0; l<3; ++l)
1229 : 0 : sum += sigma_dev[i][l]*ginv[3*l+j];
1230 : 0 : aux_mat[i][j] = sum;
1231 : : }
1232 : :
1233 : : // 4. Compute g*(dev(sigma)*inv(g))
1234 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1235 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1236 : : {
1237 : : tk::real sum = 0.0;
1238 [ - - ]: 0 : for (std::size_t l=0; l<3; ++l)
1239 : 0 : sum += g[i][l]*aux_mat[l][j];
1240 : 0 : Lp[i][j] = sum;
1241 : : }
1242 : :
1243 : : // 5. Divide by 2*mu*tau
1244 : : // 'Perfect' plasticity
1245 [ - - ]: 0 : tk::real yield_stress = getmatprop< tag::yield_stress >(k);
1246 : : tk::real equiv_stress = 0.0;
1247 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1248 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1249 : 0 : equiv_stress += sigma_dev[i][j]*sigma_dev[i][j];
1250 : 0 : equiv_stress = std::sqrt(3.0*equiv_stress/2.0);
1251 : : // rel_factor = 1/tau <- Perfect plasticity for now.
1252 : : tk::real rel_factor = 0.0;
1253 [ - - ]: 0 : if (equiv_stress >= yield_stress)
1254 : : rel_factor = 1.0e07;
1255 [ - - ]: 0 : tk::real mu = getmatprop< tag::mu >(k);
1256 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1257 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1258 : 0 : Lp[i][j] *= rel_factor/(2.0*mu);
1259 : :
1260 : : // Compute the source terms
1261 [ - - ]: 0 : std::vector< tk::real > s(9*ndof, 0.0);
1262 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1263 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1264 [ - - ]: 0 : for (std::size_t idof=0; idof<ndof; ++idof)
1265 : : {
1266 : 0 : s[(i*3+j)*ndof+idof] = B[idof] * (Lp[i][0]*g[0][j]
1267 : 0 : +Lp[i][1]*g[1][j]
1268 : 0 : +Lp[i][2]*g[2][j]);
1269 : : }
1270 : :
1271 : 0 : auto wt = wgp[igp] * geoElem(e, 0);
1272 : :
1273 : : // Contribute to the right-hand-side
1274 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1275 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1276 [ - - ]: 0 : for (std::size_t idof=0; idof<ndof; ++idof)
1277 : : {
1278 : 0 : std::size_t srcId = (i*3+j)*ndof+idof;
1279 : 0 : std::size_t dofId = solidTensorIdx(ksld,i,j)*ndof+idof;
1280 : 0 : R(e, dofId) += wt * s[srcId];
1281 : : }
1282 : :
1283 [ - - ]: 0 : ksld++;
1284 : : }
1285 : : }
1286 : :
1287 : : }
1288 : 0 : }
1289 : :
1290 : : //! Extract the velocity field at cell nodes. Currently unused.
1291 : : //! \param[in] U Solution vector at recent time step
1292 : : //! \param[in] N Element node indices
1293 : : //! \return Array of the four values of the velocity field
1294 : : std::array< std::array< tk::real, 4 >, 3 >
1295 : : velocity( const tk::Fields& U,
1296 : : const std::array< std::vector< tk::real >, 3 >&,
1297 : : const std::array< std::size_t, 4 >& N ) const
1298 : : {
1299 : : const auto rdof = g_inputdeck.get< tag::rdof >();
1300 : : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1301 : :
1302 : : std::array< std::array< tk::real, 4 >, 3 > v;
1303 : : v[0] = U.extract( momentumDofIdx(nmat, 0, rdof, 0), N );
1304 : : v[1] = U.extract( momentumDofIdx(nmat, 1, rdof, 0), N );
1305 : : v[2] = U.extract( momentumDofIdx(nmat, 2, rdof, 0), N );
1306 : :
1307 : : std::vector< std::array< tk::real, 4 > > ar;
1308 : : ar.resize(nmat);
1309 : : for (std::size_t k=0; k<nmat; ++k)
1310 : : ar[k] = U.extract( densityDofIdx(nmat, k, rdof, 0), N );
1311 : :
1312 : : std::array< tk::real, 4 > r{{ 0.0, 0.0, 0.0, 0.0 }};
1313 : : for (std::size_t i=0; i<r.size(); ++i) {
1314 : : for (std::size_t k=0; k<nmat; ++k)
1315 : : r[i] += ar[k][i];
1316 : : }
1317 : :
1318 : : std::transform( r.begin(), r.end(), v[0].begin(), v[0].begin(),
1319 : : []( tk::real s, tk::real& d ){ return d /= s; } );
1320 : : std::transform( r.begin(), r.end(), v[1].begin(), v[1].begin(),
1321 : : []( tk::real s, tk::real& d ){ return d /= s; } );
1322 : : std::transform( r.begin(), r.end(), v[2].begin(), v[2].begin(),
1323 : : []( tk::real s, tk::real& d ){ return d /= s; } );
1324 : : return v;
1325 : : }
1326 : :
1327 : : //! Return a map that associates user-specified strings to functions
1328 : : //! \return Map that associates user-specified strings to functions that
1329 : : //! compute relevant quantities to be output to file
1330 : : std::map< std::string, tk::GetVarFn > OutVarFn() const
1331 : 310 : { return MultiMatOutVarFn(); }
1332 : :
1333 : : //! Return analytic field names to be output to file
1334 : : //! \return Vector of strings labelling analytic fields output in file
1335 : : std::vector< std::string > analyticFieldNames() const {
1336 : 0 : auto nmat = g_inputdeck.get< eq, tag::nmat >();
1337 : :
1338 : 0 : return MultiMatFieldNames(nmat);
1339 : : }
1340 : :
1341 : : //! Return time history field names to be output to file
1342 : : //! \return Vector of strings labelling time history fields output in file
1343 : : std::vector< std::string > histNames() const {
1344 : 0 : return MultiMatHistNames();
1345 : : }
1346 : :
1347 : : //! Return surface field output going to file
1348 : : std::vector< std::vector< tk::real > >
1349 : : surfOutput( const std::map< int, std::vector< std::size_t > >&,
1350 : : tk::Fields& ) const
1351 : : {
1352 : : std::vector< std::vector< tk::real > > s; // punt for now
1353 : : return s;
1354 : : }
1355 : :
1356 : : //! Return time history field output evaluated at time history points
1357 : : //! \param[in] h History point data
1358 : : //! \param[in] inpoel Element-node connectivity
1359 : : //! \param[in] coord Array of nodal coordinates
1360 : : //! \param[in] U Array of unknowns
1361 : : //! \param[in] P Array of primitive quantities
1362 : : //! \return Vector of time history output of bulk flow quantities (density,
1363 : : //! velocity, total energy, and pressure) evaluated at time history points
1364 : : std::vector< std::vector< tk::real > >
1365 : 0 : histOutput( const std::vector< HistData >& h,
1366 : : const std::vector< std::size_t >& inpoel,
1367 : : const tk::UnsMesh::Coords& coord,
1368 : : const tk::Fields& U,
1369 : : const tk::Fields& P ) const
1370 : : {
1371 : 0 : const auto rdof = g_inputdeck.get< tag::rdof >();
1372 : 0 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1373 : :
1374 : : const auto& x = coord[0];
1375 : : const auto& y = coord[1];
1376 : : const auto& z = coord[2];
1377 : :
1378 : 0 : std::vector< std::vector< tk::real > > Up(h.size());
1379 : :
1380 : : std::size_t j = 0;
1381 [ - - ]: 0 : for (const auto& p : h) {
1382 : 0 : auto e = p.get< tag::elem >();
1383 : 0 : auto chp = p.get< tag::coord >();
1384 : :
1385 : : // Evaluate inverse Jacobian
1386 : 0 : std::array< std::array< tk::real, 3>, 4 > cp{{
1387 : : {{ x[inpoel[4*e ]], y[inpoel[4*e ]], z[inpoel[4*e ]] }},
1388 : : {{ x[inpoel[4*e+1]], y[inpoel[4*e+1]], z[inpoel[4*e+1]] }},
1389 : : {{ x[inpoel[4*e+2]], y[inpoel[4*e+2]], z[inpoel[4*e+2]] }},
1390 : : {{ x[inpoel[4*e+3]], y[inpoel[4*e+3]], z[inpoel[4*e+3]] }} }};
1391 : 0 : auto J = tk::inverseJacobian( cp[0], cp[1], cp[2], cp[3] );
1392 : :
1393 : : // evaluate solution at history-point
1394 : 0 : std::array< tk::real, 3 > dc{{chp[0]-cp[0][0], chp[1]-cp[0][1],
1395 [ - - ]: 0 : chp[2]-cp[0][2]}};
1396 [ - - ]: 0 : auto B = tk::eval_basis(rdof, tk::dot(J[0],dc), tk::dot(J[1],dc),
1397 : : tk::dot(J[2],dc));
1398 [ - - ]: 0 : auto uhp = eval_state(m_ncomp, rdof, rdof, e, U, B);
1399 [ - - ]: 0 : auto php = eval_state(m_nprim, rdof, rdof, e, P, B);
1400 : :
1401 : : // store solution in history output vector
1402 [ - - ][ - - ]: 0 : Up[j].resize(6, 0.0);
1403 [ - - ]: 0 : for (std::size_t k=0; k<nmat; ++k) {
1404 : 0 : Up[j][0] += uhp[densityIdx(nmat,k)];
1405 : 0 : Up[j][4] += uhp[energyIdx(nmat,k)];
1406 : 0 : Up[j][5] += php[pressureIdx(nmat,k)];
1407 : : }
1408 [ - - ]: 0 : Up[j][1] = php[velocityIdx(nmat,0)];
1409 : 0 : Up[j][2] = php[velocityIdx(nmat,1)];
1410 : 0 : Up[j][3] = php[velocityIdx(nmat,2)];
1411 [ - - ]: 0 : ++j;
1412 : : }
1413 : :
1414 : 0 : return Up;
1415 : : }
1416 : :
1417 : : //! Return names of integral variables to be output to diagnostics file
1418 : : //! \return Vector of strings labelling integral variables output
1419 : : std::vector< std::string > names() const
1420 : : {
1421 : 12 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1422 : 12 : return MultiMatDiagNames(nmat);
1423 : : }
1424 : :
1425 : : //! Return analytic solution (if defined by Problem) at xi, yi, zi, t
1426 : : //! \param[in] xi X-coordinate at which to evaluate the analytic solution
1427 : : //! \param[in] yi Y-coordinate at which to evaluate the analytic solution
1428 : : //! \param[in] zi Z-coordinate at which to evaluate the analytic solution
1429 : : //! \param[in] t Physical time at which to evaluate the analytic solution
1430 : : //! \return Vector of analytic solution at given location and time
1431 : : std::vector< tk::real >
1432 : : analyticSolution( tk::real xi, tk::real yi, tk::real zi, tk::real t ) const
1433 : 0 : { return Problem::analyticSolution( m_ncomp, m_mat_blk, xi, yi, zi, t ); }
1434 : :
1435 : : //! Return analytic solution for conserved variables
1436 : : //! \param[in] xi X-coordinate at which to evaluate the analytic solution
1437 : : //! \param[in] yi Y-coordinate at which to evaluate the analytic solution
1438 : : //! \param[in] zi Z-coordinate at which to evaluate the analytic solution
1439 : : //! \param[in] t Physical time at which to evaluate the analytic solution
1440 : : //! \return Vector of analytic solution at given location and time
1441 : : std::vector< tk::real >
1442 : : solution( tk::real xi, tk::real yi, tk::real zi, tk::real t ) const
1443 : 1322920 : { return Problem::initialize( m_ncomp, m_mat_blk, xi, yi, zi, t ); }
1444 : :
1445 : : //! Return cell-averaged specific total energy for an element
1446 : : //! \param[in] e Element id for which total energy is required
1447 : : //! \param[in] unk Vector of conserved quantities
1448 : : //! \return Cell-averaged specific total energy for given element
1449 : : tk::real sp_totalenergy(std::size_t e, const tk::Fields& unk) const
1450 : : {
1451 : 822640 : const auto rdof = g_inputdeck.get< tag::rdof >();
1452 : 822640 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1453 : :
1454 : : tk::real sp_te(0.0);
1455 : : // sum each material total energy
1456 [ - - ][ + + ]: 2858500 : for (std::size_t k=0; k<nmat; ++k) {
[ + + ][ + + ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
1457 : 2035860 : sp_te += unk(e, energyDofIdx(nmat,k,rdof,0));
1458 : : }
1459 : : return sp_te;
1460 : : }
1461 : :
1462 : : private:
1463 : : //! Number of components in this PDE system
1464 : : const ncomp_t m_ncomp;
1465 : : //! Number of primitive quantities stored in this PDE system
1466 : : const ncomp_t m_nprim;
1467 : : //! Riemann solver
1468 : : tk::RiemannFluxFn m_riemann;
1469 : : //! BC configuration
1470 : : BCStateFn m_bc;
1471 : : //! EOS material block
1472 : : std::vector< EOS > m_mat_blk;
1473 : :
1474 : : //! Evaluate conservative part of physical flux function for this PDE system
1475 : : //! \param[in] ncomp Number of scalar components in this PDE system
1476 : : //! \param[in] ugp Numerical solution at the Gauss point at which to
1477 : : //! evaluate the flux
1478 : : //! \return Flux vectors for all components in this PDE system
1479 : : //! \note The function signature must follow tk::FluxFn
1480 : : static tk::FluxFn::result_type
1481 : 7374000 : flux( ncomp_t ncomp,
1482 : : const std::vector< EOS >& mat_blk,
1483 : : const std::vector< tk::real >& ugp,
1484 : : const std::vector< std::array< tk::real, 3 > >& )
1485 : : {
1486 : 7374000 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1487 : :
1488 : 7374000 : return tk::fluxTerms(ncomp, nmat, mat_blk, ugp);
1489 : : }
1490 : :
1491 : : //! \brief Boundary state function providing the left and right state of a
1492 : : //! face at Dirichlet boundaries
1493 : : //! \param[in] ncomp Number of scalar components in this PDE system
1494 : : //! \param[in] mat_blk EOS material block
1495 : : //! \param[in] ul Left (domain-internal) state
1496 : : //! \param[in] x X-coordinate at which to compute the states
1497 : : //! \param[in] y Y-coordinate at which to compute the states
1498 : : //! \param[in] z Z-coordinate at which to compute the states
1499 : : //! \param[in] t Physical time
1500 : : //! \return Left and right states for all scalar components in this PDE
1501 : : //! system
1502 : : //! \note The function signature must follow tk::StateFn. For multimat, the
1503 : : //! left or right state is the vector of conserved quantities, followed by
1504 : : //! the vector of primitive quantities appended to it.
1505 : : static tk::StateFn::result_type
1506 : 48000 : dirichlet( ncomp_t ncomp,
1507 : : const std::vector< EOS >& mat_blk,
1508 : : const std::vector< tk::real >& ul, tk::real x, tk::real y,
1509 : : tk::real z, tk::real t, const std::array< tk::real, 3 >& )
1510 : : {
1511 : 48000 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1512 : : const auto& solidx = g_inputdeck.get<
1513 : : tag::matidxmap, tag::solidx >();
1514 : :
1515 : 48000 : [[maybe_unused]] auto nsld = numSolids(nmat, solidx);
1516 : :
1517 : 48000 : auto ur = Problem::initialize( ncomp, mat_blk, x, y, z, t );
1518 : : Assert( ur.size() == ncomp, "Incorrect size for boundary state vector" );
1519 : :
1520 [ + - ]: 48000 : ur.resize(ul.size());
1521 : :
1522 : : tk::real rho(0.0);
1523 [ + + ]: 192000 : for (std::size_t k=0; k<nmat; ++k)
1524 : 144000 : rho += ur[densityIdx(nmat, k)];
1525 : :
1526 : : // get primitives in boundary state
1527 : :
1528 : : // velocity
1529 : 48000 : ur[ncomp+velocityIdx(nmat, 0)] = ur[momentumIdx(nmat, 0)] / rho;
1530 : 48000 : ur[ncomp+velocityIdx(nmat, 1)] = ur[momentumIdx(nmat, 1)] / rho;
1531 : 48000 : ur[ncomp+velocityIdx(nmat, 2)] = ur[momentumIdx(nmat, 2)] / rho;
1532 : :
1533 : : // material pressures
1534 [ + + ]: 192000 : for (std::size_t k=0; k<nmat; ++k)
1535 : : {
1536 [ + - ]: 144000 : auto gk = getDeformGrad(nmat, k, ur);
1537 [ + - ]: 144000 : ur[ncomp+pressureIdx(nmat, k)] = mat_blk[k].compute< EOS::pressure >(
1538 : : ur[densityIdx(nmat, k)], ur[ncomp+velocityIdx(nmat, 0)],
1539 : : ur[ncomp+velocityIdx(nmat, 1)], ur[ncomp+velocityIdx(nmat, 2)],
1540 [ + - ]: 144000 : ur[energyIdx(nmat, k)], ur[volfracIdx(nmat, k)], k, gk );
1541 : : }
1542 : :
1543 : : Assert( ur.size() == ncomp+nmat+3+nsld*6, "Incorrect size for appended "
1544 : : "boundary state vector" );
1545 : :
1546 [ + - ]: 48000 : return {{ std::move(ul), std::move(ur) }};
1547 : : }
1548 : :
1549 : : // Other boundary condition types that do not depend on "Problem" should be
1550 : : // added in BCFunctions.hpp
1551 : : };
1552 : :
1553 : : } // dg::
1554 : :
1555 : : } // inciter::
1556 : :
1557 : : #endif // DGMultiMat_h
|