Line data Source code
1 1 : // Distributed under the MIT License. 2 : // See LICENSE.txt for details. 3 : 4 : /// \file 5 : /// Defines functions logical_coordinates and interface_logical_coordinates 6 : 7 : #pragma once 8 : 9 : #include <cstddef> 10 : 11 : #include "DataStructures/DataBox/Tag.hpp" 12 : #include "DataStructures/Tensor/TypeAliases.hpp" 13 : #include "Utilities/TMPL.hpp" 14 : 15 : /// \cond 16 : namespace domain { 17 : namespace Tags { 18 : template<size_t Dim> 19 : struct Mesh; 20 : template <size_t, typename> 21 : struct Coordinates; 22 : } // namespace Tags 23 : } // namespace domain 24 : 25 : template <size_t Dim> 26 : class Mesh; 27 : class DataVector; 28 : template <size_t Dim> 29 : class Direction; 30 : 31 : namespace gsl { 32 : template <typename> 33 : struct not_null; 34 : } // namespace gsl 35 : /// \endcond 36 : 37 : /// @{ 38 : /*! 39 : * \ingroup ComputationalDomainGroup 40 : * \brief Compute the logical coordinates of a Mesh in an Element. 41 : * 42 : * \details The logical coordinates are the collocation points 43 : * associated with the Spectral::Basis and Spectral::Quadrature of the \p mesh. 44 : * The Spectral::Basis determines the domain of the logical coordinates, while 45 : * the Spectral::Quadrature determines their distribution. For Legendre or 46 : * Chebyshev bases, the logical coordinates are in the interval \f$[-1, 1]\f$. 47 : * These bases may have either GaussLobatto or Gauss quadrature, which are not 48 : * uniformly distributed, and either include (GaussLobatto) or do not include 49 : * (Gauss) the end points. For the FiniteDifference basis, the logical 50 : * coordinates are again in the interval \f$[-1, 1]\f$. This basis may have 51 : * either FaceCentered or CellCentered quadrature, which are uniformly 52 : * distributed, and either include (FaceCentered) or do not include 53 : * (CellCentered) the end points. The SphericalHarmonic basis corresponds to 54 : * the spherical coordinates \f$(\theta, \phi)\f$ where the polar angle 55 : * \f$\theta\f$ is in the interval \f$[0, \pi]\f$ and the azimuth \f$\phi\f$ is 56 : * in the interval \f$[0, 2 \pi]\f$. The polar angle has Gauss quadrature 57 : * corresponding to the Legendre Gauss points of \f$\cos \theta\f$ (and thus 58 : * have no points at the poles), while the azimuth has Equiangular quadrature 59 : * which are distributed uniformly including the left endpoint, but not the 60 : * right. For the Cartoon basis, only one logical coordinate is allowed and it 61 : * is never used for computations. 62 : * 63 : * \example 64 : * \snippet Test_LogicalCoordinates.cpp logical_coordinates_example 65 : */ 66 : template <size_t VolumeDim> 67 1 : void logical_coordinates( 68 : gsl::not_null<tnsr::I<DataVector, VolumeDim, Frame::ElementLogical>*> 69 : logical_coords, 70 : const Mesh<VolumeDim>& mesh); 71 : 72 : template <size_t VolumeDim> 73 1 : tnsr::I<DataVector, VolumeDim, Frame::ElementLogical> logical_coordinates( 74 : const Mesh<VolumeDim>& mesh); 75 : /// @} 76 : 77 : namespace domain { 78 : namespace Tags { 79 : /// \ingroup DataBoxTagsGroup 80 : /// \ingroup ComputationalDomainGroup 81 : /// The logical coordinates in the Element 82 : template <size_t VolumeDim> 83 1 : struct LogicalCoordinates : Coordinates<VolumeDim, Frame::ElementLogical>, 84 : db::ComputeTag { 85 0 : using base = Coordinates<VolumeDim, Frame::ElementLogical>; 86 0 : using return_type = typename base::type; 87 0 : using argument_tags = tmpl::list<Mesh<VolumeDim>>; 88 0 : static constexpr auto function = static_cast<void (*)( 89 : gsl::not_null<return_type*>, const ::Mesh<VolumeDim>&)>( 90 : &logical_coordinates<VolumeDim>); 91 : }; 92 : } // namespace Tags 93 : } // namespace domain