DDC 0.0.0

a discrete domain computation library

ddc Namespace Reference

Namespaces

namespace  reducer
 

Classes

class  AlignedAllocator
 
struct  cartesian_prod
 
struct  cartesian_prod< ddc::DiscreteDomain< DDim1... >, ddc::DiscreteDomain< DDim2... > >
 
struct  cartesian_prod< DDom >
 
struct  cartesian_prod< DDom1, DDom2, Tail... >
 
class  Chunk
 
class  Chunk< ElementType, DiscreteDomain< DDims... >, Allocator >
 
struct  chunk_traits
 
class  ChunkCommon
 
class  ChunkCommon< ElementType, DiscreteDomain< DDims... >, LayoutStridedPolicy >
 
class  ChunkSpan
 
class  ChunkSpan< ElementType, DiscreteDomain< DDims... >, LayoutStridedPolicy, MemorySpace >
 
struct  ConstantExtrapolationRule
 
struct  ConstantExtrapolationRule< DimI >
 A functor for describing a spline boundary value by a constant extrapolation for 1D evaluator. More...
 
struct  ConstantExtrapolationRule< DimI, DimNI >
 A functor for describing a spline boundary value by a constant extrapolation for 2D evaluator. More...
 
struct  coordinate_of
 
struct  Deriv
 
class  DiscreteDomain
 
class  DiscreteDomain<>
 
struct  DiscreteDomainIterator
 
class  DiscreteElement
 A DiscreteElement identifies an element of the discrete dimension. More...
 
class  DiscreteVector
 A DiscreteVector is a vector in the discrete dimension. More...
 
struct  Fourier
 
class  GrevilleInterpolationPoints
 A class which provides helper functions to initialise the Greville points from a B-Spline definition. More...
 
struct  is_discrete_domain
 
struct  is_discrete_domain< DiscreteDomain< Tags... > >
 
struct  is_discrete_element
 
struct  is_discrete_element< DiscreteElement< Tags... > >
 
struct  is_discrete_vector
 
struct  is_discrete_vector< DiscreteVector< Tags... > >
 
struct  is_non_uniform_sampling
 
struct  is_non_uniform_sampling< NonUniformPointSampling< CDim > >
 
struct  is_periodic_domain
 
struct  is_periodic_domain< DiscreteDomain< DDims... > >
 
struct  is_periodic_sampling
 
struct  is_periodic_sampling< PeriodicSampling< CDim > >
 
struct  is_rectilinear_domain
 
struct  is_rectilinear_domain< DiscreteDomain< DDims... > >
 
struct  is_uniform_domain
 
struct  is_uniform_domain< DiscreteDomain< DDims... > >
 
struct  is_uniform_sampling
 
struct  is_uniform_sampling< UniformPointSampling< CDim > >
 
class  KnotsAsInterpolationPoints
 Helper class for the initialisation of the mesh of interpolation points. More...
 
class  KokkosAllocator
 
struct  kwArgs_fft
 
class  NonUniformBSplines
 NonUniformPointSampling specialization of BSplines. More...
 
class  NonUniformPointSampling
 NonUniformPointSampling models a non-uniform discretization of the CDim segment $[a, b]$. More...
 
struct  NullExtrapolationRule
 
class  PdiEvent
 
struct  PeriodicExtrapolationRule
 
class  PeriodicSampling
 PeriodicSampling models a periodic discretization of the provided continuous dimension. More...
 
class  ScopeGuard
 
class  SplineBuilder
 A class for creating a spline approximation of a function. More...
 
class  SplineBuilder2D
 A class for creating a 2D spline approximation of a function. More...
 
class  SplineEvaluator
 
class  SplineEvaluator2D
 Define an evaluator 2D on B-splines. More...
 
class  UniformBSplines
 
class  UniformPointSampling
 UniformPointSampling models a uniform discretization of the provided continuous dimension. More...
 

Typedefs

template<class ElementType , class SupportType , class LayoutStridedPolicy = std::experimental::layout_right, class MemorySpace = Kokkos::HostSpace>
using ChunkView = ChunkSpan< ElementType const, SupportType, LayoutStridedPolicy, MemorySpace >
 
template<class T >
using chunk_value_t = typename chunk_traits< T >::value_type
 
template<class T >
using chunk_pointer_t = typename chunk_traits< T >::pointer_type
 
template<class T >
using chunk_reference_t = typename chunk_traits< T >::reference_type
 
using CoordinateElement = Real
 A CoordinateElement the type of the scalar used to represent elements of coordinates in the continuous space. More...
 
template<class... CDims>
using Coordinate = detail::TaggedVector< CoordinateElement, CDims... >
 A Coordinate represents a coordinate in the continuous space. More...
 
template<class T >
using coordinate_of_t = typename coordinate_of< T >::type
 Helper type of ddc::coordinate_of. More...
 
template<typename... DDom>
using cartesian_prod_t = typename cartesian_prod< DDom... >::type
 
using DiscreteElementType = std::size_t
 A DiscreteCoordElement is a scalar that identifies an element of the discrete dimension. More...
 
using DiscreteVectorElement = std::ptrdiff_t
 A DiscreteVectorElement is a scalar that represents the difference between two coordinates. More...
 
template<std::size_t N, class ElementType >
using SpanND = std::experimental::mdspan< ElementType, std::experimental::dextents< std::size_t, N > >
 
template<std::size_t N, class ElementType >
using ViewND = SpanND< N, ElementType const >
 
template<class ElementType >
using Span1D = SpanND< 1, ElementType >
 
template<class ElementType >
using Span2D = SpanND< 2, ElementType >
 
template<class ElementType >
using View1D = ViewND< 1, ElementType >
 
template<class ElementType >
using View2D = ViewND< 2, ElementType >
 
using DSpan1D = ddc::Span1D< double >
 
using DSpan2D = ddc::Span2D< double >
 
using CDSpan1D = ddc::Span1D< double const >
 
using CDSpan2D = ddc::Span2D< double const >
 
using DView1D = View1D< double >
 
using DView2D = View2D< double >
 
template<class T >
using DeviceAllocator = KokkosAllocator< T, Kokkos::DefaultExecutionSpace::memory_space >
 
template<class T >
using HostAllocator = KokkosAllocator< T, Kokkos::HostSpace >
 
using Real = float
 

Enumerations

enum class  FFT_Direction {
  FORWARD ,
  BACKWARD
}
 
enum class  FFT_Normalization {
  OFF ,
  FORWARD ,
  BACKWARD ,
  ORTHO ,
  FULL
}
 
enum class  BoundCond {
  PERIODIC ,
  HERMITE ,
  GREVILLE ,
  NATURAL
}
 
enum class  SplineSolver { GINKGO }
 

Functions

template<class T , std::size_t NT, class U , std::size_t NU>
constexpr bool operator== (AlignedAllocator< T, NT > const &, AlignedAllocator< U, NU > const &) noexcept
 
template<class T , std::size_t NT, class U , std::size_t NU>
constexpr bool operator!= (AlignedAllocator< T, NT > const &, AlignedAllocator< U, NU > const &) noexcept
 
template<class... DDims, class Allocator >
 Chunk (std::string const &, DiscreteDomain< DDims... > const &, Allocator) -> Chunk< typename Allocator::value_type, DiscreteDomain< DDims... >, Allocator >
 
template<class... DDims, class Allocator >
 Chunk (DiscreteDomain< DDims... > const &, Allocator) -> Chunk< typename Allocator::value_type, DiscreteDomain< DDims... >, Allocator >
 
template<class... QueryDDims, class ChunkType >
KOKKOS_FUNCTION auto get_domain (ChunkType const &chunk) noexcept
 Access the domain (or subdomain) of a view. More...
 
template<class KokkosView , class... DDims, class = std::enable_if_t<Kokkos::is_view<KokkosView>::value>>
 ChunkSpan (KokkosView const &view, DiscreteDomain< DDims... > domain) -> ChunkSpan< detail::kokkos_to_mdspan_element_t< typename KokkosView::data_type >, DiscreteDomain< DDims... >, detail::kokkos_to_mdspan_layout_t< typename KokkosView::array_layout >, typename KokkosView::memory_space >
 
template<class... DDim, std::enable_if_t<(sizeof...(DDim) > 1), int > = 0>
KOKKOS_FUNCTION Coordinate< typename DDim::continuous_dimension_type... > coordinate (DiscreteElement< DDim... > const &c)
 
template<class... QueryDDims, class... DDims>
KOKKOS_FUNCTION constexpr DiscreteDomain< QueryDDims... > select (DiscreteDomain< DDims... > const &domain)
 
template<class... DDimsA, class... DDimsB>
KOKKOS_FUNCTION constexpr auto remove_dims_of (DiscreteDomain< DDimsA... > const &DDom_a, DiscreteDomain< DDimsB... > const &DDom_b) noexcept
 
template<typename DDim1 , typename DDim2 , typename DDimA , typename... DDimsB>
KOKKOS_FUNCTION constexpr std::conditional_t< std::is_same_v< DDimA, DDim1 >, ddc::DiscreteDomain< DDim2 >, ddc::DiscreteDomain< DDimA > > replace_dim_of_1d (DiscreteDomain< DDimA > const &DDom_a, DiscreteDomain< DDimsB... > const &DDom_b) noexcept
 
template<typename DDim1 , typename DDim2 , typename... DDimsA, typename... DDimsB>
KOKKOS_FUNCTION constexpr auto replace_dim_of (DiscreteDomain< DDimsA... > const &DDom_a, DiscreteDomain< DDimsB... > const &DDom_b) noexcept
 
