Branch data Line data Source code
1 : : // ***************************************************************************** 2 : : /*! 3 : : \file src/PDE/Riemann/HLLC.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 Harten-Lax-van Leer-Contact (HLLC) Riemann flux function 9 : : \details This file implements the Harten-Lax-van Leer-Contact (HLLC) Riemann 10 : : solver. 11 : : */ 12 : : // ***************************************************************************** 13 : : #ifndef HLLC_h 14 : : #define HLLC_h 15 : : 16 : : #include <vector> 17 : : 18 : : #include "Types.hpp" 19 : : #include "Fields.hpp" 20 : : #include "Tags.hpp" 21 : : #include "FunctionPrototypes.hpp" 22 : : #include "Inciter/Options/Flux.hpp" 23 : : #include "EoS/EoS.hpp" 24 : : 25 : : namespace inciter { 26 : : 27 : : //! HLLC approximate Riemann solver 28 : : //! \details This class can be used polymorphically with inciter::RiemannSolver 29 : : struct HLLC { 30 : : 31 : : //! HLLC approximate Riemann solver flux function 32 : : //! \param[in] fn Face/Surface normal 33 : : //! \param[in] u Left and right unknown/state vector 34 : : //! \return Riemann solution according to Harten-Lax-van Leer-Contact 35 : : //! \note The function signature must follow tk::RiemannFluxFn 36 : : static tk::RiemannFluxFn::result_type 37 : 30633342 : flux( const std::array< tk::real, 3 >& fn, 38 : : const std::array< std::vector< tk::real >, 2 >& u, 39 : : const std::vector< std::array< tk::real, 3 > >& = {} ) 40 : : { 41 [ + - ]: 30633342 : std::vector< tk::real > flx( u[0].size(), 0 ); 42 : : 43 : : // Primitive variables 44 : 30633342 : auto rhol = u[0][0]; 45 : 30633342 : auto rhor = u[1][0]; 46 : : 47 : 30633342 : auto ul = u[0][1]/rhol; 48 : 30633342 : auto vl = u[0][2]/rhol; 49 : 30633342 : auto wl = u[0][3]/rhol; 50 : : 51 : 30633342 : auto ur = u[1][1]/rhor; 52 : 30633342 : auto vr = u[1][2]/rhor; 53 : 30633342 : auto wr = u[1][3]/rhor; 54 : : 55 [ + - ]: 30633342 : auto pl = eos_pressure< tag::compflow >( 0, rhol, ul, vl, wl, u[0][4] ); 56 [ + - ]: 30633342 : auto pr = eos_pressure< tag::compflow >( 0, rhor, ur, vr, wr, u[1][4] ); 57 : : 58 [ + - ]: 30633342 : auto al = eos_soundspeed< tag::compflow >( 0, rhol, pl ); 59 [ + - ]: 30633342 : auto ar = eos_soundspeed< tag::compflow >( 0, rhor, pr ); 60 : : 61 : : // Face-normal velocities 62 : 30633342 : tk::real vnl = ul*fn[0] + vl*fn[1] + wl*fn[2]; 63 : 30633342 : tk::real vnr = ur*fn[0] + vr*fn[1] + wr*fn[2]; 64 : : 65 : : // Roe-averaged variables 66 : 30633342 : auto rlr = sqrt(rhor/rhol); 67 : 30633342 : auto rlr1 = 1.0 + rlr; 68 : : 69 : 30633342 : auto vnroe = (vnr*rlr + vnl)/rlr1 ; 70 : 30633342 : auto aroe = (ar*rlr + al)/rlr1 ; 71 : : 72 : : // Signal velocities 73 : 30633342 : auto Sl = fmin(vnl-al, vnroe-aroe); 74 : 30633342 : auto Sr = fmax(vnr+ar, vnroe+aroe); 75 : 30633342 : auto Sm = ( rhor*vnr*(Sr-vnr) - rhol*vnl*(Sl-vnl) + pl-pr ) 76 : 30633342 : /( rhor*(Sr-vnr) - rhol*(Sl-vnl) ); 77 : : 78 : : // Middle-zone (star) variables 79 : 30633342 : auto pStar = rhol*(vnl-Sl)*(vnl-Sm) + pl; 80 [ + - ]: 61266684 : auto uStar = u; 81 : : 82 : 30633342 : uStar[0][0] = (Sl-vnl) * rhol/ (Sl-Sm); 83 : 30633342 : uStar[0][1] = ((Sl-vnl) * u[0][1] + (pStar-pl)*fn[0]) / (Sl-Sm); 84 : 30633342 : uStar[0][2] = ((Sl-vnl) * u[0][2] + (pStar-pl)*fn[1]) / (Sl-Sm); 85 : 30633342 : uStar[0][3] = ((Sl-vnl) * u[0][3] + (pStar-pl)*fn[2]) / (Sl-Sm); 86 : 30633342 : uStar[0][4] = ((Sl-vnl) * u[0][4] - pl*vnl + pStar*Sm) / (Sl-Sm); 87 : : 88 : 30633342 : uStar[1][0] = (Sr-vnr) * rhor/ (Sr-Sm); 89 : 30633342 : uStar[1][1] = ((Sr-vnr) * u[1][1] + (pStar-pr)*fn[0]) / (Sr-Sm); 90 : 30633342 : uStar[1][2] = ((Sr-vnr) * u[1][2] + (pStar-pr)*fn[1]) / (Sr-Sm); 91 : 30633342 : uStar[1][3] = ((Sr-vnr) * u[1][3] + (pStar-pr)*fn[2]) / (Sr-Sm); 92 : 30633342 : uStar[1][4] = ((Sr-vnr) * u[1][4] - pr*vnr + pStar*Sm) / (Sr-Sm); 93 : : 94 : : // Numerical fluxes 95 [ + + ]: 30633342 : if (Sl > 0.0) { 96 : 215529 : flx[0] = u[0][0] * vnl; 97 : 215529 : flx[1] = u[0][1] * vnl + pl*fn[0]; 98 : 215529 : flx[2] = u[0][2] * vnl + pl*fn[1]; 99 : 215529 : flx[3] = u[0][3] * vnl + pl*fn[2]; 100 : 215529 : flx[4] = ( u[0][4] + pl ) * vnl; 101 : : } 102 [ + - ][ + + ]: 30417813 : else if (Sl <= 0.0 && Sm > 0.0) { 103 : 11220200 : flx[0] = uStar[0][0] * Sm; 104 : 11220200 : flx[1] = uStar[0][1] * Sm + pStar*fn[0]; 105 : 11220200 : flx[2] = uStar[0][2] * Sm + pStar*fn[1]; 106 : 11220200 : flx[3] = uStar[0][3] * Sm + pStar*fn[2]; 107 : 11220200 : flx[4] = ( uStar[0][4] + pStar ) * Sm; 108 : : } 109 [ + - ][ + + ]: 19197613 : else if (Sm <= 0.0 && Sr >= 0.0) { 110 : 18966923 : flx[0] = uStar[1][0] * Sm; 111 : 18966923 : flx[1] = uStar[1][1] * Sm + pStar*fn[0]; 112 : 18966923 : flx[2] = uStar[1][2] * Sm + pStar*fn[1]; 113 : 18966923 : flx[3] = uStar[1][3] * Sm + pStar*fn[2]; 114 : 18966923 : flx[4] = ( uStar[1][4] + pStar ) * Sm; 115 : : } 116 : : else { 117 : 230690 : flx[0] = u[1][0] * vnr; 118 : 230690 : flx[1] = u[1][1] * vnr + pr*fn[0]; 119 : 230690 : flx[2] = u[1][2] * vnr + pr*fn[1]; 120 : 230690 : flx[3] = u[1][3] * vnr + pr*fn[2]; 121 : 230690 : flx[4] = ( u[1][4] + pr ) * vnr; 122 : : } 123 : : 124 : 61266684 : return flx; 125 : : } 126 : : 127 : : //! Flux type accessor 128 : : //! \return Flux type 129 : 80 : static ctr::FluxType type() noexcept { return ctr::FluxType::HLLC; } 130 : : }; 131 : : 132 : : } // inciter:: 133 : : 134 : : #endif // HLLC_h