Gromacs  2024.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members | Public Types | Public Member Functions | Static Public Member Functions
gmx::basic_mdspan< ElementType, Extents, LayoutPolicy, AccessorPolicy > Class Template Reference

#include <gromacs/mdspan/mdspan.h>

Description

template<class ElementType, class Extents, class LayoutPolicy = layout_right, class AccessorPolicy = accessor_basic<ElementType>>
class gmx::basic_mdspan< ElementType, Extents, LayoutPolicy, AccessorPolicy >

Multidimensional array indexing and memory access with flexible mapping and access model.

Template Parameters
ElementTypeType of elemnt to be viewed
ExtentsThe dimensions of the multidimenisonal array to view.
LayoutPolicyDescribes is the memory layout of the multidimensional array; right by default.
AccessorPolicyDescribes memory access model.

Public Types

using extents_type = Extents
 Expose type used to define the extents of the data.
 
using layout_type = LayoutPolicy
 Expose type used to define the layout of the data.
 
using accessor_type = AccessorPolicy
 Expose type used to define the memory access model of the data.
 
using mapping_type = typename layout_type::template mapping< extents_type >
 Expose type used to map multidimensional indices to one-dimensioal indices.
 
using element_type = typename accessor_type::element_type
 Exposes the type of stored element.
 
using value_type = std::remove_cv_t< element_type >
 Expose the underlying type of the stored elements.
 
using index_type = ptrdiff_t
 Expose the type used for indexing.
 
using difference_type = ptrdiff_t
 Expose type for index differences.
 
using pointer = typename accessor_type::pointer
 Expose underlying pointer to data type.
 
using reference = typename accessor_type::reference
 Expose reference to data type.
 

Public Member Functions

constexpr basic_mdspan () noexcept
 Trivial constructor.
 
constexpr basic_mdspan (basic_mdspan &&other) noexcept=default
 Move constructor.
 
constexpr basic_mdspan (const basic_mdspan &other) noexcept=default
 copy constructor
 
basic_mdspanoperator= (const basic_mdspan &other) noexcept=default
 Copy assignment.
 
basic_mdspanoperator= (basic_mdspan &&other) noexcept=default
 Move assignment.
 
template<class OtherElementType , class OtherExtents , class OtherLayoutPolicy , class OtherAccessor >
constexpr basic_mdspan (const basic_mdspan< OtherElementType, OtherExtents, OtherLayoutPolicy, OtherAccessor > &rhs) noexcept
 Copy constructor.
 
template<class OtherElementType , class OtherExtents , class OtherLayoutPolicy , class OtherAccessor >
basic_mdspanoperator= (const basic_mdspan< OtherElementType, OtherExtents, OtherLayoutPolicy, OtherAccessor > &rhs) noexcept
 Copy assignment constructor.
 
template<class... IndexType>
constexpr basic_mdspan (pointer ptr, IndexType...DynamicExtents) noexcept
 Construct mdspan by setting the dynamic extents and pointer to data. More...
 
constexpr basic_mdspan (pointer ptr, const std::array< ptrdiff_t, extents_type::rank_dynamic()> &dynamic_extents)
 Construct from array describing dynamic extents. More...
 
constexpr basic_mdspan (pointer ptr, const mapping_type &m) noexcept
 Construct from pointer and mapping. More...
 
constexpr basic_mdspan (pointer ptr, const mapping_type &m, const accessor_type &a) noexcept
 Construct with pointer, mapping and accessor. More...
 
template<typename U , typename = std::enable_if_t<std::is_same_v<typename std::remove_reference_t<U>::view_type::element_type, ElementType>>>
constexpr basic_mdspan (U &&other)
 Construct mdspan from multidimensional arrays implemented with mdspan. More...
 
template<typename U , typename = std::enable_if_t<std::is_same_v<typename std::remove_reference_t<U>::const_view_type::element_type, ElementType>>>
constexpr basic_mdspan (const U &other)
 Construct mdspan of const Elements from multidimensional arrays implemented with mdspan. More...
 
template<class... IndexType>
constexpr std::enable_if_t
< sizeof...(IndexType)==extents_type::rank(),
reference
operator() (IndexType...indices) const noexcept
 Brace operator to access multidimensional array element. More...
 
template<class IndexType >
constexpr std::enable_if_t
< std::is_integral_v
< IndexType >
&&extents_type::rank()==1,
reference
operator[] (const IndexType &i) const noexcept
 Canonical bracket operator for one-dimensional arrays. Allows mdspan to act like array in one-dimension. Enabled only when rank==1. More...
 