template<class... QueryDDims, class... DDims>
KOKKOS_FUNCTION constexpr DiscreteVector< QueryDDims... > extents (DiscreteDomain< DDims... > const &domain) noexcept
 
template<class... QueryDDims, class... DDims>
KOKKOS_FUNCTION constexpr DiscreteElement< QueryDDims... > front (DiscreteDomain< DDims... > const &domain) noexcept
 
template<class... QueryDDims, class... DDims>
KOKKOS_FUNCTION constexpr DiscreteElement< QueryDDims... > back (DiscreteDomain< DDims... > const &domain) noexcept
 
template<class QueryDDimSeq , class... DDims>
KOKKOS_FUNCTION constexpr auto select_by_type_seq (DiscreteDomain< DDims... > const &domain)
 
template<class Tag >
KOKKOS_FUNCTION constexpr DiscreteElementType const & uid (DiscreteElement< Tag > const &tuple) noexcept
 
template<class Tag >
KOKKOS_FUNCTION constexpr DiscreteElementTypeuid (DiscreteElement< Tag > &tuple) noexcept
 
template<class QueryTag , class... Tags>
KOKKOS_FUNCTION constexpr DiscreteElementType const & uid (DiscreteElement< Tags... > const &tuple) noexcept
 
template<class QueryTag , class... Tags>
KOKKOS_FUNCTION constexpr DiscreteElementTypeuid (DiscreteElement< Tags... > &tuple) noexcept
 
template<class QueryTag , class... Tags>
KOKKOS_FUNCTION constexpr DiscreteElementType const & uid_or (DiscreteElement< Tags... > const &tuple, DiscreteElementType const &default_value) noexcept
 
template<class... QueryTags, class... Tags>
KOKKOS_FUNCTION constexpr DiscreteElement< QueryTags... > select (DiscreteElement< Tags... > const &arr) noexcept
 
template<class... QueryTags, class... Tags>
KOKKOS_FUNCTION constexpr DiscreteElement< QueryTags... > select (DiscreteElement< Tags... > &&arr) noexcept
 
template<class QueryTag , class HeadDElem , class... TailDElems, std::enable_if_t< is_discrete_element_v< HeadDElem > &&(is_discrete_element_v< TailDElems > &&...), int > = 1>
KOKKOS_FUNCTION constexpr auto const & take (HeadDElem const &head, TailDElems const &... tail)
 Returns a reference towards the DiscreteElement that contains the QueryTag. More...
 
std::ostream & operator<< (std::ostream &out, DiscreteElement<> const &)
 
template<class Head , class... Tags>
std::ostream & operator<< (std::ostream &out, DiscreteElement< Head, Tags... > const &arr)
 
template<class... Tags, class... OTags>
KOKKOS_FUNCTION constexpr bool operator== (DiscreteElement< Tags... > const &lhs, DiscreteElement< OTags... > const &rhs) noexcept
 
template<class... Tags, class... OTags>
KOKKOS_FUNCTION constexpr bool operator!= (DiscreteElement< Tags... > const &lhs, DiscreteElement< OTags... > const &rhs) noexcept
 
template<class Tag >
KOKKOS_FUNCTION constexpr bool operator< (DiscreteElement< Tag > const &lhs, DiscreteElement< Tag > const &rhs)
 
template<class Tag >
KOKKOS_FUNCTION constexpr bool operator<= (DiscreteElement< Tag > const &lhs, DiscreteElement< Tag > const &rhs)
 
template<class Tag >
KOKKOS_FUNCTION constexpr bool operator> (DiscreteElement< Tag > const &lhs, DiscreteElement< Tag > const &rhs)
 
template<class Tag >
KOKKOS_FUNCTION constexpr bool operator>= (DiscreteElement< Tag > const &lhs, DiscreteElement< Tag > const &rhs)
 
template<class... Tags, class... OTags>
KOKKOS_FUNCTION constexpr DiscreteElement< Tags... > operator+ (DiscreteElement< Tags... > const &lhs, DiscreteVector< OTags... > const &rhs)
 right external binary operators: +, - More...
 
template<class Tag , class IntegralType , class = std::enable_if_t<std::is_integral_v<IntegralType>>, class = std::enable_if_t<!is_discrete_vector_v<IntegralType>>>
KOKKOS_FUNCTION constexpr DiscreteElement< Tag > operator+ (DiscreteElement< Tag > const &lhs, IntegralType const &rhs)
 
template<class... Tags, class... OTags>
KOKKOS_FUNCTION constexpr DiscreteElement< Tags... > operator- (DiscreteElement< Tags... > const &lhs, DiscreteVector< OTags... > const &rhs)
 
template<class Tag , class IntegralType , class = std::enable_if_t<std::is_integral_v<IntegralType>>, class = std::enable_if_t<!is_discrete_vector_v<IntegralType>>>
KOKKOS_FUNCTION constexpr DiscreteElement< Tag > operator- (DiscreteElement< Tag > const &lhs, IntegralType const &rhs)
 
template<class... Tags, class... OTags>
KOKKOS_FUNCTION constexpr DiscreteVector< Tags... > operator- (DiscreteElement< Tags... > const &lhs, DiscreteElement< OTags... > const &rhs)
 binary operator: - More...
 
template<class DDim , class... Args>
void init_discrete_space (Args &&... args)
 Initialize (emplace) a global singleton discrete space. More...
 
template<class DDimImpl , class Arg >
Arg init_discrete_space (std::tuple< DDimImpl, Arg > &&a)
 Move construct a global singleton discrete space and pass through the other argument. More...
 
template<class DDimImpl , class... Args>
std::enable_if_t< 2<=sizeof...(Args), std::tuple< Args... > > init_discrete_space(std::tuple< DDimImpl, Args... > &&a){ using DDim=typename DDimImpl::discrete_dimension_type;init_discrete_space< DDim >(std::move(std::get< 0 >(a)));return detail::extract_after(std::move(a), std::index_sequence_for< Args... >());}template< class DDim, class MemorySpace=DDC_CURRENT_KOKKOS_SPACE > KOKKOS_FUNCTION detail::ddim_impl_t< DDim, MemorySpace > const & discrete_space ()
 Move construct a global singleton discrete space and pass through remaining arguments. More...
 
template<class DDim >
bool is_discrete_space_initialized () noexcept
 
template<class DDim >
detail::ddim_impl_t< DDim, Kokkos::HostSpace > const & host_discrete_space ()
 
template<class QueryTag , class... Tags>
KOKKOS_FUNCTION constexpr DiscreteVectorElement const & get (DiscreteVector< Tags... > const &tuple) noexcept
 
template<class QueryTag , class... Tags>
KOKKOS_FUNCTION constexpr DiscreteVectorElementget (DiscreteVector< Tags... > &tuple) noexcept
 
template<class QueryTag , class... Tags>
KOKKOS_FUNCTION constexpr DiscreteVectorElement const & get_or (DiscreteVector< Tags... > const &tuple, DiscreteVectorElement const &default_value) noexcept
 
template<class... Tags>
KOKKOS_FUNCTION constexpr DiscreteVector< Tags... > operator+ (DiscreteVector< Tags... > const &x)
 Unary operators: +, -. More...
 
template<class... Tags>
KOKKOS_FUNCTION constexpr DiscreteVector< Tags... > operator- (DiscreteVector< Tags... > const &x)
 
template<class... Tags, class... OTags>
KOKKOS_FUNCTION constexpr auto operator+ (DiscreteVector< Tags... > const &lhs, DiscreteVector< OTags... > const &rhs)
 Internal binary operators: +, -. More...
 
template<class... Tags, class... OTags>
KOKKOS_FUNCTION constexpr auto operator- (DiscreteVector< Tags... > const &lhs, DiscreteVector< OTags... > const &rhs)
 
template<class Tag , class IntegralType , class = std::enable_if_t<std::is_integral_v<IntegralType>>>
KOKKOS_FUNCTION constexpr DiscreteVector< Tag > operator+ (DiscreteVector< Tag > const &lhs, IntegralType const &rhs)
 
template<class IntegralType , class Tag , class = std::enable_if_t<std::is_integral_v<IntegralType>>>
KOKKOS_FUNCTION constexpr DiscreteVector< Tag > operator+ (IntegralType const &lhs, DiscreteVector< Tag > const &rhs)
 
template<class Tag , class IntegralType , class = std::enable_if_t<std::is_integral_v<IntegralType>>>
KOKKOS_FUNCTION constexpr DiscreteVector< Tag > operator- (DiscreteVector< Tag > const &lhs, IntegralType const &rhs)
 
template<class IntegralType , class Tag , class = std::enable_if_t<std::is_integral_v<IntegralType>>>
KOKKOS_FUNCTION constexpr DiscreteVector< Tag > operator- (IntegralType const &lhs, DiscreteVector< Tag > const &rhs)
 
template<class IntegralType , class... Tags, class = std::enable_if_t<std::is_integral_v<IntegralType>>>
KOKKOS_FUNCTION constexpr auto operator* (IntegralType const &lhs, DiscreteVector< Tags... > const &rhs)
 external left binary operator: * More...
 
template<class... QueryTags, class... Tags>
KOKKOS_FUNCTION constexpr DiscreteVector< QueryTags... > select (DiscreteVector< Tags... > const &arr) noexcept
 
template<class... QueryTags, class... Tags>
KOKKOS_FUNCTION constexpr DiscreteVector< QueryTags... > select (DiscreteVector< Tags... > &&arr) noexcept
 
