10#ifndef XTENSOR_CONTAINER_HPP
11#define XTENSOR_CONTAINER_HPP
17#include <xtl/xmeta_utils.hpp>
18#include <xtl/xsequence.hpp>
20#include "../core/xaccessible.hpp"
21#include "../core/xiterable.hpp"
22#include "../core/xiterator.hpp"
23#include "../core/xmath.hpp"
24#include "../core/xoperation.hpp"
25#include "../core/xstrides.hpp"
26#include "../core/xtensor_config.hpp"
27#include "../core/xtensor_forward.hpp"
42 struct allocator_type_impl
44 using type =
typename T::allocator_type;
47 template <
class T, std::
size_t N>
48 struct allocator_type_impl<std::array<T, N>>
50 using type = std::allocator<T>;
55 using allocator_type_t =
typename detail::allocator_type_impl<T>::type;
70 private xaccessible<D>
74 using derived_type = D;
77 using storage_type =
typename inner_types::storage_type;
78 using allocator_type = allocator_type_t<std::decay_t<storage_type>>;
79 using value_type =
typename storage_type::value_type;
80 using reference =
typename inner_types::reference;
81 using const_reference =
typename inner_types::const_reference;
82 using pointer =
typename storage_type::pointer;
83 using const_pointer =
typename storage_type::const_pointer;
84 using size_type =
typename inner_types::size_type;
85 using difference_type =
typename storage_type::difference_type;
86 using simd_value_type = xt_simd::simd_type<value_type>;
87 using bool_load_type = xt::bool_load_type<value_type>;
89 using shape_type =
typename inner_types::shape_type;
90 using strides_type =
typename inner_types::strides_type;
91 using backstrides_type =
typename inner_types::backstrides_type;
93 using inner_shape_type =
typename inner_types::inner_shape_type;
94 using inner_strides_type =
typename inner_types::inner_strides_type;
95 using inner_backstrides_type =
typename inner_types::inner_backstrides_type;
98 using stepper =
typename iterable_base::stepper;
99 using const_stepper =
typename iterable_base::const_stepper;
101 using accessible_base = xaccessible<D>;
103 static constexpr layout_type static_layout = inner_types::layout;
105 using data_alignment = xt_simd::container_alignment_t<storage_type>;
106 using simd_type = xt_simd::simd_type<value_type>;
108 using linear_iterator =
typename iterable_base::linear_iterator;
109 using const_linear_iterator =
typename iterable_base::const_linear_iterator;
110 using reverse_linear_iterator =
typename iterable_base::reverse_linear_iterator;
111 using const_reverse_linear_iterator =
typename iterable_base::const_reverse_linear_iterator;
113 static_assert(static_layout !=
layout_type::any,
"Container layout can never be layout_type::any!");
115 size_type
size() const noexcept;
119 constexpr const inner_shape_type&
shape() const noexcept;
120 constexpr const inner_strides_type&
strides() const noexcept;
121 constexpr const inner_backstrides_type&
backstrides() const noexcept;
126 template <class... Args>
127 reference operator()(Args... args);
129 template <class... Args>
130 const_reference operator()(Args... args) const;
132 template <class... Args>
133 reference unchecked(Args... args);
135 template <class... Args>
136 const_reference unchecked(Args... args) const;
138 using accessible_base::at;
139 using accessible_base::
shape;
140 using accessible_base::operator[];
141 using accessible_base::
back;
142 using accessible_base::
front;
144 using accessible_base::periodic;
147 reference element(It first, It last);
149 const_reference element(It first, It last) const;
155 const_pointer
data() const noexcept;
164 stepper stepper_begin(const S&
shape) noexcept;
169 const_stepper stepper_begin(const S&
shape) const noexcept;
173 reference data_element(size_type i);
174 const_reference data_element(size_type i) const;
177 const_reference
flat(size_type i) const;
179 template <class requested_type>
180 using simd_return_type = xt_simd::simd_return_type<value_type, requested_type>;
182 template <class align, class simd>
183 void store_simd(size_type i, const simd& e);
184 template <class align, class requested_type = value_type, std::
size_t N = xt_simd::simd_traits<requested_type>::
size>
185 container_simd_return_type_t<storage_type, value_type, requested_type>
186 load_simd(size_type i) const;
188 linear_iterator linear_begin() noexcept;
189 linear_iterator linear_end() noexcept;
191 const_linear_iterator linear_begin() const noexcept;
192 const_linear_iterator linear_end() const noexcept;
193 const_linear_iterator linear_cbegin() const noexcept;
194 const_linear_iterator linear_cend() const noexcept;
196 reverse_linear_iterator linear_rbegin() noexcept;
197 reverse_linear_iterator linear_rend() noexcept;
199 const_reverse_linear_iterator linear_rbegin() const noexcept;
200 const_reverse_linear_iterator linear_rend() const noexcept;
201 const_reverse_linear_iterator linear_crbegin() const noexcept;
202 const_reverse_linear_iterator linear_crend() const noexcept;
204 using container_iterator = linear_iterator;
205 using const_container_iterator = const_linear_iterator;
209 xcontainer() = default;
210 ~xcontainer() = default;
212 xcontainer(const xcontainer&) = default;
213 xcontainer& operator=(const xcontainer&) = default;
215 xcontainer(xcontainer&&) = default;
216 xcontainer& operator=(xcontainer&&) = default;
218 container_iterator data_xbegin() noexcept;
219 const_container_iterator data_xbegin() const noexcept;
220 container_iterator data_xend(
layout_type l, size_type offset) noexcept;
221 const_container_iterator data_xend(
layout_type l, size_type offset) const noexcept;
225 derived_type& derived_cast() & noexcept;
226 const derived_type& derived_cast() const& noexcept;
227 derived_type derived_cast() && noexcept;
232 It data_xend_impl(It begin,
layout_type l, size_type offset) const noexcept;
234 inner_shape_type& mutable_shape();
235 inner_strides_type& mutable_strides();
236 inner_backstrides_type& mutable_backstrides();
239 friend class xstepper;
241 friend class xaccessible<D>;
258 class xstrided_container : public xcontainer<D>
262 using base_type = xcontainer<D>;
263 using storage_type =
typename base_type::storage_type;
264 using value_type =
typename base_type::value_type;
265 using reference =
typename base_type::reference;
266 using const_reference =
typename base_type::const_reference;
267 using pointer =
typename base_type::pointer;
268 using const_pointer =
typename base_type::const_pointer;
269 using size_type =
typename base_type::size_type;
270 using shape_type =
typename base_type::shape_type;
271 using strides_type =
typename base_type::strides_type;
272 using inner_shape_type =
typename base_type::inner_shape_type;
273 using inner_strides_type =
typename base_type::inner_strides_type;
274 using inner_backstrides_type =
typename base_type::inner_backstrides_type;
276 template <
class S = shape_type>
278 template <
class S = shape_type>
280 template <
class S = shape_type>
283 template <
class S = shape_type>
290 bool is_contiguous() const noexcept;
294 xstrided_container() noexcept;
295 ~xstrided_container() = default;
297 xstrided_container(const xstrided_container&) = default;
298 xstrided_container& operator=(const xstrided_container&) = default;
300 xstrided_container(xstrided_container&&) = default;
301 xstrided_container& operator=(xstrided_container&&) = default;
303 explicit xstrided_container(inner_shape_type&&, inner_strides_type&&) noexcept;
304 explicit xstrided_container(inner_shape_type&&, inner_strides_type&&, inner_backstrides_type&&,
layout_type&&) noexcept;
306 inner_shape_type& shape_impl() noexcept;
307 const inner_shape_type& shape_impl() const noexcept;
309 inner_strides_type& strides_impl() noexcept;
310 const inner_strides_type& strides_impl() const noexcept;
312 inner_backstrides_type& backstrides_impl() noexcept;
313 const inner_backstrides_type& backstrides_impl() const noexcept;
315 template <class S = shape_type>
317 template <class S = shape_type>
324 inner_shape_type m_shape;
325 inner_strides_type m_strides;
326 inner_backstrides_type m_backstrides;
336 inline It xcontainer<D>::data_xend_impl(It begin,
layout_type l, size_type offset) const noexcept
338 return strided_data_end(*
this, begin, l, offset);
342 inline auto xcontainer<D>::mutable_shape() -> inner_shape_type&
344 return derived_cast().shape_impl();
348 inline auto xcontainer<D>::mutable_strides() -> inner_strides_type&
350 return derived_cast().strides_impl();
354 inline auto xcontainer<D>::mutable_backstrides() -> inner_backstrides_type&
356 return derived_cast().backstrides_impl();
369 return contiguous_layout ?
storage().size() : compute_size(
shape());
378 return shape().size();
387 return derived_cast().shape_impl();
396 return derived_cast().strides_impl();
405 return derived_cast().backstrides_impl();
423 if (contiguous_layout)
425 std::fill(this->linear_begin(), this->linear_end(), value);
429 std::fill(this->begin(), this->end(), value);
440 template <
class... Args>
441 inline auto xcontainer<D>::operator()(Args... args) -> reference
443 XTENSOR_TRY(check_index(
shape(), args...));
444 XTENSOR_CHECK_DIMENSION(
shape(), args...);
445 size_type index = xt::data_offset<size_type>(
strides(), args...);
456 template <
class... Args>
457 inline auto xcontainer<D>::operator()(Args... args)
const -> const_reference
459 XTENSOR_TRY(check_index(
shape(), args...));
460 XTENSOR_CHECK_DIMENSION(
shape(), args...);
461 size_type index = xt::data_offset<size_type>(
strides(), args...);
485 template <
class... Args>
486 inline auto xcontainer<D>::unchecked(Args... args) -> reference
488 size_type index = xt::unchecked_data_offset<size_type, static_layout>(
490 static_cast<std::ptrdiff_t
>(args)...
515 template <
class... Args>
516 inline auto xcontainer<D>::unchecked(Args... args)
const -> const_reference
518 size_type index = xt::unchecked_data_offset<size_type, static_layout>(
520 static_cast<std::ptrdiff_t
>(args)...
534 inline auto xcontainer<D>::element(It first, It last) -> reference
536 XTENSOR_TRY(check_element_index(
shape(), first, last));
537 return storage()[element_offset<size_type>(
strides(), first, last)];
549 inline auto xcontainer<D>::element(It first, It last)
const -> const_reference
551 XTENSOR_TRY(check_element_index(
shape(), first, last));
552 return storage()[element_offset<size_type>(
strides(), first, last)];
561 return derived_cast().storage_impl();
571 return derived_cast().storage_impl();
621 return xt::broadcast_shape(this->
shape(), shape);
633 return str.size() ==
strides().size() && std::equal(str.cbegin(), str.cend(),
strides().begin());
639 inline auto xcontainer<D>::derived_cast() const& noexcept -> const derived_type&
641 return *
static_cast<const derived_type*
>(
this);
645 inline auto xcontainer<D>::derived_cast() &&
noexcept -> derived_type
647 return *
static_cast<derived_type*
>(
this);
651 inline auto xcontainer<D>::data_element(size_type i) -> reference
657 inline auto xcontainer<D>::data_element(size_type i)
const -> const_reference
671 XTENSOR_ASSERT(i <
size());
684 XTENSOR_ASSERT(i <
size());
694 inline auto xcontainer<D>::stepper_begin(
const S& shape)
noexcept -> stepper
696 size_type offset = shape.size() -
dimension();
697 return stepper(
static_cast<derived_type*
>(
this), data_xbegin(), offset);
702 inline auto xcontainer<D>::stepper_end(
const S&
shape,
layout_type l)
noexcept -> stepper
705 return stepper(
static_cast<derived_type*
>(
this), data_xend(l, offset), offset);
710 inline auto xcontainer<D>::stepper_begin(
const S&
shape)
const noexcept -> const_stepper
713 return const_stepper(
static_cast<const derived_type*
>(
this), data_xbegin(), offset);
718 inline auto xcontainer<D>::stepper_end(
const S&
shape,
layout_type l)
const noexcept -> const_stepper
721 return const_stepper(
static_cast<const derived_type*
>(
this), data_xend(l, offset), offset);
725 inline auto xcontainer<D>::data_xbegin() noexcept -> container_iterator
731 inline auto xcontainer<D>::data_xbegin() const noexcept -> const_container_iterator
737 inline auto xcontainer<D>::data_xend(
layout_type l, size_type offset)
noexcept -> container_iterator
739 return data_xend_impl(
storage().begin(), l, offset);
743 inline auto xcontainer<D>::data_xend(
layout_type l, size_type offset)
const noexcept
744 -> const_container_iterator
746 return data_xend_impl(
storage().cbegin(), l, offset);
750 template <
class alignment,
class simd>
751 inline void xcontainer<D>::store_simd(size_type i,
const simd& e)
753 using align_mode = driven_align_mode_t<alignment, data_alignment>;
754 xt_simd::store_as(std::addressof(
storage()[i]), e, align_mode());
758 template <
class alignment,
class requested_type, std::
size_t N>
759 inline auto xcontainer<D>::load_simd(size_type i)
const
760 -> container_simd_return_type_t<storage_type, value_type, requested_type>
762 using align_mode = driven_align_mode_t<alignment, data_alignment>;
763 return xt_simd::load_as<requested_type>(std::addressof(
storage()[i]), align_mode());
767 inline auto xcontainer<D>::linear_begin() noexcept -> linear_iterator
773 inline auto xcontainer<D>::linear_end() noexcept -> linear_iterator
779 inline auto xcontainer<D>::linear_begin() const noexcept -> const_linear_iterator
785 inline auto xcontainer<D>::linear_end() const noexcept -> const_linear_iterator
791 inline auto xcontainer<D>::linear_cbegin() const noexcept -> const_linear_iterator
797 inline auto xcontainer<D>::linear_cend() const noexcept -> const_linear_iterator
803 inline auto xcontainer<D>::linear_rbegin() noexcept -> reverse_linear_iterator
809 inline auto xcontainer<D>::linear_rend() noexcept -> reverse_linear_iterator
815 inline auto xcontainer<D>::linear_rbegin() const noexcept -> const_reverse_linear_iterator
821 inline auto xcontainer<D>::linear_rend() const noexcept -> const_reverse_linear_iterator
827 inline auto xcontainer<D>::linear_crbegin() const noexcept -> const_reverse_linear_iterator
833 inline auto xcontainer<D>::linear_crend() const noexcept -> const_reverse_linear_iterator
839 inline auto xcontainer<D>::derived_cast() &
noexcept -> derived_type&
841 return *
static_cast<derived_type*
>(
this);
849 inline xstrided_container<D>::xstrided_container() noexcept
858 inline xstrided_container<D>::xstrided_container(inner_shape_type&&
shape, inner_strides_type&&
strides) noexcept
860 , m_shape(std::move(
shape))
861 , m_strides(std::move(
strides))
863 m_backstrides = xtl::make_sequence<inner_backstrides_type>(m_shape.size(), 0);
864 adapt_strides(m_shape, m_strides, m_backstrides);
868 inline xstrided_container<D>::xstrided_container(
869 inner_shape_type&&
shape,
875 , m_shape(std::move(
shape))
876 , m_strides(std::move(
strides))
878 , m_layout(std::move(layout))
883 inline auto xstrided_container<D>::shape_impl() noexcept -> inner_shape_type&
889 inline auto xstrided_container<D>::shape_impl() const noexcept -> const inner_shape_type&
895 inline auto xstrided_container<D>::strides_impl() noexcept -> inner_strides_type&
901 inline auto xstrided_container<D>::strides_impl() const noexcept -> const inner_strides_type&
907 inline auto xstrided_container<D>::backstrides_impl() noexcept -> inner_backstrides_type&
909 return m_backstrides;
913 inline auto xstrided_container<D>::backstrides_impl() const noexcept -> const inner_backstrides_type&
915 return m_backstrides;
929 inline bool xstrided_container<D>::is_contiguous() const noexcept
931 using str_type =
typename inner_strides_type::value_type;
932 auto is_zero = [](
auto i)
936 if (!is_contiguous_container<storage_type>::value)
944 auto it = std::find_if_not(m_strides.rbegin(), m_strides.rend(), is_zero);
946 return it == m_strides.rend() || *it == str_type(1);
950 auto it = std::find_if_not(m_strides.begin(), m_strides.end(), is_zero);
952 return it == m_strides.end() || *it == str_type(1);
956 return m_strides.empty();
962 template <
class C,
class S>
963 inline void resize_data_container(C& c, S size)
965 xt::resize_container(c, size);
968 template <
class C,
class S>
969 inline void resize_data_container(
const C& c, S size)
973 XTENSOR_ASSERT_MSG(c.size() == size,
"Trying to resize const data container with wrong size.");
976 template <
class S,
class T>
977 constexpr bool check_resize_dimension(
const S&,
const T&)
982 template <
class T,
size_t N,
class S>
983 constexpr bool check_resize_dimension(
const std::array<T, N>&,
const S& s)
985 return N == s.size();
1001 detail::check_resize_dimension(m_shape,
shape),
1002 "cannot change the number of dimensions of xtensor"
1004 std::size_t dim =
shape.size();
1005 if (m_shape.size() != dim || !std::equal(std::begin(
shape), std::end(
shape), std::begin(m_shape))
1010 m_layout = XTENSOR_DEFAULT_LAYOUT;
1012 m_shape = xtl::forward_sequence<shape_type, S>(
shape);
1014 resize_container(m_strides, dim);
1015 resize_container(m_backstrides, dim);
1017 detail::resize_data_container(this->
storage(), data_size);
1033 detail::check_resize_dimension(m_shape,
shape),
1034 "cannot change the number of dimensions of xtensor"
1040 "Cannot change layout_type if template parameter not layout_type::dynamic."
1059 detail::check_resize_dimension(m_shape,
shape),
1060 "cannot change the number of dimensions of xtensor"
1066 "Cannot resize with custom strides when layout() is != layout_type::dynamic."
1069 m_shape = xtl::forward_sequence<shape_type, S>(
shape);
1071 resize_container(m_backstrides, m_strides.size());
1072 adapt_strides(m_shape, m_strides, m_backstrides);
1074 detail::resize_data_container(this->
storage(), compute_size(m_shape));
1095 std::forward<S>(
shape),
1096 xtl::is_signed<std::decay_t<
typename std::decay_t<S>::value_type>>(),
1097 std::forward<layout_type>(
layout)
1099 return this->derived_cast();
1106 using sh_type = rebind_container_t<T, shape_type>;
1107 sh_type sh = xtl::make_sequence<sh_type>(shape.size());
1108 std::copy(shape.begin(), shape.end(), sh.begin());
1109 reshape_impl(std::move(sh), xtl::is_signed<T>(), std::forward<layout_type>(layout));
1110 return this->derived_cast();
1116 xstrided_container<D>::reshape_impl(S&&
shape, std::false_type ,
layout_type layout)
1118 if (compute_size(
shape) != this->
size())
1122 "Cannot reshape with incorrect number of elements. Do you mean to resize?"
1127 layout = XTENSOR_DEFAULT_LAYOUT;
1131 XTENSOR_THROW(std::runtime_error,
"Cannot reshape with different layout if static layout != dynamic.");
1134 m_shape = xtl::forward_sequence<shape_type, S>(
shape);
1135 resize_container(m_strides, m_shape.size());
1136 resize_container(m_backstrides, m_shape.size());
1143 xstrided_container<D>::reshape_impl(S&& _shape, std::true_type ,
layout_type layout)
1145 using tmp_value_type =
typename std::decay_t<S>::value_type;
1146 auto new_size = compute_size(_shape);
1147 if (this->
size() % new_size)
1149 XTENSOR_THROW(std::runtime_error,
"Negative axis size cannot be inferred. Shape mismatch.");
1151 std::decay_t<S>
shape = _shape;
1152 tmp_value_type accumulator = 1;
1153 std::size_t neg_idx = 0;
1155 for (
auto it =
shape.begin(); it !=
shape.end(); ++it, i++)
1160 XTENSOR_ASSERT(dim == -1 && !neg_idx);
1165 if (accumulator < 0)
1167 shape[neg_idx] =
static_cast<tmp_value_type
>(this->
size()) / std::abs(accumulator);
1169 else if (this->
size() != new_size)
1173 "Cannot reshape with incorrect number of elements. Do you mean to resize?"
1177 m_shape = xtl::forward_sequence<shape_type, S>(
shape);
1178 resize_container(m_strides, m_shape.size());
1179 resize_container(m_backstrides, m_shape.size());
1184 inline auto xstrided_container<D>::mutable_layout() noexcept ->
layout_type&
Base class for implementation of common expression constant access methods.
size_type dimension() const noexcept
Returns the number of dimensions of the expression.
bool in_bounds(Args... args) const
storage_type & storage() noexcept
size_type size() const noexcept
Returns the number of element in the container.
constexpr const inner_backstrides_type & backstrides() const noexcept
constexpr size_type dimension() const noexcept
reference flat(size_type i)
constexpr const inner_strides_type & strides() const noexcept
bool has_linear_assign(const S &strides) const noexcept
bool broadcast_shape(S &shape, bool reuse_cache=false) const
constexpr const inner_shape_type & shape() const noexcept
const size_type data_offset() const noexcept
void fill(const T &value)
Base class for multidimensional iterable expressions with contiguous storage.
layout_type layout() const noexcept
void resize(S &&shape, const strides_type &strides)
Resizes the container.
void resize(S &&shape, layout_type l)
Resizes the container.
auto & reshape(S &&shape, layout_type layout=base_type::static_layout) &
Reshapes the container and keeps old elements.
void resize(S &&shape, bool force=false)
Resizes the container.
std::size_t compute_strides(const shape_type &shape, layout_type l, strides_type &strides)
Compute the strides given the shape and the layout of an array.
standard mathematical functions for xexpressions