Branch data Line data Source code
1 : : // *****************************************************************************
2 : : /*!
3 : : \file src/PDE/CompFlow/DGCompFlow.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 single-material flow using discontinuous Galerkin
9 : : finite elements
10 : : \details This file implements calls to the physics operators governing
11 : : compressible single-material flow using discontinuous Galerkin
12 : : discretizations.
13 : : */
14 : : // *****************************************************************************
15 : : #ifndef DGCompFlow_h
16 : : #define DGCompFlow_h
17 : :
18 : : #include <cmath>
19 : : #include <algorithm>
20 : : #include <unordered_set>
21 : : #include <map>
22 : :
23 : : #include <brigand/algorithms/for_each.hpp>
24 : :
25 : : #include "Macro.hpp"
26 : : #include "Exception.hpp"
27 : : #include "Vector.hpp"
28 : : #include "ContainerUtil.hpp"
29 : : #include "UnsMesh.hpp"
30 : : #include "Inciter/InputDeck/InputDeck.hpp"
31 : : #include "Integrate/Basis.hpp"
32 : : #include "Integrate/Quadrature.hpp"
33 : : #include "Integrate/Initialize.hpp"
34 : : #include "Integrate/Mass.hpp"
35 : : #include "Integrate/Surface.hpp"
36 : : #include "Integrate/Boundary.hpp"
37 : : #include "Integrate/Volume.hpp"
38 : : #include "Integrate/Source.hpp"
39 : : #include "RiemannChoice.hpp"
40 : : #include "EoS/EOS.hpp"
41 : : #include "Reconstruction.hpp"
42 : : #include "Limiter.hpp"
43 : : #include "PrefIndicator.hpp"
44 : :
45 : : namespace inciter {
46 : :
47 : : extern ctr::InputDeck g_inputdeck;
48 : :
49 : : namespace dg {
50 : :
51 : : //! \brief CompFlow used polymorphically with tk::DGPDE
52 : : //! \details The template arguments specify policies and are used to configure
53 : : //! the behavior of the class. The policies are:
54 : : //! - Physics - physics configuration, see PDE/CompFlow/Physics.h
55 : : //! - Problem - problem configuration, see PDE/CompFlow/Problem.h
56 : : //! \note The default physics is Euler, set in inciter::deck::check_compflow()
57 : : template< class Physics, class Problem >
58 : : class CompFlow {
59 : :
60 : : private:
61 : : using eq = tag::compflow;
62 : :
63 : : public:
64 : : //! Constructor
65 : 106 : explicit CompFlow() :
66 : : m_physics(),
67 : : m_problem(),
68 : : m_ncomp( g_inputdeck.get< tag::ncomp >() ),
69 : : m_riemann( compflowRiemannSolver(
70 [ + - ]: 106 : g_inputdeck.get< tag::flux >() ) )
71 : : {
72 : : // associate boundary condition configurations with state functions, the
73 : : // order in which the state functions listed matters, see ctr::bc::Keys
74 [ + - ][ + - ]: 1378 : brigand::for_each< ctr::bclist::Keys >( ConfigBC( m_bc,
[ + - ][ + + ]
[ + - ][ + + ]
[ + - ][ - - ]
[ - - ][ - - ]
[ - - ]
75 : : // BC State functions
76 : : { dirichlet
77 : : , symmetry
78 : : , invalidBC // Outlet BC not implemented
79 : : , farfield
80 : : , extrapolate
81 : : , invalidBC }, // No slip wall BC not implemented
82 : : // BC Gradient functions
83 : : { noOpGrad
84 : : , noOpGrad
85 : : , noOpGrad
86 : : , noOpGrad
87 : : , noOpGrad
88 : : , noOpGrad }
89 : : ) );
90 : :
91 : : // EoS initialization
92 : : const auto& matprop =
93 : : g_inputdeck.get< tag::material >();
94 : : const auto& matidxmap =
95 : : g_inputdeck.get< tag::matidxmap >();
96 [ + - ]: 106 : auto mateos = matprop[matidxmap.get< tag::eosidx >()[0]].get<tag::eos>();
97 [ + - ]: 106 : m_mat_blk.emplace_back(mateos, EqType::compflow, 0);
98 : :
99 : 106 : }
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 : : // compflow does not need/store any primitive quantities currently
107 : : return 0;
108 : : }
109 : :
110 : : //! Find the number of materials set up for this PDE system
111 : : //! \return The number of materials set up for this PDE system
112 : : std::size_t nmat() const
113 : : {
114 : : // compflow does not need nmat
115 : : return 0;
116 : : }
117 : :
118 : : //! Assign number of DOFs per equation in the PDE system
119 : : //! \param[in,out] numEqDof Array storing number of Dofs for each PDE
120 : : //! equation
121 : : void numEquationDofs(std::vector< std::size_t >& numEqDof) const
122 : : {
123 : : // all equation-dofs initialized to ndof
124 [ - - ][ + + ]: 1236 : for (std::size_t i=0; i<m_ncomp; ++i) {
[ + + ][ - - ]
[ + + ][ + + ]
[ - - ][ + + ]
[ + + ][ + + ]
[ - - ]
125 : 1030 : numEqDof.push_back(g_inputdeck.get< tag::ndof >());
126 : : }
127 : : }
128 : :
129 : : //! Determine elements that lie inside the user-defined IC box
130 : : //! \param[in] geoElem Element geometry array
131 : : //! \param[in] nielem Number of internal elements
132 : : //! \param[in,out] inbox List of nodes at which box user ICs are set for
133 : : //! each IC box
134 : : void IcBoxElems( const tk::Fields& geoElem,
135 : : std::size_t nielem,
136 : : std::vector< std::unordered_set< std::size_t > >& inbox ) const
137 : : {
138 : 206 : tk::BoxElems< eq >(geoElem, nielem, inbox);
139 : : }
140 : :
141 : : //! Find how 'stiff equations', which we currently
142 : : //! have none for Compflow
143 : : //! \return number of stiff equations
144 : : std::size_t nstiffeq() const
145 : : { return 0; }
146 : :
147 : : //! Find how 'nonstiff equations', which we currently
148 : : //! don't use for Compflow
149 : : //! \return number of non-stiff equations
150 : : std::size_t nnonstiffeq() const
151 : : { return 0; }
152 : :
153 : : //! Locate the stiff equations. Unused for compflow.
154 : : //! \param[out] stiffEqIdx list
155 : : void setStiffEqIdx( std::vector< std::size_t >& stiffEqIdx ) const
156 : : {
157 : 0 : stiffEqIdx.resize(0);
158 : : }
159 : :
160 : : //! Locate the nonstiff equations. Unused for compflow.
161 : : //! \param[out] nonStiffEqIdx list
162 : : void setNonStiffEqIdx( std::vector< std::size_t >& nonStiffEqIdx ) const
163 : : {
164 : 0 : nonStiffEqIdx.resize(0);
165 : : }
166 : :
167 : : //! Initalize the compressible flow equations, prepare for time integration
168 : : //! \param[in] L Block diagonal mass matrix
169 : : //! \param[in] inpoel Element-node connectivity
170 : : //! \param[in] coord Array of nodal coordinates
171 : : //! \param[in] inbox List of elements at which box user ICs are set for
172 : : //! each IC box
173 : : //! \param[in,out] unk Array of unknowns
174 : : //! \param[in] t Physical time
175 : : //! \param[in] nielem Number of internal elements
176 : : void
177 [ + - ]: 214 : initialize( const tk::Fields& L,
178 : : const std::vector< std::size_t >& inpoel,
179 : : const tk::UnsMesh::Coords& coord,
180 : : const std::vector< std::unordered_set< std::size_t > >& inbox,
181 : : const std::unordered_map< std::size_t,
182 : : std::set< std::size_t > >&,
183 : : tk::Fields& unk,
184 : : tk::real t,
185 : : const std::size_t nielem ) const
186 : : {
187 [ + - ]: 214 : tk::initialize( m_ncomp, m_mat_blk, L, inpoel, coord,
188 : : Problem::initialize, unk, t, nielem );
189 : :
190 : 214 : const auto rdof = g_inputdeck.get< tag::rdof >();
191 : : const auto& ic = g_inputdeck.get< tag::ic >();
192 : : const auto& icbox = ic.get< tag::box >();
193 : : const auto& bgpreic = ic.get< tag::pressure >();
194 : 214 : auto c_v = getmatprop< tag::cv >();
195 : :
196 : : // Set initial conditions inside user-defined IC box
197 : 214 : std::vector< tk::real > s(m_ncomp, 0.0);
198 [ + + ]: 93920 : for (std::size_t e=0; e<nielem; ++e) {
199 [ - + ]: 93706 : if (icbox.size() > 0) {
200 : : std::size_t bcnt = 0;
201 [ - - ]: 0 : for (const auto& b : icbox) { // for all boxes
202 [ - - ][ - - ]: 0 : if (inbox.size() > bcnt && inbox[bcnt].find(e) != inbox[bcnt].end())
203 : : {
204 [ - - ]: 0 : std::vector< tk::real > box
205 : : { b.template get< tag::xmin >(), b.template get< tag::xmax >(),
206 : : b.template get< tag::ymin >(), b.template get< tag::ymax >(),
207 : : b.template get< tag::zmin >(), b.template get< tag::zmax >() };
208 : 0 : auto V_ex = (box[1]-box[0]) * (box[3]-box[2]) * (box[5]-box[4]);
209 [ - - ]: 0 : for (std::size_t c=0; c<m_ncomp; ++c) {
210 : 0 : auto mark = c*rdof;
211 : 0 : s[c] = unk(e,mark);
212 : : // set high-order DOFs to zero
213 [ - - ]: 0 : for (std::size_t i=1; i<rdof; ++i)
214 : 0 : unk(e,mark+i) = 0.0;
215 : : }
216 [ - - ]: 0 : initializeBox<ctr::boxList>( m_mat_blk, 1.0, V_ex,
217 : : t, b, bgpreic, c_v, s );
218 : : // store box-initialization in solution vector
219 [ - - ]: 0 : for (std::size_t c=0; c<m_ncomp; ++c) {
220 : 0 : auto mark = c*rdof;
221 : 0 : unk(e,mark) = s[c];
222 : : }
223 : : }
224 : 0 : ++bcnt;
225 : : }
226 : : }
227 : : }
228 : 214 : }
229 : :
230 : : //! Compute density constraint for a given material
231 : : // //! \param[in] nelem Number of elements
232 : : // //! \param[in] unk Array of unknowns
233 : : //! \param[out] densityConstr Density Constraint: rho/(rho0*det(g))
234 : : void computeDensityConstr( std::size_t /*nelem*/,
235 : : tk::Fields& /*unk*/,
236 : : std::vector< tk::real >& densityConstr) const
237 : : {
238 : 874 : densityConstr.resize(0);
239 : : }
240 : :
241 : : //! Compute the left hand side block-diagonal mass matrix
242 : : //! \param[in] geoElem Element geometry array
243 : : //! \param[in,out] l Block diagonal mass matrix
244 : : void lhs( const tk::Fields& geoElem, tk::Fields& l ) const {
245 : 214 : const auto ndof = g_inputdeck.get< tag::ndof >();
246 : 214 : tk::mass( m_ncomp, ndof, geoElem, l );
247 : : }
248 : :
249 : : //! Update the interface cells to first order dofs
250 : : //! \details This function resets the high-order terms in interface cells,
251 : : //! and is currently not used in compflow.
252 : : void updateInterfaceCells( tk::Fields&,
253 : : std::size_t,
254 : : std::vector< std::size_t >&,
255 : : std::vector< std::size_t >& ) const {}
256 : :
257 : : //! Update the primitives for this PDE system
258 : : //! \details This function computes and stores the dofs for primitive
259 : : //! quantities, which is currently unused for compflow. But if a limiter
260 : : //! requires primitive variables for example, this would be the place to
261 : : //! add the computation of the primitive variables.
262 : : void updatePrimitives( const tk::Fields&,
263 : : const tk::Fields&,
264 : : const tk::Fields&,
265 : : tk::Fields&,
266 : : std::size_t,
267 : : std::vector< std::size_t >& ) const {}
268 : :
269 : : //! Clean up the state of trace materials for this PDE system
270 : : //! \details This function cleans up the state of materials present in trace
271 : : //! quantities in each cell. This is unused for compflow.
272 : : void cleanTraceMaterial( tk::real,
273 : : const tk::Fields&,
274 : : tk::Fields&,
275 : : tk::Fields&,
276 : : std::size_t ) const {}
277 : :
278 : : //! Reconstruct second-order solution from first-order using least-squares
279 : : // //! \param[in] t Physical time
280 : : // //! \param[in] geoFace Face geometry array
281 : : // //! \param[in] geoElem Element geometry array
282 : : // //! \param[in] fd Face connectivity and boundary conditions object
283 : : // //! \param[in] inpoel Element-node connectivity
284 : : // //! \param[in] coord Array of nodal coordinates
285 : : // //! \param[in,out] U Solution vector at recent time step
286 : : // //! \param[in,out] P Primitive vector at recent time step
287 : 14160 : void reconstruct( tk::real,
288 : : const tk::Fields&,
289 : : const tk::Fields&,
290 : : const inciter::FaceData&,
291 : : const std::map< std::size_t, std::vector< std::size_t > >&,
292 : : const std::vector< std::size_t >&,
293 : : const tk::UnsMesh::Coords&,
294 : : tk::Fields&,
295 : : tk::Fields&,
296 : : const bool,
297 : : const std::vector< std::size_t >& ) const
298 : : {
299 : : // do reconstruction only if P0P1
300 [ + + ]: 14160 : if (g_inputdeck.get< tag::rdof >() == 4 &&
301 [ - + ]: 7050 : g_inputdeck.get< tag::ndof >() == 1)
302 [ - - ][ - - ]: 0 : Throw("P0P1 not supported for CompFlow.");
[ - - ][ - - ]
[ - - ][ - - ]
303 : 14160 : }
304 : :
305 : : //! Limit second-order solution
306 : : //! \param[in] t Physical time
307 : : //! \param[in] geoFace Face geometry array
308 : : //! \param[in] geoElem Element geometry array
309 : : //! \param[in] fd Face connectivity and boundary conditions object
310 : : //! \param[in] esup Elements surrounding points
311 : : //! \param[in] inpoel Element-node connectivity
312 : : //! \param[in] coord Array of nodal coordinates
313 : : //! \param[in] ndofel Vector of local number of degrees of freedome
314 : : //! \param[in] gid Local->global node id map
315 : : //! \param[in] bid Local chare-boundary node ids (value) associated to
316 : : //! global node ids (key)
317 : : //! \param[in] uNodalExtrm Chare-boundary nodal extrema for conservative
318 : : //! variables
319 : : //! \param[in] mtInv Inverse of Taylor mass matrix
320 : : //! \param[in,out] U Solution vector at recent time step
321 : : //! \param[in,out] shockmarker Vector of shock-marker values
322 : 14160 : void limit( [[maybe_unused]] tk::real t,
323 : : [[maybe_unused]] const bool pref,
324 : : [[maybe_unused]] const tk::Fields& geoFace,
325 : : const tk::Fields& geoElem,
326 : : const inciter::FaceData& fd,
327 : : const std::map< std::size_t, std::vector< std::size_t > >& esup,
328 : : const std::vector< std::size_t >& inpoel,
329 : : const tk::UnsMesh::Coords& coord,
330 : : const std::vector< std::size_t >& ndofel,
331 : : const std::vector< std::size_t >& gid,
332 : : const std::unordered_map< std::size_t, std::size_t >& bid,
333 : : const std::vector< std::vector<tk::real> >& uNodalExtrm,
334 : : const std::vector< std::vector<tk::real> >&,
335 : : const std::vector< std::vector<tk::real> >& mtInv,
336 : : tk::Fields& U,
337 : : tk::Fields&,
338 : : std::vector< std::size_t >& shockmarker) const
339 : : {
340 : 14160 : const auto limiter = g_inputdeck.get< tag::limiter >();
341 : 14160 : const auto rdof = g_inputdeck.get< tag::rdof >();
342 : : const auto& solidx = g_inputdeck.get< tag::matidxmap, tag::solidx >();
343 : :
344 [ - + ]: 14160 : if (limiter == ctr::LimiterType::WENOP1)
345 : 0 : WENO_P1( fd.Esuel(), U );
346 [ + + ]: 14160 : else if (limiter == ctr::LimiterType::SUPERBEEP1)
347 : 1650 : Superbee_P1( fd.Esuel(), inpoel, ndofel, coord, U );
348 [ - + ]: 12510 : else if (limiter == ctr::LimiterType::VERTEXBASEDP1 && rdof == 4)
349 : 0 : VertexBasedCompflow_P1( esup, inpoel, ndofel, fd.Esuel().size()/4,
350 [ - - ]: 0 : m_mat_blk, fd, geoFace, geoElem, coord, flux, solidx, U,
351 : : shockmarker);
352 [ + + ]: 12510 : else if (limiter == ctr::LimiterType::VERTEXBASEDP1 && rdof == 10)
353 : 2580 : VertexBasedCompflow_P2( esup, inpoel, ndofel, fd.Esuel().size()/4,
354 [ + - ]: 2580 : m_mat_blk, fd, geoFace, geoElem, coord, gid, bid,
355 : : uNodalExtrm, mtInv, flux, solidx, U, shockmarker);
356 : 14160 : }
357 : :
358 : : //! Update the conservative variable solution for this PDE system
359 : : //! \details This function computes the updated dofs for conservative
360 : : //! quantities based on the limited solution and is currently not used in
361 : : //! compflow.
362 : : void CPL( const tk::Fields&,
363 : : const tk::Fields&,
364 : : const std::vector< std::size_t >&,
365 : : const tk::UnsMesh::Coords&,
366 : : tk::Fields&,
367 : : std::size_t ) const {}
368 : :
369 : : //! Return cell-average deformation gradient tensor (no-op for compflow)
370 : : //! \details This function is a no-op in compflow.
371 : : std::array< std::vector< tk::real >, 9 > cellAvgDeformGrad(
372 : : const tk::Fields&,
373 : : std::size_t ) const
374 : : {
375 : : return {};
376 : : }
377 : :
378 : : //! Reset the high order solution for p-adaptive scheme
379 : : //! \param[in] fd Face connectivity and boundary conditions object
380 : : //! \param[in,out] unk Solution vector at recent time step
381 : : //! \param[in] ndofel Vector of local number of degrees of freedome
382 : : //! \details This function reset the high order coefficient for p-adaptive
383 : : //! solution polynomials.
384 : 1770 : void resetAdapSol( const inciter::FaceData& fd,
385 : : tk::Fields& unk,
386 : : tk::Fields&,
387 : : const std::vector< std::size_t >& ndofel ) const
388 : : {
389 : 1770 : const auto rdof = g_inputdeck.get< tag::rdof >();
390 : 1770 : const auto ncomp = unk.nprop() / rdof;
391 : :
392 [ + + ]: 414490 : for(std::size_t e=0; e<fd.Esuel().size()/4; e++)
393 : : {
394 [ + + ]: 412720 : if(ndofel[e] == 1)
395 : : {
396 [ + + ]: 1671966 : for (std::size_t c=0; c<ncomp; ++c)
397 : : {
398 : 1393305 : auto mark = c*rdof;
399 : 1393305 : unk(e, mark+1) = 0.0;
400 : 1393305 : unk(e, mark+2) = 0.0;
401 : 1393305 : unk(e, mark+3) = 0.0;
402 : : }
403 [ + + ]: 134059 : } else if(ndofel[e] == 4)
404 : : {
405 [ + + ]: 274422 : for (std::size_t c=0; c<ncomp; ++c)
406 : : {
407 : 228685 : auto mark = c*rdof;
408 : 228685 : unk(e, mark+4) = 0.0;
409 : 228685 : unk(e, mark+5) = 0.0;
410 : 228685 : unk(e, mark+6) = 0.0;
411 : 228685 : unk(e, mark+7) = 0.0;
412 : 228685 : unk(e, mark+8) = 0.0;
413 : 228685 : unk(e, mark+9) = 0.0;
414 : : }
415 : : }
416 : : }
417 : 1770 : }
418 : :
419 : : //! Compute right hand side
420 : : //! \param[in] t Physical time
421 : : //! \param[in] pref Indicator for p-adaptive algorithm
422 : : //! \param[in] geoFace Face geometry array
423 : : //! \param[in] geoElem Element geometry array
424 : : //! \param[in] fd Face connectivity and boundary conditions object
425 : : //! \param[in] inpoel Element-node connectivity
426 : : //! \param[in] boxelems Mesh node ids within user-defined IC boxes
427 : : //! \param[in] coord Array of nodal coordinates
428 : : //! \param[in] U Solution vector at recent time step
429 : : //! \param[in] P Primitive vector at recent time step
430 : : //! \param[in] ndofel Vector of local number of degrees of freedom
431 : : //! \param[in] dt Delta time
432 : : //! \param[in,out] R Right-hand side vector computed
433 : 15720 : void rhs( tk::real t,
434 : : const bool pref,
435 : : const tk::Fields& geoFace,
436 : : const tk::Fields& geoElem,
437 : : const inciter::FaceData& fd,
438 : : const std::vector< std::size_t >& inpoel,
439 : : const std::vector< std::unordered_set< std::size_t > >& boxelems,
440 : : const tk::UnsMesh::Coords& coord,
441 : : const tk::Fields& U,
442 : : const tk::Fields& P,
443 : : const std::vector< std::size_t >& ndofel,
444 : : const tk::real dt,
445 : : tk::Fields& R ) const
446 : : {
447 : 15720 : const auto ndof = g_inputdeck.get< tag::ndof >();
448 : 15720 : const auto rdof = g_inputdeck.get< tag::rdof >();
449 : :
450 : : const auto& solidx = g_inputdeck.get< tag::matidxmap, tag::solidx >();
451 : :
452 : : Assert( U.nunk() == P.nunk(), "Number of unknowns in solution "
453 : : "vector and primitive vector at recent time step incorrect" );
454 : : Assert( U.nunk() == R.nunk(), "Number of unknowns in solution "
455 : : "vector and right-hand side at recent time step incorrect" );
456 : : Assert( U.nprop() == rdof*5, "Number of components in solution "
457 : : "vector must equal "+ std::to_string(rdof*5) );
458 : : Assert( P.nprop() == 0, "Number of components in primitive "
459 : : "vector must equal "+ std::to_string(0) );
460 : : Assert( R.nprop() == ndof*5, "Number of components in right-hand "
461 : : "side vector must equal "+ std::to_string(ndof*5) );
462 : : Assert( fd.Inpofa().size()/3 == fd.Esuf().size()/2,
463 : : "Mismatch in inpofa size" );
464 : :
465 : : // set rhs to zero
466 : : R.fill(0.0);
467 : :
468 : : // empty vector for non-conservative terms. This vector is unused for
469 : : // single-material hydrodynamics since, there are no non-conservative
470 : : // terms in the system of PDEs.
471 : 15720 : std::vector< std::vector < tk::real > > riemannDeriv;
472 : :
473 : : // configure a no-op lambda for prescribed velocity
474 : : auto velfn = []( ncomp_t, tk::real, tk::real, tk::real, tk::real ){
475 : : return tk::VelFn::result_type(); };
476 : :
477 : : // compute internal surface flux integrals
478 [ + - ]: 15720 : tk::surfInt( pref, 1, m_mat_blk, t, ndof, rdof, inpoel, solidx,
479 [ + - ]: 15720 : coord, fd, geoFace, geoElem, m_riemann, velfn, U, P, ndofel,
480 : : dt, R, riemannDeriv );
481 : :
482 : : // compute optional source term
483 [ + - ]: 15720 : tk::srcInt( m_mat_blk, t, ndof, fd.Esuel().size()/4,
484 : : inpoel, coord, geoElem, Problem::src, ndofel, R );
485 : :
486 [ + + ]: 15720 : if(ndof > 1)
487 : : // compute volume integrals
488 [ + - ][ + - ]: 42480 : tk::volInt( 1, t, m_mat_blk, ndof, rdof,
[ - - ]
489 [ + - ]: 14160 : fd.Esuel().size()/4, inpoel, coord, geoElem, flux, velfn,
490 : : U, P, ndofel, R );
491 : :
492 : : // compute boundary surface flux integrals
493 [ + + ]: 110040 : for (const auto& b : m_bc)
494 [ + - ]: 188640 : tk::bndSurfInt( pref, 1, m_mat_blk, ndof, rdof, std::get<0>(b),
495 : : fd, geoFace, geoElem, inpoel, coord, t, m_riemann,
496 : : velfn, std::get<1>(b), U, P, ndofel, R, riemannDeriv );
497 : :
498 : : // compute external (energy) sources
499 : : const auto& ic = g_inputdeck.get< tag::ic >();
500 : : const auto& icbox = ic.get< tag::box >();
501 : :
502 [ - + ][ - - ]: 15720 : if (!icbox.empty() && !boxelems.empty()) {
503 : : std::size_t bcnt = 0;
504 [ - - ]: 0 : for (const auto& b : icbox) { // for all boxes for this eq
505 [ - - ]: 0 : std::vector< tk::real > box
506 : : { b.template get< tag::xmin >(), b.template get< tag::xmax >(),
507 : : b.template get< tag::ymin >(), b.template get< tag::ymax >(),
508 : : b.template get< tag::zmin >(), b.template get< tag::zmax >() };
509 : :
510 : : const auto& initiate = b.template get< tag::initiate >();
511 [ - - ]: 0 : if (initiate == ctr::InitiateType::LINEAR) {
512 [ - - ]: 0 : boxSrc( t, inpoel, boxelems[bcnt], coord, geoElem, ndofel, R );
513 : : }
514 [ - - ]: 0 : ++bcnt;
515 : : }
516 : : }
517 : 15720 : }
518 : :
519 : : //! Evaluate the adaptive indicator and mark the ndof for each element
520 : : //! \param[in] nunk Number of unknowns
521 : : //! \param[in] coord Array of nodal coordinates
522 : : //! \param[in] inpoel Element-node connectivity
523 : : //! \param[in] fd Face connectivity and boundary conditions object
524 : : //! \param[in] unk Array of unknowns
525 : : //! \param[in] prim Array of primitive quantities
526 : : //! \param[in] indicator p-refinement indicator type
527 : : //! \param[in] ndof Number of degrees of freedom in the solution
528 : : //! \param[in] ndofmax Max number of degrees of freedom for p-refinement
529 : : //! \param[in] tolref Tolerance for p-refinement
530 : : //! \param[in,out] ndofel Vector of local number of degrees of freedome
531 [ + - ]: 1636 : void eval_ndof( std::size_t nunk,
532 : : const tk::UnsMesh::Coords& coord,
533 : : const std::vector< std::size_t >& inpoel,
534 : : const inciter::FaceData& fd,
535 : : const tk::Fields& unk,
536 : : [[maybe_unused]] const tk::Fields& prim,
537 : : inciter::ctr::PrefIndicatorType indicator,
538 : : std::size_t ndof,
539 : : std::size_t ndofmax,
540 : : tk::real tolref,
541 : : std::vector< std::size_t >& ndofel ) const
542 : : {
543 : : const auto& esuel = fd.Esuel();
544 : :
545 [ + - ]: 1636 : if(indicator == inciter::ctr::PrefIndicatorType::SPECTRAL_DECAY)
546 : 1636 : spectral_decay( 1, nunk, esuel, unk, ndof, ndofmax, tolref, ndofel );
547 [ - - ]: 0 : else if(indicator == inciter::ctr::PrefIndicatorType::NON_CONFORMITY)
548 : 0 : non_conformity( nunk, fd.Nbfac(), inpoel, coord, esuel, fd.Esuf(),
549 : : fd.Inpofa(), unk, ndof, ndofmax, ndofel );
550 : : else
551 [ - - ][ - - ]: 0 : Throw( "No such adaptive indicator type" );
[ - - ][ - - ]
[ - - ][ - - ]
552 : 1636 : }
553 : :
554 : : //! Compute the minimum time step size
555 : : //! \param[in] coord Mesh node coordinates
556 : : //! \param[in] inpoel Mesh element connectivity
557 : : //! \param[in] fd Face connectivity and boundary conditions object
558 : : //! \param[in] geoFace Face geometry array
559 : : //! \param[in] geoElem Element geometry array
560 : : //! \param[in] ndofel Vector of local number of degrees of freedom
561 : : //! \param[in] U Solution vector at recent time step
562 : : //! \return Minimum time step size
563 : 2070 : tk::real dt( const std::array< std::vector< tk::real >, 3 >& coord,
564 : : const std::vector< std::size_t >& inpoel,
565 : : const inciter::FaceData& fd,
566 : : const tk::Fields& geoFace,
567 : : const tk::Fields& geoElem,
568 : : const std::vector< std::size_t >& ndofel,
569 : : const tk::Fields& U,
570 : : const tk::Fields&,
571 : : const std::size_t /*nielem*/ ) const
572 : : {
573 : 2070 : const auto rdof = g_inputdeck.get< tag::rdof >();
574 : :
575 : : const auto& esuf = fd.Esuf();
576 : : const auto& inpofa = fd.Inpofa();
577 : :
578 : : tk::real rho, u, v, w, rhoE, p, a, vn, dSV_l, dSV_r;
579 : 2070 : std::vector< tk::real > delt( U.nunk(), 0.0 );
580 : :
581 : : const auto& cx = coord[0];
582 : : const auto& cy = coord[1];
583 : : const auto& cz = coord[2];
584 : :
585 : : // compute internal surface maximum characteristic speed
586 [ + + ]: 1487410 : for (std::size_t f=0; f<esuf.size()/2; ++f)
587 : : {
588 : :
589 [ + + ]: 1485340 : std::size_t el = static_cast< std::size_t >(esuf[2*f]);
590 : 1485340 : auto er = esuf[2*f+1];
591 : :
592 : : // Number of quadrature points for face integration
593 : : std::size_t ng;
594 : :
595 [ + + ]: 1485340 : if(er > -1)
596 : : {
597 : 1101260 : auto eR = static_cast< std::size_t >( er );
598 : :
599 [ + - ]: 1101260 : auto ng_l = tk::NGfa(ndofel[el]);
600 [ + - ][ + + ]: 1101260 : auto ng_r = tk::NGfa(ndofel[eR]);
601 : :
602 : : // When the number of gauss points for the left and right element are
603 : : // different, choose the larger ng
604 : 1101260 : ng = std::max( ng_l, ng_r );
605 : : }
606 : : else
607 : : {
608 [ + - ]: 384080 : ng = tk::NGfa(ndofel[el]);
609 : : }
610 : :
611 : : // arrays for quadrature points
612 : : std::array< std::vector< tk::real >, 2 > coordgp;
613 : : std::vector< tk::real > wgp;
614 : :
615 [ + - ]: 1485340 : coordgp[0].resize( ng );
616 [ + - ]: 1485340 : coordgp[1].resize( ng );
617 [ + - ]: 1485340 : wgp.resize( ng );
618 : :
619 : : // get quadrature point weights and coordinates for triangle
620 [ + - ]: 1485340 : tk::GaussQuadratureTri( ng, coordgp, wgp );
621 : :
622 : : // Extract the left element coordinates
623 : 1485340 : std::array< std::array< tk::real, 3>, 4 > coordel_l {{
624 : : {{ cx[inpoel[4*el ]], cy[inpoel[4*el ]], cz[inpoel[4*el ]] }},
625 : : {{ cx[inpoel[4*el+1]], cy[inpoel[4*el+1]], cz[inpoel[4*el+1]] }},
626 : : {{ cx[inpoel[4*el+2]], cy[inpoel[4*el+2]], cz[inpoel[4*el+2]] }},
627 : : {{ cx[inpoel[4*el+3]], cy[inpoel[4*el+3]], cz[inpoel[4*el+3]] }} }};
628 : :
629 : : // Compute the determinant of Jacobian matrix
630 : : auto detT_l =
631 : 1485340 : tk::Jacobian(coordel_l[0], coordel_l[1], coordel_l[2], coordel_l[3]);
632 : :
633 : : // Extract the face coordinates
634 : 1485340 : std::array< std::array< tk::real, 3>, 3 > coordfa {{
635 : : {{ cx[ inpofa[3*f ] ], cy[ inpofa[3*f ] ], cz[ inpofa[3*f ] ] }},
636 : : {{ cx[ inpofa[3*f+1] ], cy[ inpofa[3*f+1] ], cz[ inpofa[3*f+1] ] }},
637 : : {{ cx[ inpofa[3*f+2] ], cy[ inpofa[3*f+2] ], cz[ inpofa[3*f+2] ] }}
638 : : }};
639 : :
640 : 1485340 : dSV_l = 0.0;
641 : 1485340 : dSV_r = 0.0;
642 : :
643 : : // Gaussian quadrature
644 [ + + ]: 5637109 : for (std::size_t igp=0; igp<ng; ++igp)
645 : : {
646 : : // Compute the coordinates of quadrature point at physical domain
647 [ + - ]: 4151769 : auto gp = tk::eval_gp( igp, coordfa, coordgp );
648 : :
649 : : // Compute the basis function for the left element
650 : 4151769 : auto B_l = tk::eval_basis( ndofel[el],
651 : 4151769 : tk::Jacobian(coordel_l[0], gp, coordel_l[2], coordel_l[3])/detT_l,
652 : 4151769 : tk::Jacobian(coordel_l[0], coordel_l[1], gp, coordel_l[3])/detT_l,
653 [ + - ]: 4151769 : tk::Jacobian(coordel_l[0], coordel_l[1], coordel_l[2], gp)/detT_l );
654 : :
655 : 4151769 : auto wt = wgp[igp] * geoFace(f,0);
656 : :
657 : : std::array< std::vector< tk::real >, 2 > ugp;
658 : :
659 : : // left element
660 [ + + ]: 24910614 : for (ncomp_t c=0; c<5; ++c)
661 : : {
662 [ + - ]: 20758845 : auto mark = c*rdof;
663 [ + - ]: 20758845 : ugp[0].push_back( U(el, mark) );
664 : :
665 [ + + ]: 20758845 : if(ndofel[el] > 1) //DG(P1)
666 : 17318295 : ugp[0][c] += U(el, mark+1) * B_l[1]
667 : 17318295 : + U(el, mark+2) * B_l[2]
668 : 17318295 : + U(el, mark+3) * B_l[3];
669 : :
670 [ + + ]: 20758845 : if(ndofel[el] > 4) //DG(P2)
671 : 10487970 : ugp[0][c] += U(el, mark+4) * B_l[4]
672 : 10487970 : + U(el, mark+5) * B_l[5]
673 : 10487970 : + U(el, mark+6) * B_l[6]
674 : 10487970 : + U(el, mark+7) * B_l[7]
675 : 10487970 : + U(el, mark+8) * B_l[8]
676 : 10487970 : + U(el, mark+9) * B_l[9];
677 : : }
678 : :
679 : 4151769 : rho = ugp[0][0];
680 : 4151769 : u = ugp[0][1]/rho;
681 : 4151769 : v = ugp[0][2]/rho;
682 : 4151769 : w = ugp[0][3]/rho;
683 : 4151769 : rhoE = ugp[0][4];
684 [ + - ]: 4151769 : p = m_mat_blk[0].compute< EOS::pressure >( rho, u, v, w, rhoE );
685 : :
686 [ + - ]: 4151769 : a = m_mat_blk[0].compute< EOS::soundspeed >( rho, p );
687 : :
688 : 4151769 : vn = u*geoFace(f,1) + v*geoFace(f,2) + w*geoFace(f,3);
689 : :
690 : 4151769 : dSV_l = wt * (std::fabs(vn) + a);
691 : :
692 : : // right element
693 [ + + ]: 4151769 : if (er > -1) {
694 : :
695 : : // nodal coordinates of the right element
696 : 3082564 : std::size_t eR = static_cast< std::size_t >( er );
697 : :
698 : : // Extract the left element coordinates
699 [ + - ]: 3082564 : std::array< std::array< tk::real, 3>, 4 > coordel_r {{
700 : : {{ cx[inpoel[4*eR ]], cy[inpoel[4*eR ]], cz[inpoel[4*eR ]] }},
701 : : {{ cx[inpoel[4*eR+1]], cy[inpoel[4*eR+1]], cz[inpoel[4*eR+1]] }},
702 : : {{ cx[inpoel[4*eR+2]], cy[inpoel[4*eR+2]], cz[inpoel[4*eR+2]] }},
703 : : {{ cx[inpoel[4*eR+3]], cy[inpoel[4*eR+3]], cz[inpoel[4*eR+3]] }}
704 : : }};
705 : :
706 : : // Compute the determinant of Jacobian matrix
707 : : auto detT_r =
708 : 3082564 : tk::Jacobian(coordel_r[0],coordel_r[1],coordel_r[2],coordel_r[3]);
709 : :
710 : : // Compute the coordinates of quadrature point at physical domain
711 [ + - ]: 3082564 : gp = tk::eval_gp( igp, coordfa, coordgp );
712 : :
713 : : // Compute the basis function for the right element
714 : 3082564 : auto B_r = tk::eval_basis( ndofel[eR],
715 : 3082564 : tk::Jacobian(coordel_r[0],gp,coordel_r[2],coordel_r[3])/detT_r,
716 : 3082564 : tk::Jacobian(coordel_r[0],coordel_r[1],gp,coordel_r[3])/detT_r,
717 [ + - ]: 3082564 : tk::Jacobian(coordel_r[0],coordel_r[1],coordel_r[2],gp)/detT_r );
718 : :
719 [ + + ]: 18495384 : for (ncomp_t c=0; c<5; ++c)
720 : : {
721 [ + - ]: 15412820 : auto mark = c*rdof;
722 [ + - ]: 15412820 : ugp[1].push_back( U(eR, mark) );
723 : :
724 [ + + ]: 15412820 : if(ndofel[eR] > 1) //DG(P1)
725 : 12829620 : ugp[1][c] += U(eR, mark+1) * B_r[1]
726 : 12829620 : + U(eR, mark+2) * B_r[2]
727 : 12829620 : + U(eR, mark+3) * B_r[3];
728 : :
729 [ + + ]: 15412820 : if(ndofel[eR] > 4) //DG(P2)
730 : 7864410 : ugp[1][c] += U(eR, mark+4) * B_r[4]
731 : 7864410 : + U(eR, mark+5) * B_r[5]
732 : 7864410 : + U(eR, mark+6) * B_r[6]
733 : 7864410 : + U(eR, mark+7) * B_r[7]
734 : 7864410 : + U(eR, mark+8) * B_r[8]
735 : 7864410 : + U(eR, mark+9) * B_r[9];
736 : : }
737 : :
738 : 3082564 : rho = ugp[1][0];
739 : 3082564 : u = ugp[1][1]/rho;
740 : 3082564 : v = ugp[1][2]/rho;
741 : 3082564 : w = ugp[1][3]/rho;
742 : 3082564 : rhoE = ugp[1][4];
743 [ + - ]: 3082564 : p = m_mat_blk[0].compute< EOS::pressure >( rho, u, v, w, rhoE );
744 [ + - ]: 3082564 : a = m_mat_blk[0].compute< EOS::soundspeed >( rho, p );
745 : :
746 : 3082564 : vn = u*geoFace(f,1) + v*geoFace(f,2) + w*geoFace(f,3);
747 : :
748 [ + + ]: 3082564 : dSV_r = wt * (std::fabs(vn) + a);
749 [ + - ]: 3082564 : delt[eR] += std::max( dSV_l, dSV_r );
750 : : }
751 : :
752 : 4151769 : delt[el] += std::max( dSV_l, dSV_r );
753 : : }
754 : : }
755 : :
756 : 2070 : tk::real mindt = std::numeric_limits< tk::real >::max();
757 : : tk::real dgp = 0.0;
758 : :
759 : : // compute allowable dt
760 [ + + ]: 618910 : for (std::size_t e=0; e<fd.Esuel().size()/4; ++e)
761 : : {
762 : : dgp = 0.0;
763 [ + + ]: 616840 : if (ndofel[e] == 4)
764 : : {
765 : : dgp = 1.0;
766 : : }
767 [ + + ]: 425383 : else if (ndofel[e] == 10)
768 : : {
769 : : dgp = 2.0;
770 : : }
771 : :
772 : : // Scale smallest dt with CFL coefficient and the CFL is scaled by (2*p+1)
773 : : // where p is the order of the DG polynomial by linear stability theory.
774 [ + + ]: 627915 : mindt = std::min( mindt, geoElem(e,0)/ (delt[e] * (2.0*dgp + 1.0)) );
775 : : }
776 : :
777 [ + - ]: 4140 : return mindt;
778 : : }
779 : :
780 : : //! Compute stiff terms for a single element, not implemented here
781 : : // //! \param[in] e Element number
782 : : // //! \param[in] geoElem Element geometry array
783 : : // //! \param[in] inpoel Element-node connectivity
784 : : // //! \param[in] coord Array of nodal coordinates
785 : : // //! \param[in] U Solution vector at recent time step
786 : : // //! \param[in] P Primitive vector at recent time step
787 : : // //! \param[in] ndofel Vector of local number of degrees of freedom
788 : : // //! \param[in,out] R Right-hand side vector computed
789 : : void stiff_rhs( std::size_t /*e*/,
790 : : const tk::Fields& /*geoElem*/,
791 : : const std::vector< std::size_t >& /*inpoel*/,
792 : : const tk::UnsMesh::Coords& /*coord*/,
793 : : const tk::Fields& /*U*/,
794 : : const tk::Fields& /*P*/,
795 : : const std::vector< std::size_t >& /*ndofel*/,
796 : : tk::Fields& /*R*/ ) const
797 : : {}
798 : :
799 : : //! Extract the velocity field at cell nodes. Currently unused.
800 : : //! \param[in] U Solution vector at recent time step
801 : : //! \param[in] N Element node indices
802 : : //! \return Array of the four values of the velocity field
803 : : std::array< std::array< tk::real, 4 >, 3 >
804 : : velocity( const tk::Fields& U,
805 : : const std::array< std::vector< tk::real >, 3 >&,
806 : : const std::array< std::size_t, 4 >& N ) const
807 : : {
808 : : std::array< std::array< tk::real, 4 >, 3 > v;
809 : : v[0] = U.extract( 1, N );
810 : : v[1] = U.extract( 2, N );
811 : : v[2] = U.extract( 3, N );
812 : : auto r = U.extract( 0, N );
813 : : std::transform( r.begin(), r.end(), v[0].begin(), v[0].begin(),
814 : : []( tk::real s, tk::real& d ){ return d /= s; } );
815 : : std::transform( r.begin(), r.end(), v[1].begin(), v[1].begin(),
816 : : []( tk::real s, tk::real& d ){ return d /= s; } );
817 : : std::transform( r.begin(), r.end(), v[2].begin(), v[2].begin(),
818 : : []( tk::real s, tk::real& d ){ return d /= s; } );
819 : : return v;
820 : : }
821 : :
822 : : //! Return a map that associates user-specified strings to functions
823 : : //! \return Map that associates user-specified strings to functions that
824 : : //! compute relevant quantities to be output to file
825 : : std::map< std::string, tk::GetVarFn > OutVarFn() const
826 : 1748 : { return CompFlowOutVarFn(); }
827 : :
828 : : //! Return analytic field names to be output to file
829 : : //! \return Vector of strings labelling analytic fields output in file
830 : : std::vector< std::string > analyticFieldNames() const
831 : 588 : { return m_problem.analyticFieldNames( m_ncomp ); }
832 : :
833 : : //! Return time history field names to be output to file
834 : : //! \return Vector of strings labeling time history fields output in file
835 : : std::vector< std::string > histNames() const
836 : 0 : { return CompFlowHistNames(); }
837 : :
838 : : //! Return surface field output going to file
839 : : std::vector< std::vector< tk::real > >
840 : : surfOutput( const std::map< int, std::vector< std::size_t > >&,
841 : : tk::Fields& ) const
842 : : {
843 : : std::vector< std::vector< tk::real > > s; // punt for now
844 : : return s;
845 : : }
846 : :
847 : : //! Return time history field output evaluated at time history points
848 : : //! \param[in] h History point data
849 : : //! \param[in] inpoel Element-node connectivity
850 : : //! \param[in] coord Array of nodal coordinates
851 : : //! \param[in] U Array of unknowns
852 : : std::vector< std::vector< tk::real > >
853 : 0 : histOutput( const std::vector< HistData >& h,
854 : : const std::vector< std::size_t >& inpoel,
855 : : const tk::UnsMesh::Coords& coord,
856 : : const tk::Fields& U,
857 : : const tk::Fields& ) const
858 : : {
859 : 0 : const auto rdof = g_inputdeck.get< tag::rdof >();
860 : :
861 : : const auto& x = coord[0];
862 : : const auto& y = coord[1];
863 : : const auto& z = coord[2];
864 : :
865 : 0 : std::vector< std::vector< tk::real > > Up(h.size());
866 : :
867 : : std::size_t j = 0;
868 [ - - ]: 0 : for (const auto& p : h) {
869 : 0 : auto e = p.get< tag::elem >();
870 : 0 : auto chp = p.get< tag::coord >();
871 : :
872 : : // Evaluate inverse Jacobian
873 : 0 : std::array< std::array< tk::real, 3>, 4 > cp{{
874 : : {{ x[inpoel[4*e ]], y[inpoel[4*e ]], z[inpoel[4*e ]] }},
875 : : {{ x[inpoel[4*e+1]], y[inpoel[4*e+1]], z[inpoel[4*e+1]] }},
876 : : {{ x[inpoel[4*e+2]], y[inpoel[4*e+2]], z[inpoel[4*e+2]] }},
877 : : {{ x[inpoel[4*e+3]], y[inpoel[4*e+3]], z[inpoel[4*e+3]] }} }};
878 : 0 : auto J = tk::inverseJacobian( cp[0], cp[1], cp[2], cp[3] );
879 : :
880 : : // evaluate solution at history-point
881 : 0 : std::array< tk::real, 3 > dc{{chp[0]-cp[0][0], chp[1]-cp[0][1],
882 [ - - ]: 0 : chp[2]-cp[0][2]}};
883 [ - - ]: 0 : auto B = tk::eval_basis(rdof, tk::dot(J[0],dc), tk::dot(J[1],dc),
884 : : tk::dot(J[2],dc));
885 [ - - ]: 0 : auto uhp = eval_state(m_ncomp, rdof, rdof, e, U, B);
886 : :
887 : : // store solution in history output vector
888 [ - - ]: 0 : Up[j].resize(6, 0.0);
889 [ - - ]: 0 : Up[j][0] = uhp[0];
890 : 0 : Up[j][1] = uhp[1]/uhp[0];
891 : 0 : Up[j][2] = uhp[2]/uhp[0];
892 : 0 : Up[j][3] = uhp[3]/uhp[0];
893 : 0 : Up[j][4] = uhp[4]/uhp[0];
894 [ - - ]: 0 : Up[j][5] = m_mat_blk[0].compute< EOS::pressure >( uhp[0], uhp[1]/uhp[0],
895 [ - - ][ - - ]: 0 : uhp[2]/uhp[0], uhp[3]/uhp[0], uhp[4] );
896 [ - - ]: 0 : ++j;
897 : : }
898 : :
899 : 0 : return Up;
900 : : }
901 : :
902 : : //! Return names of integral variables to be output to diagnostics file
903 : : //! \return Vector of strings labelling integral variables output
904 : : std::vector< std::string > names() const
905 : 32 : { return m_problem.names( m_ncomp ); }
906 : :
907 : : //! Return analytic solution (if defined by Problem) at xi, yi, zi, t
908 : : //! \param[in] xi X-coordinate at which to evaluate the analytic solution
909 : : //! \param[in] yi Y-coordinate at which to evaluate the analytic solution
910 : : //! \param[in] zi Z-coordinate at which to evaluate the analytic solution
911 : : //! \param[in] t Physical time at which to evaluate the analytic solution
912 : : //! \return Vector of analytic solution at given location and time
913 : : std::vector< tk::real >
914 : : analyticSolution( tk::real xi, tk::real yi, tk::real zi, tk::real t ) const
915 : 348546 : { return Problem::analyticSolution( m_ncomp, m_mat_blk, xi, yi,
916 : 348546 : zi, t ); }
917 : :
918 : : //! Return analytic solution for conserved variables
919 : : //! \param[in] xi X-coordinate at which to evaluate the analytic solution
920 : : //! \param[in] yi Y-coordinate at which to evaluate the analytic solution
921 : : //! \param[in] zi Z-coordinate at which to evaluate the analytic solution
922 : : //! \param[in] t Physical time at which to evaluate the analytic solution
923 : : //! \return Vector of analytic solution at given location and time
924 : : std::vector< tk::real >
925 : : solution( tk::real xi, tk::real yi, tk::real zi, tk::real t ) const
926 : 2469998 : { return Problem::initialize( m_ncomp, m_mat_blk, xi, yi, zi, t ); }
927 : :
928 : : //! Return cell-averaged specific total energy for an element
929 : : //! \param[in] e Element id for which total energy is required
930 : : //! \param[in] unk Vector of conserved quantities
931 : : //! \return Cell-averaged specific total energy for given element
932 : : tk::real sp_totalenergy(std::size_t e, const tk::Fields& unk) const
933 : : {
934 : 1015102 : const auto rdof = g_inputdeck.get< tag::rdof >();
935 : :
936 : 1015102 : return unk(e,4*rdof);
937 : : }
938 : :
939 : : private:
940 : : //! Physics policy
941 : : const Physics m_physics;
942 : : //! Problem policy
943 : : const Problem m_problem;
944 : : //! Number of components in this PDE system
945 : : const ncomp_t m_ncomp;
946 : : //! Riemann solver
947 : : tk::RiemannFluxFn m_riemann;
948 : : //! BC configuration
949 : : BCStateFn m_bc;
950 : : //! EOS material block
951 : : std::vector< EOS > m_mat_blk;
952 : :
953 : : //! Evaluate physical flux function for this PDE system
954 : : //! \param[in] ncomp Number of scalar components in this PDE system
955 : : //! \param[in] mat_blk EOS material block
956 : : //! \param[in] ugp Numerical solution at the Gauss point at which to
957 : : //! evaluate the flux
958 : : //! \return Flux vectors for all components in this PDE system
959 : : //! \note The function signature must follow tk::FluxFn
960 : : static tk::FluxFn::result_type
961 : 17110593 : flux( [[maybe_unused]] ncomp_t ncomp,
962 : : const std::vector< EOS >& mat_blk,
963 : : const std::vector< tk::real >& ugp,
964 : : const std::vector< std::array< tk::real, 3 > >& )
965 : : {
966 : : Assert( ugp.size() == ncomp, "Size mismatch" );
967 : :
968 : 17110593 : auto u = ugp[1] / ugp[0];
969 : 17110593 : auto v = ugp[2] / ugp[0];
970 : 17110593 : auto w = ugp[3] / ugp[0];
971 : 17110593 : auto p = mat_blk[0].compute< EOS::pressure >( ugp[0], u, v, w, ugp[4] );
972 : :
973 : 17110593 : std::vector< std::array< tk::real, 3 > > fl( ugp.size() );
974 : :
975 : 17110593 : fl[0][0] = ugp[1];
976 : 17110593 : fl[1][0] = ugp[1] * u + p;
977 : 17110593 : fl[2][0] = ugp[1] * v;
978 : 17110593 : fl[3][0] = ugp[1] * w;
979 : 17110593 : fl[4][0] = u * (ugp[4] + p);
980 : :
981 : 17110593 : fl[0][1] = ugp[2];
982 : 17110593 : fl[1][1] = ugp[2] * u;
983 : 17110593 : fl[2][1] = ugp[2] * v + p;
984 : 17110593 : fl[3][1] = ugp[2] * w;
985 : 17110593 : fl[4][1] = v * (ugp[4] + p);
986 : :
987 : 17110593 : fl[0][2] = ugp[3];
988 : 17110593 : fl[1][2] = ugp[3] * u;
989 : 17110593 : fl[2][2] = ugp[3] * v;
990 : 17110593 : fl[3][2] = ugp[3] * w + p;
991 : 17110593 : fl[4][2] = w * (ugp[4] + p);
992 : :
993 : 17110593 : return fl;
994 : : }
995 : :
996 : : //! \brief Boundary state function providing the left and right state of a
997 : : //! face at Dirichlet boundaries
998 : : //! \param[in] ncomp Number of scalar components in this PDE system
999 : : //! \param[in] mat_blk EOS material block
1000 : : //! \param[in] ul Left (domain-internal) state
1001 : : //! \param[in] x X-coordinate at which to compute the states
1002 : : //! \param[in] y Y-coordinate at which to compute the states
1003 : : //! \param[in] z Z-coordinate at which to compute the states
1004 : : //! \param[in] t Physical time
1005 : : //! \return Left and right states for all scalar components in this PDE
1006 : : //! system
1007 : : //! \note The function signature must follow tk::StateFn
1008 : : static tk::StateFn::result_type
1009 : 3241770 : dirichlet( ncomp_t ncomp,
1010 : : const std::vector< EOS >& mat_blk,
1011 : : const std::vector< tk::real >& ul, tk::real x, tk::real y,
1012 : : tk::real z, tk::real t, const std::array< tk::real, 3 >& )
1013 : : {
1014 : 3241770 : return {{ ul, Problem::initialize( ncomp, mat_blk, x, y, z, t ) }};
1015 : : }
1016 : :
1017 : : //! \brief Boundary state function providing the left and right state of a
1018 : : //! face at symmetry boundaries
1019 : : //! \param[in] ul Left (domain-internal) state
1020 : : //! \param[in] fn Unit face normal
1021 : : //! \return Left and right states for all scalar components in this PDE
1022 : : //! system
1023 : : //! \note The function signature must follow tk::StateFn
1024 : : static tk::StateFn::result_type
1025 : 2867445 : symmetry( ncomp_t, const std::vector< EOS >&,
1026 : : const std::vector< tk::real >& ul, tk::real, tk::real, tk::real,
1027 : : tk::real, const std::array< tk::real, 3 >& fn )
1028 : : {
1029 : 2867445 : std::vector< tk::real > ur(5);
1030 : : // Internal cell velocity components
1031 : 2867445 : auto v1l = ul[1]/ul[0];
1032 : 2867445 : auto v2l = ul[2]/ul[0];
1033 : 2867445 : auto v3l = ul[3]/ul[0];
1034 : : // Normal component of velocity
1035 : 2867445 : auto vnl = v1l*fn[0] + v2l*fn[1] + v3l*fn[2];
1036 : : // Ghost state velocity components
1037 : 2867445 : auto v1r = v1l - 2.0*vnl*fn[0];
1038 : 2867445 : auto v2r = v2l - 2.0*vnl*fn[1];
1039 : 2867445 : auto v3r = v3l - 2.0*vnl*fn[2];
1040 : : // Boundary condition
1041 : 2867445 : ur[0] = ul[0];
1042 : 2867445 : ur[1] = ur[0] * v1r;
1043 : 2867445 : ur[2] = ur[0] * v2r;
1044 : 2867445 : ur[3] = ur[0] * v3r;
1045 : 2867445 : ur[4] = ul[4];
1046 [ + - ]: 2867445 : return {{ std::move(ul), std::move(ur) }};
1047 : : }
1048 : :
1049 : : //! \brief Boundary state function providing the left and right state of a
1050 : : //! face at farfield boundaries
1051 : : //! \param[in] mat_blk EOS material block
1052 : : //! \param[in] ul Left (domain-internal) state
1053 : : //! \param[in] fn Unit face normal
1054 : : //! \return Left and right states for all scalar components in this PDE
1055 : : //! system
1056 : : //! \note The function signature must follow tk::StateFn
1057 : : static tk::StateFn::result_type
1058 : 40800 : farfield( ncomp_t, const std::vector< EOS >& mat_blk,
1059 : : const std::vector< tk::real >& ul, tk::real, tk::real, tk::real,
1060 : : tk::real, const std::array< tk::real, 3 >& fn )
1061 : : {
1062 : : // Primitive variables from farfield
1063 : 40800 : const auto& bc = g_inputdeck.get< tag::bc >()[0];
1064 : 40800 : auto frho = bc.get< tag::density >();
1065 : 40800 : auto fp = bc.get< tag::pressure >();
1066 : : const auto& fu = bc.get< tag::velocity >();
1067 : :
1068 : : // Speed of sound from farfield
1069 : 40800 : auto fa = mat_blk[0].compute< EOS::soundspeed >( frho, fp );
1070 : :
1071 : : // Normal component from farfield
1072 : 40800 : auto fvn = fu[0]*fn[0] + fu[1]*fn[1] + fu[2]*fn[2];
1073 : :
1074 : : // Mach number from farfield
1075 : 40800 : auto fM = fvn / fa;
1076 : :
1077 : : // Specific total energy from farfield
1078 : 40800 : auto frhoE = mat_blk[0].compute< EOS::totalenergy >( frho, fu[0], fu[1],
1079 : : fu[2], fp );
1080 : :
1081 : : // Pressure from internal cell
1082 : 40800 : auto p = mat_blk[0].compute< EOS::pressure >( ul[0], ul[1]/ul[0],
1083 : 40800 : ul[2]/ul[0], ul[3]/ul[0], ul[4] );
1084 : :
1085 : 40800 : auto ur = ul;
1086 : :
1087 [ - + ]: 40800 : if(fM <= -1) // Supersonic inflow
1088 : : {
1089 : : // For supersonic inflow, all the characteristics are from outside.
1090 : : // Therefore, we calculate the ghost cell state using the primitive
1091 : : // variables from outside.
1092 : 0 : ur[0] = frho;
1093 : 0 : ur[1] = frho * fu[0];
1094 : 0 : ur[2] = frho * fu[1];
1095 : 0 : ur[3] = frho * fu[2];
1096 : 0 : ur[4] = frhoE;
1097 [ + - ][ - + ]: 40800 : } else if(fM > -1 && fM < 0) // Subsonic inflow
1098 : : {
1099 : : // For subsonic inflow, there are 1 outgoing characteristcs and 4
1100 : : // incoming characteristic. Therefore, we calculate the ghost cell state
1101 : : // by taking pressure from the internal cell and other quantities from
1102 : : // the outside.
1103 : 0 : ur[0] = frho;
1104 : 0 : ur[1] = frho * fu[0];
1105 : 0 : ur[2] = frho * fu[1];
1106 : 0 : ur[3] = frho * fu[2];
1107 [ - - ]: 0 : ur[4] = mat_blk[0].compute< EOS::totalenergy >( frho, fu[0], fu[1],
1108 : : fu[2], p );
1109 [ + - ][ + - ]: 40800 : } else if(fM >= 0 && fM < 1) // Subsonic outflow
1110 : : {
1111 : : // For subsonic outflow, there are 1 incoming characteristcs and 4
1112 : : // outgoing characteristic. Therefore, we calculate the ghost cell state
1113 : : // by taking pressure from the outside and other quantities from the
1114 : : // internal cell.
1115 : 81600 : ur[4] = mat_blk[0].compute< EOS::totalenergy >( ul[0], ul[1]/ul[0],
1116 [ + - ][ - - ]: 40800 : ul[2]/ul[0], ul[3]/ul[0], fp );
1117 : : }
1118 : : // Otherwise, for supersonic outflow, all the characteristics are from
1119 : : // internal cell. Therefore, we calculate the ghost cell state using the
1120 : : // conservative variables from outside.
1121 : :
1122 [ + - ][ + - ]: 81600 : return {{ ul, ur }};
1123 : : }
1124 : :
1125 : : //! \brief Boundary state function providing the left and right state of a
1126 : : //! face at extrapolation boundaries
1127 : : //! \param[in] ul Left (domain-internal) state
1128 : : //! \return Left and right states for all scalar components in this PDE
1129 : : //! system
1130 : : //! \note The function signature must follow tk::StateFn
1131 : : static tk::StateFn::result_type
1132 : 93360 : extrapolate( ncomp_t, const std::vector< EOS >&,
1133 : : const std::vector< tk::real >& ul, tk::real, tk::real,
1134 : : tk::real, tk::real, const std::array< tk::real, 3 >& )
1135 : : {
1136 : 93360 : return {{ ul, ul }};
1137 : : }
1138 : :
1139 : : //! \brief Boundary gradient function copying the left gradient to the right
1140 : : //! gradient at a face
1141 : : //! \param[in] dul Left (domain-internal) state
1142 : : //! \return Left and right states for all scalar components in this PDE
1143 : : //! system
1144 : : //! \note The function signature must follow tk::StateFn.
1145 : : static tk::StateFn::result_type
1146 : 0 : noOpGrad( ncomp_t,
1147 : : const std::vector< EOS >&,
1148 : : const std::vector< tk::real >& dul,
1149 : : tk::real, tk::real, tk::real, tk::real,
1150 : : const std::array< tk::real, 3 >& )
1151 : : {
1152 : 0 : return {{ dul, dul }};
1153 : : }
1154 : :
1155 : : //! Compute sources corresponding to a propagating front in user-defined box
1156 : : //! \param[in] t Physical time
1157 : : //! \param[in] inpoel Element point connectivity
1158 : : //! \param[in] boxelems Mesh node ids within user-defined box
1159 : : //! \param[in] coord Mesh node coordinates
1160 : : //! \param[in] geoElem Element geometry array
1161 : : //! \param[in] ndofel Vector of local number of degrees of freedome
1162 : : //! \param[in] R Right-hand side vector
1163 : : //! \details This function add the energy source corresponding to a planar
1164 : : //! wave-front propagating along the z-direction with a user-specified
1165 : : //! velocity, within a box initial condition, configured by the user.
1166 : : //! Example (SI) units of the quantities involved:
1167 : : //! * internal energy content (energy per unit volume): J/m^3
1168 : : //! * specific energy (internal energy per unit mass): J/kg
1169 : 0 : void boxSrc( tk::real t,
1170 : : const std::vector< std::size_t >& inpoel,
1171 : : const std::unordered_set< std::size_t >& boxelems,
1172 : : const tk::UnsMesh::Coords& coord,
1173 : : const tk::Fields& geoElem,
1174 : : const std::vector< std::size_t >& ndofel,
1175 : : tk::Fields& R ) const
1176 : : {
1177 : 0 : const auto ndof = g_inputdeck.get< tag::ndof >();
1178 : : const auto& ic = g_inputdeck.get< tag::ic >();
1179 : : const auto& icbox = ic.get< tag::box >();
1180 : :
1181 [ - - ]: 0 : for (const auto& b : icbox) { // for all boxes for this eq
1182 [ - - ]: 0 : std::vector< tk::real > box
1183 : : { b.template get< tag::xmin >(), b.template get< tag::xmax >(),
1184 : : b.template get< tag::ymin >(), b.template get< tag::ymax >(),
1185 : : b.template get< tag::zmin >(), b.template get< tag::zmax >() };
1186 : :
1187 : 0 : auto boxenc = b.template get< tag::energy_content >();
1188 : : Assert( boxenc > 0.0, "Box energy content must be nonzero" );
1189 : :
1190 : 0 : auto V_ex = (box[1]-box[0]) * (box[3]-box[2]) * (box[5]-box[4]);
1191 : :
1192 : : // determine times at which sourcing is initialized and terminated
1193 : 0 : auto iv = b.template get< tag::front_speed >();
1194 : : auto wFront = 0.1;
1195 : : auto tInit = 0.0;
1196 : 0 : auto tFinal = tInit + (box[5] - box[4] - 2.0*wFront) / std::fabs(iv);
1197 : : auto aBox = (box[1]-box[0]) * (box[3]-box[2]);
1198 : :
1199 : : const auto& cx = coord[0];
1200 : : const auto& cy = coord[1];
1201 : : const auto& cz = coord[2];
1202 : :
1203 [ - - ][ - - ]: 0 : if (t >= tInit && t <= tFinal) {
1204 : : // The energy front is assumed to have a half-sine-wave shape. The
1205 : : // half wave-length is the width of the front. At t=0, the center of
1206 : : // this front (i.e. the peak of the partial-sine-wave) is at X_0 +
1207 : : // W_0. W_0 is calculated based on the width of the front and the
1208 : : // direction of propagation (which is assumed to be along the
1209 : : // z-direction). If the front propagation velocity is positive, it
1210 : : // is assumed that the initial position of the energy source is the
1211 : : // minimum z-coordinate of the box; whereas if this velocity is
1212 : : // negative, the initial position is the maximum z-coordinate of the
1213 : : // box.
1214 : :
1215 : : // Orientation of box
1216 : 0 : std::array< tk::real, 3 > b_orientn{{
1217 : : b.template get< tag::orientation >()[0],
1218 : : b.template get< tag::orientation >()[1],
1219 : : b.template get< tag::orientation >()[2] }};
1220 : 0 : std::array< tk::real, 3 > b_centroid{{ 0.5*(box[0]+box[1]),
1221 : 0 : 0.5*(box[2]+box[3]), 0.5*(box[4]+box[5]) }};
1222 : : // Transform box to reference space
1223 : 0 : std::array< tk::real, 3 > b_min{{box[0], box[2], box[4]}};
1224 : 0 : std::array< tk::real, 3 > b_max{{box[1], box[3], box[5]}};
1225 : : tk::movePoint(b_centroid, b_min);
1226 : : tk::movePoint(b_centroid, b_max);
1227 : :
1228 : : // initial center of front
1229 : 0 : tk::real zInit(b_min[2]);
1230 [ - - ]: 0 : if (iv < 0.0) zInit = b_max[2];
1231 : : // current location of front
1232 : 0 : auto z0 = zInit + iv*t;
1233 : 0 : auto z1 = z0 + std::copysign(wFront, iv);
1234 : : tk::real s0(z0), s1(z1);
1235 : : // if velocity of propagation is negative, initial position is z1
1236 [ - - ]: 0 : if (iv < 0.0) {
1237 : : s0 = z1;
1238 : : s1 = z0;
1239 : : }
1240 : : // Sine-wave (positive part of the wave) source term amplitude
1241 : : auto pi = 4.0 * std::atan(1.0);
1242 : 0 : auto amplE = boxenc * V_ex * pi
1243 : 0 : / (aBox * wFront * 2.0 * (tFinal-tInit));
1244 : : //// Square wave (constant) source term amplitude
1245 : : //auto amplE = boxenc * V_ex
1246 : : // / (aBox * wFront * (tFinal-tInit));
1247 : :
1248 : : // add source
1249 [ - - ]: 0 : for (auto e : boxelems) {
1250 : 0 : std::array< tk::real, 3 > node{{ geoElem(e,1), geoElem(e,2),
1251 : : geoElem(e,3) }};
1252 : : // Transform node to reference space of box
1253 : : tk::movePoint(b_centroid, node);
1254 : 0 : tk::rotatePoint({{-b_orientn[0], -b_orientn[1], -b_orientn[2]}},
1255 : : node);
1256 : :
1257 [ - - ][ - - ]: 0 : if (node[2] >= s0 && node[2] <= s1) {
1258 [ - - ]: 0 : auto ng = tk::NGvol(ndofel[e]);
1259 : :
1260 : : // arrays for quadrature points
1261 : : std::array< std::vector< tk::real >, 3 > coordgp;
1262 : : std::vector< tk::real > wgp;
1263 : :
1264 [ - - ]: 0 : coordgp[0].resize( ng );
1265 [ - - ]: 0 : coordgp[1].resize( ng );
1266 [ - - ]: 0 : coordgp[2].resize( ng );
1267 [ - - ]: 0 : wgp.resize( ng );
1268 : :
1269 [ - - ]: 0 : tk::GaussQuadratureTet( ng, coordgp, wgp );
1270 : :
1271 : : // Extract the element coordinates
1272 : 0 : std::array< std::array< tk::real, 3>, 4 > coordel{{
1273 : : {{ cx[inpoel[4*e ]], cy[inpoel[4*e ]], cz[inpoel[4*e ]] }},
1274 : : {{ cx[inpoel[4*e+1]], cy[inpoel[4*e+1]], cz[inpoel[4*e+1]] }},
1275 : : {{ cx[inpoel[4*e+2]], cy[inpoel[4*e+2]], cz[inpoel[4*e+2]] }},
1276 : : {{ cx[inpoel[4*e+3]], cy[inpoel[4*e+3]], cz[inpoel[4*e+3]] }}}};
1277 : :
1278 [ - - ]: 0 : for (std::size_t igp=0; igp<ng; ++igp) {
1279 : : // Compute the coordinates of quadrature point at physical
1280 : : // domain
1281 [ - - ]: 0 : auto gp = tk::eval_gp( igp, coordel, coordgp );
1282 : :
1283 : : // Transform quadrature point to reference space of box
1284 : : tk::movePoint(b_centroid, gp);
1285 : 0 : tk::rotatePoint({{-b_orientn[0], -b_orientn[1], -b_orientn[2]}},
1286 : : gp);
1287 : :
1288 : : // Compute the basis function
1289 [ - - ]: 0 : auto B = tk::eval_basis( ndofel[e], coordgp[0][igp],
1290 : : coordgp[1][igp], coordgp[2][igp] );
1291 : :
1292 : : // Compute the source term variable
1293 [ - - ][ - - ]: 0 : std::vector< tk::real > s(5, 0.0);
1294 : 0 : s[4] = amplE * std::sin(pi*(gp[2]-s0)/wFront);
1295 : :
1296 [ - - ]: 0 : auto wt = wgp[igp] * geoElem(e, 0);
1297 : :
1298 [ - - ]: 0 : tk::update_rhs( ndof, ndofel[e], wt, e, B, s, R );
1299 : : }
1300 : : }
1301 : : }
1302 : : }
1303 : : }
1304 : 0 : }
1305 : : };
1306 : :
1307 : : } // dg::
1308 : :
1309 : : } // inciter::
1310 : :
1311 : : #endif // DGCompFlow_h
|