template<class QueryTag , class HeadDVect , class... TailDVects, std::enable_if_t< is_discrete_vector_v< HeadDVect > &&(is_discrete_vector_v< TailDVects > &&...), int > = 1>
KOKKOS_FUNCTION constexpr auto const & take (HeadDVect const &head, TailDVects const &... tail)
 Returns a reference towards the DiscreteVector that contains the QueryTag. More...
 
template<class Tag >
KOKKOS_FUNCTION constexpr bool operator< (DiscreteVector< Tag > const &lhs, DiscreteVector< Tag > const &rhs)
 
template<class Tag , class IntegralType >
KOKKOS_FUNCTION constexpr bool operator< (DiscreteVector< Tag > const &lhs, IntegralType const &rhs)
 
std::ostream & operator<< (std::ostream &out, DiscreteVector<> const &)
 
template<class Head , class... Tags>
std::ostream & operator<< (std::ostream &out, DiscreteVector< Head, Tags... > const &arr)
 
template<class... DDims, class Functor >
void for_each (DiscreteDomain< DDims... > const &domain, Functor &&f) noexcept
 iterates over a nD domain in serial More...
 
template<typename... X>
void init_fourier_space (ddc::DiscreteDomain< ddc::UniformPointSampling< X >... > x_mesh)
 
template<typename... X>
ddc::DiscreteDomain< ddc::PeriodicSampling< ddc::Fourier< X > >... > FourierMesh (ddc::DiscreteDomain< ddc::UniformPointSampling< X >... > x_mesh, bool C2C)
 
template<typename Tin , typename Tout , typename... X, typename ExecSpace , typename MemorySpace , typename layout_in , typename layout_out >
void fft (ExecSpace const &execSpace, ddc::ChunkSpan< Tout, ddc::DiscreteDomain< ddc::PeriodicSampling< ddc::Fourier< X > >... >, layout_out, MemorySpace > out, ddc::ChunkSpan< Tin, ddc::DiscreteDomain< ddc::UniformPointSampling< X >... >, layout_in, MemorySpace > in, ddc::kwArgs_fft kwargs={ddc::FFT_Normalization::OFF})
 
template<typename Tin , typename Tout , typename... X, typename ExecSpace , typename MemorySpace , typename layout_in , typename layout_out >
void ifft (ExecSpace const &execSpace, ddc::ChunkSpan< Tout, ddc::DiscreteDomain< ddc::UniformPointSampling< X >... >, layout_out, MemorySpace > out, ddc::ChunkSpan< Tin, ddc::DiscreteDomain< ddc::PeriodicSampling< ddc::Fourier< X > >... >, layout_in, MemorySpace > in, ddc::kwArgs_fft kwargs={ddc::FFT_Normalization::OFF})
 
constexpr int n_boundary_equations (ddc::BoundCond const bc, std::size_t const degree)
 
constexpr int n_user_input (ddc::BoundCond const bc, std::size_t const degree)
 
constexpr bool is_spline_interpolation_mesh_uniform (bool const is_uniform, ddc::BoundCond const BcXmin, ddc::BoundCond const BcXmax, int degree)
 
template<class ElementType , std::size_t N>
Span1D< ElementType > as_span (std::array< ElementType, N > &arr) noexcept
 
template<class ElementType , std::size_t N>
Span1D< const ElementType > as_span (std::array< ElementType, N > const &arr) noexcept
 
template<class ElementType , class Extents , class Layout , class Accessor >
std::ostream & operator<< (std::ostream &os, std::experimental::mdspan< ElementType, Extents, Layout, Accessor > const &s)
 Convenient function to dump a mdspan, it recursively prints all dimensions. More...
 
template<class T , class MST , class U , class MSU >
constexpr bool operator== (KokkosAllocator< T, MST > const &, KokkosAllocator< U, MSU > const &) noexcept
 
template<class T , class MST , class U , class MSU >
constexpr bool operator!= (KokkosAllocator< T, MST > const &, KokkosAllocator< U, MSU > const &) noexcept
 
template<class Space , class ElementType , class Support , class Layout , class MemorySpace >
auto create_mirror (Space const &space, ChunkSpan< ElementType, Support, Layout, MemorySpace > const &src)
 Returns a new Chunk with the same layout as src allocated on the memory space Space::memory_space. More...
 
template<class ElementType , class Support , class Layout , class MemorySpace >
auto create_mirror (ChunkSpan< ElementType, Support, Layout, MemorySpace > const &src)
 Returns a new host Chunk with the same layout as src. More...
 
template<class Space , class ElementType , class Support , class Layout , class MemorySpace >
auto create_mirror_and_copy (Space const &space, ChunkSpan< ElementType, Support, Layout, MemorySpace > const &src)
 Returns a new Chunk with the same layout as src allocated on the memory space Space::memory_space and operates a deep copy between the two. More...
 
template<class ElementType , class Support , class Layout , class MemorySpace >
auto create_mirror_and_copy (ChunkSpan< ElementType, Support, Layout, MemorySpace > const &src)
 Returns a new host Chunk with the same layout as src and operates a deep copy between the two. More...
 
template<class Space , class ElementType , class Support , class Layout , class MemorySpace >
auto create_mirror_view (Space const &space, ChunkSpan< ElementType, Support, Layout, MemorySpace > const &src)
 If src is accessible from space then returns a copy of src, otherwise returns a new Chunk with the same layout as src allocated on the memory space Space::memory_space. More...
 
template<class ElementType , class Support , class Layout , class MemorySpace >
auto create_mirror_view (ChunkSpan< ElementType, Support, Layout, MemorySpace > const &src)
 If src is host accessible then returns a copy of src, otherwise returns a new host Chunk with the same layout. More...
 
template<class Space , class ElementType , class Support , class Layout , class MemorySpace >
auto create_mirror_view_and_copy (Space const &space, ChunkSpan< ElementType, Support, Layout, MemorySpace > const &src)
 If src is accessible from space then returns a copy of src, otherwise returns a new Chunk with the same layout as src allocated on the memory space Space::memory_space and operates a deep copy between the two. More...
 
template<class ElementType , class Support , class Layout , class MemorySpace >
auto create_mirror_view_and_copy (ChunkSpan< ElementType, Support, Layout, MemorySpace > const &src)
 If src is host accessible then returns a copy of src, otherwise returns a new host Chunk with the same layout as src and operates a deep copy between the two. More...
 
template<class DDimImpl , std::enable_if_t< is_non_uniform_sampling_v< typename DDimImpl::discrete_dimension_type >, int > = 0>
std::ostream & operator<< (std::ostream &out, DDimImpl const &mesh)
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > coordinate (DiscreteElement< NonUniformPointSampling< CDim > > const &c)
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > distance_at_left (DiscreteElement< NonUniformPointSampling< CDim > > i)
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > distance_at_right (DiscreteElement< NonUniformPointSampling< CDim > > i)
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > rmin (DiscreteDomain< NonUniformPointSampling< CDim > > const &d)
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > rmax (DiscreteDomain< NonUniformPointSampling< CDim > > const &d)
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > rlength (DiscreteDomain< NonUniformPointSampling< CDim > > const &d)
 
template<class ChunkDst , class ChunkSrc >
auto parallel_deepcopy (ChunkDst &&dst, ChunkSrc &&src)
 Copy the content of a borrowed chunk into another. More...
 
template<class ExecSpace , class ChunkDst , class ChunkSrc >
auto parallel_deepcopy (ExecSpace const &execution_space, ChunkDst &&dst, ChunkSrc &&src)
 Copy the content of a borrowed chunk into another. More...
 
template<class ChunkDst , class T >
auto parallel_fill (ChunkDst &&dst, T const &value)
 Fill a borrowed chunk with a given value. More...
 
template<class ExecSpace , class ChunkDst , class T >
auto parallel_fill (ExecSpace const &execution_space, ChunkDst &&dst, T const &value)
 Fill a borrowed chunk with a given value. More...
 
template<class ExecSpace , class... DDims, class Functor >
void parallel_for_each (ExecSpace const &execution_space, DiscreteDomain< DDims... > const &domain, Functor &&f) noexcept
 iterates over a nD domain using a given Kokkos execution space More...
 
template<class... DDims, class Functor >
void parallel_for_each (DiscreteDomain< DDims... > const &domain, Functor &&f) noexcept
 iterates over a nD domain using the Kokkos default execution space More...
 
template<class ExecSpace , class... DDims, class T , class BinaryReductionOp , class UnaryTransformOp >
parallel_transform_reduce (ExecSpace const &execution_space, DiscreteDomain< DDims... > const &domain, T neutral, BinaryReductionOp &&reduce, UnaryTransformOp &&transform) noexcept
 A reduction over a nD domain using a given Kokkos execution space. More...
 
template<class... DDims, class T , class BinaryReductionOp , class UnaryTransformOp >
parallel_transform_reduce (DiscreteDomain< DDims... > const &domain, T neutral, BinaryReductionOp &&reduce, UnaryTransformOp &&transform) noexcept
 A reduction over a nD domain using the Kokkos default execution space. More...
 
template<PDI_inout_t access, class DataType >
void expose_to_pdi (std::string const &name, DataType &&data)
 
template<class DataType >
void expose_to_pdi (std::string const &name, DataType &&data)
 
template<class DDim >
KOKKOS_FUNCTION std::enable_if_t< is_periodic_sampling_v< DDim >, typename DDim::continuous_element_type > origin () noexcept
 Lower bound index of the mesh. More...
 
template<class DDim >
KOKKOS_FUNCTION std::enable_if_t< is_periodic_sampling_v< DDim >, typename DDim::discrete_element_type > front () noexcept
 Lower bound index of the mesh. More...
 
