10#ifndef XTENSOR_CONTAINER_HPP 
   11#define XTENSOR_CONTAINER_HPP 
   19#include <xtl/xmeta_utils.hpp> 
   20#include <xtl/xsequence.hpp> 
   22#include "../core/xaccessible.hpp" 
   23#include "../core/xiterable.hpp" 
   24#include "../core/xiterator.hpp" 
   25#include "../core/xmath.hpp" 
   26#include "../core/xoperation.hpp" 
   27#include "../core/xstrides.hpp" 
   28#include "../core/xtensor_config.hpp" 
   29#include "../core/xtensor_forward.hpp" 
   44        struct allocator_type_impl
 
   46            using type = 
typename T::allocator_type;
 
   49        template <
class T, std::
size_t N>
 
   50        struct allocator_type_impl<std::array<T, N>>
 
   52            using type = std::allocator<T>;  
 
   57    using allocator_type_t = 
typename detail::allocator_type_impl<T>::type;
 
   72                       private xaccessible<D>
 
   76        using derived_type = D;
 
   79        using storage_type = 
typename inner_types::storage_type;
 
   80        using allocator_type = allocator_type_t<std::decay_t<storage_type>>;
 
   81        using value_type = 
typename storage_type::value_type;
 
   82        using reference = 
typename inner_types::reference;
 
   83        using const_reference = 
typename inner_types::const_reference;
 
   84        using pointer = 
typename storage_type::pointer;
 
   85        using const_pointer = 
typename storage_type::const_pointer;
 
   86        using size_type = 
typename inner_types::size_type;
 
   87        using difference_type = 
typename storage_type::difference_type;
 
   88        using simd_value_type = xt_simd::simd_type<value_type>;
 
   89        using bool_load_type = xt::bool_load_type<value_type>;
 
   91        using shape_type = 
typename inner_types::shape_type;
 
   92        using strides_type = 
typename inner_types::strides_type;
 
   93        using backstrides_type = 
typename inner_types::backstrides_type;
 
   95        using inner_shape_type = 
typename inner_types::inner_shape_type;
 
   96        using inner_strides_type = 
typename inner_types::inner_strides_type;
 
   97        using inner_backstrides_type = 
typename inner_types::inner_backstrides_type;
 
  100        using stepper = 
typename iterable_base::stepper;
 
  101        using const_stepper = 
typename iterable_base::const_stepper;
 
  103        using accessible_base = xaccessible<D>;
 
  105        static constexpr layout_type static_layout = inner_types::layout;
 
  107        using data_alignment = xt_simd::container_alignment_t<storage_type>;
 
  108        using simd_type = xt_simd::simd_type<value_type>;
 
  110        using linear_iterator = 
typename iterable_base::linear_iterator;
 
  111        using const_linear_iterator = 
typename iterable_base::const_linear_iterator;
 
  112        using reverse_linear_iterator = 
typename iterable_base::reverse_linear_iterator;
 
  113        using const_reverse_linear_iterator = 
typename iterable_base::const_reverse_linear_iterator;
 
  115        static_assert(static_layout != 
layout_type::any, 
"Container layout can never be layout_type::any!");
 
  117        size_type 
size() const noexcept;
 
  119        XTENSOR_CONSTEXPR_RETURN size_type 
dimension() const noexcept;
 
  121        XTENSOR_CONSTEXPR_RETURN const inner_shape_type& 
shape() const noexcept;
 
  122        XTENSOR_CONSTEXPR_RETURN const inner_strides_type& 
strides() const noexcept;
 
  123        XTENSOR_CONSTEXPR_RETURN const inner_backstrides_type& 
backstrides() const noexcept;
 
  128        template <class... Args>
 
  129        reference operator()(Args... args);
 
  131        template <class... Args>
 
  132        const_reference operator()(Args... args) const;
 
  134        template <class... Args>
 
  135        reference unchecked(Args... args);
 
  137        template <class... Args>
 
  138        const_reference unchecked(Args... args) const;
 
  140        using accessible_base::at;
 
  141        using accessible_base::
shape;
 
  142        using accessible_base::operator[];
 
  143        using accessible_base::
back;
 
  144        using accessible_base::
front;
 
  146        using accessible_base::periodic;
 
  149        reference element(It first, It last);
 
  151        const_reference element(It first, It last) const;
 
  157        const_pointer 
data() const noexcept;
 
  166        stepper stepper_begin(const S& 
shape) noexcept;
 
  171        const_stepper stepper_begin(const S& 
shape) const noexcept;
 
  175        reference data_element(size_type i);
 
  176        const_reference data_element(size_type i) const;
 
  179        const_reference 
flat(size_type i) const;
 
  181        template <class requested_type>
 
  182        using simd_return_type = xt_simd::simd_return_type<value_type, requested_type>;
 
  184        template <class align, class simd>
 
  185        void store_simd(size_type i, const simd& e);
 
  186        template <class align, class requested_type = value_type, std::
size_t N = xt_simd::simd_traits<requested_type>::
size>
 
  187        container_simd_return_type_t<storage_type, value_type, requested_type>
 
  188         load_simd(size_type i) const;
 
  190        linear_iterator linear_begin() noexcept;
 
  191        linear_iterator linear_end() noexcept;
 
  193        const_linear_iterator linear_begin() const noexcept;
 
  194        const_linear_iterator linear_end() const noexcept;
 
  195        const_linear_iterator linear_cbegin() const noexcept;
 
  196        const_linear_iterator linear_cend() const noexcept;
 
  198        reverse_linear_iterator linear_rbegin() noexcept;
 
  199        reverse_linear_iterator linear_rend() noexcept;
 
  201        const_reverse_linear_iterator linear_rbegin() const noexcept;
 
  202        const_reverse_linear_iterator linear_rend() const noexcept;
 
  203        const_reverse_linear_iterator linear_crbegin() const noexcept;
 
  204        const_reverse_linear_iterator linear_crend() const noexcept;
 
  206        using container_iterator = linear_iterator;
 
  207        using const_container_iterator = const_linear_iterator;
 
  211        xcontainer() = default;
 
  212        ~xcontainer() = default;
 
  214        xcontainer(const xcontainer&) = default;
 
  215        xcontainer& operator=(const xcontainer&) = default;
 
  217        xcontainer(xcontainer&&) = default;
 
  218        xcontainer& operator=(xcontainer&&) = default;
 
  220        container_iterator data_xbegin() noexcept;
 
  221        const_container_iterator data_xbegin() const noexcept;
 
  222        container_iterator data_xend(
layout_type l, size_type offset) noexcept;
 
  223        const_container_iterator data_xend(
layout_type l, size_type offset) const noexcept;
 
  227        derived_type& derived_cast() & noexcept;
 
  228        const derived_type& derived_cast() const& noexcept;
 
  229        derived_type derived_cast() && noexcept;
 
  234        It data_xend_impl(It begin, 
layout_type l, size_type offset) const noexcept;
 
  236        inner_shape_type& mutable_shape();
 
  237        inner_strides_type& mutable_strides();
 
  238        inner_backstrides_type& mutable_backstrides();
 
  241        friend class xstepper;
 
  243        friend class xaccessible<D>;
 
 
  260    class xstrided_container : public xcontainer<D>
 
  264        using base_type = xcontainer<D>;
 
  265        using storage_type = 
typename base_type::storage_type;
 
  266        using value_type = 
typename base_type::value_type;
 
  267        using reference = 
typename base_type::reference;
 
  268        using const_reference = 
typename base_type::const_reference;
 
  269        using pointer = 
typename base_type::pointer;
 
  270        using const_pointer = 
typename base_type::const_pointer;
 
  271        using size_type = 
typename base_type::size_type;
 
  272        using shape_type = 
typename base_type::shape_type;
 
  273        using strides_type = 
typename base_type::strides_type;
 
  274        using inner_shape_type = 
typename base_type::inner_shape_type;
 
  275        using inner_strides_type = 
typename base_type::inner_strides_type;
 
  276        using inner_backstrides_type = 
typename base_type::inner_backstrides_type;
 
  278        template <
class S = shape_type>
 
  280        template <
class S = shape_type>
 
  282        template <
class S = shape_type>
 
  285        template <
class S = shape_type>
 
  292        bool is_contiguous() const noexcept;
 
  296        xstrided_container() noexcept;
 
  297        ~xstrided_container() = default;
 
  299        xstrided_container(const xstrided_container&) = default;
 
  300        xstrided_container& operator=(const xstrided_container&) = default;
 
  302        xstrided_container(xstrided_container&&) = default;
 
  303        xstrided_container& operator=(xstrided_container&&) = default;
 
  305        explicit xstrided_container(inner_shape_type&&, inner_strides_type&&) noexcept;
 
  306        explicit xstrided_container(inner_shape_type&&, inner_strides_type&&, inner_backstrides_type&&, 
layout_type&&) noexcept;
 
  308        inner_shape_type& shape_impl() noexcept;
 
  309        const inner_shape_type& shape_impl() const noexcept;
 
  311        inner_strides_type& strides_impl() noexcept;
 
  312        const inner_strides_type& strides_impl() const noexcept;
 
  314        inner_backstrides_type& backstrides_impl() noexcept;
 
  315        const inner_backstrides_type& backstrides_impl() const noexcept;
 
  317        template <class S = shape_type>
 
  319        template <class S = shape_type>
 
  326        inner_shape_type m_shape;
 
  327        inner_strides_type m_strides;
 
  328        inner_backstrides_type m_backstrides;
 
 
  338    inline It xcontainer<D>::data_xend_impl(It begin, 
layout_type l, size_type offset) const noexcept
 
  340        return strided_data_end(*
this, begin, l, offset);
 
  344    inline auto xcontainer<D>::mutable_shape() -> inner_shape_type&
 
  346        return derived_cast().shape_impl();
 
  350    inline auto xcontainer<D>::mutable_strides() -> inner_strides_type&
 
  352        return derived_cast().strides_impl();
 
  356    inline auto xcontainer<D>::mutable_backstrides() -> inner_backstrides_type&
 
  358        return derived_cast().backstrides_impl();
 
  371        return contiguous_layout ? 
storage().size() : compute_size(
shape());
 
 
  380        return shape().size();
 
 
  389        return derived_cast().shape_impl();
 
 
  398        return derived_cast().strides_impl();
 
 
  407        return derived_cast().backstrides_impl();
 
 
  425        if (contiguous_layout)
 
  427            std::fill(this->linear_begin(), this->linear_end(), value);
 
  431            std::fill(this->begin(), this->end(), value);
 
 
  442    template <
class... Args>
 
  443    inline auto xcontainer<D>::operator()(Args... args) -> reference
 
  445        XTENSOR_TRY(check_index(
shape(), args...));
 
  446        XTENSOR_CHECK_DIMENSION(
shape(), args...);
 
  447        size_type index = xt::data_offset<size_type>(
strides(), args...);
 
 
  458    template <
class... Args>
 
  459    inline auto xcontainer<D>::operator()(Args... args) 
const -> const_reference
 
  461        XTENSOR_TRY(check_index(
shape(), args...));
 
  462        XTENSOR_CHECK_DIMENSION(
shape(), args...);
 
  463        size_type index = xt::data_offset<size_type>(
strides(), args...);
 
 
  487    template <
class... Args>
 
  488    inline auto xcontainer<D>::unchecked(Args... args) -> reference
 
  490        size_type index = xt::unchecked_data_offset<size_type, static_layout>(
 
  492            static_cast<std::ptrdiff_t
>(args)...
 
 
  517    template <
class... Args>
 
  518    inline auto xcontainer<D>::unchecked(Args... args) 
const -> const_reference
 
  520        size_type index = xt::unchecked_data_offset<size_type, static_layout>(
 
  522            static_cast<std::ptrdiff_t
>(args)...
 
 
  536    inline auto xcontainer<D>::element(It first, It last) -> reference
 
  538        XTENSOR_TRY(check_element_index(
shape(), first, last));
 
  539        return storage()[element_offset<size_type>(
strides(), first, last)];
 
 
  551    inline auto xcontainer<D>::element(It first, It last) 
const -> const_reference
 
  553        XTENSOR_TRY(check_element_index(
shape(), first, last));
 
  554        return storage()[element_offset<size_type>(
strides(), first, last)];
 
 
  563        return derived_cast().storage_impl();
 
 
  573        return derived_cast().storage_impl();
 
 
  623        return xt::broadcast_shape(this->
shape(), shape);
 
 
  635        return str.size() == 
strides().size() && std::equal(str.cbegin(), str.cend(), 
strides().begin());
 
 
  641    inline auto xcontainer<D>::derived_cast() const& noexcept -> const derived_type&
 
  643        return *
static_cast<const derived_type*
>(
this);
 
  647    inline auto xcontainer<D>::derived_cast() && 
noexcept -> derived_type
 
  649        return *
static_cast<derived_type*
>(
this);
 
  653    inline auto xcontainer<D>::data_element(size_type i) -> reference
 
  659    inline auto xcontainer<D>::data_element(size_type i) 
const -> const_reference
 
  673        XTENSOR_ASSERT(i < 
size());
 
 
  686        XTENSOR_ASSERT(i < 
size());
 
 
  696    inline auto xcontainer<D>::stepper_begin(
const S& shape) 
noexcept -> stepper
 
  698        size_type offset = shape.size() - 
dimension();
 
  699        return stepper(
static_cast<derived_type*
>(
this), data_xbegin(), offset);
 
  704    inline auto xcontainer<D>::stepper_end(
const S& 
shape, 
layout_type l) 
noexcept -> stepper
 
  707        return stepper(
static_cast<derived_type*
>(
this), data_xend(l, offset), offset);
 
  712    inline auto xcontainer<D>::stepper_begin(
const S& 
shape) 
const noexcept -> const_stepper
 
  715        return const_stepper(
static_cast<const derived_type*
>(
this), data_xbegin(), offset);
 
  720    inline auto xcontainer<D>::stepper_end(
const S& 
shape, 
layout_type l) 
const noexcept -> const_stepper
 
  723        return const_stepper(
static_cast<const derived_type*
>(
this), data_xend(l, offset), offset);
 
  727    inline auto xcontainer<D>::data_xbegin() noexcept -> container_iterator
 
  733    inline auto xcontainer<D>::data_xbegin() const noexcept -> const_container_iterator
 
  739    inline auto xcontainer<D>::data_xend(
layout_type l, size_type offset) 
noexcept -> container_iterator
 
  741        return data_xend_impl(
storage().begin(), l, offset);
 
  745    inline auto xcontainer<D>::data_xend(
layout_type l, size_type offset) 
const noexcept 
  746        -> const_container_iterator
 
  748        return data_xend_impl(
storage().cbegin(), l, offset);
 
  752    template <
class alignment, 
class simd>
 
  753    inline void xcontainer<D>::store_simd(size_type i, 
const simd& e)
 
  755        using align_mode = driven_align_mode_t<alignment, data_alignment>;
 
  756        xt_simd::store_as(std::addressof(
storage()[i]), e, align_mode());
 
  760    template <
class alignment, 
class requested_type, std::
size_t N>
 
  761    inline auto xcontainer<D>::load_simd(size_type i) 
const 
  762        -> container_simd_return_type_t<storage_type, value_type, requested_type>
 
  764        using align_mode = driven_align_mode_t<alignment, data_alignment>;
 
  765        return xt_simd::load_as<requested_type>(std::addressof(
storage()[i]), align_mode());
 
  769    inline auto xcontainer<D>::linear_begin() noexcept -> linear_iterator
 
  775    inline auto xcontainer<D>::linear_end() noexcept -> linear_iterator
 
  781    inline auto xcontainer<D>::linear_begin() const noexcept -> const_linear_iterator
 
  787    inline auto xcontainer<D>::linear_end() const noexcept -> const_linear_iterator
 
  793    inline auto xcontainer<D>::linear_cbegin() const noexcept -> const_linear_iterator
 
  799    inline auto xcontainer<D>::linear_cend() const noexcept -> const_linear_iterator
 
  805    inline auto xcontainer<D>::linear_rbegin() noexcept -> reverse_linear_iterator
 
  811    inline auto xcontainer<D>::linear_rend() noexcept -> reverse_linear_iterator
 
  817    inline auto xcontainer<D>::linear_rbegin() const noexcept -> const_reverse_linear_iterator
 
  823    inline auto xcontainer<D>::linear_rend() const noexcept -> const_reverse_linear_iterator
 
  829    inline auto xcontainer<D>::linear_crbegin() const noexcept -> const_reverse_linear_iterator
 
  835    inline auto xcontainer<D>::linear_crend() const noexcept -> const_reverse_linear_iterator
 
  841    inline auto xcontainer<D>::derived_cast() & 
noexcept -> derived_type&
 
  843        return *
static_cast<derived_type*
>(
this);
 
  851    inline xstrided_container<D>::xstrided_container() noexcept
 
  860    inline xstrided_container<D>::xstrided_container(inner_shape_type&& 
shape, inner_strides_type&& 
strides) noexcept
 
  862        , m_shape(std::move(
shape))
 
  863        , m_strides(std::move(
strides))
 
  865        m_backstrides = xtl::make_sequence<inner_backstrides_type>(m_shape.size(), 0);
 
  866        adapt_strides(m_shape, m_strides, m_backstrides);
 
  870    inline xstrided_container<D>::xstrided_container(
 
  871        inner_shape_type&& 
shape,
 
  877        , m_shape(std::move(
shape))
 
  878        , m_strides(std::move(
strides))
 
  880        , m_layout(std::move(layout))
 
  885    inline auto xstrided_container<D>::shape_impl() noexcept -> inner_shape_type&
 
  891    inline auto xstrided_container<D>::shape_impl() const noexcept -> const inner_shape_type&
 
  897    inline auto xstrided_container<D>::strides_impl() noexcept -> inner_strides_type&
 
  903    inline auto xstrided_container<D>::strides_impl() const noexcept -> const inner_strides_type&
 
  909    inline auto xstrided_container<D>::backstrides_impl() noexcept -> inner_backstrides_type&
 
  911        return m_backstrides;
 
  915    inline auto xstrided_container<D>::backstrides_impl() const noexcept -> const inner_backstrides_type&
 
  917        return m_backstrides;
 
  931    inline bool xstrided_container<D>::is_contiguous() const noexcept
 
  933        using str_type = 
typename inner_strides_type::value_type;
 
  934        auto is_zero = [](
auto i)
 
  938        if (!is_contiguous_container<storage_type>::value)
 
  946            auto it = std::find_if_not(m_strides.rbegin(), m_strides.rend(), is_zero);
 
  948            return it == m_strides.rend() || *it == str_type(1);
 
  952            auto it = std::find_if_not(m_strides.begin(), m_strides.end(), is_zero);
 
  954            return it == m_strides.end() || *it == str_type(1);
 
  958            return m_strides.empty();
 
  964        template <
class C, 
class S>
 
  965        inline void resize_data_container(C& c, S size)
 
  967            xt::resize_container(c, size);
 
  970        template <
class C, 
class S>
 
  971        inline void resize_data_container(
const C& c, S size)
 
  975            XTENSOR_ASSERT_MSG(c.size() == size, 
"Trying to resize const data container with wrong size.");
 
  978        template <
class S, 
class T>
 
  979        constexpr bool check_resize_dimension(
const S&, 
const T&)
 
  984        template <
class T, 
size_t N, 
class S>
 
  985        constexpr bool check_resize_dimension(
const std::array<T, N>&, 
const S& s)
 
  987            return N == s.size();
 
 1003            detail::check_resize_dimension(m_shape, 
shape),
 
 1004            "cannot change the number of dimensions of xtensor" 
 1006        std::size_t dim = 
shape.size();
 
 1007        if (m_shape.size() != dim || !std::equal(std::begin(
shape), std::end(
shape), std::begin(m_shape))
 
 1012                m_layout = XTENSOR_DEFAULT_LAYOUT;  
 
 1014            m_shape = xtl::forward_sequence<shape_type, S>(
shape);
 
 1016            resize_container(m_strides, dim);
 
 1017            resize_container(m_backstrides, dim);
 
 1019            detail::resize_data_container(this->
storage(), data_size);
 
 
 1035            detail::check_resize_dimension(m_shape, 
shape),
 
 1036            "cannot change the number of dimensions of xtensor" 
 1042                "Cannot change layout_type if template parameter not layout_type::dynamic." 
 
 1061            detail::check_resize_dimension(m_shape, 
shape),
 
 1062            "cannot change the number of dimensions of xtensor" 
 1068                "Cannot resize with custom strides when layout() is != layout_type::dynamic." 
 1071        m_shape = xtl::forward_sequence<shape_type, S>(
shape);
 
 1073        resize_container(m_backstrides, m_strides.size());
 
 1074        adapt_strides(m_shape, m_strides, m_backstrides);
 
 1076        detail::resize_data_container(this->
storage(), compute_size(m_shape));
 
 
 1097            std::forward<S>(
shape),
 
 1098            xtl::is_signed<std::decay_t<
typename std::decay_t<S>::value_type>>(),
 
 1099            std::forward<layout_type>(
layout)
 
 1101        return this->derived_cast();
 
 
 1108        using sh_type = rebind_container_t<T, shape_type>;
 
 1109        sh_type sh = xtl::make_sequence<sh_type>(shape.size());
 
 1110        std::copy(shape.begin(), shape.end(), sh.begin());
 
 1111        reshape_impl(std::move(sh), xtl::is_signed<T>(), std::forward<layout_type>(layout));
 
 1112        return this->derived_cast();
 
 1118    xstrided_container<D>::reshape_impl(S&& 
shape, std::false_type , 
layout_type layout)
 
 1120        if (compute_size(
shape) != this->
size())
 
 1124                "Cannot reshape with incorrect number of elements. Do you mean to resize?" 
 1129            layout = XTENSOR_DEFAULT_LAYOUT;  
 
 1133            XTENSOR_THROW(std::runtime_error, 
"Cannot reshape with different layout if static layout != dynamic.");
 
 1136        m_shape = xtl::forward_sequence<shape_type, S>(
shape);
 
 1137        resize_container(m_strides, m_shape.size());
 
 1138        resize_container(m_backstrides, m_shape.size());
 
 1145    xstrided_container<D>::reshape_impl(S&& _shape, std::true_type , 
layout_type layout)
 
 1147        using tmp_value_type = 
typename std::decay_t<S>::value_type;
 
 1148        auto new_size = compute_size(_shape);
 
 1149        if (this->
size() % new_size)
 
 1151            XTENSOR_THROW(std::runtime_error, 
"Negative axis size cannot be inferred. Shape mismatch.");
 
 1153        std::decay_t<S> 
shape = _shape;
 
 1154        tmp_value_type accumulator = 1;
 
 1155        std::size_t neg_idx = 0;
 
 1157        for (
auto it = 
shape.begin(); it != 
shape.end(); ++it, i++)
 
 1162                XTENSOR_ASSERT(dim == -1 && !neg_idx);
 
 1167        if (accumulator < 0)
 
 1169            shape[neg_idx] = 
static_cast<tmp_value_type
>(this->
size()) / std::abs(accumulator);
 
 1171        else if (this->
size() != new_size)
 
 1175                "Cannot reshape with incorrect number of elements. Do you mean to resize?" 
 1179        m_shape = xtl::forward_sequence<shape_type, S>(
shape);
 
 1180        resize_container(m_strides, m_shape.size());
 
 1181        resize_container(m_backstrides, m_shape.size());
 
 1186    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