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 : 142 : explicit MultiMat() :
72 : : m_ncomp( g_inputdeck.get< tag::ncomp >() ),
73 : : m_nprim(nprim()),
74 : : m_riemann( multimatRiemannSolver(
75 [ + - ]: 142 : g_inputdeck.get< tag::flux >() ) )
76 : : {
77 : : // associate boundary condition configurations with state functions
78 [ + - ][ + - ]: 2130 : 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 [ + - ][ + - ]: 426 : 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 [ + - ][ + - ]: 284 : ConfigBackPressureBC(m_bc, back_pressure, noOpGrad);
[ - - ]
104 : :
105 : : // EoS initialization
106 [ + - ]: 142 : initializeMaterialEoS( m_mat_blk );
107 : 142 : }
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 : 211 : 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 : 211 : std::size_t np(nmat+3);
120 : :
121 [ + + ][ + + ]: 717 : for (std::size_t k=0; k<nmat; ++k) {
[ + + ][ + + ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ + + ][ + + ]
[ + + ][ + + ]
122 [ - + ][ - + ]: 506 : 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 : 69 : 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 : 414 : std::size_t nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
175 : 414 : 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 : 138 : 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 [ + - ]: 69 : 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 [ + - ]: 69 : tk::initialize( m_ncomp, m_mat_blk, L, inpoel, coord,
236 : : Problem::initialize, unk, t, nielem );
237 : :
238 : 69 : 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 : 69 : std::vector< tk::real > s(m_ncomp, 0.0);
248 [ + + ]: 23411 : for (std::size_t e=0; e<nielem; ++e) {
249 : : // inside user-defined box
250 [ + + ]: 23342 : if (!icbox.empty()) {
251 : : std::size_t bcnt = 0;
252 [ + + ]: 1792 : for (const auto& b : icbox) { // for all boxes
253 [ + - ][ + + ]: 1792 : if (inbox.size() > bcnt && inbox[bcnt].find(e) != inbox[bcnt].end())
254 : : {
255 [ + - ][ - - ]: 448 : 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 : 448 : auto V_ex = (box[1]-box[0]) * (box[3]-box[2]) * (box[5]-box[4]);
260 [ + + ]: 4480 : for (std::size_t c=0; c<m_ncomp; ++c) {
261 : 4032 : auto mark = c*rdof;
262 : 4032 : s[c] = unk(e,mark);
263 : : // set high-order DOFs to zero
264 [ - + ]: 4032 : for (std::size_t i=1; i<rdof; ++i)
265 : 0 : unk(e,mark+i) = 0.0;
266 : : }
267 [ + - ]: 448 : initializeBox<ctr::boxList>( m_mat_blk, V_ex, t, b, bgpre,
268 : : bgtemp, s );
269 : : // store box-initialization in solution vector
270 [ + + ]: 4480 : for (std::size_t c=0; c<m_ncomp; ++c) {
271 : 4032 : auto mark = c*rdof;
272 : 4032 : unk(e,mark) = s[c];
273 : : }
274 : : }
275 : 896 : ++bcnt;
276 : : }
277 : : }
278 : :
279 : : // inside user-specified mesh blocks
280 [ - + ]: 23342 : 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 : 69 : }
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 : 179 : 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 : 179 : std::size_t rdof = g_inputdeck.get< tag::rdof >();
311 : 179 : std::size_t nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
312 [ + + ]: 69861 : for (std::size_t e=0; e<nelem; ++e)
313 : 69682 : densityConstr[e] = 0.0;
314 [ + + ]: 552 : for (std::size_t imat=0; imat<nmat; ++imat)
315 [ - + ]: 373 : 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 [ + + ]: 161595 : for (std::size_t e=0; e<nelem; ++e)
334 : : {
335 : : // Retrieve alpha and add it to the constraint measure
336 : 161222 : tk::real alpha = unk(e, volfracDofIdx(nmat, imat, rdof, 0));
337 : 161222 : densityConstr[e] += alpha;
338 : : }
339 : : }
340 : 179 : }
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 : 69 : 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 : 69 : 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 : 6705 : 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 : 6705 : auto intsharp =
365 : : g_inputdeck.get< tag::multimat, tag::intsharp >();
366 : : // If this cell is not material interface, return this function
367 [ + + ]: 6705 : 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 : 6774 : 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 : : const std::vector< std::size_t >& ndofel ) const
449 : : {
450 : 6774 : const auto rdof = g_inputdeck.get< tag::rdof >();
451 : 6774 : const auto ndof = g_inputdeck.get< tag::ndof >();
452 : 6774 : 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 [ + + ]: 2527796 : for (std::size_t e=0; e<nielem; ++e)
463 : : {
464 : 2521022 : std::vector< tk::real > R(m_nprim*ndof, 0.0);
465 : :
466 [ + - ]: 2521022 : 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 [ + - ]: 2521022 : coordgp[0].resize( ng );
473 [ + - ]: 2521022 : coordgp[1].resize( ng );
474 [ + - ]: 2521022 : coordgp[2].resize( ng );
475 [ + - ]: 2521022 : wgp.resize( ng );
476 : :
477 [ + - ]: 2521022 : tk::GaussQuadratureTet( ng, coordgp, wgp );
478 : :
479 : : // Local degree of freedom
480 : 2521022 : auto dof_el = ndofel[e];
481 : :
482 : : // Loop over quadrature points in element e
483 [ + + ]: 7073484 : for (std::size_t igp=0; igp<ng; ++igp)
484 : : {
485 : : // Compute the basis function
486 [ + - ]: 4552462 : auto B =
487 [ + - ]: 4552462 : tk::eval_basis( dof_el, coordgp[0][igp], coordgp[1][igp], coordgp[2][igp] );
488 : :
489 [ + - ]: 4552462 : auto w = wgp[igp] * geoElem(e, 0);
490 : :
491 [ + - ]: 4552462 : 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 [ + + ]: 14757572 : for (std::size_t k=0; k<nmat; ++k)
496 : 10205110 : rhob += state[densityIdx(nmat, k)];
497 : :
498 : : // velocity vector at quadrature point
499 : : std::array< tk::real, 3 >
500 [ + - ]: 4552462 : vel{ state[momentumIdx(nmat, 0)]/rhob,
501 : 4552462 : state[momentumIdx(nmat, 1)]/rhob,
502 : 4552462 : state[momentumIdx(nmat, 2)]/rhob };
503 : :
504 [ + - ][ - - ]: 4552462 : std::vector< tk::real > pri(m_nprim, 0.0);
505 : :
506 : : // Evaluate material pressure at quadrature point
507 [ + + ]: 14757572 : for(std::size_t imat = 0; imat < nmat; imat++)
508 : : {
509 [ + - ]: 10205110 : auto alphamat = state[volfracIdx(nmat, imat)];
510 [ + - ]: 10205110 : auto arhomat = state[densityIdx(nmat, imat)];
511 : 10205110 : auto arhoemat = state[energyIdx(nmat, imat)];
512 [ + - ]: 10205110 : auto gmat = getDeformGrad(nmat, imat, state);
513 [ + - ]: 10205110 : pri[pressureIdx(nmat,imat)] = m_mat_blk[imat].compute<
514 [ + - ]: 10205110 : EOS::pressure >( arhomat, vel[0], vel[1], vel[2], arhoemat,
515 : : alphamat, imat, gmat );
516 : :
517 [ + - ][ - + ]: 10205110 : pri[pressureIdx(nmat,imat)] = constrain_pressure( m_mat_blk,
518 : : pri[pressureIdx(nmat,imat)], arhomat, alphamat, imat);
519 : :
520 [ - + ]: 10205110 : 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 [ + + ]: 18209848 : for (std::size_t idir=0; idir<3; ++idir) {
536 : 13657386 : pri[velocityIdx(nmat,idir)] = vel[idir];
537 : : }
538 : :
539 [ + + ]: 28414958 : for(std::size_t k = 0; k < m_nprim; k++)
540 : : {
541 : 23862496 : auto mark = k * ndof;
542 [ + + ]: 85814492 : for(std::size_t idof = 0; idof < dof_el; idof++)
543 : 61951996 : R[mark+idof] += w * pri[k] * B[idof];
544 : : }
545 : : }
546 : :
547 : : // Update the DG solution of primitive variables
548 [ + + ]: 16226318 : for(std::size_t k = 0; k < m_nprim; k++)
549 : : {
550 : 13705296 : auto mark = k * ndof;
551 : 13705296 : auto rmark = k * rdof;
552 [ + + ]: 35028492 : for(std::size_t idof = 0; idof < dof_el; idof++)
553 : : {
554 [ + + ]: 21323196 : prim(e, rmark+idof) = R[mark+idof] / L(e, mark+idof);
555 [ + + ]: 21323196 : if(fabs(prim(e, rmark+idof)) < 1e-16)
556 : 5296062 : prim(e, rmark+idof) = 0;
557 : : }
558 : : }
559 : : }
560 : 6774 : }
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 : 6705 : 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 : 6705 : 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 : 6705 : auto neg_density = cleanTraceMultiMat(t, nielem, m_mat_blk, geoElem, nmat,
593 : : unk, prim);
594 : :
595 [ - + ][ - - ]: 6705 : if (neg_density) Throw("Negative partial density.");
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ]
596 : 6705 : }
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 : : Assert( P.nprop() == rdof*m_nprim, "Number of components in primitive "
634 : : "vector must equal "+ std::to_string(rdof*m_nprim) );
635 : :
636 : : //----- reconstruction of conserved quantities -----
637 : : //--------------------------------------------------
638 : :
639 [ + + ]: 845460 : for (std::size_t e=0; e<nelem; ++e)
640 : : {
641 : : // 1. specify how many variables need to be reconstructed
642 : : std::vector< std::size_t > vars;
643 : : // for p-adaptive DG
644 [ - + ]: 841380 : if (pref) {
645 : : // If DG is applied, reconstruct only volume fractions
646 [ - - ]: 0 : if(ndofel[e] > 1) {
647 [ - - ][ - - ]: 0 : for (std::size_t k=0; k<nmat; ++k) vars.push_back(volfracIdx(nmat, k));
648 : : }
649 : : else // If P0P1 is applied for this element
650 [ - - ][ - - ]: 0 : for (std::size_t c=0; c<m_ncomp; ++c) vars.push_back(c);
651 : : }
652 : : else {
653 : : // for P0P1, reconstruct all variables
654 [ + + ]: 841380 : if (is_p0p1)
655 [ + + ][ + - ]: 3411000 : for (std::size_t c=0; c<m_ncomp; ++c) vars.push_back(c);
656 : : // for high-order DG, reconstruct only volume fractions
657 [ + - ]: 500280 : else if (ndof > 1)
658 [ + + ][ + - ]: 1500840 : for (std::size_t k=0; k<nmat; ++k) vars.push_back(volfracIdx(nmat, k));
[ - - ]
659 : : }
660 : :
661 : : // 2. solve 3x3 least-squares system
662 : : // Reconstruct second-order dofs of volume-fractions in Taylor space
663 : : // using nodal-stencils, for a good interface-normal estimate
664 [ + - ]: 841380 : tk::recoLeastSqExtStencil( rdof, e, esup, inpoel, geoElem, U, vars );
665 : :
666 : : // 3. transform reconstructed derivatives to Dubiner dofs
667 [ + - ]: 841380 : tk::transform_P0P1( rdof, e, inpoel, coord, U, vars );
668 : : }
669 : :
670 : : //----- reconstruction of primitive quantities -----
671 : : //--------------------------------------------------
672 : : // For multimat, conserved and primitive quantities are reconstructed
673 : : // separately.
674 : :
675 [ + + ]: 845460 : for (std::size_t e=0; e<nelem; ++e)
676 : : {
677 : : // There are two conditions that requires the reconstruction of the
678 : : // primitive variables:
679 : : // 1. p-adaptive is triggered and P0P1 scheme is applied to specific
680 : : // elements
681 : : // 2. p-adaptive is not triggered and P0P1 scheme is applied to the
682 : : // whole computation domain
683 [ - + ][ - - ]: 841380 : if ((pref && ndofel[e] == 1) || (!pref && is_p0p1)) {
[ + + ]
684 : : std::vector< std::size_t > vars;
685 [ + + ][ + - ]: 2046600 : for (std::size_t c=0; c<m_nprim; ++c) vars.push_back(c);
686 : :
687 : : // 1.
688 : : // Reconstruct second-order dofs of volume-fractions in Taylor space
689 : : // using nodal-stencils, for a good interface-normal estimate
690 [ + - ]: 341100 : tk::recoLeastSqExtStencil( rdof, e, esup, inpoel, geoElem, P, vars );
691 : :
692 : : // 2.
693 [ + - ]: 341100 : tk::transform_P0P1(rdof, e, inpoel, coord, P, vars);
694 : : }
695 : : }
696 : :
697 : 4080 : }
698 : :
699 : : //! Limit second-order solution, and primitive quantities separately
700 : : //! \param[in] t Physical time
701 : : //! \param[in] pref Indicator for p-adaptive algorithm
702 : : //! \param[in] geoFace Face geometry array
703 : : //! \param[in] geoElem Element geometry array
704 : : //! \param[in] fd Face connectivity and boundary conditions object
705 : : //! \param[in] esup Elements-surrounding-nodes connectivity
706 : : //! \param[in] inpoel Element-node connectivity
707 : : //! \param[in] coord Array of nodal coordinates
708 : : //! \param[in] ndofel Vector of local number of degrees of freedome
709 : : //! \param[in] gid Local->global node id map
710 : : //! \param[in] bid Local chare-boundary node ids (value) associated to
711 : : //! global node ids (key)
712 : : //! \param[in] uNodalExtrm Chare-boundary nodal extrema for conservative
713 : : //! variables
714 : : //! \param[in] pNodalExtrm Chare-boundary nodal extrema for primitive
715 : : //! variables
716 : : //! \param[in] mtInv Inverse of Taylor mass matrix
717 : : //! \param[in,out] U Solution vector at recent time step
718 : : //! \param[in,out] P Vector of primitives at recent time step
719 : : //! \param[in,out] shockmarker Vector of shock-marker values
720 : 4080 : void limit( [[maybe_unused]] tk::real t,
721 : : const bool pref,
722 : : const tk::Fields& geoFace,
723 : : const tk::Fields& geoElem,
724 : : const inciter::FaceData& fd,
725 : : const std::map< std::size_t, std::vector< std::size_t > >& esup,
726 : : const std::vector< std::size_t >& inpoel,
727 : : const tk::UnsMesh::Coords& coord,
728 : : const std::vector< std::size_t >& ndofel,
729 : : const std::vector< std::size_t >& gid,
730 : : const std::unordered_map< std::size_t, std::size_t >& bid,
731 : : const std::vector< std::vector<tk::real> >& uNodalExtrm,
732 : : const std::vector< std::vector<tk::real> >& pNodalExtrm,
733 : : const std::vector< std::vector<tk::real> >& mtInv,
734 : : tk::Fields& U,
735 : : tk::Fields& P,
736 : : std::vector< std::size_t >& shockmarker ) const
737 : : {
738 : : Assert( U.nunk() == P.nunk(), "Number of unknowns in solution "
739 : : "vector and primitive vector at recent time step incorrect" );
740 : :
741 : 4080 : const auto limiter = g_inputdeck.get< tag::limiter >();
742 : 4080 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
743 : 4080 : const auto rdof = g_inputdeck.get< tag::rdof >();
744 : : const auto& solidx = g_inputdeck.get<
745 : : tag::matidxmap, tag::solidx >();
746 : :
747 : : // limit vectors of conserved and primitive quantities
748 [ - + ]: 4080 : if (limiter == ctr::LimiterType::SUPERBEEP1)
749 : : {
750 : 0 : SuperbeeMultiMat_P1( fd.Esuel(), inpoel, ndofel,
751 : : coord, solidx, U, P, nmat );
752 : : }
753 [ + - ]: 4080 : else if (limiter == ctr::LimiterType::VERTEXBASEDP1 && rdof == 4)
754 : : {
755 : 4080 : VertexBasedMultiMat_P1( esup, inpoel, ndofel, fd.Esuel().size()/4,
756 [ + - ]: 4080 : m_mat_blk, fd, geoFace, geoElem, coord, flux, solidx, U, P,
757 : : nmat, shockmarker );
758 : : }
759 [ - - ]: 0 : else if (limiter == ctr::LimiterType::VERTEXBASEDP1 && rdof == 10)
760 : : {
761 : 0 : VertexBasedMultiMat_P2( pref, esup, inpoel, ndofel, fd.Esuel().size()/4,
762 [ - - ]: 0 : m_mat_blk, fd, geoFace, geoElem, coord, gid, bid,
763 : : uNodalExtrm, pNodalExtrm, mtInv, flux, solidx, U, P, nmat,
764 : : shockmarker );
765 : : }
766 [ - - ]: 0 : else if (limiter != ctr::LimiterType::NOLIMITER)
767 : : {
768 [ - - ][ - - ]: 0 : Throw("Limiter type not configured for multimat.");
[ - - ][ - - ]
[ - - ][ - - ]
769 : : }
770 : 4080 : }
771 : :
772 : : //! Apply CPL to the conservative variable solution for this PDE system
773 : : //! \param[in] prim Array of primitive variables
774 : : //! \param[in] geoElem Element geometry array
775 : : //! \param[in] inpoel Element-node connectivity
776 : : //! \param[in] coord Array of nodal coordinates
777 : : //! \param[in,out] unk Array of conservative variables
778 : : //! \param[in] nielem Number of internal elements
779 : : //! \details This function applies CPL to obtain consistent dofs for
780 : : //! conservative quantities based on the limited primitive quantities.
781 : : //! See Pandare et al. (2023). On the Design of Stable,
782 : : //! Consistent, and Conservative High-Order Methods for Multi-Material
783 : : //! Hydrodynamics. J Comp Phys, 112313.
784 : : void CPL( const tk::Fields& prim,
785 : : const tk::Fields& geoElem,
786 : : const std::vector< std::size_t >& inpoel,
787 : : const tk::UnsMesh::Coords& coord,
788 : : tk::Fields& unk,
789 : : std::size_t nielem ) const
790 : : {
791 : : [[maybe_unused]] const auto rdof = g_inputdeck.get< tag::rdof >();
792 : 30 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
793 : :
794 : : Assert( unk.nunk() == prim.nunk(), "Number of unknowns in solution "
795 : : "vector and primitive vector at recent time step incorrect" );
796 : : Assert( unk.nprop() == rdof*m_ncomp, "Number of components in solution "
797 : : "vector must equal "+ std::to_string(rdof*m_ncomp) );
798 : : Assert( prim.nprop() == rdof*m_nprim, "Number of components in vector of "
799 : : "primitive quantities must equal "+ std::to_string(rdof*m_nprim) );
800 : :
801 : 30 : correctLimConservMultiMat(nielem, m_mat_blk, nmat, inpoel,
802 : : coord, geoElem, prim, unk);
803 : : }
804 : :
805 : : //! Return cell-average deformation gradient tensor
806 : : //! \param[in] unk Solution vector at recent time step
807 : : //! \param[in] nielem Number of internal elements
808 : : //! \details This function returns the bulk cell-average inverse
809 : : //! deformation gradient tensor
810 : 0 : std::array< std::vector< tk::real >, 9 > cellAvgDeformGrad(
811 : : const tk::Fields& unk,
812 : : std::size_t nielem ) const
813 : : {
814 : 0 : const auto rdof = g_inputdeck.get< tag::rdof >();
815 : 0 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
816 : : const auto& solidx = g_inputdeck.get< tag::matidxmap, tag::solidx >();
817 : :
818 : : std::array< std::vector< tk::real >, 9 > gb;
819 [ - - ][ - - ]: 0 : if (inciter::haveSolid(nmat, solidx)) {
820 [ - - ]: 0 : for (auto& gij : gb)
821 [ - - ]: 0 : gij.resize(nielem, 0.0);
822 [ - - ]: 0 : for (std::size_t e=0; e<nielem; ++e) {
823 [ - - ]: 0 : for (std::size_t k=0; k<nmat; ++k) {
824 [ - - ]: 0 : if (solidx[k] > 0) {
825 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
826 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
827 : 0 : gb[3*i+j][e] += unk(e, volfracDofIdx(nmat,k,rdof,0)) *
828 : : unk(e,deformDofIdx(nmat,solidx[k],i,j,rdof,0));
829 : : }
830 : : }
831 : : }
832 : : }
833 : :
834 : 0 : return gb;
835 : : }
836 : :
837 : :
838 : : //! Reset the high order solution for p-adaptive scheme
839 : : //! \param[in] fd Face connectivity and boundary conditions object
840 : : //! \param[in,out] unk Solution vector at recent time step
841 : : //! \param[in,out] prim Primitive vector at recent time step
842 : : //! \param[in] ndofel Vector of local number of degrees of freedome
843 : : //! \details This function reset the high order coefficient for p-adaptive
844 : : //! solution polynomials. Unlike compflow class, the high order of fv
845 : : //! solution will not be reset since p0p1 is the base scheme for
846 : : //! multi-material p-adaptive DG method.
847 : 0 : void resetAdapSol( const inciter::FaceData& fd,
848 : : tk::Fields& unk,
849 : : tk::Fields& prim,
850 : : const std::vector< std::size_t >& ndofel ) const
851 : : {
852 : 0 : const auto rdof = g_inputdeck.get< tag::rdof >();
853 : 0 : const auto ncomp = unk.nprop() / rdof;
854 : 0 : const auto nprim = prim.nprop() / rdof;
855 : :
856 [ - - ]: 0 : for(std::size_t e = 0; e < fd.Esuel().size()/4; e++)
857 : : {
858 [ - - ]: 0 : if(ndofel[e] < 10)
859 : : {
860 [ - - ]: 0 : for (std::size_t c=0; c<ncomp; ++c)
861 : : {
862 : 0 : auto mark = c*rdof;
863 : 0 : unk(e, mark+4) = 0.0;
864 : 0 : unk(e, mark+5) = 0.0;
865 : 0 : unk(e, mark+6) = 0.0;
866 : 0 : unk(e, mark+7) = 0.0;
867 : 0 : unk(e, mark+8) = 0.0;
868 : 0 : unk(e, mark+9) = 0.0;
869 : : }
870 [ - - ]: 0 : for (std::size_t c=0; c<nprim; ++c)
871 : : {
872 : 0 : auto mark = c*rdof;
873 : 0 : prim(e, mark+4) = 0.0;
874 : 0 : prim(e, mark+5) = 0.0;
875 : 0 : prim(e, mark+6) = 0.0;
876 : 0 : prim(e, mark+7) = 0.0;
877 : 0 : prim(e, mark+8) = 0.0;
878 : 0 : prim(e, mark+9) = 0.0;
879 : : }
880 : : }
881 : : }
882 : 0 : }
883 : :
884 : : //! Compute right hand side
885 : : //! \param[in] t Physical time
886 : : //! \param[in] pref Indicator for p-adaptive algorithm
887 : : //! \param[in] geoFace Face geometry array
888 : : //! \param[in] geoElem Element geometry array
889 : : //! \param[in] fd Face connectivity and boundary conditions object
890 : : //! \param[in] inpoel Element-node connectivity
891 : : //! \param[in] coord Array of nodal coordinates
892 : : //! \param[in] U Solution vector at recent time step
893 : : //! \param[in] P Primitive vector at recent time step
894 : : //! \param[in] ndofel Vector of local number of degrees of freedom
895 : : //! \param[in] dt Delta time
896 : : //! \param[in,out] R Right-hand side vector computed
897 : 6705 : void rhs( tk::real t,
898 : : const bool pref,
899 : : const tk::Fields& geoFace,
900 : : const tk::Fields& geoElem,
901 : : const inciter::FaceData& fd,
902 : : const std::vector< std::size_t >& inpoel,
903 : : const std::vector< std::unordered_set< std::size_t > >&,
904 : : const tk::UnsMesh::Coords& coord,
905 : : const tk::Fields& U,
906 : : const tk::Fields& P,
907 : : const std::vector< std::size_t >& ndofel,
908 : : const tk::real dt,
909 : : tk::Fields& R ) const
910 : : {
911 : 6705 : const auto ndof = g_inputdeck.get< tag::ndof >();
912 : 6705 : const auto rdof = g_inputdeck.get< tag::rdof >();
913 : 6705 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
914 : 6705 : const auto intsharp =
915 : : g_inputdeck.get< tag::multimat, tag::intsharp >();
916 : : const auto& solidx = inciter::g_inputdeck.get<
917 : : tag::matidxmap, tag::solidx >();
918 : 6705 : auto nsld = numSolids(nmat, solidx);
919 : :
920 : 6705 : const auto nelem = fd.Esuel().size()/4;
921 : :
922 : : Assert( U.nunk() == P.nunk(), "Number of unknowns in solution "
923 : : "vector and primitive vector at recent time step incorrect" );
924 : : Assert( U.nunk() == R.nunk(), "Number of unknowns in solution "
925 : : "vector and right-hand side at recent time step incorrect" );
926 : : Assert( U.nprop() == rdof*m_ncomp, "Number of components in solution "
927 : : "vector must equal "+ std::to_string(rdof*m_ncomp) );
928 : : Assert( P.nprop() == rdof*m_nprim, "Number of components in primitive "
929 : : "vector must equal "+ std::to_string(rdof*m_nprim) );
930 : : Assert( R.nprop() == ndof*m_ncomp, "Number of components in right-hand "
931 : : "side vector must equal "+ std::to_string(ndof*m_ncomp) );
932 : : Assert( fd.Inpofa().size()/3 == fd.Esuf().size()/2,
933 : : "Mismatch in inpofa size" );
934 : :
935 : : // set rhs to zero
936 : : R.fill(0.0);
937 : :
938 : : // Allocate space for Riemann derivatives used in non-conservative terms.
939 : : // The following Riemann derivatives are stored, in order:
940 : : // 1) 3*nmat terms: derivatives of partial pressure of each material,
941 : : // for the energy equations.
942 : : // 2) ndof terms: derivatives of Riemann velocity times the basis
943 : : // function, for the volume fraction equations.
944 : : // 3) nmat*3*3*9 terms: 3 derivatives of u_l*g_ij for each material, for
945 : : // the deformation gradient equations.
946 : : // 4) 3*nsld terms: 3 derivatives of \alpha \sigma_ij for each solid
947 : : // material, for the energy equations.
948 : : // 5) 27*nsld terms: all combinations of d(g_il)/d(x_j) - d(g_ij)/d(x_l)
949 : : // for each solid material, for the deformation equations.
950 : : std::vector< std::vector< tk::real > >
951 [ + - ][ + - ]: 20115 : riemannDeriv(3*nmat+ndof+3*nsld+27*nsld, std::vector<tk::real>(U.nunk(),0.0));
[ + - ]
952 : :
953 : : // configure a no-op lambda for prescribed velocity
954 : : auto velfn = []( ncomp_t, tk::real, tk::real, tk::real, tk::real ){
955 : : return tk::VelFn::result_type(); };
956 : :
957 : : // compute internal surface flux integrals
958 [ + - ]: 6705 : tk::surfInt( pref, nmat, m_mat_blk, t, ndof, rdof, inpoel, solidx,
959 [ + - ]: 6705 : coord, fd, geoFace, geoElem, m_riemann, velfn, U, P, ndofel,
960 : : dt, R, riemannDeriv, intsharp );
961 : :
962 : : // compute optional source term
963 [ + - ]: 6705 : tk::srcInt( m_mat_blk, t, ndof, fd.Esuel().size()/4, inpoel,
964 : : coord, geoElem, Problem::src, ndofel, R, nmat );
965 : :
966 [ + + ]: 6705 : if(ndof > 1)
967 : : // compute volume integrals
968 [ + - ][ + - ]: 2340 : tk::volInt( nmat, t, m_mat_blk, ndof, rdof, nelem,
[ - - ]
969 : : inpoel, coord, geoElem, flux, velfn, U, P, ndofel, R,
970 : : intsharp );
971 : :
972 : : // compute boundary surface flux integrals
973 [ + + ]: 60345 : for (const auto& b : m_bc)
974 [ + - ]: 107280 : tk::bndSurfInt( pref, nmat, m_mat_blk, ndof, rdof,
975 : : std::get<0>(b), fd, geoFace, geoElem, inpoel, coord, t,
976 : : m_riemann, velfn, std::get<1>(b), U, P, ndofel, R,
977 : : riemannDeriv, intsharp );
978 : :
979 : : Assert( riemannDeriv.size() == 3*nmat+ndof+3*nsld+27*nsld, "Size of "
980 : : "Riemann derivative vector incorrect" );
981 : :
982 : : // get derivatives from riemannDeriv
983 [ + + ]: 58230 : for (std::size_t k=0; k<riemannDeriv.size(); ++k)
984 : : {
985 : : Assert( riemannDeriv[k].size() == U.nunk(), "Riemann derivative vector "
986 : : "for non-conservative terms has incorrect size" );
987 [ + + ]: 33037725 : for (std::size_t e=0; e<U.nunk(); ++e)
988 : 32986200 : riemannDeriv[k][e] /= geoElem(e, 0);
989 : : }
990 : :
991 : : // compute volume integrals of non-conservative terms
992 [ + - ]: 6705 : tk::nonConservativeInt( pref, nmat, m_mat_blk, ndof, rdof, nelem,
993 : : inpoel, coord, geoElem, U, P, riemannDeriv,
994 : : ndofel, R, intsharp );
995 : :
996 : : // Compute integrals for inverse deformation correction in solid materials
997 [ + - ][ - + ]: 6705 : if (inciter::haveSolid(nmat, solidx) &&
[ - - ]
998 : : g_inputdeck.get< tag::multimat, tag::rho0constraint >())
999 [ - - ]: 0 : tk::solidTermsVolInt( nmat, m_mat_blk, ndof, rdof, nelem,
1000 : : inpoel, coord, geoElem, U, P, ndofel,
1001 : : dt, R);
1002 : :
1003 : : // compute finite pressure relaxation terms
1004 [ + + ]: 6705 : if (g_inputdeck.get< tag::multimat, tag::prelax >())
1005 : : {
1006 : 375 : const auto ct = g_inputdeck.get< tag::multimat,
1007 : : tag::prelax_timescale >();
1008 [ + - ]: 375 : tk::pressureRelaxationInt( pref, nmat, m_mat_blk, ndof,
1009 : : rdof, nelem, inpoel, coord, geoElem, U, P,
1010 : : ndofel, ct, R, intsharp );
1011 : : }
1012 : 6705 : }
1013 : :
1014 : : //! Evaluate the adaptive indicator and mark the ndof for each element
1015 : : //! \param[in] nunk Number of unknowns
1016 : : //! \param[in] coord Array of nodal coordinates
1017 : : //! \param[in] inpoel Element-node connectivity
1018 : : //! \param[in] fd Face connectivity and boundary conditions object
1019 : : //! \param[in] unk Array of unknowns
1020 : : //! \param[in] prim Array of primitive quantities
1021 : : //! \param[in] indicator p-refinement indicator type
1022 : : //! \param[in] ndof Number of degrees of freedom in the solution
1023 : : //! \param[in] ndofmax Max number of degrees of freedom for p-refinement
1024 : : //! \param[in] tolref Tolerance for p-refinement
1025 : : //! \param[in,out] ndofel Vector of local number of degrees of freedome
1026 [ - - ]: 0 : void eval_ndof( std::size_t nunk,
1027 : : [[maybe_unused]] const tk::UnsMesh::Coords& coord,
1028 : : [[maybe_unused]] const std::vector< std::size_t >& inpoel,
1029 : : const inciter::FaceData& fd,
1030 : : const tk::Fields& unk,
1031 : : [[maybe_unused]] const tk::Fields& prim,
1032 : : inciter::ctr::PrefIndicatorType indicator,
1033 : : std::size_t ndof,
1034 : : std::size_t ndofmax,
1035 : : tk::real tolref,
1036 : : std::vector< std::size_t >& ndofel ) const
1037 : : {
1038 : : const auto& esuel = fd.Esuel();
1039 : 0 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1040 : :
1041 [ - - ]: 0 : if(indicator == inciter::ctr::PrefIndicatorType::SPECTRAL_DECAY)
1042 : 0 : spectral_decay(nmat, nunk, esuel, unk, ndof, ndofmax, tolref, ndofel);
1043 : : else
1044 [ - - ][ - - ]: 0 : Throw( "No such adaptive indicator type" );
[ - - ][ - - ]
[ - - ][ - - ]
1045 : 0 : }
1046 : :
1047 : : //! Compute the minimum time step size
1048 : : //! \param[in] fd Face connectivity and boundary conditions object
1049 : : //! \param[in] geoFace Face geometry array
1050 : : //! \param[in] geoElem Element geometry array
1051 : : // //! \param[in] ndofel Vector of local number of degrees of freedom
1052 : : //! \param[in] U Solution vector at recent time step
1053 : : //! \param[in] P Vector of primitive quantities at recent time step
1054 : : //! \param[in] nielem Number of internal elements
1055 : : //! \return Minimum time step size
1056 : : //! \details The allowable dt is calculated by looking at the maximum
1057 : : //! wave-speed in elements surrounding each face, times the area of that
1058 : : //! face. Once the maximum of this quantity over the mesh is determined,
1059 : : //! the volume of each cell is divided by this quantity. A minimum of this
1060 : : //! ratio is found over the entire mesh, which gives the allowable dt.
1061 : : tk::real dt( const std::array< std::vector< tk::real >, 3 >&,
1062 : : const std::vector< std::size_t >&,
1063 : : const inciter::FaceData& fd,
1064 : : const tk::Fields& geoFace,
1065 : : const tk::Fields& geoElem,
1066 : : const std::vector< std::size_t >& /*ndofel*/,
1067 : : const tk::Fields& U,
1068 : : const tk::Fields& P,
1069 : : const std::size_t nielem ) const
1070 : : {
1071 : : const auto ndof = g_inputdeck.get< tag::ndof >();
1072 : : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1073 : :
1074 : : auto mindt = timeStepSizeMultiMat( m_mat_blk, fd.Esuf(), geoFace, geoElem,
1075 : : nielem, nmat, U, P);
1076 : :
1077 : : tk::real dgp = 0.0;
1078 : : if (ndof == 4)
1079 : : {
1080 : : dgp = 1.0;
1081 : : }
1082 : : else if (ndof == 10)
1083 : : {
1084 : : dgp = 2.0;
1085 : : }
1086 : :
1087 : : // Scale smallest dt with CFL coefficient and the CFL is scaled by (2*p+1)
1088 : : // where p is the order of the DG polynomial by linear stability theory.
1089 : : mindt /= (2.0*dgp + 1.0);
1090 : : return mindt;
1091 : : }
1092 : :
1093 : : //! Compute stiff terms for a single element
1094 : : //! \param[in] e Element number
1095 : : //! \param[in] geoElem Element geometry array
1096 : : //! \param[in] inpoel Element-node connectivity
1097 : : //! \param[in] coord Array of nodal coordinates
1098 : : //! \param[in] U Solution vector at recent time step
1099 : : //! \param[in] P Primitive vector at recent time step
1100 : : //! \param[in] ndofel Vector of local number of degrees of freedom
1101 : : //! \param[in,out] R Right-hand side vector computed
1102 : 0 : void stiff_rhs( std::size_t e,
1103 : : const tk::Fields& geoElem,
1104 : : const std::vector< std::size_t >& inpoel,
1105 : : const tk::UnsMesh::Coords& coord,
1106 : : const tk::Fields& U,
1107 : : const tk::Fields& P,
1108 : : const std::vector< std::size_t >& ndofel,
1109 : : tk::Fields& R ) const
1110 : : {
1111 : 0 : const auto ndof = g_inputdeck.get< tag::ndof >();
1112 : 0 : const auto rdof = g_inputdeck.get< tag::rdof >();
1113 : 0 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1114 : 0 : const auto intsharp =
1115 : : g_inputdeck.get< tag::multimat, tag::intsharp >();
1116 : : const auto& solidx = inciter::g_inputdeck.get<
1117 : : tag::matidxmap, tag::solidx >();
1118 : :
1119 : : Assert( U.nunk() == P.nunk(), "Number of unknowns in solution "
1120 : : "vector and primitive vector at recent time step incorrect" );
1121 : : Assert( U.nprop() == rdof*m_ncomp, "Number of components in solution "
1122 : : "vector must equal "+ std::to_string(rdof*m_ncomp) );
1123 : : Assert( P.nprop() == rdof*m_nprim, "Number of components in primitive "
1124 : : "vector must equal "+ std::to_string(rdof*m_nprim) );
1125 : : Assert( R.nprop() == ndof*nstiffeq(), "Number of components in "
1126 : : "right-hand side must equal "+ std::to_string(ndof*nstiffeq()) );
1127 : :
1128 : : // set rhs to zero for element e
1129 [ - - ]: 0 : for (std::size_t i=0; i<ndof*nstiffeq(); ++i)
1130 : 0 : R(e, i) = 0.0;
1131 : :
1132 : : const auto& cx = coord[0];
1133 : : const auto& cy = coord[1];
1134 : : const auto& cz = coord[2];
1135 : :
1136 : 0 : auto ncomp = U.nprop()/rdof;
1137 : 0 : auto nprim = P.nprop()/rdof;
1138 : :
1139 : 0 : auto ng = tk::NGvol(ndofel[e]);
1140 : :
1141 : : // arrays for quadrature points
1142 : : std::array< std::vector< tk::real >, 3 > coordgp;
1143 : : std::vector< tk::real > wgp;
1144 : :
1145 [ - - ]: 0 : coordgp[0].resize( ng );
1146 [ - - ]: 0 : coordgp[1].resize( ng );
1147 [ - - ]: 0 : coordgp[2].resize( ng );
1148 [ - - ]: 0 : wgp.resize( ng );
1149 : :
1150 [ - - ]: 0 : tk::GaussQuadratureTet( ng, coordgp, wgp );
1151 : :
1152 : : // Extract the element coordinates
1153 : 0 : std::array< std::array< tk::real, 3>, 4 > coordel {{
1154 : : {{ cx[ inpoel[4*e ] ], cy[ inpoel[4*e ] ], cz[ inpoel[4*e ] ] }},
1155 : : {{ cx[ inpoel[4*e+1] ], cy[ inpoel[4*e+1] ], cz[ inpoel[4*e+1] ] }},
1156 : : {{ cx[ inpoel[4*e+2] ], cy[ inpoel[4*e+2] ], cz[ inpoel[4*e+2] ] }},
1157 : : {{ cx[ inpoel[4*e+3] ], cy[ inpoel[4*e+3] ], cz[ inpoel[4*e+3] ] }}
1158 : : }};
1159 : :
1160 : : // Gaussian quadrature
1161 [ - - ]: 0 : for (std::size_t igp=0; igp<ng; ++igp)
1162 : : {
1163 : : // Compute the coordinates of quadrature point at physical domain
1164 [ - - ]: 0 : auto gp = tk::eval_gp( igp, coordel, coordgp );
1165 : :
1166 : : // Compute the basis function
1167 [ - - ]: 0 : auto B = tk::eval_basis( ndofel[e], coordgp[0][igp], coordgp[1][igp],
1168 : : coordgp[2][igp] );
1169 : :
1170 [ - - ]: 0 : auto state = tk::evalPolynomialSol(m_mat_blk, intsharp, ncomp, nprim,
1171 : : rdof, nmat, e, ndofel[e], inpoel, coord, geoElem, gp, B, U, P);
1172 : :
1173 : : // compute source
1174 : : // Loop through materials
1175 : : std::size_t ksld = 0;
1176 [ - - ]: 0 : for (std::size_t k=0; k<nmat; ++k)
1177 : : {
1178 [ - - ]: 0 : if (solidx[k] > 0)
1179 : : {
1180 : 0 : tk::real alpha = state[inciter::volfracIdx(nmat, k)];
1181 : : std::array< std::array< tk::real, 3 >, 3 > g;
1182 : : // Compute the source terms
1183 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1184 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1185 : 0 : g[i][j] = state[inciter::deformIdx(nmat,solidx[k],i,j)];
1186 : :
1187 : : // Compute Lp
1188 : : // Reference: Ortega, A. L., Lombardini, M., Pullin, D. I., &
1189 : : // Meiron, D. I. (2014). Numerical simulation of elastic–plastic
1190 : : // solid mechanics using an Eulerian stretch tensor approach and
1191 : : // HLLD Riemann solver. Journal of Computational Physics, 257,
1192 : : // 414-441
1193 : : std::array< std::array< tk::real, 3 >, 3 > Lp;
1194 : :
1195 : : // 1. Compute dev(sigma)
1196 [ - - ]: 0 : auto sigma_dev = m_mat_blk[k].computeTensor< EOS::CauchyStress >(
1197 [ - - ]: 0 : state[inciter::densityIdx(nmat, k)], 0.0, 0.0, 0.0, 0.0, alpha, k,
1198 : : g );
1199 : 0 : tk::real apr = state[ncomp+inciter::pressureIdx(nmat, k)];
1200 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i) sigma_dev[i][i] -= apr;
1201 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1202 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1203 : 0 : sigma_dev[i][j] /= alpha;
1204 : 0 : tk::real sigma_trace =
1205 : 0 : sigma_dev[0][0]+sigma_dev[1][1]+sigma_dev[2][2];
1206 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1207 : 0 : sigma_dev[i][i] -= sigma_trace/3.0;
1208 : :
1209 : : // 2. Compute inv(g)
1210 : : double ginv[9];
1211 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1212 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1213 : 0 : ginv[3*i+j] = g[i][j];
1214 : : lapack_int ipiv[3];
1215 : : #ifndef NDEBUG
1216 : : lapack_int ierr =
1217 : : #endif
1218 [ - - ]: 0 : LAPACKE_dgetrf(LAPACK_ROW_MAJOR, 3, 3, ginv, 3, ipiv);
1219 : : Assert(ierr==0, "Lapack error in LU factorization of g");
1220 : : #ifndef NDEBUG
1221 : : lapack_int jerr =
1222 : : #endif
1223 [ - - ]: 0 : LAPACKE_dgetri(LAPACK_ROW_MAJOR, 3, ginv, 3, ipiv);
1224 : : Assert(jerr==0, "Lapack error in inverting g");
1225 : :
1226 : : // 3. Compute dev(sigma)*inv(g)
1227 : : std::array< std::array< tk::real, 3 >, 3 > aux_mat;
1228 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1229 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1230 : : {
1231 : : tk::real sum = 0.0;
1232 [ - - ]: 0 : for (std::size_t l=0; l<3; ++l)
1233 : 0 : sum += sigma_dev[i][l]*ginv[3*l+j];
1234 : 0 : aux_mat[i][j] = sum;
1235 : : }
1236 : :
1237 : : // 4. Compute g*(dev(sigma)*inv(g))
1238 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1239 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1240 : : {
1241 : : tk::real sum = 0.0;
1242 [ - - ]: 0 : for (std::size_t l=0; l<3; ++l)
1243 : 0 : sum += g[i][l]*aux_mat[l][j];
1244 : 0 : Lp[i][j] = sum;
1245 : : }
1246 : :
1247 : : // 5. Divide by 2*mu*tau
1248 : : // 'Perfect' plasticity
1249 [ - - ]: 0 : tk::real yield_stress = getmatprop< tag::yield_stress >(k);
1250 : : tk::real equiv_stress = 0.0;
1251 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1252 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1253 : 0 : equiv_stress += sigma_dev[i][j]*sigma_dev[i][j];
1254 : 0 : equiv_stress = std::sqrt(3.0*equiv_stress/2.0);
1255 : : // rel_factor = 1/tau <- Perfect plasticity for now.
1256 : : tk::real rel_factor = 0.0;
1257 [ - - ]: 0 : if (equiv_stress >= yield_stress)
1258 : : rel_factor = 1.0e07;
1259 [ - - ]: 0 : tk::real mu = getmatprop< tag::mu >(k);
1260 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1261 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1262 : 0 : Lp[i][j] *= rel_factor/(2.0*mu);
1263 : :
1264 : : // Compute the source terms
1265 [ - - ]: 0 : std::vector< tk::real > s(9*ndof, 0.0);
1266 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1267 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1268 [ - - ]: 0 : for (std::size_t idof=0; idof<ndof; ++idof)
1269 : : {
1270 : 0 : s[(i*3+j)*ndof+idof] = B[idof] * (Lp[i][0]*g[0][j]
1271 : 0 : +Lp[i][1]*g[1][j]
1272 : 0 : +Lp[i][2]*g[2][j]);
1273 : : }
1274 : :
1275 : 0 : auto wt = wgp[igp] * geoElem(e, 0);
1276 : :
1277 : : // Contribute to the right-hand-side
1278 [ - - ]: 0 : for (std::size_t i=0; i<3; ++i)
1279 [ - - ]: 0 : for (std::size_t j=0; j<3; ++j)
1280 [ - - ]: 0 : for (std::size_t idof=0; idof<ndof; ++idof)
1281 : : {
1282 : 0 : std::size_t srcId = (i*3+j)*ndof+idof;
1283 : 0 : std::size_t dofId = solidTensorIdx(ksld,i,j)*ndof+idof;
1284 : 0 : R(e, dofId) += wt * s[srcId];
1285 : : }
1286 : :
1287 [ - - ]: 0 : ksld++;
1288 : : }
1289 : : }
1290 : :
1291 : : }
1292 : 0 : }
1293 : :
1294 : : //! Extract the velocity field at cell nodes. Currently unused.
1295 : : //! \param[in] U Solution vector at recent time step
1296 : : //! \param[in] N Element node indices
1297 : : //! \return Array of the four values of the velocity field
1298 : : std::array< std::array< tk::real, 4 >, 3 >
1299 : : velocity( const tk::Fields& U,
1300 : : const std::array< std::vector< tk::real >, 3 >&,
1301 : : const std::array< std::size_t, 4 >& N ) const
1302 : : {
1303 : : const auto rdof = g_inputdeck.get< tag::rdof >();
1304 : : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1305 : :
1306 : : std::array< std::array< tk::real, 4 >, 3 > v;
1307 : : v[0] = U.extract( momentumDofIdx(nmat, 0, rdof, 0), N );
1308 : : v[1] = U.extract( momentumDofIdx(nmat, 1, rdof, 0), N );
1309 : : v[2] = U.extract( momentumDofIdx(nmat, 2, rdof, 0), N );
1310 : :
1311 : : std::vector< std::array< tk::real, 4 > > ar;
1312 : : ar.resize(nmat);
1313 : : for (std::size_t k=0; k<nmat; ++k)
1314 : : ar[k] = U.extract( densityDofIdx(nmat, k, rdof, 0), N );
1315 : :
1316 : : std::array< tk::real, 4 > r{{ 0.0, 0.0, 0.0, 0.0 }};
1317 : : for (std::size_t i=0; i<r.size(); ++i) {
1318 : : for (std::size_t k=0; k<nmat; ++k)
1319 : : r[i] += ar[k][i];
1320 : : }
1321 : :
1322 : : std::transform( r.begin(), r.end(), v[0].begin(), v[0].begin(),
1323 : : []( tk::real s, tk::real& d ){ return d /= s; } );
1324 : : std::transform( r.begin(), r.end(), v[1].begin(), v[1].begin(),
1325 : : []( tk::real s, tk::real& d ){ return d /= s; } );
1326 : : std::transform( r.begin(), r.end(), v[2].begin(), v[2].begin(),
1327 : : []( tk::real s, tk::real& d ){ return d /= s; } );
1328 : : return v;
1329 : : }
1330 : :
1331 : : //! Return a map that associates user-specified strings to functions
1332 : : //! \return Map that associates user-specified strings to functions that
1333 : : //! compute relevant quantities to be output to file
1334 : : std::map< std::string, tk::GetVarFn > OutVarFn() const
1335 : 358 : { return MultiMatOutVarFn(); }
1336 : :
1337 : : //! Return analytic field names to be output to file
1338 : : //! \return Vector of strings labelling analytic fields output in file
1339 : : std::vector< std::string > analyticFieldNames() const {
1340 : 0 : auto nmat = g_inputdeck.get< eq, tag::nmat >();
1341 : :
1342 : 0 : return MultiMatFieldNames(nmat);
1343 : : }
1344 : :
1345 : : //! Return time history field names to be output to file
1346 : : //! \return Vector of strings labelling time history fields output in file
1347 : : std::vector< std::string > histNames() const {
1348 : 0 : return MultiMatHistNames();
1349 : : }
1350 : :
1351 : : //! Return surface field output going to file
1352 : : std::vector< std::vector< tk::real > >
1353 : : surfOutput( const std::map< int, std::vector< std::size_t > >&,
1354 : : tk::Fields& ) const
1355 : : {
1356 : : std::vector< std::vector< tk::real > > s; // punt for now
1357 : : return s;
1358 : : }
1359 : :
1360 : : //! Return time history field output evaluated at time history points
1361 : : //! \param[in] h History point data
1362 : : //! \param[in] inpoel Element-node connectivity
1363 : : //! \param[in] coord Array of nodal coordinates
1364 : : //! \param[in] U Array of unknowns
1365 : : //! \param[in] P Array of primitive quantities
1366 : : //! \return Vector of time history output of bulk flow quantities (density,
1367 : : //! velocity, total energy, and pressure) evaluated at time history points
1368 : : std::vector< std::vector< tk::real > >
1369 : 0 : histOutput( const std::vector< HistData >& h,
1370 : : const std::vector< std::size_t >& inpoel,
1371 : : const tk::UnsMesh::Coords& coord,
1372 : : const tk::Fields& U,
1373 : : const tk::Fields& P ) const
1374 : : {
1375 : 0 : const auto rdof = g_inputdeck.get< tag::rdof >();
1376 : 0 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1377 : :
1378 : : const auto& x = coord[0];
1379 : : const auto& y = coord[1];
1380 : : const auto& z = coord[2];
1381 : :
1382 : 0 : std::vector< std::vector< tk::real > > Up(h.size());
1383 : :
1384 : : std::size_t j = 0;
1385 [ - - ]: 0 : for (const auto& p : h) {
1386 : 0 : auto e = p.get< tag::elem >();
1387 : 0 : auto chp = p.get< tag::coord >();
1388 : :
1389 : : // Evaluate inverse Jacobian
1390 : 0 : std::array< std::array< tk::real, 3>, 4 > cp{{
1391 : : {{ x[inpoel[4*e ]], y[inpoel[4*e ]], z[inpoel[4*e ]] }},
1392 : : {{ x[inpoel[4*e+1]], y[inpoel[4*e+1]], z[inpoel[4*e+1]] }},
1393 : : {{ x[inpoel[4*e+2]], y[inpoel[4*e+2]], z[inpoel[4*e+2]] }},
1394 : : {{ x[inpoel[4*e+3]], y[inpoel[4*e+3]], z[inpoel[4*e+3]] }} }};
1395 : 0 : auto J = tk::inverseJacobian( cp[0], cp[1], cp[2], cp[3] );
1396 : :
1397 : : // evaluate solution at history-point
1398 : 0 : std::array< tk::real, 3 > dc{{chp[0]-cp[0][0], chp[1]-cp[0][1],
1399 [ - - ]: 0 : chp[2]-cp[0][2]}};
1400 [ - - ]: 0 : auto B = tk::eval_basis(rdof, tk::dot(J[0],dc), tk::dot(J[1],dc),
1401 : : tk::dot(J[2],dc));
1402 [ - - ]: 0 : auto uhp = eval_state(m_ncomp, rdof, rdof, e, U, B);
1403 [ - - ]: 0 : auto php = eval_state(m_nprim, rdof, rdof, e, P, B);
1404 : :
1405 : : // store solution in history output vector
1406 [ - - ][ - - ]: 0 : Up[j].resize(6, 0.0);
1407 [ - - ]: 0 : for (std::size_t k=0; k<nmat; ++k) {
1408 : 0 : Up[j][0] += uhp[densityIdx(nmat,k)];
1409 : 0 : Up[j][4] += uhp[energyIdx(nmat,k)];
1410 : 0 : Up[j][5] += php[pressureIdx(nmat,k)];
1411 : : }
1412 [ - - ]: 0 : Up[j][1] = php[velocityIdx(nmat,0)];
1413 : 0 : Up[j][2] = php[velocityIdx(nmat,1)];
1414 : 0 : Up[j][3] = php[velocityIdx(nmat,2)];
1415 [ - - ]: 0 : ++j;
1416 : : }
1417 : :
1418 : 0 : return Up;
1419 : : }
1420 : :
1421 : : //! Return names of integral variables to be output to diagnostics file
1422 : : //! \return Vector of strings labelling integral variables output
1423 : : std::vector< std::string > names() const
1424 : : {
1425 : 13 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1426 : 13 : return MultiMatDiagNames(nmat);
1427 : : }
1428 : :
1429 : : //! Return analytic solution (if defined by Problem) at xi, yi, zi, t
1430 : : //! \param[in] xi X-coordinate at which to evaluate the analytic solution
1431 : : //! \param[in] yi Y-coordinate at which to evaluate the analytic solution
1432 : : //! \param[in] zi Z-coordinate at which to evaluate the analytic solution
1433 : : //! \param[in] t Physical time at which to evaluate the analytic solution
1434 : : //! \return Vector of analytic solution at given location and time
1435 : : std::vector< tk::real >
1436 : : analyticSolution( tk::real xi, tk::real yi, tk::real zi, tk::real t ) const
1437 : 0 : { return Problem::analyticSolution( m_ncomp, m_mat_blk, xi, yi, zi, t ); }
1438 : :
1439 : : //! Return analytic solution for conserved variables
1440 : : //! \param[in] xi X-coordinate at which to evaluate the analytic solution
1441 : : //! \param[in] yi Y-coordinate at which to evaluate the analytic solution
1442 : : //! \param[in] zi Z-coordinate at which to evaluate the analytic solution
1443 : : //! \param[in] t Physical time at which to evaluate the analytic solution
1444 : : //! \return Vector of analytic solution at given location and time
1445 : : std::vector< tk::real >
1446 : : solution( tk::real xi, tk::real yi, tk::real zi, tk::real t ) const
1447 : 1358616 : { return Problem::initialize( m_ncomp, m_mat_blk, xi, yi, zi, t ); }
1448 : :
1449 : : //! Return cell-averaged specific total energy for an element
1450 : : //! \param[in] e Element id for which total energy is required
1451 : : //! \param[in] unk Vector of conserved quantities
1452 : : //! \return Cell-averaged specific total energy for given element
1453 : : tk::real sp_totalenergy(std::size_t e, const tk::Fields& unk) const
1454 : : {
1455 : 858336 : const auto rdof = g_inputdeck.get< tag::rdof >();
1456 : 858336 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1457 : :
1458 : : tk::real sp_te(0.0);
1459 : : // sum each material total energy
1460 [ + + ][ + + ]: 2965588 : for (std::size_t k=0; k<nmat; ++k) {
[ + + ][ + + ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
[ - - ][ - - ]
1461 : 2107252 : sp_te += unk(e, energyDofIdx(nmat,k,rdof,0));
1462 : : }
1463 : : return sp_te;
1464 : : }
1465 : :
1466 : : private:
1467 : : //! Number of components in this PDE system
1468 : : const ncomp_t m_ncomp;
1469 : : //! Number of primitive quantities stored in this PDE system
1470 : : const ncomp_t m_nprim;
1471 : : //! Riemann solver
1472 : : tk::RiemannFluxFn m_riemann;
1473 : : //! BC configuration
1474 : : BCStateFn m_bc;
1475 : : //! EOS material block
1476 : : std::vector< EOS > m_mat_blk;
1477 : :
1478 : : //! Evaluate conservative part of physical flux function for this PDE system
1479 : : //! \param[in] ncomp Number of scalar components in this PDE system
1480 : : //! \param[in] ugp Numerical solution at the Gauss point at which to
1481 : : //! evaluate the flux
1482 : : //! \return Flux vectors for all components in this PDE system
1483 : : //! \note The function signature must follow tk::FluxFn
1484 : : static tk::FluxFn::result_type
1485 : 7374000 : flux( ncomp_t ncomp,
1486 : : const std::vector< EOS >& mat_blk,
1487 : : const std::vector< tk::real >& ugp,
1488 : : const std::vector< std::array< tk::real, 3 > >& )
1489 : : {
1490 : 7374000 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1491 : :
1492 : 7374000 : return tk::fluxTerms(ncomp, nmat, mat_blk, ugp);
1493 : : }
1494 : :
1495 : : //! \brief Boundary state function providing the left and right state of a
1496 : : //! face at Dirichlet boundaries
1497 : : //! \param[in] ncomp Number of scalar components in this PDE system
1498 : : //! \param[in] mat_blk EOS material block
1499 : : //! \param[in] ul Left (domain-internal) state
1500 : : //! \param[in] x X-coordinate at which to compute the states
1501 : : //! \param[in] y Y-coordinate at which to compute the states
1502 : : //! \param[in] z Z-coordinate at which to compute the states
1503 : : //! \param[in] t Physical time
1504 : : //! \return Left and right states for all scalar components in this PDE
1505 : : //! system
1506 : : //! \note The function signature must follow tk::StateFn. For multimat, the
1507 : : //! left or right state is the vector of conserved quantities, followed by
1508 : : //! the vector of primitive quantities appended to it.
1509 : : static tk::StateFn::result_type
1510 : 48000 : dirichlet( ncomp_t ncomp,
1511 : : const std::vector< EOS >& mat_blk,
1512 : : const std::vector< tk::real >& ul, tk::real x, tk::real y,
1513 : : tk::real z, tk::real t, const std::array< tk::real, 3 >& )
1514 : : {
1515 : 48000 : auto nmat = g_inputdeck.get< tag::multimat, tag::nmat >();
1516 : : const auto& solidx = g_inputdeck.get<
1517 : : tag::matidxmap, tag::solidx >();
1518 : :
1519 : 48000 : [[maybe_unused]] auto nsld = numSolids(nmat, solidx);
1520 : :
1521 : 48000 : auto ur = Problem::initialize( ncomp, mat_blk, x, y, z, t );
1522 : : Assert( ur.size() == ncomp, "Incorrect size for boundary state vector" );
1523 : :
1524 [ + - ]: 48000 : ur.resize(ul.size());
1525 : :
1526 : : tk::real rho(0.0);
1527 [ + + ]: 192000 : for (std::size_t k=0; k<nmat; ++k)
1528 : 144000 : rho += ur[densityIdx(nmat, k)];
1529 : :
1530 : : // get primitives in boundary state
1531 : :
1532 : : // velocity
1533 : 48000 : ur[ncomp+velocityIdx(nmat, 0)] = ur[momentumIdx(nmat, 0)] / rho;
1534 : 48000 : ur[ncomp+velocityIdx(nmat, 1)] = ur[momentumIdx(nmat, 1)] / rho;
1535 : 48000 : ur[ncomp+velocityIdx(nmat, 2)] = ur[momentumIdx(nmat, 2)] / rho;
1536 : :
1537 : : // material pressures
1538 [ + + ]: 192000 : for (std::size_t k=0; k<nmat; ++k)
1539 : : {
1540 [ + - ]: 144000 : auto gk = getDeformGrad(nmat, k, ur);
1541 [ + - ]: 144000 : ur[ncomp+pressureIdx(nmat, k)] = mat_blk[k].compute< EOS::pressure >(
1542 : : ur[densityIdx(nmat, k)], ur[ncomp+velocityIdx(nmat, 0)],
1543 : : ur[ncomp+velocityIdx(nmat, 1)], ur[ncomp+velocityIdx(nmat, 2)],
1544 [ + - ]: 144000 : ur[energyIdx(nmat, k)], ur[volfracIdx(nmat, k)], k, gk );
1545 : : }
1546 : :
1547 : : Assert( ur.size() == ncomp+nmat+3+nsld*6, "Incorrect size for appended "
1548 : : "boundary state vector" );
1549 : :
1550 [ + - ]: 48000 : return {{ std::move(ul), std::move(ur) }};
1551 : : }
1552 : :
1553 : : // Other boundary condition types that do not depend on "Problem" should be
1554 : : // added in BCFunctions.hpp
1555 : : };
1556 : :
1557 : : } // dg::
1558 : :
1559 : : } // inciter::
1560 : :
1561 : : #endif // DGMultiMat_h
|