template<class DDim >
KOKKOS_FUNCTION std::enable_if_t< is_periodic_sampling_v< DDim >, Realstep () noexcept
 Spacing step of the mesh. More...
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > coordinate (DiscreteElement< PeriodicSampling< CDim > > const &c)
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > distance_at_left (DiscreteElement< PeriodicSampling< CDim > >)
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > distance_at_right (DiscreteElement< PeriodicSampling< CDim > >)
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > rmin (DiscreteDomain< PeriodicSampling< CDim > > const &d)
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > rmax (DiscreteDomain< PeriodicSampling< CDim > > const &d)
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > rlength (DiscreteDomain< PeriodicSampling< CDim > > const &d)
 
template<class... DDims, class T , class BinaryReductionOp , class UnaryTransformOp >
transform_reduce (DiscreteDomain< DDims... > const &domain, T neutral, BinaryReductionOp &&reduce, UnaryTransformOp &&transform) noexcept
 A reduction over a nD domain in serial. More...
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > coordinate (DiscreteElement< UniformPointSampling< CDim > > const &c)
 
template<class DDim >
KOKKOS_FUNCTION std::enable_if_t< is_uniform_sampling_v< DDim >, typename DDim::continuous_element_type > origin () noexcept
 Lower bound index of the mesh. More...
 
template<class DDim >
KOKKOS_FUNCTION std::enable_if_t< is_uniform_sampling_v< DDim >, typename DDim::discrete_element_type > front () noexcept
 Lower bound index of the mesh. More...
 
template<class DDim >
KOKKOS_FUNCTION std::enable_if_t< is_uniform_sampling_v< DDim >, Realstep () noexcept
 Spacing step of the mesh. More...
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > distance_at_left (DiscreteElement< UniformPointSampling< CDim > >)
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > distance_at_right (DiscreteElement< UniformPointSampling< CDim > >)
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > rmin (DiscreteDomain< UniformPointSampling< CDim > > const &d)
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > rmax (DiscreteDomain< UniformPointSampling< CDim > > const &d)
 
template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > rlength (DiscreteDomain< UniformPointSampling< CDim > > const &d)
 

Variables

template<class ElementType , class SupportType , class Allocator >
constexpr bool enable_chunk< Chunk< ElementType, SupportType, Allocator > > = true
 
template<class ElementType , class SupportType , class LayoutStridedPolicy , class MemorySpace >
constexpr bool enable_chunk< ChunkSpan< ElementType, SupportType, LayoutStridedPolicy, MemorySpace > > = true
 
template<class ElementType , class SupportType , class LayoutStridedPolicy , class MemorySpace >
constexpr bool enable_borrowed_chunk< ChunkSpan< ElementType, SupportType, LayoutStridedPolicy, MemorySpace > > = true
 
template<class T >
constexpr bool enable_borrowed_chunk = false
 
template<class T >
constexpr bool enable_chunk = false
 
template<class T >
constexpr bool is_chunk_v = enable_chunk<std::remove_const_t<std::remove_reference_t<T>>>
 
template<class T >
constexpr bool is_borrowed_chunk_v
 
template<class T >
constexpr bool is_writable_chunk_v = !std::is_const_v<std::remove_pointer_t<chunk_pointer_t<T>>>
 
template<class T >
constexpr bool is_discrete_domain_v = is_discrete_domain<T>::value
 
template<class T >
constexpr bool is_discrete_element_v = is_discrete_element<T>::value
 
template<class T >
constexpr bool is_discrete_vector_v = is_discrete_vector<T>::value
 
template<class DDim >
constexpr bool is_non_uniform_sampling_v = is_non_uniform_sampling<DDim>::value
 
template<class DDim >
constexpr bool is_periodic_sampling_v = is_periodic_sampling<DDim>::value
 
template<class T >
constexpr bool is_periodic_domain_v = is_periodic_domain<T>::value
 
template<class T >
constexpr bool is_rectilinear_domain_v = is_rectilinear_domain<T>::value
 
template<class DDim >
constexpr bool is_uniform_sampling_v = is_uniform_sampling<DDim>::value
 
template<class T >
constexpr bool is_uniform_domain_v = is_uniform_domain<T>::value
 

Class Documentation

◆ ddc::cartesian_prod

struct ddc::cartesian_prod
template<typename... DDom>
struct ddc::cartesian_prod< DDom >

◆ ddc::cartesian_prod< ddc::DiscreteDomain< DDim1... >, ddc::DiscreteDomain< DDim2... > >

struct ddc::cartesian_prod< ddc::DiscreteDomain< DDim1... >, ddc::DiscreteDomain< DDim2... > >
template<typename... DDim1, typename... DDim2>
struct ddc::cartesian_prod< ddc::DiscreteDomain< DDim1... >, ddc::DiscreteDomain< DDim2... > >
Class Members
typedef DiscreteDomain< DDim1..., DDim2... > type

◆ ddc::cartesian_prod< DDom >

struct ddc::cartesian_prod< DDom >
template<typename DDom>
struct ddc::cartesian_prod< DDom >
Class Members
typedef DDom type

◆ ddc::cartesian_prod< DDom1, DDom2, Tail... >

struct ddc::cartesian_prod< DDom1, DDom2, Tail... >
template<typename DDom1, typename DDom2, typename... Tail>
struct ddc::cartesian_prod< DDom1, DDom2, Tail... >
Class Members
typedef typename type, Tail... >::type type

◆ ddc::Chunk

class ddc::Chunk
template<class ElementType, class, class Allocator = HostAllocator<ElementType>>
class ddc::Chunk< ElementType, class, Allocator >

◆ ddc::chunk_traits

struct ddc::chunk_traits
template<class T>
struct ddc::chunk_traits< T >
Class Members
typedef remove_cv_t< remove_pointer_t< decltype(declval< T >().data_handle())> > value_type
typedef decltype(declval< T >().data_handle()) pointer_type
typedef decltype(*declval< T >().data_handle()) reference_type

◆ ddc::ChunkCommon

class ddc::ChunkCommon
template<class ElementType, class SupportType, class LayoutStridedPolicy>
class ddc::ChunkCommon< ElementType, SupportType, LayoutStridedPolicy >

◆ ddc::ChunkSpan

class ddc::ChunkSpan
template<class ElementType, class SupportType, class LayoutStridedPolicy = std::experimental::layout_right, class MemorySpace = Kokkos::DefaultHostExecutionSpace::memory_space>
class ddc::ChunkSpan< ElementType, SupportType, LayoutStridedPolicy, MemorySpace >

◆ ddc::ConstantExtrapolationRule

struct ddc::ConstantExtrapolationRule
template<class DimI, class... Dim>
struct ddc::ConstantExtrapolationRule< DimI, Dim >

◆ ddc::coordinate_of

struct ddc::coordinate_of
template<class T>
struct ddc::coordinate_of< T >
Class Members
typedef decltype(coordinate(declval< T >())) type

◆ ddc::Deriv

struct ddc::Deriv
template<class Tag>
struct ddc::Deriv< Tag >

◆ ddc::Fourier

struct ddc::Fourier
template<typename Dim>
struct ddc::Fourier< Dim >

◆ ddc::kwArgs_fft

struct ddc::kwArgs_fft
Class Members
FFT_Normalization normalization

◆ ddc::NonUniformPointSampling

class ddc::NonUniformPointSampling
template<class CDim>
class ddc::NonUniformPointSampling< CDim >

NonUniformPointSampling models a non-uniform discretization of the CDim segment $[a, b]$.

Class Members
typedef CDim continuous_dimension_type
typedef Coordinate< CDim > continuous_element_type
typedef NonUniformPointSampling discrete_dimension_type
typedef DiscreteDomain< NonUniformPointSampling > discrete_domain_type
typedef DiscreteElement< NonUniformPointSampling > discrete_element_type
typedef DiscreteVector< NonUniformPointSampling > discrete_vector_type

Typedef Documentation

◆ ChunkView

template<class ElementType , class SupportType , class LayoutStridedPolicy = std::experimental::layout_right, class MemorySpace = Kokkos::HostSpace>
using ddc::ChunkView = typedef ChunkSpan<ElementType const, SupportType, LayoutStridedPolicy, MemorySpace>

◆ chunk_value_t

template<class T >
using ddc::chunk_value_t = typedef typename chunk_traits<T>::value_type

◆ chunk_pointer_t

template<class T >
using ddc::chunk_pointer_t = typedef typename chunk_traits<T>::pointer_type

◆ chunk_reference_t

template<class T >
using ddc::chunk_reference_t = typedef typename chunk_traits<T>::reference_type

◆ CoordinateElement

using ddc::CoordinateElement = typedef Real

A CoordinateElement the type of the scalar used to represent elements of coordinates in the continuous space.

◆ Coordinate

template<class... CDims>
using ddc::Coordinate = typedef detail::TaggedVector<CoordinateElement, CDims...>

A Coordinate represents a coordinate in the continuous space.

It is tagged by its dimensions.

◆ coordinate_of_t

template<class T >
using ddc::coordinate_of_t = typedef typename coordinate_of<T>::type

Helper type of ddc::coordinate_of.

◆ cartesian_prod_t

template<typename... DDom>
using ddc::cartesian_prod_t = typedef typename cartesian_prod<DDom...>::type

◆ DiscreteElementType

using ddc::DiscreteElementType = typedef std::size_t

A DiscreteCoordElement is a scalar that identifies an element of the discrete dimension.

◆ DiscreteVectorElement

using ddc::DiscreteVectorElement = typedef std::ptrdiff_t

A DiscreteVectorElement is a scalar that represents the difference between two coordinates.

◆ SpanND