template<class IndexType , typename sliced_mdspan_type = basic_mdspan<element_type, decltype(extents_type().sliced_extents()), LayoutPolicy, AccessorPolicy>>
constexpr std::enable_if_t
< std::is_integral_v
< IndexType >
&&(extents_type::rank() >
1)&&std::is_same_v
< LayoutPolicy, layout_right >
, sliced_mdspan_type > 
operator[] (const IndexType index) const noexcept
 Bracket operator for multi-dimensional arrays. More...
 
constexpr index_type static_extent (size_t k) const noexcept
 Return the static extent. More...
 
constexpr index_type extent (int k) const noexcept
 Return the extent. More...
 
constexpr const extents_typeextents () const noexcept
 Return all extents.
 
constexpr bool is_unique () const noexcept
 Report if the currently applied map is unique.
 
constexpr bool is_strided () const noexcept
 Report if the currently applied map is strided.
 
constexpr bool is_contiguous () const noexcept
 Report if the currently applied map is contiguous.
 
constexpr index_type stride (size_t r) const noexcept
 Report stride along a specific rank.
 
constexpr mapping_type mapping () const noexcept
 Return the currently applied mapping.
 
constexpr accessor_type accessor () const noexcept
 Return the memory access model.
 
constexpr pointer data () const noexcept
 Return pointer to underlying data.
 

Static Public Member Functions

static constexpr int rank () noexcept
 Report the rank.
 
static constexpr int rank_dynamic () noexcept
 Report the dynamic rank.
 
static constexpr bool is_always_unique () noexcept
 Report if mappings for this basic_span is always unique.
 
static constexpr bool is_always_strided () noexcept
 Report if mapping for this basic_span is always strided.
 
static constexpr bool is_always_contiguous () noexcept
 Report if mapping for this basic_span is always is_contiguous.
 

Constructor & Destructor Documentation

template<class ElementType, class Extents, class LayoutPolicy = layout_right, class AccessorPolicy = accessor_basic<ElementType>>
template<class... IndexType>
constexpr gmx::basic_mdspan< ElementType, Extents, LayoutPolicy, AccessorPolicy >::basic_mdspan ( pointer  ptr,
IndexType...  DynamicExtents 
)
inlineexplicitnoexcept

Construct mdspan by setting the dynamic extents and pointer to data.

Parameters
[in]ptrPointer to data to be accessed by this span
[in]DynamicExtents
Template Parameters
IndexTypeindex type to describe dynamic extents
template<class ElementType, class Extents, class LayoutPolicy = layout_right, class AccessorPolicy = accessor_basic<ElementType>>
constexpr gmx::basic_mdspan< ElementType, Extents, LayoutPolicy, AccessorPolicy >::basic_mdspan ( pointer  ptr,
const std::array< ptrdiff_t, extents_type::rank_dynamic()> &  dynamic_extents 
)
inline

Construct from array describing dynamic extents.

Parameters
[in]ptrPointer to data to be accessed by this span
[in]dynamic_extentsArray the size of dynamic extents.
template<class ElementType, class Extents, class LayoutPolicy = layout_right, class AccessorPolicy = accessor_basic<ElementType>>
constexpr gmx::basic_mdspan< ElementType, Extents, LayoutPolicy, AccessorPolicy >::basic_mdspan ( pointer  ptr,
const mapping_type m 
)
inlinenoexcept

Construct from pointer and mapping.

Parameters
[in]ptrPointer to data to be accessed by this span
[in]mMapping from multidimenisonal indices to one-dimensional offset.
template<class ElementType, class Extents, class LayoutPolicy = layout_right, class AccessorPolicy = accessor_basic<ElementType>>
constexpr gmx::basic_mdspan< ElementType, Extents, LayoutPolicy, AccessorPolicy >::basic_mdspan ( pointer  ptr,
const mapping_type m,
const accessor_type a 
)
inlinenoexcept

Construct with pointer, mapping and accessor.

Parameters
[in]ptrPointer to data to be accessed by this span
[in]mMapping from multidimenisonal indices to one-dimensional offset.
[in]aAccessor implementing memory access model.
template<class ElementType, class Extents, class LayoutPolicy = layout_right, class AccessorPolicy = accessor_basic<ElementType>>
template<typename U , typename = std::enable_if_t<std::is_same_v<typename std::remove_reference_t<U>::view_type::element_type, ElementType>>>
constexpr gmx::basic_mdspan< ElementType, Extents, LayoutPolicy, AccessorPolicy >::basic_mdspan ( U &&  other)
inline

