Line data Source code
1 0 : // Distributed under the MIT License. 2 : // See LICENSE.txt for details. 3 : 4 : #pragma once 5 : 6 : #include <memory> 7 : #include <optional> 8 : #include <pup.h> 9 : 10 : #include "DataStructures/DataBox/Prefixes.hpp" 11 : #include "DataStructures/DataVector.hpp" 12 : #include "DataStructures/Tensor/TypeAliases.hpp" 13 : #include "Evolution/Systems/ScalarAdvection/BoundaryCorrections/BoundaryCorrection.hpp" 14 : #include "Evolution/Systems/ScalarAdvection/Tags.hpp" 15 : #include "NumericalAlgorithms/DiscontinuousGalerkin/Formulation.hpp" 16 : #include "Options/String.hpp" 17 : #include "Utilities/Gsl.hpp" 18 : #include "Utilities/Serialization/CharmPupable.hpp" 19 : #include "Utilities/TMPL.hpp" 20 : 21 : namespace ScalarAdvection { 22 : namespace BoundaryCorrections { 23 : /*! 24 : * \brief A Rusanov/local Lax-Friedrichs Riemann solver 25 : * 26 : * Let \f$U\f$ be the evolved scalar variable, \f$F^i\f$ the corresponding 27 : * fluxes, and \f$n_i\f$ be the outward directed unit normal to the interface. 28 : * Denoting \f$F := n_i F^i\f$, the %Rusanov boundary correction is 29 : * 30 : * \f{align*} 31 : * G_\text{Rusanov} = \frac{F_\text{int} - F_\text{ext}}{2} - 32 : * \frac{\text{max}\left(|\lambda_\text{int}|, 33 : * |\lambda_\text{ext}|\right)}{2} \left(U_\text{ext} - U_\text{int}\right), 34 : * \f} 35 : * 36 : * where "int" and "ext" stand for interior and exterior, and 37 : * \f$\lambda\f$ is the characteristic/signal speed. The minus sign in 38 : * front of the \f$F_{\text{ext}}\f$ is necessary because the outward directed 39 : * normal of the neighboring element has the opposite sign, i.e. 40 : * \f$n_i^{\text{ext}}=-n_i^{\text{int}}\f$. 41 : * 42 : * For the ScalarAdvection system, \f$\lambda\f$ is given as 43 : * 44 : * \f{align*} 45 : * \lambda = |v| = \sqrt{v^iv_i} 46 : * \f} 47 : * 48 : * where \f$v^i\f$ is the advection velocity field. 49 : * 50 : * \note In the strong form the `dg_boundary_terms` function returns 51 : * \f$G - F_\text{int}\f$ 52 : */ 53 : template <size_t Dim> 54 1 : class Rusanov final : public BoundaryCorrection<Dim> { 55 : private: 56 0 : struct AbsCharSpeed : db::SimpleTag { 57 0 : using type = Scalar<DataVector>; 58 : }; 59 : 60 : public: 61 0 : using options = tmpl::list<>; 62 0 : static constexpr Options::String help = { 63 : "Computes the Rusanov or local Lax-Friedrichs boundary correction term " 64 : "for the ScalarAdvection system."}; 65 : 66 0 : Rusanov() = default; 67 0 : Rusanov(const Rusanov&) = default; 68 0 : Rusanov& operator=(const Rusanov&) = default; 69 0 : Rusanov(Rusanov&&) = default; 70 0 : Rusanov& operator=(Rusanov&&) = default; 71 0 : ~Rusanov() override = default; 72 : 73 : /// \cond 74 : explicit Rusanov(CkMigrateMessage* msg); 75 : using PUP::able::register_constructor; 76 : WRAPPED_PUPable_decl_template(Rusanov); // NOLINT 77 : /// \endcond 78 0 : void pup(PUP::er& p) override; // NOLINT 79 : 80 0 : std::unique_ptr<BoundaryCorrection<Dim>> get_clone() const override; 81 : 82 0 : using dg_package_field_tags = 83 : tmpl::list<Tags::U, ::Tags::NormalDotFlux<Tags::U>, AbsCharSpeed>; 84 0 : using dg_package_data_temporary_tags = tmpl::list<Tags::VelocityField<Dim>>; 85 0 : using dg_package_data_volume_tags = tmpl::list<>; 86 0 : using dg_boundary_terms_volume_tags = tmpl::list<>; 87 : 88 0 : static double dg_package_data( 89 : gsl::not_null<Scalar<DataVector>*> packaged_u, 90 : gsl::not_null<Scalar<DataVector>*> packaged_normal_dot_flux_u, 91 : gsl::not_null<Scalar<DataVector>*> packaged_abs_char_speed, 92 : const Scalar<DataVector>& u, 93 : const tnsr::I<DataVector, Dim, Frame::Inertial>& flux_u, 94 : const tnsr::I<DataVector, Dim, Frame::Inertial>& velocity_field, 95 : const tnsr::i<DataVector, Dim, Frame::Inertial>& normal_covector, 96 : const std::optional<tnsr::I<DataVector, Dim, Frame::Inertial>>& 97 : mesh_velocity, 98 : const std::optional<Scalar<DataVector>>& normal_dot_mesh_velocity); 99 : 100 0 : static void dg_boundary_terms( 101 : gsl::not_null<Scalar<DataVector>*> boundary_correction_u, 102 : const Scalar<DataVector>& u_int, 103 : const Scalar<DataVector>& normal_dot_flux_u_int, 104 : const Scalar<DataVector>& abs_char_speed_int, 105 : const Scalar<DataVector>& u_ext, 106 : const Scalar<DataVector>& normal_dot_flux_u_ext, 107 : const Scalar<DataVector>& abs_char_speed_ext, 108 : dg::Formulation dg_formulation); 109 : }; 110 : } // namespace BoundaryCorrections 111 : } // namespace ScalarAdvection