template<std::size_t N, class ElementType >
using ddc::SpanND = typedef std::experimental::mdspan<ElementType, std::experimental::dextents<std::size_t, N> >

◆ ViewND

template<std::size_t N, class ElementType >
using ddc::ViewND = typedef SpanND<N, ElementType const>

◆ Span1D

template<class ElementType >
using ddc::Span1D = typedef SpanND<1, ElementType>

◆ Span2D

template<class ElementType >
using ddc::Span2D = typedef SpanND<2, ElementType>

◆ View1D

template<class ElementType >
using ddc::View1D = typedef ViewND<1, ElementType>

◆ View2D

template<class ElementType >
using ddc::View2D = typedef ViewND<2, ElementType>

◆ DSpan1D

using ddc::DSpan1D = typedef ddc::Span1D<double>

◆ DSpan2D

using ddc::DSpan2D = typedef ddc::Span2D<double>

◆ CDSpan1D

using ddc::CDSpan1D = typedef ddc::Span1D<double const>

◆ CDSpan2D

using ddc::CDSpan2D = typedef ddc::Span2D<double const>

◆ DView1D

using ddc::DView1D = typedef View1D<double>

◆ DView2D

using ddc::DView2D = typedef View2D<double>

◆ DeviceAllocator

template<class T >
using ddc::DeviceAllocator = typedef KokkosAllocator<T, Kokkos::DefaultExecutionSpace::memory_space>

◆ HostAllocator

template<class T >
using ddc::HostAllocator = typedef KokkosAllocator<T, Kokkos::HostSpace>

◆ Real

using ddc::Real = typedef float

Enumeration Type Documentation

◆ FFT_Direction

enum class ddc::FFT_Direction
strong
Enumerator
FORWARD 
BACKWARD 

◆ FFT_Normalization

enum class ddc::FFT_Normalization
strong
Enumerator
OFF 
FORWARD 
BACKWARD 
ORTHO 
FULL 

◆ BoundCond

enum class ddc::BoundCond
strong
Enumerator
PERIODIC 
HERMITE 
GREVILLE 
NATURAL 

◆ SplineSolver

enum class ddc::SplineSolver
strong
Enumerator
GINKGO 

Function Documentation

◆ operator==() [1/3]

template<class T , std::size_t NT, class U , std::size_t NU>
constexpr bool ddc::operator== ( AlignedAllocator< T, NT > const &  ,
AlignedAllocator< U, NU > const &   
)
constexprnoexcept

◆ operator!=() [1/3]

template<class T , std::size_t NT, class U , std::size_t NU>
constexpr bool ddc::operator!= ( AlignedAllocator< T, NT > const &  ,
AlignedAllocator< U, NU > const &   
)
constexprnoexcept

◆ Chunk() [1/2]

template<class... DDims, class Allocator >
ddc::Chunk ( std::string const &  ,
DiscreteDomain< DDims... > const &  ,
Allocator   
) -> Chunk< typename Allocator::value_type, DiscreteDomain< DDims... >, Allocator >

◆ Chunk() [2/2]

template<class... DDims, class Allocator >
ddc::Chunk ( DiscreteDomain< DDims... > const &  ,
Allocator   
) -> Chunk< typename Allocator::value_type, DiscreteDomain< DDims... >, Allocator >

◆ get_domain()

template<class... QueryDDims, class ChunkType >
KOKKOS_FUNCTION auto ddc::get_domain ( ChunkType const &  chunk)
noexcept

Access the domain (or subdomain) of a view.

Parameters
[in]chunkthe view whose domain to access
Returns
the domain of view in the queried dimensions

◆ ChunkSpan()

template<class KokkosView , class... DDims, class = std::enable_if_t<Kokkos::is_view<KokkosView>::value>>
ddc::ChunkSpan ( KokkosView const &  view,
DiscreteDomain< DDims... >  domain 
) -> ChunkSpan< detail::kokkos_to_mdspan_element_t< typename KokkosView::data_type >, DiscreteDomain< DDims... >, detail::kokkos_to_mdspan_layout_t< typename KokkosView::array_layout >, typename KokkosView::memory_space >

◆ coordinate() [1/4]

template<class... DDim, std::enable_if_t<(sizeof...(DDim) > 1), int > = 0>
KOKKOS_FUNCTION Coordinate< typename DDim::continuous_dimension_type... > ddc::coordinate ( DiscreteElement< DDim... > const &  c)

◆ select() [1/5]

template<class... QueryDDims, class... DDims>
KOKKOS_FUNCTION constexpr DiscreteDomain< QueryDDims... > ddc::select ( DiscreteDomain< DDims... > const &  domain)
constexpr

◆ remove_dims_of()

template<class... DDimsA, class... DDimsB>
KOKKOS_FUNCTION constexpr auto ddc::remove_dims_of ( DiscreteDomain< DDimsA... > const &  DDom_a,
DiscreteDomain< DDimsB... > const &  DDom_b 
)
constexprnoexcept

◆ replace_dim_of_1d()

template<typename DDim1 , typename DDim2 , typename DDimA , typename... DDimsB>
KOKKOS_FUNCTION constexpr std::conditional_t< std::is_same_v< DDimA, DDim1 >, ddc::DiscreteDomain< DDim2 >, ddc::DiscreteDomain< DDimA > > ddc::replace_dim_of_1d ( DiscreteDomain< DDimA > const &  DDom_a,
DiscreteDomain< DDimsB... > const &  DDom_b 
)
constexprnoexcept

◆ replace_dim_of()

template<typename DDim1 , typename DDim2 , typename... DDimsA, typename... DDimsB>
KOKKOS_FUNCTION constexpr auto ddc::replace_dim_of ( DiscreteDomain< DDimsA... > const &  DDom_a,
DiscreteDomain< DDimsB... > const &  DDom_b 
)
constexprnoexcept

◆ extents()

template<class... QueryDDims, class... DDims>
KOKKOS_FUNCTION constexpr DiscreteVector< QueryDDims... > ddc::extents ( DiscreteDomain< DDims... > const &  domain)
constexprnoexcept

◆ front() [1/3]

template<class... QueryDDims, class... DDims>
KOKKOS_FUNCTION constexpr DiscreteElement< QueryDDims... > ddc::front ( DiscreteDomain< DDims... > const &  domain)
constexprnoexcept

◆ back()

template<class... QueryDDims, class... DDims>
KOKKOS_FUNCTION constexpr DiscreteElement< QueryDDims... > ddc::back ( DiscreteDomain< DDims... > const &  domain)
constexprnoexcept

◆ select_by_type_seq()

template<class QueryDDimSeq , class... DDims>
KOKKOS_FUNCTION constexpr auto ddc::select_by_type_seq ( DiscreteDomain< DDims... > const &  domain)
constexpr

◆ uid() [1/4]

template<class Tag >
KOKKOS_FUNCTION constexpr DiscreteElementType const & ddc::uid ( DiscreteElement< Tag > const &  tuple)
constexprnoexcept

◆ uid() [2/4]

template<class Tag >
KOKKOS_FUNCTION constexpr DiscreteElementType & ddc::uid ( DiscreteElement< Tag > &  tuple)
constexprnoexcept

◆ uid() [3/4]

template<class QueryTag , class... Tags>
KOKKOS_FUNCTION constexpr DiscreteElementType const & ddc::uid ( DiscreteElement< Tags... > const &  tuple)
constexprnoexcept

◆ uid() [4/4]

template<class QueryTag , class... Tags>
KOKKOS_FUNCTION constexpr DiscreteElementType & ddc::uid ( DiscreteElement< Tags... > &  tuple)
constexprnoexcept

◆ uid_or()

template<class QueryTag , class... Tags>
KOKKOS_FUNCTION constexpr DiscreteElementType const & ddc::uid_or ( DiscreteElement< Tags... > const &  tuple,
DiscreteElementType const &  default_value 
)
constexprnoexcept

◆ select() [2/5]

template<class... QueryTags, class... Tags>
KOKKOS_FUNCTION constexpr DiscreteElement< QueryTags... > ddc::select ( DiscreteElement< Tags... > const &  arr)
constexprnoexcept

◆ select() [3/5]

template<class... QueryTags, class... Tags>
KOKKOS_FUNCTION constexpr DiscreteElement< QueryTags... > ddc::select ( DiscreteElement< Tags... > &&  arr)
constexprnoexcept

◆ take() [1/2]

template<class QueryTag , class HeadDElem , class... TailDElems, std::enable_if_t< is_discrete_element_v< HeadDElem > &&(is_discrete_element_v< TailDElems > &&...), int > = 1>
KOKKOS_FUNCTION constexpr auto const & ddc::take ( HeadDElem const &  head,
TailDElems const &...  tail 
)
constexpr

Returns a reference towards the DiscreteElement that contains the QueryTag.

◆ operator<<() [1/6]

std::ostream & ddc::operator<< ( std::ostream &  out,
DiscreteElement<> const &   
)
inline

◆ operator<<() [2/6]

template<class Head , class... Tags>
std::ostream & ddc::operator<< ( std::ostream &  out,
DiscreteElement< Head, Tags... > const &  arr 
)

◆ operator==() [2/3]

template<class... Tags, class... OTags>
KOKKOS_FUNCTION constexpr bool ddc::operator== ( DiscreteElement< Tags... > const &  lhs,
DiscreteElement< OTags... > const &  rhs 
)
constexprnoexcept

◆ operator!=() [2/3]

template<class... Tags, class... OTags>
KOKKOS_FUNCTION constexpr bool ddc::operator!= ( DiscreteElement< Tags... > const &  lhs,
DiscreteElement< OTags... > const &  rhs 
)
constexprnoexcept