Construct mdspan from multidimensional arrays implemented with mdspan.

Requires the container to have a view_type describing the mdspan, which is accessible through an asView() call

This allows functions to declare mdspans as arguments, but take e.g. multidimensional arrays implicitly during the function call

Template Parameters
Ucontainer type
Parameters
[in]othermdspan-implementing container
template<class ElementType, class Extents, class LayoutPolicy = layout_right, class AccessorPolicy = accessor_basic<ElementType>>
template<typename U , typename = std::enable_if_t<std::is_same_v<typename std::remove_reference_t<U>::const_view_type::element_type, ElementType>>>
constexpr gmx::basic_mdspan< ElementType, Extents, LayoutPolicy, AccessorPolicy >::basic_mdspan ( const U &  other)
inline

Construct mdspan of const Elements from multidimensional arrays implemented with mdspan.

Requires the container to have a const_view_type describing the mdspan, which is accessible through an asConstView() call

This allows functions to declare mdspans as arguments, but take e.g. multidimensional arrays implicitly during the function call

Template Parameters
Ucontainer type
Parameters
[in]othermdspan-implementing container

Member Function Documentation

template<class ElementType, class Extents, class LayoutPolicy = layout_right, class AccessorPolicy = accessor_basic<ElementType>>
constexpr index_type gmx::basic_mdspan< ElementType, Extents, LayoutPolicy, AccessorPolicy >::extent ( int  k) const
inlinenoexcept

Return the extent.

Parameters
[in]kdimension to query for extent
Returns
extent along specified dimension
template<class ElementType, class Extents, class LayoutPolicy = layout_right, class AccessorPolicy = accessor_basic<ElementType>>
template<class... IndexType>
constexpr std::enable_if_t<sizeof...(IndexType) == extents_type::rank(), reference> gmx::basic_mdspan< ElementType, Extents, LayoutPolicy, AccessorPolicy >::operator() ( IndexType...  indices) const
inlinenoexcept

Brace operator to access multidimensional array element.

Parameters
[in]indicesThe multidimensional indices of the object. Requires rank() == sizeof...(IndexType). Slicing is implemented via sub_span.
Returns
reference to element at indices.
template<class ElementType, class Extents, class LayoutPolicy = layout_right, class AccessorPolicy = accessor_basic<ElementType>>
template<class IndexType >
constexpr std::enable_if_t<std::is_integral_v<IndexType> && extents_type::rank() == 1, reference> gmx::basic_mdspan< ElementType, Extents, LayoutPolicy, AccessorPolicy >::operator[] ( const IndexType &  i) const
inlinenoexcept

Canonical bracket operator for one-dimensional arrays. Allows mdspan to act like array in one-dimension. Enabled only when rank==1.

Parameters
[in]ione-dimensional index
Returns
reference to element stored at position i
template<class ElementType, class Extents, class LayoutPolicy = layout_right, class AccessorPolicy = accessor_basic<ElementType>>
template<class IndexType , typename sliced_mdspan_type = basic_mdspan<element_type, decltype(extents_type().sliced_extents()), LayoutPolicy, AccessorPolicy>>
constexpr std::enable_if_t<std::is_integral_v<IndexType> && (extents_type::rank() > 1) && std::is_same_v<LayoutPolicy, layout_right>, sliced_mdspan_type> gmx::basic_mdspan< ElementType, Extents, LayoutPolicy, AccessorPolicy >::operator[] ( const IndexType  index) const
inlinenoexcept

Bracket operator for multi-dimensional arrays.

Note
Prefer operator() for better compile-time and run-time performance

Slices two- and higher-dimensional arrays along a given slice by returning a new basic_mdspan that drops the first extent and indexes the remaining extents

Note
Currently only implemented for layout_right
For layout_right this implementation has significant performance benefits over implementing a more general slicing operator with a strided layout
Enabled only when rank() > 1
Template Parameters
IndexTypeintegral tyoe for the index that enables indexing with, e.g., int or size_t
Parameters
[in]indexone-dimensional index of the slice to be indexed
Returns
basic_mdspan that is sliced at the given index
template<class ElementType, class Extents, class LayoutPolicy = layout_right, class AccessorPolicy = accessor_basic<ElementType>>
constexpr index_type gmx::basic_mdspan< ElementType, Extents, LayoutPolicy, AccessorPolicy >::static_extent ( size_t  k) const
inlinenoexcept

Return the static extent.

Parameters
[in]kdimension to query for static extent
Returns
static extent along specified dimension

The documentation for this class was generated from the following file: