Quinoa all test code coverage report
Current view: top level - PDE/CompFlow - DGCompFlow.hpp (source / functions) Hit Total Coverage
Commit: -128-NOTFOUND Lines: 230 334 68.9 %
Date: 2024-12-12 08:36:05 Functions: 61 176 34.7 %
Legend: Lines: hit not hit | Branches: + taken - not taken # not executed Branches: 130 328 39.6 %

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

Generated by: LCOV version 1.14