◆ operator<() [1/3]

template<class Tag >
KOKKOS_FUNCTION constexpr bool ddc::operator< ( DiscreteElement< Tag > const &  lhs,
DiscreteElement< Tag > const &  rhs 
)
constexpr

◆ operator<=()

template<class Tag >
KOKKOS_FUNCTION constexpr bool ddc::operator<= ( DiscreteElement< Tag > const &  lhs,
DiscreteElement< Tag > const &  rhs 
)
constexpr

◆ operator>()

template<class Tag >
KOKKOS_FUNCTION constexpr bool ddc::operator> ( DiscreteElement< Tag > const &  lhs,
DiscreteElement< Tag > const &  rhs 
)
constexpr

◆ operator>=()

template<class Tag >
KOKKOS_FUNCTION constexpr bool ddc::operator>= ( DiscreteElement< Tag > const &  lhs,
DiscreteElement< Tag > const &  rhs 
)
constexpr

◆ operator+() [1/6]

template<class... Tags, class... OTags>
KOKKOS_FUNCTION constexpr DiscreteElement< Tags... > ddc::operator+ ( DiscreteElement< Tags... > const &  lhs,
DiscreteVector< OTags... > const &  rhs 
)
constexpr

right external binary operators: +, -

◆ operator+() [2/6]

template<class Tag , class IntegralType , class = std::enable_if_t<std::is_integral_v<IntegralType>>, class = std::enable_if_t<!is_discrete_vector_v<IntegralType>>>
KOKKOS_FUNCTION constexpr DiscreteElement< Tag > ddc::operator+ ( DiscreteElement< Tag > const &  lhs,
IntegralType const &  rhs 
)
constexpr

◆ operator-() [1/7]

template<class... Tags, class... OTags>
KOKKOS_FUNCTION constexpr DiscreteElement< Tags... > ddc::operator- ( DiscreteElement< Tags... > const &  lhs,
DiscreteVector< OTags... > const &  rhs 
)
constexpr

◆ operator-() [2/7]

template<class Tag , class IntegralType , class = std::enable_if_t<std::is_integral_v<IntegralType>>, class = std::enable_if_t<!is_discrete_vector_v<IntegralType>>>
KOKKOS_FUNCTION constexpr DiscreteElement< Tag > ddc::operator- ( DiscreteElement< Tag > const &  lhs,
IntegralType const &  rhs 
)
constexpr

◆ operator-() [3/7]

template<class... Tags, class... OTags>
KOKKOS_FUNCTION constexpr DiscreteVector< Tags... > ddc::operator- ( DiscreteElement< Tags... > const &  lhs,
DiscreteElement< OTags... > const &  rhs 
)
constexpr

binary operator: -

◆ init_discrete_space() [1/2]

template<class DDim , class... Args>
void ddc::init_discrete_space ( Args &&...  args)

Initialize (emplace) a global singleton discrete space.

Parameters
argsthe constructor arguments

◆ init_discrete_space() [2/2]

template<class DDimImpl , class Arg >
Arg ddc::init_discrete_space ( std::tuple< DDimImpl, Arg > &&  a)

Move construct a global singleton discrete space and pass through the other argument.

Parameters
a- the discrete space to move at index 0
  • the arguments to pass through at index 1

◆ discrete_space()

template<class DDimImpl , class... Args>
std::enable_if_t< 2<=sizeof...(Args), std::tuple< Args... > > init_discrete_space( std::tuple< DDimImpl, Args... > &&a){ using DDim=typename DDimImpl::discrete_dimension_type; init_discrete_space< DDim >(std::move(std::get< 0 >(a))); return detail::extract_after(std::move(a), std::index_sequence_for< Args... >());}template< class DDim, class MemorySpace=DDC_CURRENT_KOKKOS_SPACE > KOKKOS_FUNCTION detail::ddim_impl_t< DDim, MemorySpace > const & ddc::discrete_space ( )

Move construct a global singleton discrete space and pass through remaining arguments.

Parameters
a- the discrete space to move at index 0
  • the (2+) arguments to pass through in other indices

◆ is_discrete_space_initialized()

template<class DDim >
bool ddc::is_discrete_space_initialized ( )
noexcept

◆ host_discrete_space()

template<class DDim >
detail::ddim_impl_t< DDim, Kokkos::HostSpace > const & ddc::host_discrete_space ( )

◆ get() [1/2]

template<class QueryTag , class... Tags>
KOKKOS_FUNCTION constexpr DiscreteVectorElement const & ddc::get ( DiscreteVector< Tags... > const &  tuple)
constexprnoexcept

◆ get() [2/2]

template<class QueryTag , class... Tags>
KOKKOS_FUNCTION constexpr DiscreteVectorElement & ddc::get ( DiscreteVector< Tags... > &  tuple)
constexprnoexcept

◆ get_or()

template<class QueryTag , class... Tags>
KOKKOS_FUNCTION constexpr DiscreteVectorElement const & ddc::get_or ( DiscreteVector< Tags... > const &  tuple,
DiscreteVectorElement const &  default_value 
)
constexprnoexcept

◆ operator+() [3/6]

template<class... Tags>
KOKKOS_FUNCTION constexpr DiscreteVector< Tags... > ddc::operator+ ( DiscreteVector< Tags... > const &  x)
constexpr

Unary operators: +, -.

◆ operator-() [4/7]

template<class... Tags>
KOKKOS_FUNCTION constexpr DiscreteVector< Tags... > ddc::operator- ( DiscreteVector< Tags... > const &  x)
constexpr

◆ operator+() [4/6]

template<class... Tags, class... OTags>
KOKKOS_FUNCTION constexpr auto ddc::operator+ ( DiscreteVector< Tags... > const &  lhs,
DiscreteVector< OTags... > const &  rhs 
)
constexpr

Internal binary operators: +, -.

◆ operator-() [5/7]

template<class... Tags, class... OTags>
KOKKOS_FUNCTION constexpr auto ddc::operator- ( DiscreteVector< Tags... > const &  lhs,
DiscreteVector< OTags... > const &  rhs 
)
constexpr

◆ operator+() [5/6]

template<class Tag , class IntegralType , class = std::enable_if_t<std::is_integral_v<IntegralType>>>
KOKKOS_FUNCTION constexpr DiscreteVector< Tag > ddc::operator+ ( DiscreteVector< Tag > const &  lhs,
IntegralType const &  rhs 
)
constexpr

◆ operator+() [6/6]

template<class IntegralType , class Tag , class = std::enable_if_t<std::is_integral_v<IntegralType>>>
KOKKOS_FUNCTION constexpr DiscreteVector< Tag > ddc::operator+ ( IntegralType const &  lhs,
DiscreteVector< Tag > const &  rhs 
)
constexpr

◆ operator-() [6/7]

template<class Tag , class IntegralType , class = std::enable_if_t<std::is_integral_v<IntegralType>>>
KOKKOS_FUNCTION constexpr DiscreteVector< Tag > ddc::operator- ( DiscreteVector< Tag > const &  lhs,
IntegralType const &  rhs 
)
constexpr

◆ operator-() [7/7]

template<class IntegralType , class Tag , class = std::enable_if_t<std::is_integral_v<IntegralType>>>
KOKKOS_FUNCTION constexpr DiscreteVector< Tag > ddc::operator- ( IntegralType const &  lhs,
DiscreteVector< Tag > const &  rhs 
)
constexpr

◆ operator*()

template<class IntegralType , class... Tags, class = std::enable_if_t<std::is_integral_v<IntegralType>>>
KOKKOS_FUNCTION constexpr auto ddc::operator* ( IntegralType const &  lhs,
DiscreteVector< Tags... > const &  rhs 
)
constexpr

external left binary operator: *

◆ select() [4/5]

template<class... QueryTags, class... Tags>
KOKKOS_FUNCTION constexpr DiscreteVector< QueryTags... > ddc::select ( DiscreteVector< Tags... > const &  arr)
constexprnoexcept

◆ select() [5/5]

template<class... QueryTags, class... Tags>
KOKKOS_FUNCTION constexpr DiscreteVector< QueryTags... > ddc::select ( DiscreteVector< Tags... > &&  arr)
constexprnoexcept

◆ take() [2/2]

template<class QueryTag , class HeadDVect , class... TailDVects, std::enable_if_t< is_discrete_vector_v< HeadDVect > &&(is_discrete_vector_v< TailDVects > &&...), int > = 1>
KOKKOS_FUNCTION constexpr auto const & ddc::take ( HeadDVect const &  head,
TailDVects const &...  tail 
)
constexpr

Returns a reference towards the DiscreteVector that contains the QueryTag.

◆ operator<() [2/3]

template<class Tag >
KOKKOS_FUNCTION constexpr bool ddc::operator< ( DiscreteVector< Tag > const &  lhs,
DiscreteVector< Tag > const &  rhs 
)
constexpr

◆ operator<() [3/3]

template<class Tag , class IntegralType >
KOKKOS_FUNCTION constexpr bool ddc::operator< ( DiscreteVector< Tag > const &  lhs,
IntegralType const &  rhs 
)
constexpr

◆ operator<<() [3/6]

std::ostream & ddc::operator<< ( std::ostream &  out,
DiscreteVector<> const &   
)
inline

◆ operator<<() [4/6]

template<class Head , class... Tags>
std::ostream & ddc::operator<< ( std::ostream &  out,
DiscreteVector< Head, Tags... > const &  arr 
)

◆ for_each()

template<class... DDims, class Functor >
void ddc::for_each ( DiscreteDomain< DDims... > const &  domain,
Functor &&  f 
)
noexcept

