xtensor
Loading...
Searching...
No Matches
xcontainer.hpp
1/***************************************************************************
2 * Copyright (c) Johan Mabille, Sylvain Corlay and Wolf Vollprecht *
3 * Copyright (c) QuantStack *
4 * *
5 * Distributed under the terms of the BSD 3-Clause License. *
6 * *
7 * The full license is in the file LICENSE, distributed with this software. *
8 ****************************************************************************/
9
10#ifndef XTENSOR_CONTAINER_HPP
11#define XTENSOR_CONTAINER_HPP
12
13#include <algorithm>
14#include <functional>
15#include <memory>
16#include <numeric>
17#include <stdexcept>
18
19#include <xtl/xmeta_utils.hpp>
20#include <xtl/xsequence.hpp>
21
22#include "xaccessible.hpp"
23#include "xiterable.hpp"
24#include "xiterator.hpp"
25#include "xmath.hpp"
26#include "xoperation.hpp"
27#include "xstrides.hpp"
28#include "xtensor_config.hpp"
29#include "xtensor_forward.hpp"
30
31namespace xt
32{
33 template <class D>
35 {
36 using inner_shape_type = typename xcontainer_inner_types<D>::inner_shape_type;
37 using stepper = xstepper<D>;
39 };
40
41 namespace detail
42 {
43 template <class T>
44 struct allocator_type_impl
45 {
46 using type = typename T::allocator_type;
47 };
48
49 template <class T, std::size_t N>
50 struct allocator_type_impl<std::array<T, N>>
51 {
52 using type = std::allocator<T>; // fake allocator for testing
53 };
54 }
55
56 template <class T>
57 using allocator_type_t = typename detail::allocator_type_impl<T>::type;
58
70 template <class D>
72 private xaccessible<D>
73 {
74 public:
75
76 using derived_type = D;
77
79 using storage_type = typename inner_types::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>;
90
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;
94
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;
98
100 using stepper = typename iterable_base::stepper;
101 using const_stepper = typename iterable_base::const_stepper;
102
104
105 static constexpr layout_type static_layout = inner_types::layout;
106 static constexpr bool contiguous_layout = static_layout != layout_type::dynamic;
107 using data_alignment = xt_simd::container_alignment_t<storage_type>;
108 using simd_type = xt_simd::simd_type<value_type>;
109
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;
114
115 static_assert(static_layout != layout_type::any, "Container layout can never be layout_type::any!");
116
117 size_type size() const noexcept;
118
119 XTENSOR_CONSTEXPR_RETURN size_type dimension() const noexcept;
120
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;
124
125 template <class T>
126 void fill(const T& value);
127
128 template <class... Args>
129 reference operator()(Args... args);
130
131 template <class... Args>
132 const_reference operator()(Args... args) const;
133
134 template <class... Args>
135 reference unchecked(Args... args);
136
137 template <class... Args>
138 const_reference unchecked(Args... args) const;
139
146 using accessible_base::periodic;
147
149 reference element(It first, It last);
151 const_reference element(It first, It last) const;
152
153 storage_type& storage() noexcept;
154 const storage_type& storage() const noexcept;
155
156 pointer data() noexcept;
157 const_pointer data() const noexcept;
159
162
166 stepper stepper_begin(const S& shape) noexcept;
168 stepper stepper_end(const S& shape, layout_type l) noexcept;
169
171 const_stepper stepper_begin(const S& shape) const noexcept;
173 const_stepper stepper_end(const S& shape, layout_type l) const noexcept;
174
175 reference data_element(size_type i);
176 const_reference data_element(size_type i) const;
177
178 reference flat(size_type i);
179 const_reference flat(size_type i) const;
180
182 using simd_return_type = xt_simd::simd_return_type<value_type, requested_type>;
183
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 /*simd_return_type<requested_type>*/ load_simd(size_type i) const;
189
190 linear_iterator linear_begin() noexcept;
191 linear_iterator linear_end() noexcept;
192
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;
197
198 reverse_linear_iterator linear_rbegin() noexcept;
199 reverse_linear_iterator linear_rend() noexcept;
200
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;
205
206 using container_iterator = linear_iterator;
207 using const_container_iterator = const_linear_iterator;
208
209 protected:
210
212 ~xcontainer() = default;
213
216
219
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;
224
225 protected:
226
227 derived_type& derived_cast() & noexcept;
228 const derived_type& derived_cast() const& noexcept;
229 derived_type derived_cast() && noexcept;
230
231 private:
232
234 It data_xend_impl(It begin, layout_type l, size_type offset) const noexcept;
235
236 inner_shape_type& mutable_shape();
237 inner_strides_type& mutable_strides();
238 inner_backstrides_type& mutable_backstrides();
239
242
245 };
246
261 {
262 public:
263
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;
277
278 template <class S = shape_type>
279 void resize(S&& shape, bool force = false);
280 template <class S = shape_type>
282 template <class S = shape_type>
283 void resize(S&& shape, const strides_type& strides);
284
285 template <class S = shape_type>
286 auto& reshape(S&& shape, layout_type layout = base_type::static_layout) &;
287
288 template <class T>
289 auto& reshape(std::initializer_list<T> shape, layout_type layout = base_type::static_layout) &;
290
292 bool is_contiguous() const noexcept;
293
294 protected:
295
298
301
304
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;
307
308 inner_shape_type& shape_impl() noexcept;
309 const inner_shape_type& shape_impl() const noexcept;
310
311 inner_strides_type& strides_impl() noexcept;
312 const inner_strides_type& strides_impl() const noexcept;
313
314 inner_backstrides_type& backstrides_impl() noexcept;
315 const inner_backstrides_type& backstrides_impl() const noexcept;
316
317 template <class S = shape_type>
318 void reshape_impl(S&& shape, std::true_type, layout_type layout = base_type::static_layout);
319 template <class S = shape_type>
320 void reshape_impl(S&& shape, std::false_type, layout_type layout = base_type::static_layout);
321
322 layout_type& mutable_layout() noexcept;
323
324 private:
325
326 inner_shape_type m_shape;
327 inner_strides_type m_strides;
328 inner_backstrides_type m_backstrides;
329 layout_type m_layout = base_type::static_layout;
330 };
331
332 /******************************
333 * xcontainer implementation *
334 ******************************/
335
338 inline It xcontainer<D>::data_xend_impl(It begin, layout_type l, size_type offset) const noexcept
339 {
340 return strided_data_end(*this, begin, l, offset);
341 }
342
343 template <class D>
344 inline auto xcontainer<D>::mutable_shape() -> inner_shape_type&
345 {
346 return derived_cast().shape_impl();
347 }
348
349 template <class D>
350 inline auto xcontainer<D>::mutable_strides() -> inner_strides_type&
351 {
352 return derived_cast().strides_impl();
353 }
354
355 template <class D>
356 inline auto xcontainer<D>::mutable_backstrides() -> inner_backstrides_type&
357 {
358 return derived_cast().backstrides_impl();
359 }
360
368 template <class D>
369 inline auto xcontainer<D>::size() const noexcept -> size_type
370 {
371 return contiguous_layout ? storage().size() : compute_size(shape());
372 }
373
377 template <class D>
378 XTENSOR_CONSTEXPR_RETURN auto xcontainer<D>::dimension() const noexcept -> size_type
379 {
380 return shape().size();
381 }
382
386 template <class D>
387 XTENSOR_CONSTEXPR_RETURN auto xcontainer<D>::shape() const noexcept -> const inner_shape_type&
388 {
389 return derived_cast().shape_impl();
390 }
391
395 template <class D>
396 XTENSOR_CONSTEXPR_RETURN auto xcontainer<D>::strides() const noexcept -> const inner_strides_type&
397 {
398 return derived_cast().strides_impl();
399 }
400
404 template <class D>
405 XTENSOR_CONSTEXPR_RETURN auto xcontainer<D>::backstrides() const noexcept -> const inner_backstrides_type&
406 {
407 return derived_cast().backstrides_impl();
408 }
409
411
416
421 template <class D>
422 template <class T>
423 inline void xcontainer<D>::fill(const T& value)
424 {
425 if (contiguous_layout)
426 {
427 std::fill(this->linear_begin(), this->linear_end(), value);
428 }
429 else
430 {
431 std::fill(this->begin(), this->end(), value);
432 }
433 }
434
441 template <class D>
442 template <class... Args>
443 inline auto xcontainer<D>::operator()(Args... args) -> reference
444 {
445 XTENSOR_TRY(check_index(shape(), args...));
446 XTENSOR_CHECK_DIMENSION(shape(), args...);
447 size_type index = xt::data_offset<size_type>(strides(), args...);
448 return storage()[index];
449 }
450
457 template <class D>
458 template <class... Args>
459 inline auto xcontainer<D>::operator()(Args... args) const -> const_reference
460 {
461 XTENSOR_TRY(check_index(shape(), args...));
462 XTENSOR_CHECK_DIMENSION(shape(), args...);
463 size_type index = xt::data_offset<size_type>(strides(), args...);
464 return storage()[index];
465 }
466
486 template <class D>
487 template <class... Args>
488 inline auto xcontainer<D>::unchecked(Args... args) -> reference
489 {
491 strides(),
492 static_cast<std::ptrdiff_t>(args)...
493 );
494 return storage()[index];
495 }
496
516 template <class D>
517 template <class... Args>
518 inline auto xcontainer<D>::unchecked(Args... args) const -> const_reference
519 {
521 strides(),
522 static_cast<std::ptrdiff_t>(args)...
523 );
524 return storage()[index];
525 }
526
534 template <class D>
535 template <class It>
536 inline auto xcontainer<D>::element(It first, It last) -> reference
537 {
538 XTENSOR_TRY(check_element_index(shape(), first, last));
540 }
541
549 template <class D>
550 template <class It>
551 inline auto xcontainer<D>::element(It first, It last) const -> const_reference
552 {
553 XTENSOR_TRY(check_element_index(shape(), first, last));
555 }
556
560 template <class D>
561 inline auto xcontainer<D>::storage() noexcept -> storage_type&
562 {
563 return derived_cast().storage_impl();
564 }
565
570 template <class D>
571 inline auto xcontainer<D>::storage() const noexcept -> const storage_type&
572 {
573 return derived_cast().storage_impl();
574 }
575
581 template <class D>
582 inline auto xcontainer<D>::data() noexcept -> pointer
583 {
584 return storage().data();
585 }
586
592 template <class D>
593 inline auto xcontainer<D>::data() const noexcept -> const_pointer
594 {
595 return storage().data();
596 }
597
601 template <class D>
603 {
604 return size_type(0);
605 }
606
608
619 template <class D>
620 template <class S>
621 inline bool xcontainer<D>::broadcast_shape(S& shape, bool) const
622 {
623 return xt::broadcast_shape(this->shape(), shape);
624 }
625
631 template <class D>
632 template <class S>
633 inline bool xcontainer<D>::has_linear_assign(const S& str) const noexcept
634 {
635 return str.size() == strides().size() && std::equal(str.cbegin(), str.cend(), strides().begin());
636 }
637
639
640 template <class D>
641 inline auto xcontainer<D>::derived_cast() const& noexcept -> const derived_type&
642 {
643 return *static_cast<const derived_type*>(this);
644 }
645
646 template <class D>
647 inline auto xcontainer<D>::derived_cast() && noexcept -> derived_type
648 {
649 return *static_cast<derived_type*>(this);
650 }
651
652 template <class D>
653 inline auto xcontainer<D>::data_element(size_type i) -> reference
654 {
655 return storage()[i];
656 }
657
658 template <class D>
659 inline auto xcontainer<D>::data_element(size_type i) const -> const_reference
660 {
661 return storage()[i];
662 }
663
670 template <class D>
671 inline auto xcontainer<D>::flat(size_type i) -> reference
672 {
673 XTENSOR_ASSERT(i < size());
674 return storage()[i];
675 }
676
683 template <class D>
684 inline auto xcontainer<D>::flat(size_type i) const -> const_reference
685 {
686 XTENSOR_ASSERT(i < size());
687 return storage()[i];
688 }
689
690 /***************
691 * stepper api *
692 ***************/
693
694 template <class D>
695 template <class S>
696 inline auto xcontainer<D>::stepper_begin(const S& shape) noexcept -> stepper
697 {
698 size_type offset = shape.size() - dimension();
699 return stepper(static_cast<derived_type*>(this), data_xbegin(), offset);
700 }
701
702 template <class D>
703 template <class S>
704 inline auto xcontainer<D>::stepper_end(const S& shape, layout_type l) noexcept -> stepper
705 {
706 size_type offset = shape.size() - dimension();
707 return stepper(static_cast<derived_type*>(this), data_xend(l, offset), offset);
708 }
709
710 template <class D>
711 template <class S>
712 inline auto xcontainer<D>::stepper_begin(const S& shape) const noexcept -> const_stepper
713 {
714 size_type offset = shape.size() - dimension();
715 return const_stepper(static_cast<const derived_type*>(this), data_xbegin(), offset);
716 }
717
718 template <class D>
719 template <class S>
720 inline auto xcontainer<D>::stepper_end(const S& shape, layout_type l) const noexcept -> const_stepper
721 {
722 size_type offset = shape.size() - dimension();
723 return const_stepper(static_cast<const derived_type*>(this), data_xend(l, offset), offset);
724 }
725
726 template <class D>
727 inline auto xcontainer<D>::data_xbegin() noexcept -> container_iterator
728 {
729 return storage().begin();
730 }
731
732 template <class D>
733 inline auto xcontainer<D>::data_xbegin() const noexcept -> const_container_iterator
734 {
735 return storage().cbegin();
736 }
737
738 template <class D>
739 inline auto xcontainer<D>::data_xend(layout_type l, size_type offset) noexcept -> container_iterator
740 {
741 return data_xend_impl(storage().begin(), l, offset);
742 }
743
744 template <class D>
745 inline auto xcontainer<D>::data_xend(layout_type l, size_type offset) const noexcept
746 -> const_container_iterator
747 {
748 return data_xend_impl(storage().cbegin(), l, offset);
749 }
750
751 template <class D>
752 template <class alignment, class simd>
753 inline void xcontainer<D>::store_simd(size_type i, const simd& e)
754 {
755 using align_mode = driven_align_mode_t<alignment, data_alignment>;
756 xt_simd::store_as(std::addressof(storage()[i]), e, align_mode());
757 }
758
759 template <class D>
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>
763 {
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());
766 }
767
768 template <class D>
769 inline auto xcontainer<D>::linear_begin() noexcept -> linear_iterator
770 {
771 return storage().begin();
772 }
773
774 template <class D>
775 inline auto xcontainer<D>::linear_end() noexcept -> linear_iterator
776 {
777 return storage().end();
778 }
779
780 template <class D>
781 inline auto xcontainer<D>::linear_begin() const noexcept -> const_linear_iterator
782 {
783 return storage().begin();
784 }
785
786 template <class D>
787 inline auto xcontainer<D>::linear_end() const noexcept -> const_linear_iterator
788 {
789 return storage().cend();
790 }
791
792 template <class D>
793 inline auto xcontainer<D>::linear_cbegin() const noexcept -> const_linear_iterator
794 {
795 return storage().cbegin();
796 }
797
798 template <class D>
799 inline auto xcontainer<D>::linear_cend() const noexcept -> const_linear_iterator
800 {
801 return storage().cend();
802 }
803
804 template <class D>
805 inline auto xcontainer<D>::linear_rbegin() noexcept -> reverse_linear_iterator
806 {
807 return storage().rbegin();
808 }
809
810 template <class D>
811 inline auto xcontainer<D>::linear_rend() noexcept -> reverse_linear_iterator
812 {
813 return storage().rend();
814 }
815
816 template <class D>
817 inline auto xcontainer<D>::linear_rbegin() const noexcept -> const_reverse_linear_iterator
818 {
819 return storage().rbegin();
820 }
821
822 template <class D>
823 inline auto xcontainer<D>::linear_rend() const noexcept -> const_reverse_linear_iterator
824 {
825 return storage().rend();
826 }
827
828 template <class D>
829 inline auto xcontainer<D>::linear_crbegin() const noexcept -> const_reverse_linear_iterator
830 {
831 return storage().crbegin();
832 }
833
834 template <class D>
835 inline auto xcontainer<D>::linear_crend() const noexcept -> const_reverse_linear_iterator
836 {
837 return storage().crend();
838 }
839
840 template <class D>
841 inline auto xcontainer<D>::derived_cast() & noexcept -> derived_type&
842 {
843 return *static_cast<derived_type*>(this);
844 }
845
846 /*************************************
847 * xstrided_container implementation *
848 *************************************/
849
850 template <class D>
851 inline xstrided_container<D>::xstrided_container() noexcept
852 : base_type()
853 {
854 m_shape = xtl::make_sequence<inner_shape_type>(base_type::dimension(), 0);
855 m_strides = xtl::make_sequence<inner_strides_type>(base_type::dimension(), 0);
856 m_backstrides = xtl::make_sequence<inner_backstrides_type>(base_type::dimension(), 0);
857 }
858
859 template <class D>
860 inline xstrided_container<D>::xstrided_container(inner_shape_type&& shape, inner_strides_type&& strides) noexcept
861 : base_type()
862 , m_shape(std::move(shape))
863 , m_strides(std::move(strides))
864 {
865 m_backstrides = xtl::make_sequence<inner_backstrides_type>(m_shape.size(), 0);
866 adapt_strides(m_shape, m_strides, m_backstrides);
867 }
868
869 template <class D>
870 inline xstrided_container<D>::xstrided_container(
871 inner_shape_type&& shape,
872 inner_strides_type&& strides,
873 inner_backstrides_type&& backstrides,
874 layout_type&& layout
875 ) noexcept
876 : base_type()
877 , m_shape(std::move(shape))
878 , m_strides(std::move(strides))
879 , m_backstrides(std::move(backstrides))
880 , m_layout(std::move(layout))
881 {
882 }
883
884 template <class D>
885 inline auto xstrided_container<D>::shape_impl() noexcept -> inner_shape_type&
886 {
887 return m_shape;
888 }
889
890 template <class D>
891 inline auto xstrided_container<D>::shape_impl() const noexcept -> const inner_shape_type&
892 {
893 return m_shape;
894 }
895
896 template <class D>
897 inline auto xstrided_container<D>::strides_impl() noexcept -> inner_strides_type&
898 {
899 return m_strides;
900 }
901
902 template <class D>
903 inline auto xstrided_container<D>::strides_impl() const noexcept -> const inner_strides_type&
904 {
905 return m_strides;
906 }
907
908 template <class D>
909 inline auto xstrided_container<D>::backstrides_impl() noexcept -> inner_backstrides_type&
910 {
911 return m_backstrides;
912 }
913
914 template <class D>
915 inline auto xstrided_container<D>::backstrides_impl() const noexcept -> const inner_backstrides_type&
916 {
917 return m_backstrides;
918 }
919
924 template <class D>
926 {
927 return m_layout;
928 }
929
930 template <class D>
932 {
933 using str_type = typename inner_strides_type::value_type;
934 auto is_zero = [](auto i)
935 {
936 return i == 0;
937 };
938 if (!is_contiguous_container<storage_type>::value)
939 {
940 return false;
941 }
942 // We need to make sure the inner-most non-zero stride is one.
943 // Trailing zero strides are ignored because they indicate bradcasted dimensions.
944 if (m_layout == layout_type::row_major)
945 {
946 auto it = std::find_if_not(m_strides.rbegin(), m_strides.rend(), is_zero);
947 // If the array has strides of zero, it is a constant, and therefore contiguous.
948 return it == m_strides.rend() || *it == str_type(1);
949 }
950 else if (m_layout == layout_type::column_major)
951 {
952 auto it = std::find_if_not(m_strides.begin(), m_strides.end(), is_zero);
953 // If the array has strides of zero, it is a constant, and therefore contiguous.
954 return it == m_strides.end() || *it == str_type(1);
955 }
956 else
957 {
958 return m_strides.empty();
959 }
960 }
961
962 namespace detail
963 {
964 template <class C, class S>
965 inline void resize_data_container(C& c, S size)
966 {
967 xt::resize_container(c, size);
968 }
969
970 template <class C, class S>
971 inline void resize_data_container(const C& c, S size)
972 {
973 (void) c; // remove unused parameter warning
974 (void) size;
975 XTENSOR_ASSERT_MSG(c.size() == size, "Trying to resize const data container with wrong size.");
976 }
977
978 template <class S, class T>
979 constexpr bool check_resize_dimension(const S&, const T&)
980 {
981 return true;
982 }
983
984 template <class T, size_t N, class S>
985 constexpr bool check_resize_dimension(const std::array<T, N>&, const S& s)
986 {
987 return N == s.size();
988 }
989 }
990
998 template <class D>
999 template <class S>
1001 {
1002 XTENSOR_ASSERT_MSG(
1003 detail::check_resize_dimension(m_shape, shape),
1004 "cannot change the number of dimensions of xtensor"
1005 )
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))
1008 || force)
1009 {
1010 if (D::static_layout == layout_type::dynamic && m_layout == layout_type::dynamic)
1011 {
1012 m_layout = XTENSOR_DEFAULT_LAYOUT; // fall back to default layout
1013 }
1014 m_shape = xtl::forward_sequence<shape_type, S>(shape);
1015
1016 resize_container(m_strides, dim);
1017 resize_container(m_backstrides, dim);
1018 size_type data_size = compute_strides<D::static_layout>(m_shape, m_layout, m_strides, m_backstrides);
1019 detail::resize_data_container(this->storage(), data_size);
1020 }
1021 }
1022
1030 template <class D>
1031 template <class S>
1033 {
1034 XTENSOR_ASSERT_MSG(
1035 detail::check_resize_dimension(m_shape, shape),
1036 "cannot change the number of dimensions of xtensor"
1037 )
1038 if (base_type::static_layout != layout_type::dynamic && l != base_type::static_layout)
1039 {
1040 XTENSOR_THROW(
1041 std::runtime_error,
1042 "Cannot change layout_type if template parameter not layout_type::dynamic."
1043 );
1044 }
1045 m_layout = l;
1046 resize(std::forward<S>(shape), true);
1047 }
1048
1056 template <class D>
1057 template <class S>
1058 inline void xstrided_container<D>::resize(S&& shape, const strides_type& strides)
1059 {
1060 XTENSOR_ASSERT_MSG(
1061 detail::check_resize_dimension(m_shape, shape),
1062 "cannot change the number of dimensions of xtensor"
1063 )
1064 if (base_type::static_layout != layout_type::dynamic)
1065 {
1066 XTENSOR_THROW(
1067 std::runtime_error,
1068 "Cannot resize with custom strides when layout() is != layout_type::dynamic."
1069 );
1070 }
1071 m_shape = xtl::forward_sequence<shape_type, S>(shape);
1072 m_strides = strides;
1073 resize_container(m_backstrides, m_strides.size());
1074 adapt_strides(m_shape, m_strides, m_backstrides);
1075 m_layout = layout_type::dynamic;
1076 detail::resize_data_container(this->storage(), compute_size(m_shape));
1077 }
1078
1092 template <class D>
1093 template <class S>
1095 {
1096 reshape_impl(
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)
1100 );
1101 return this->derived_cast();
1102 }
1103
1104 template <class D>
1105 template <class T>
1106 inline auto& xstrided_container<D>::reshape(std::initializer_list<T> shape, layout_type layout) &
1107 {
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();
1113 }
1114
1115 template <class D>
1116 template <class S>
1117 inline void
1118 xstrided_container<D>::reshape_impl(S&& shape, std::false_type /* is unsigned */, layout_type layout)
1119 {
1120 if (compute_size(shape) != this->size())
1121 {
1122 XTENSOR_THROW(
1123 std::runtime_error,
1124 "Cannot reshape with incorrect number of elements. Do you mean to resize?"
1125 );
1126 }
1127 if (D::static_layout == layout_type::dynamic && layout == layout_type::dynamic)
1128 {
1129 layout = XTENSOR_DEFAULT_LAYOUT; // fall back to default layout
1130 }
1131 if (D::static_layout != layout_type::dynamic && layout != D::static_layout)
1132 {
1133 XTENSOR_THROW(std::runtime_error, "Cannot reshape with different layout if static layout != dynamic.");
1134 }
1135 m_layout = layout;
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());
1139 compute_strides<D::static_layout>(m_shape, m_layout, m_strides, m_backstrides);
1140 }
1141
1142 template <class D>
1143 template <class S>
1144 inline void
1145 xstrided_container<D>::reshape_impl(S&& _shape, std::true_type /* is signed */, layout_type layout)
1146 {
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)
1150 {
1151 XTENSOR_THROW(std::runtime_error, "Negative axis size cannot be inferred. Shape mismatch.");
1152 }
1153 std::decay_t<S> shape = _shape;
1154 tmp_value_type accumulator = 1;
1155 std::size_t neg_idx = 0;
1156 std::size_t i = 0;
1157 for (auto it = shape.begin(); it != shape.end(); ++it, i++)
1158 {
1159 auto&& dim = *it;
1160 if (dim < 0)
1161 {
1162 XTENSOR_ASSERT(dim == -1 && !neg_idx);
1163 neg_idx = i;
1164 }
1165 accumulator *= dim;
1166 }
1167 if (accumulator < 0)
1168 {
1169 shape[neg_idx] = static_cast<tmp_value_type>(this->size()) / std::abs(accumulator);
1170 }
1171 else if (this->size() != new_size)
1172 {
1173 XTENSOR_THROW(
1174 std::runtime_error,
1175 "Cannot reshape with incorrect number of elements. Do you mean to resize?"
1176 );
1177 }
1178 m_layout = layout;
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());
1182 compute_strides<D::static_layout>(m_shape, m_layout, m_strides, m_backstrides);
1183 }
1184
1185 template <class D>
1186 inline auto xstrided_container<D>::mutable_layout() noexcept -> layout_type&
1187 {
1188 return m_layout;
1189 }
1190}
1191
1192#endif
Base class for implementation of common expression access methods.
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
Returns true only if the the specified position is a valid entry in the expression.
Base class for dense multidimensional containers.
storage_type & storage() noexcept
Returns a reference to the buffer containing the elements of the container.
size_type size() const noexcept
Returns the number of element in the container.
constexpr const inner_backstrides_type & backstrides() const noexcept
Returns the backstrides of the container.
constexpr size_type dimension() const noexcept
Returns the number of dimensions of the container.
reference flat(size_type i)
Returns a reference to the element at the specified position in the container storage (as if it was o...
constexpr const inner_strides_type & strides() const noexcept
Returns the strides of the container.
bool has_linear_assign(const S &strides) const noexcept
Checks whether the xcontainer can be linearly assigned to an expression with the specified strides.
bool broadcast_shape(S &shape, bool reuse_cache=false) const
Broadcast the shape of the container to the specified parameter.
constexpr const inner_shape_type & shape() const noexcept
Returns the shape of the container.
pointer data() noexcept
Returns a pointer to the underlying array serving as element storage.
const size_type data_offset() const noexcept
Returns the offset to the first element in the container.
reference back()
Returns a reference to the last element of the expression.
reference front()
Returns a reference to the first element of the expression.
void fill(const T &value)
Fills the container with the given value.
Base class for multidimensional iterable expressions with contiguous storage.
Partial implementation of xcontainer that embeds the strides and the shape.
layout_type layout() const noexcept
Return the layout_type of the container.
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.
standard mathematical functions for xexpressions
layout_type
Definition xlayout.hpp:24