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