iterates over a nD domain in serial

Parameters
[in]domainthe domain over which to iterate
[in]fa functor taking an index as parameter

◆ init_fourier_space()

template<typename... X>
void ddc::init_fourier_space ( ddc::DiscreteDomain< ddc::UniformPointSampling< X >... >  x_mesh)

◆ FourierMesh()

template<typename... X>
ddc::DiscreteDomain< ddc::PeriodicSampling< ddc::Fourier< X > >... > ddc::FourierMesh ( ddc::DiscreteDomain< ddc::UniformPointSampling< X >... >  x_mesh,
bool  C2C 
)

◆ fft()

template<typename Tin , typename Tout , typename... X, typename ExecSpace , typename MemorySpace , typename layout_in , typename layout_out >
void ddc::fft ( ExecSpace const &  execSpace,
ddc::ChunkSpan< Tout, ddc::DiscreteDomain< ddc::PeriodicSampling< ddc::Fourier< X > >... >, layout_out, MemorySpace >  out,
ddc::ChunkSpan< Tin, ddc::DiscreteDomain< ddc::UniformPointSampling< X >... >, layout_in, MemorySpace >  in,
ddc::kwArgs_fft  kwargs = {ddc::FFT_Normalization::OFF} 
)

◆ ifft()

template<typename Tin , typename Tout , typename... X, typename ExecSpace , typename MemorySpace , typename layout_in , typename layout_out >
void ddc::ifft ( ExecSpace const &  execSpace,
ddc::ChunkSpan< Tout, ddc::DiscreteDomain< ddc::UniformPointSampling< X >... >, layout_out, MemorySpace >  out,
ddc::ChunkSpan< Tin, ddc::DiscreteDomain< ddc::PeriodicSampling< ddc::Fourier< X > >... >, layout_in, MemorySpace >  in,
ddc::kwArgs_fft  kwargs = {ddc::FFT_Normalization::OFF} 
)

◆ n_boundary_equations()

constexpr int ddc::n_boundary_equations ( ddc::BoundCond const  bc,
std::size_t const  degree 
)
constexpr

◆ n_user_input()

constexpr int ddc::n_user_input ( ddc::BoundCond const  bc,
std::size_t const  degree 
)
constexpr

◆ is_spline_interpolation_mesh_uniform()

constexpr bool ddc::is_spline_interpolation_mesh_uniform ( bool const  is_uniform,
ddc::BoundCond const  BcXmin,
ddc::BoundCond const  BcXmax,
int  degree 
)
constexpr

◆ as_span() [1/2]

template<class ElementType , std::size_t N>
Span1D< ElementType > ddc::as_span ( std::array< ElementType, N > &  arr)
noexcept

◆ as_span() [2/2]

template<class ElementType , std::size_t N>
Span1D< const ElementType > ddc::as_span ( std::array< ElementType, N > const &  arr)
noexcept

◆ operator<<() [5/6]

template<class ElementType , class Extents , class Layout , class Accessor >
std::ostream & ddc::operator<< ( std::ostream &  os,
std::experimental::mdspan< ElementType, Extents, Layout, Accessor > const &  s 
)

Convenient function to dump a mdspan, it recursively prints all dimensions.

Disclaimer: use with caution for large arrays

◆ operator==() [3/3]

template<class T , class MST , class U , class MSU >
constexpr bool ddc::operator== ( KokkosAllocator< T, MST > const &  ,
KokkosAllocator< U, MSU > const &   
)
constexprnoexcept

◆ operator!=() [3/3]

template<class T , class MST , class U , class MSU >
constexpr bool ddc::operator!= ( KokkosAllocator< T, MST > const &  ,
KokkosAllocator< U, MSU > const &   
)
constexprnoexcept

◆ create_mirror() [1/2]

template<class Space , class ElementType , class Support , class Layout , class MemorySpace >
auto ddc::create_mirror ( Space const &  space,
ChunkSpan< ElementType, Support, Layout, MemorySpace > const &  src 
)

Returns a new Chunk with the same layout as src allocated on the memory space Space::memory_space.

◆ create_mirror() [2/2]

template<class ElementType , class Support , class Layout , class MemorySpace >
auto ddc::create_mirror ( ChunkSpan< ElementType, Support, Layout, MemorySpace > const &  src)

Returns a new host Chunk with the same layout as src.

◆ create_mirror_and_copy() [1/2]

template<class Space , class ElementType , class Support , class Layout , class MemorySpace >
auto ddc::create_mirror_and_copy ( Space const &  space,
ChunkSpan< ElementType, Support, Layout, MemorySpace > const &  src 
)

Returns a new Chunk with the same layout as src allocated on the memory space Space::memory_space and operates a deep copy between the two.

◆ create_mirror_and_copy() [2/2]

template<class ElementType , class Support , class Layout , class MemorySpace >
auto ddc::create_mirror_and_copy ( ChunkSpan< ElementType, Support, Layout, MemorySpace > const &  src)

Returns a new host Chunk with the same layout as src and operates a deep copy between the two.

◆ create_mirror_view() [1/2]

template<class Space , class ElementType , class Support , class Layout , class MemorySpace >
auto ddc::create_mirror_view ( Space const &  space,
ChunkSpan< ElementType, Support, Layout, MemorySpace > const &  src 
)

If src is accessible from space then returns a copy of src, otherwise returns a new Chunk with the same layout as src allocated on the memory space Space::memory_space.

◆ create_mirror_view() [2/2]

template<class ElementType , class Support , class Layout , class MemorySpace >
auto ddc::create_mirror_view ( ChunkSpan< ElementType, Support, Layout, MemorySpace > const &  src)

If src is host accessible then returns a copy of src, otherwise returns a new host Chunk with the same layout.

◆ create_mirror_view_and_copy() [1/2]

template<class Space , class ElementType , class Support , class Layout , class MemorySpace >
auto ddc::create_mirror_view_and_copy ( Space const &  space,
ChunkSpan< ElementType, Support, Layout, MemorySpace > const &  src 
)

If src is accessible from space then returns a copy of src, otherwise returns a new Chunk with the same layout as src allocated on the memory space Space::memory_space and operates a deep copy between the two.

◆ create_mirror_view_and_copy() [2/2]

template<class ElementType , class Support , class Layout , class MemorySpace >
auto ddc::create_mirror_view_and_copy ( ChunkSpan< ElementType, Support, Layout, MemorySpace > const &  src)

If src is host accessible then returns a copy of src, otherwise returns a new host Chunk with the same layout as src and operates a deep copy between the two.

◆ operator<<() [6/6]

template<class DDimImpl , std::enable_if_t< is_non_uniform_sampling_v< typename DDimImpl::discrete_dimension_type >, int > = 0>
std::ostream & ddc::operator<< ( std::ostream &  out,
DDimImpl const &  mesh 
)

◆ coordinate() [2/4]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::coordinate ( DiscreteElement< NonUniformPointSampling< CDim > > const &  c)

◆ distance_at_left() [1/3]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::distance_at_left ( DiscreteElement< NonUniformPointSampling< CDim > >  i)

◆ distance_at_right() [1/3]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::distance_at_right ( DiscreteElement< NonUniformPointSampling< CDim > >  i)

◆ rmin() [1/3]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::rmin ( DiscreteDomain< NonUniformPointSampling< CDim > > const &  d)

◆ rmax() [1/3]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::rmax ( DiscreteDomain< NonUniformPointSampling< CDim > > const &  d)

◆ rlength() [1/3]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::rlength ( DiscreteDomain< NonUniformPointSampling< CDim > > const &  d)

◆ parallel_deepcopy() [1/2]

template<class ChunkDst , class ChunkSrc >
auto ddc::parallel_deepcopy ( ChunkDst &&  dst,
ChunkSrc &&  src 
)

Copy the content of a borrowed chunk into another.

Parameters
[out]dstthe borrowed chunk in which to copy
[in]srcthe borrowed chunk from which to copy
Returns
dst as a ChunkSpan

◆ parallel_deepcopy() [2/2]

template<class ExecSpace , class ChunkDst , class ChunkSrc >
auto ddc::parallel_deepcopy ( ExecSpace const &  execution_space,
ChunkDst &&  dst,
ChunkSrc &&  src 
)

Copy the content of a borrowed chunk into another.

Parameters
[in]execution_spacea Kokkos execution space where the loop will be executed on
[out]dstthe borrowed chunk in which to copy
[in]srcthe borrowed chunk from which to copy
Returns
dst as a ChunkSpan

◆ parallel_fill() [1/2]

template<class ChunkDst , class T >
auto ddc::parallel_fill ( ChunkDst &&  dst,
T const &  value 
)

Fill a borrowed chunk with a given value.

Parameters
[out]dstthe borrowed chunk in which to copy
[in]valuethe value to fill dst
Returns
dst as a ChunkSpan

◆ parallel_fill() [2/2]

template<class ExecSpace , class ChunkDst , class T >
auto ddc::parallel_fill ( ExecSpace const &  execution_space,
ChunkDst &&  dst,
T const &  value 
)

Fill a borrowed chunk with a given value.

Parameters
[in]execution_spacea Kokkos execution space where the loop will be executed on
[out]dstthe borrowed chunk in which to copy
[in]valuethe value to fill dst
Returns
dst as a ChunkSpan

◆ parallel_for_each() [1/2]

template<class ExecSpace , class... DDims, class Functor >
void ddc::parallel_for_each ( ExecSpace const &  execution_space,
DiscreteDomain< DDims... > const &  domain,
Functor &&  f 
)
noexcept

iterates over a nD domain using a given Kokkos execution space

Parameters
[in]execution_spacea Kokkos execution space where the loop will be executed on
[in]domainthe domain over which to iterate
[in]fa functor taking an index as parameter

◆ parallel_for_each() [2/2]

template<class... DDims, class Functor >
void ddc::parallel_for_each ( DiscreteDomain< DDims... > const &  domain,
Functor &&  f 
)
noexcept

iterates over a nD domain using the Kokkos default execution space

Parameters
[in]domainthe domain over which to iterate
[in]fa functor taking an index as parameter

◆ parallel_transform_reduce() [1/2]

template<class ExecSpace , class... DDims, class T , class BinaryReductionOp , class UnaryTransformOp >
T ddc::parallel_transform_reduce ( ExecSpace const &  execution_space,
DiscreteDomain< DDims... > const &  domain,
neutral,
BinaryReductionOp &&  reduce,
UnaryTransformOp &&  transform 
)
noexcept

A reduction over a nD domain using a given Kokkos execution space.

Parameters
[in]execution_spacea Kokkos execution space where the loop will be executed on
[in]domainthe range over which to apply the algorithm
[in]neutralthe neutral element of the reduction operation
[in]reducea binary FunctionObject that will be applied in unspecified order to the results of transform, the results of other reduce and neutral.
[in]transforma unary FunctionObject that will be applied to each element of the input range. The return type must be acceptable as input to reduce

◆ parallel_transform_reduce() [2/2]

template<class... DDims, class T , class BinaryReductionOp , class UnaryTransformOp >
T ddc::parallel_transform_reduce ( DiscreteDomain< DDims... > const &  domain,
neutral,
BinaryReductionOp &&  reduce,
UnaryTransformOp &&  transform 
)
noexcept

A reduction over a nD domain using the Kokkos default execution space.

Parameters
[in]domainthe range over which to apply the algorithm
[in]neutralthe neutral element of the reduction operation
[in]reducea binary FunctionObject that will be applied in unspecified order to the results of transform, the results of other reduce and neutral.
[in]transforma unary FunctionObject that will be applied to each element of the input range. The return type must be acceptable as input to reduce

◆ expose_to_pdi() [1/2]

template<PDI_inout_t access, class DataType >
void ddc::expose_to_pdi ( std::string const &  name,
DataType &&  data 
)

◆ expose_to_pdi() [2/2]

template<class DataType >
void ddc::expose_to_pdi ( std::string const &  name,
DataType &&  data 
)

◆ origin() [1/2]

template<class DDim >
KOKKOS_FUNCTION std:: enable_if_t< is_periodic_sampling_v< DDim >, typename DDim::continuous_element_type > ddc::origin ( )
noexcept

Lower bound index of the mesh.

◆ front() [2/3]

template<class DDim >
KOKKOS_FUNCTION std::enable_if_t< is_periodic_sampling_v< DDim >, typename DDim::discrete_element_type > ddc::front ( )
noexcept

Lower bound index of the mesh.

◆ step() [1/2]

template<class DDim >
KOKKOS_FUNCTION std::enable_if_t< is_periodic_sampling_v< DDim >, Real > ddc::step ( )
noexcept

Spacing step of the mesh.

◆ coordinate() [3/4]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::coordinate ( DiscreteElement< PeriodicSampling< CDim > > const &  c)

◆ distance_at_left() [2/3]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::distance_at_left ( DiscreteElement< PeriodicSampling< CDim > >  )

◆ distance_at_right() [2/3]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::distance_at_right ( DiscreteElement< PeriodicSampling< CDim > >  )

◆ rmin() [2/3]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::rmin ( DiscreteDomain< PeriodicSampling< CDim > > const &  d)

◆ rmax() [2/3]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::rmax ( DiscreteDomain< PeriodicSampling< CDim > > const &  d)

◆ rlength() [2/3]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::rlength ( DiscreteDomain< PeriodicSampling< CDim > > const &  d)

◆ transform_reduce()

template<class... DDims, class T , class BinaryReductionOp , class UnaryTransformOp >
T ddc::transform_reduce ( DiscreteDomain< DDims... > const &  domain,
neutral,
BinaryReductionOp &&  reduce,
UnaryTransformOp &&  transform 
)
noexcept

A reduction over a nD domain in serial.

Parameters
[in]domainthe range over which to apply the algorithm
[in]neutralthe neutral element of the reduction operation
[in]reducea binary FunctionObject that will be applied in unspecified order to the results of transform, the results of other reduce and neutral.
[in]transforma unary FunctionObject that will be applied to each element of the input range. The return type must be acceptable as input to reduce

◆ coordinate() [4/4]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::coordinate ( DiscreteElement< UniformPointSampling< CDim > > const &  c)

◆ origin() [2/2]

template<class DDim >
KOKKOS_FUNCTION std:: enable_if_t< is_uniform_sampling_v< DDim >, typename DDim::continuous_element_type > ddc::origin ( )
noexcept

Lower bound index of the mesh.

◆ front() [3/3]

template<class DDim >
KOKKOS_FUNCTION std::enable_if_t< is_uniform_sampling_v< DDim >, typename DDim::discrete_element_type > ddc::front ( )
noexcept

Lower bound index of the mesh.

◆ step() [2/2]

template<class DDim >
KOKKOS_FUNCTION std::enable_if_t< is_uniform_sampling_v< DDim >, Real > ddc::step ( )
noexcept

Spacing step of the mesh.

◆ distance_at_left() [3/3]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::distance_at_left ( DiscreteElement< UniformPointSampling< CDim > >  )

◆ distance_at_right() [3/3]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::distance_at_right ( DiscreteElement< UniformPointSampling< CDim > >  )

◆ rmin() [3/3]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::rmin ( DiscreteDomain< UniformPointSampling< CDim > > const &  d)

◆ rmax() [3/3]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::rmax ( DiscreteDomain< UniformPointSampling< CDim > > const &  d)

◆ rlength() [3/3]

template<class CDim >
KOKKOS_FUNCTION Coordinate< CDim > ddc::rlength ( DiscreteDomain< UniformPointSampling< CDim > > const &  d)

Variable Documentation

◆ enable_chunk< Chunk< ElementType, SupportType, Allocator > >

template<class ElementType , class SupportType , class Allocator >
constexpr bool ddc::enable_chunk< Chunk< ElementType, SupportType, Allocator > > = true
inlineconstexpr

◆ enable_chunk< ChunkSpan< ElementType, SupportType, LayoutStridedPolicy, MemorySpace > >

template<class ElementType , class SupportType , class LayoutStridedPolicy , class MemorySpace >
constexpr bool ddc::enable_chunk< ChunkSpan< ElementType, SupportType, LayoutStridedPolicy, MemorySpace > > = true
inlineconstexpr

◆ enable_borrowed_chunk< ChunkSpan< ElementType, SupportType, LayoutStridedPolicy, MemorySpace > >

template<class ElementType , class SupportType , class LayoutStridedPolicy , class MemorySpace >
constexpr bool ddc::enable_borrowed_chunk< ChunkSpan< ElementType, SupportType, LayoutStridedPolicy, MemorySpace > > = true
inlineconstexpr

◆ enable_borrowed_chunk

template<class T >
constexpr bool ddc::enable_borrowed_chunk = false
inlineconstexpr

◆ enable_chunk

template<class T >
constexpr bool ddc::enable_chunk = false
inlineconstexpr

◆ is_chunk_v

template<class T >
constexpr bool ddc::is_chunk_v = enable_chunk<std::remove_const_t<std::remove_reference_t<T>>>
inlineconstexpr

◆ is_borrowed_chunk_v

template<class T >
constexpr bool ddc::is_borrowed_chunk_v
inlineconstexpr
Initial value:
T> && (std::is_lvalue_reference_v<T> || enable_borrowed_chunk<std::remove_cv_t<std::remove_reference_t<T>>>)
constexpr bool is_chunk_v
Definition: chunk_traits.hpp:19

◆ is_writable_chunk_v

template<class T >
constexpr bool ddc::is_writable_chunk_v = !std::is_const_v<std::remove_pointer_t<chunk_pointer_t<T>>>
inlineconstexpr

◆ is_discrete_domain_v

template<class T >
constexpr bool ddc::is_discrete_domain_v = is_discrete_domain<T>::value
inlineconstexpr

◆ is_discrete_element_v

template<class T >
constexpr bool ddc::is_discrete_element_v = is_discrete_element<T>::value
inlineconstexpr

◆ is_discrete_vector_v

template<class T >
constexpr bool ddc::is_discrete_vector_v = is_discrete_vector<T>::value
inlineconstexpr

◆ is_non_uniform_sampling_v

template<class DDim >
constexpr bool ddc::is_non_uniform_sampling_v = is_non_uniform_sampling<DDim>::value
constexpr

◆ is_periodic_sampling_v

template<class DDim >
constexpr bool ddc::is_periodic_sampling_v = is_periodic_sampling<DDim>::value
constexpr

◆ is_periodic_domain_v

template<class T >
constexpr bool ddc::is_periodic_domain_v = is_periodic_domain<T>::value
constexpr

◆ is_rectilinear_domain_v

template<class T >
constexpr bool ddc::is_rectilinear_domain_v = is_rectilinear_domain<T>::value
constexpr

◆ is_uniform_sampling_v

template<class DDim >
constexpr bool ddc::is_uniform_sampling_v = is_uniform_sampling<DDim>::value
constexpr

◆ is_uniform_domain_v

template<class T >
constexpr bool ddc::is_uniform_domain_v = is_uniform_domain<T>::value
constexpr