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 <memory>
15#include <stdexcept>
16
17#include <xtl/xmeta_utils.hpp>
18#include <xtl/xsequence.hpp>
19
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"
28
29namespace xt
30{
31 template <class D>
33 {
34 using inner_shape_type = typename xcontainer_inner_types<D>::inner_shape_type;
35 using stepper = xstepper<D>;
36 using const_stepper = xstepper<const D>;
37 };
38
39 namespace detail
40 {
41 template <class T>
42 struct allocator_type_impl
43 {
44 using type = typename T::allocator_type;
45 };
46
47 template <class T, std::size_t N>
48 struct allocator_type_impl<std::array<T, N>>
49 {
50 using type = std::allocator<T>; // fake allocator for testing
51 };
52 }
53
54 template <class T>
55 using allocator_type_t = typename detail::allocator_type_impl<T>::type;
56
68 template <class D>
69 class xcontainer : public xcontiguous_iterable<D>,
70 private xaccessible<D>
71 {
72 public:
73
74 using derived_type = D;
75
76 using inner_types = xcontainer_inner_types<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>;
88
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;
92
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;
96
97 using iterable_base = xcontiguous_iterable<D>;
98 using stepper = typename iterable_base::stepper;
99 using const_stepper = typename iterable_base::const_stepper;
100
101 using accessible_base = xaccessible<D>;
102
103 static constexpr layout_type static_layout = inner_types::layout;
104 static constexpr bool contiguous_layout = static_layout != layout_type::dynamic;
105 using data_alignment = xt_simd::container_alignment_t<storage_type>;
106 using simd_type = xt_simd::simd_type<value_type>;
107
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;
112
113 static_assert(static_layout != layout_type::any, "Container layout can never be layout_type::any!");
114
115 size_type size() const noexcept;
116
117 constexpr size_type dimension() const noexcept;
118
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;
122
123 template <class T>
124 void fill(const T& value);
125
126 template <class... Args>
127 reference operator()(Args... args);
128
129 template <class... Args>
130 const_reference operator()(Args... args) const;
131
132 template <class... Args>
133 reference unchecked(Args... args);
134
135 template <class... Args>
136 const_reference unchecked(Args... args) const;
137
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;
143 using accessible_base::in_bounds;
144 using accessible_base::periodic;
145
146 template <class It>
147 reference element(It first, It last);
148 template <class It>
149 const_reference element(It first, It last) const;
150
151 storage_type& storage() noexcept;
152 const storage_type& storage() const noexcept;
153
154 pointer data() noexcept;
155 const_pointer data() const noexcept;
156 const size_type data_offset() const noexcept;
157
158 template <class S>
159 bool broadcast_shape(S& shape, bool reuse_cache = false) const;
160
161 template <class S>
162 bool has_linear_assign(const S& strides) const noexcept;
163 template <class S>
164 stepper stepper_begin(const S& shape) noexcept;
165 template <class S>
166 stepper stepper_end(const S& shape, layout_type l) noexcept;
167
168 template <class S>
169 const_stepper stepper_begin(const S& shape) const noexcept;
170 template <class S>
171 const_stepper stepper_end(const S& shape, layout_type l) const noexcept;
172
173 reference data_element(size_type i);
174 const_reference data_element(size_type i) const;
175
176 reference flat(size_type i);
177 const_reference flat(size_type i) const;
178
179 template <class requested_type>
180 using simd_return_type = xt_simd::simd_return_type<value_type, requested_type>;
181
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 /*simd_return_type<requested_type>*/ load_simd(size_type i) const;
187
188 linear_iterator linear_begin() noexcept;
189 linear_iterator linear_end() noexcept;
190
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;
195
196 reverse_linear_iterator linear_rbegin() noexcept;
197 reverse_linear_iterator linear_rend() noexcept;
198
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;
203
204 using container_iterator = linear_iterator;
205 using const_container_iterator = const_linear_iterator;
206
207 protected:
208
209 xcontainer() = default;
210 ~xcontainer() = default;
211
212 xcontainer(const xcontainer&) = default;
213 xcontainer& operator=(const xcontainer&) = default;
214
215 xcontainer(xcontainer&&) = default;
216 xcontainer& operator=(xcontainer&&) = default;
217
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;
222
223 protected:
224
225 derived_type& derived_cast() & noexcept;
226 const derived_type& derived_cast() const& noexcept;
227 derived_type derived_cast() && noexcept;
228
229 private:
230
231 template <class It>
232 It data_xend_impl(It begin, layout_type l, size_type offset) const noexcept;
233
234 inner_shape_type& mutable_shape();
235 inner_strides_type& mutable_strides();
236 inner_backstrides_type& mutable_backstrides();
237
238 template <class C>
239 friend class xstepper;
240
241 friend class xaccessible<D>;
242 friend class xconst_accessible<D>;
243 };
244
257 template <class D>
258 class xstrided_container : public xcontainer<D>
259 {
260 public:
261
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;
275
276 template <class S = shape_type>
277 void resize(S&& shape, bool force = false);
278 template <class S = shape_type>
280 template <class S = shape_type>
281 void resize(S&& shape, const strides_type& strides);
282
283 template <class S = shape_type>
284 auto& reshape(S&& shape, layout_type layout = base_type::static_layout) &;
285
286 template <class T>
287 auto& reshape(std::initializer_list<T> shape, layout_type layout = base_type::static_layout) &;
288
289 layout_type layout() const noexcept;
290 bool is_contiguous() const noexcept;
291
292 protected:
293
294 xstrided_container() noexcept;
295 ~xstrided_container() = default;
296
297 xstrided_container(const xstrided_container&) = default;
298 xstrided_container& operator=(const xstrided_container&) = default;
299
300 xstrided_container(xstrided_container&&) = default;
301 xstrided_container& operator=(xstrided_container&&) = default;
302
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;
305
306 inner_shape_type& shape_impl() noexcept;
307 const inner_shape_type& shape_impl() const noexcept;
308
309 inner_strides_type& strides_impl() noexcept;
310 const inner_strides_type& strides_impl() const noexcept;
311
312 inner_backstrides_type& backstrides_impl() noexcept;
313 const inner_backstrides_type& backstrides_impl() const noexcept;
314
315 template <class S = shape_type>
316 void reshape_impl(S&& shape, std::true_type, layout_type layout = base_type::static_layout);
317 template <class S = shape_type>
318 void reshape_impl(S&& shape, std::false_type, layout_type layout = base_type::static_layout);
319
320 layout_type& mutable_layout() noexcept;
321
322 private:
323
324 inner_shape_type m_shape;
325 inner_strides_type m_strides;
326 inner_backstrides_type m_backstrides;
327 layout_type m_layout = base_type::static_layout;
328 };
329
330 /******************************
331 * xcontainer implementation *
332 ******************************/
333
334 template <class D>
335 template <class It>
336 inline It xcontainer<D>::data_xend_impl(It begin, layout_type l, size_type offset) const noexcept
337 {
338 return strided_data_end(*this, begin, l, offset);
339 }
340
341 template <class D>
342 inline auto xcontainer<D>::mutable_shape() -> inner_shape_type&
343 {
344 return derived_cast().shape_impl();
345 }
346
347 template <class D>
348 inline auto xcontainer<D>::mutable_strides() -> inner_strides_type&
349 {
350 return derived_cast().strides_impl();
351 }
352
353 template <class D>
354 inline auto xcontainer<D>::mutable_backstrides() -> inner_backstrides_type&
355 {
356 return derived_cast().backstrides_impl();
357 }
358
363
366 template <class D>
367 inline auto xcontainer<D>::size() const noexcept -> size_type
368 {
369 return contiguous_layout ? storage().size() : compute_size(shape());
370 }
371
375 template <class D>
376 constexpr auto xcontainer<D>::dimension() const noexcept -> size_type
377 {
378 return shape().size();
379 }
380
384 template <class D>
385 constexpr auto xcontainer<D>::shape() const noexcept -> const inner_shape_type&
386 {
387 return derived_cast().shape_impl();
388 }
389
393 template <class D>
394 constexpr auto xcontainer<D>::strides() const noexcept -> const inner_strides_type&
395 {
396 return derived_cast().strides_impl();
397 }
398
402 template <class D>
403 constexpr auto xcontainer<D>::backstrides() const noexcept -> const inner_backstrides_type&
404 {
405 return derived_cast().backstrides_impl();
406 }
407
409
414
419 template <class D>
420 template <class T>
421 inline void xcontainer<D>::fill(const T& value)
422 {
423 if (contiguous_layout)
424 {
425 std::fill(this->linear_begin(), this->linear_end(), value);
426 }
427 else
428 {
429 std::fill(this->begin(), this->end(), value);
430 }
431 }
432
439 template <class D>
440 template <class... Args>
441 inline auto xcontainer<D>::operator()(Args... args) -> reference
442 {
443 XTENSOR_TRY(check_index(shape(), args...));
444 XTENSOR_CHECK_DIMENSION(shape(), args...);
445 size_type index = xt::data_offset<size_type>(strides(), args...);
446 return storage()[index];
447 }
448
455 template <class D>
456 template <class... Args>
457 inline auto xcontainer<D>::operator()(Args... args) const -> const_reference
458 {
459 XTENSOR_TRY(check_index(shape(), args...));
460 XTENSOR_CHECK_DIMENSION(shape(), args...);
461 size_type index = xt::data_offset<size_type>(strides(), args...);
462 return storage()[index];
463 }
464
484 template <class D>
485 template <class... Args>
486 inline auto xcontainer<D>::unchecked(Args... args) -> reference
487 {
488 size_type index = xt::unchecked_data_offset<size_type, static_layout>(
489 strides(),
490 static_cast<std::ptrdiff_t>(args)...
491 );
492 return storage()[index];
493 }
494
514 template <class D>
515 template <class... Args>
516 inline auto xcontainer<D>::unchecked(Args... args) const -> const_reference
517 {
518 size_type index = xt::unchecked_data_offset<size_type, static_layout>(
519 strides(),
520 static_cast<std::ptrdiff_t>(args)...
521 );
522 return storage()[index];
523 }
524
532 template <class D>
533 template <class It>
534 inline auto xcontainer<D>::element(It first, It last) -> reference
535 {
536 XTENSOR_TRY(check_element_index(shape(), first, last));
537 return storage()[element_offset<size_type>(strides(), first, last)];
538 }
539
547 template <class D>
548 template <class It>
549 inline auto xcontainer<D>::element(It first, It last) const -> const_reference
550 {
551 XTENSOR_TRY(check_element_index(shape(), first, last));
552 return storage()[element_offset<size_type>(strides(), first, last)];
553 }
554
558 template <class D>
559 inline auto xcontainer<D>::storage() noexcept -> storage_type&
560 {
561 return derived_cast().storage_impl();
562 }
563
568 template <class D>
569 inline auto xcontainer<D>::storage() const noexcept -> const storage_type&
570 {
571 return derived_cast().storage_impl();
572 }
573
579 template <class D>
580 inline auto xcontainer<D>::data() noexcept -> pointer
581 {
582 return storage().data();
583 }
584
590 template <class D>
591 inline auto xcontainer<D>::data() const noexcept -> const_pointer
592 {
593 return storage().data();
594 }
595
599 template <class D>
600 inline auto xcontainer<D>::data_offset() const noexcept -> const size_type
601 {
602 return size_type(0);
603 }
604
606
611
617 template <class D>
618 template <class S>
619 inline bool xcontainer<D>::broadcast_shape(S& shape, bool) const
620 {
621 return xt::broadcast_shape(this->shape(), shape);
622 }
623
629 template <class D>
630 template <class S>
631 inline bool xcontainer<D>::has_linear_assign(const S& str) const noexcept
632 {
633 return str.size() == strides().size() && std::equal(str.cbegin(), str.cend(), strides().begin());
634 }
635
637
638 template <class D>
639 inline auto xcontainer<D>::derived_cast() const& noexcept -> const derived_type&
640 {
641 return *static_cast<const derived_type*>(this);
642 }
643
644 template <class D>
645 inline auto xcontainer<D>::derived_cast() && noexcept -> derived_type
646 {
647 return *static_cast<derived_type*>(this);
648 }
649
650 template <class D>
651 inline auto xcontainer<D>::data_element(size_type i) -> reference
652 {
653 return storage()[i];
654 }
655
656 template <class D>
657 inline auto xcontainer<D>::data_element(size_type i) const -> const_reference
658 {
659 return storage()[i];
660 }
661
668 template <class D>
669 inline auto xcontainer<D>::flat(size_type i) -> reference
670 {
671 XTENSOR_ASSERT(i < size());
672 return storage()[i];
673 }
674
681 template <class D>
682 inline auto xcontainer<D>::flat(size_type i) const -> const_reference
683 {
684 XTENSOR_ASSERT(i < size());
685 return storage()[i];
686 }
687
688 /***************
689 * stepper api *
690 ***************/
691
692 template <class D>
693 template <class S>
694 inline auto xcontainer<D>::stepper_begin(const S& shape) noexcept -> stepper
695 {
696 size_type offset = shape.size() - dimension();
697 return stepper(static_cast<derived_type*>(this), data_xbegin(), offset);
698 }
699
700 template <class D>
701 template <class S>
702 inline auto xcontainer<D>::stepper_end(const S& shape, layout_type l) noexcept -> stepper
703 {
704 size_type offset = shape.size() - dimension();
705 return stepper(static_cast<derived_type*>(this), data_xend(l, offset), offset);
706 }
707
708 template <class D>
709 template <class S>
710 inline auto xcontainer<D>::stepper_begin(const S& shape) const noexcept -> const_stepper
711 {
712 size_type offset = shape.size() - dimension();
713 return const_stepper(static_cast<const derived_type*>(this), data_xbegin(), offset);
714 }
715
716 template <class D>
717 template <class S>
718 inline auto xcontainer<D>::stepper_end(const S& shape, layout_type l) const noexcept -> const_stepper
719 {
720 size_type offset = shape.size() - dimension();
721 return const_stepper(static_cast<const derived_type*>(this), data_xend(l, offset), offset);
722 }
723
724 template <class D>
725 inline auto xcontainer<D>::data_xbegin() noexcept -> container_iterator
726 {
727 return storage().begin();
728 }
729
730 template <class D>
731 inline auto xcontainer<D>::data_xbegin() const noexcept -> const_container_iterator
732 {
733 return storage().cbegin();
734 }
735
736 template <class D>
737 inline auto xcontainer<D>::data_xend(layout_type l, size_type offset) noexcept -> container_iterator
738 {
739 return data_xend_impl(storage().begin(), l, offset);
740 }
741
742 template <class D>
743 inline auto xcontainer<D>::data_xend(layout_type l, size_type offset) const noexcept
744 -> const_container_iterator
745 {
746 return data_xend_impl(storage().cbegin(), l, offset);
747 }
748
749 template <class D>
750 template <class alignment, class simd>
751 inline void xcontainer<D>::store_simd(size_type i, const simd& e)
752 {
753 using align_mode = driven_align_mode_t<alignment, data_alignment>;
754 xt_simd::store_as(std::addressof(storage()[i]), e, align_mode());
755 }
756
757 template <class D>
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>
761 {
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());
764 }
765
766 template <class D>
767 inline auto xcontainer<D>::linear_begin() noexcept -> linear_iterator
768 {
769 return storage().begin();
770 }
771
772 template <class D>
773 inline auto xcontainer<D>::linear_end() noexcept -> linear_iterator
774 {
775 return storage().end();
776 }
777
778 template <class D>
779 inline auto xcontainer<D>::linear_begin() const noexcept -> const_linear_iterator
780 {
781 return storage().begin();
782 }
783
784 template <class D>
785 inline auto xcontainer<D>::linear_end() const noexcept -> const_linear_iterator
786 {
787 return storage().cend();
788 }
789
790 template <class D>
791 inline auto xcontainer<D>::linear_cbegin() const noexcept -> const_linear_iterator
792 {
793 return storage().cbegin();
794 }
795
796 template <class D>
797 inline auto xcontainer<D>::linear_cend() const noexcept -> const_linear_iterator
798 {
799 return storage().cend();
800 }
801
802 template <class D>
803 inline auto xcontainer<D>::linear_rbegin() noexcept -> reverse_linear_iterator
804 {
805 return storage().rbegin();
806 }
807
808 template <class D>
809 inline auto xcontainer<D>::linear_rend() noexcept -> reverse_linear_iterator
810 {
811 return storage().rend();
812 }
813
814 template <class D>
815 inline auto xcontainer<D>::linear_rbegin() const noexcept -> const_reverse_linear_iterator
816 {
817 return storage().rbegin();
818 }
819
820 template <class D>
821 inline auto xcontainer<D>::linear_rend() const noexcept -> const_reverse_linear_iterator
822 {
823 return storage().rend();
824 }
825
826 template <class D>
827 inline auto xcontainer<D>::linear_crbegin() const noexcept -> const_reverse_linear_iterator
828 {
829 return storage().crbegin();
830 }
831
832 template <class D>
833 inline auto xcontainer<D>::linear_crend() const noexcept -> const_reverse_linear_iterator
834 {
835 return storage().crend();
836 }
837
838 template <class D>
839 inline auto xcontainer<D>::derived_cast() & noexcept -> derived_type&
840 {
841 return *static_cast<derived_type*>(this);
842 }
843
844 /*************************************
845 * xstrided_container implementation *
846 *************************************/
847
848 template <class D>
849 inline xstrided_container<D>::xstrided_container() noexcept
850 : base_type()
851 {
852 m_shape = xtl::make_sequence<inner_shape_type>(base_type::dimension(), 0);
853 m_strides = xtl::make_sequence<inner_strides_type>(base_type::dimension(), 0);
854 m_backstrides = xtl::make_sequence<inner_backstrides_type>(base_type::dimension(), 0);
855 }
856
857 template <class D>
858 inline xstrided_container<D>::xstrided_container(inner_shape_type&& shape, inner_strides_type&& strides) noexcept
859 : base_type()
860 , m_shape(std::move(shape))
861 , m_strides(std::move(strides))
862 {
863 m_backstrides = xtl::make_sequence<inner_backstrides_type>(m_shape.size(), 0);
864 adapt_strides(m_shape, m_strides, m_backstrides);
865 }
866
867 template <class D>
868 inline xstrided_container<D>::xstrided_container(
869 inner_shape_type&& shape,
870 inner_strides_type&& strides,
871 inner_backstrides_type&& backstrides,
872 layout_type&& layout
873 ) noexcept
874 : base_type()
875 , m_shape(std::move(shape))
876 , m_strides(std::move(strides))
877 , m_backstrides(std::move(backstrides))
878 , m_layout(std::move(layout))
879 {
880 }
881
882 template <class D>
883 inline auto xstrided_container<D>::shape_impl() noexcept -> inner_shape_type&
884 {
885 return m_shape;
886 }
887
888 template <class D>
889 inline auto xstrided_container<D>::shape_impl() const noexcept -> const inner_shape_type&
890 {
891 return m_shape;
892 }
893
894 template <class D>
895 inline auto xstrided_container<D>::strides_impl() noexcept -> inner_strides_type&
896 {
897 return m_strides;
898 }
899
900 template <class D>
901 inline auto xstrided_container<D>::strides_impl() const noexcept -> const inner_strides_type&
902 {
903 return m_strides;
904 }
905
906 template <class D>
907 inline auto xstrided_container<D>::backstrides_impl() noexcept -> inner_backstrides_type&
908 {
909 return m_backstrides;
910 }
911
912 template <class D>
913 inline auto xstrided_container<D>::backstrides_impl() const noexcept -> const inner_backstrides_type&
914 {
915 return m_backstrides;
916 }
917
922 template <class D>
924 {
925 return m_layout;
926 }
927
928 template <class D>
929 inline bool xstrided_container<D>::is_contiguous() const noexcept
930 {
931 using str_type = typename inner_strides_type::value_type;
932 auto is_zero = [](auto i)
933 {
934 return i == 0;
935 };
936 if (!is_contiguous_container<storage_type>::value)
937 {
938 return false;
939 }
940 // We need to make sure the inner-most non-zero stride is one.
941 // Trailing zero strides are ignored because they indicate bradcasted dimensions.
942 if (m_layout == layout_type::row_major)
943 {
944 auto it = std::find_if_not(m_strides.rbegin(), m_strides.rend(), is_zero);
945 // If the array has strides of zero, it is a constant, and therefore contiguous.
946 return it == m_strides.rend() || *it == str_type(1);
947 }
948 else if (m_layout == layout_type::column_major)
949 {
950 auto it = std::find_if_not(m_strides.begin(), m_strides.end(), is_zero);
951 // If the array has strides of zero, it is a constant, and therefore contiguous.
952 return it == m_strides.end() || *it == str_type(1);
953 }
954 else
955 {
956 return m_strides.empty();
957 }
958 }
959
960 namespace detail
961 {
962 template <class C, class S>
963 inline void resize_data_container(C& c, S size)
964 {
965 xt::resize_container(c, size);
966 }
967
968 template <class C, class S>
969 inline void resize_data_container(const C& c, S size)
970 {
971 (void) c; // remove unused parameter warning
972 (void) size;
973 XTENSOR_ASSERT_MSG(c.size() == size, "Trying to resize const data container with wrong size.");
974 }
975
976 template <class S, class T>
977 constexpr bool check_resize_dimension(const S&, const T&)
978 {
979 return true;
980 }
981
982 template <class T, size_t N, class S>
983 constexpr bool check_resize_dimension(const std::array<T, N>&, const S& s)
984 {
985 return N == s.size();
986 }
987 }
988
996 template <class D>
997 template <class S>
998 inline void xstrided_container<D>::resize(S&& shape, bool force)
999 {
1000 XTENSOR_ASSERT_MSG(
1001 detail::check_resize_dimension(m_shape, shape),
1002 "cannot change the number of dimensions of xtensor"
1003 )
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))
1006 || force)
1007 {
1008 if (D::static_layout == layout_type::dynamic && m_layout == layout_type::dynamic)
1009 {
1010 m_layout = XTENSOR_DEFAULT_LAYOUT; // fall back to default layout
1011 }
1012 m_shape = xtl::forward_sequence<shape_type, S>(shape);
1013
1014 resize_container(m_strides, dim);
1015 resize_container(m_backstrides, dim);
1016 size_type data_size = compute_strides<D::static_layout>(m_shape, m_layout, m_strides, m_backstrides);
1017 detail::resize_data_container(this->storage(), data_size);
1018 }
1019 }
1020
1028 template <class D>
1029 template <class S>
1031 {
1032 XTENSOR_ASSERT_MSG(
1033 detail::check_resize_dimension(m_shape, shape),
1034 "cannot change the number of dimensions of xtensor"
1035 )
1036 if (base_type::static_layout != layout_type::dynamic && l != base_type::static_layout)
1037 {
1038 XTENSOR_THROW(
1039 std::runtime_error,
1040 "Cannot change layout_type if template parameter not layout_type::dynamic."
1041 );
1042 }
1043 m_layout = l;
1044 resize(std::forward<S>(shape), true);
1045 }
1046
1054 template <class D>
1055 template <class S>
1056 inline void xstrided_container<D>::resize(S&& shape, const strides_type& strides)
1057 {
1058 XTENSOR_ASSERT_MSG(
1059 detail::check_resize_dimension(m_shape, shape),
1060 "cannot change the number of dimensions of xtensor"
1061 )
1062 if (base_type::static_layout != layout_type::dynamic)
1063 {
1064 XTENSOR_THROW(
1065 std::runtime_error,
1066 "Cannot resize with custom strides when layout() is != layout_type::dynamic."
1067 );
1068 }
1069 m_shape = xtl::forward_sequence<shape_type, S>(shape);
1070 m_strides = strides;
1071 resize_container(m_backstrides, m_strides.size());
1072 adapt_strides(m_shape, m_strides, m_backstrides);
1073 m_layout = layout_type::dynamic;
1074 detail::resize_data_container(this->storage(), compute_size(m_shape));
1075 }
1076
1090 template <class D>
1091 template <class S>
1093 {
1094 reshape_impl(
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)
1098 );
1099 return this->derived_cast();
1100 }
1101
1102 template <class D>
1103 template <class T>
1104 inline auto& xstrided_container<D>::reshape(std::initializer_list<T> shape, layout_type layout) &
1105 {
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();
1111 }
1112
1113 template <class D>
1114 template <class S>
1115 inline void
1116 xstrided_container<D>::reshape_impl(S&& shape, std::false_type /* is unsigned */, layout_type layout)
1117 {
1118 if (compute_size(shape) != this->size())
1119 {
1120 XTENSOR_THROW(
1121 std::runtime_error,
1122 "Cannot reshape with incorrect number of elements. Do you mean to resize?"
1123 );
1124 }
1125 if (D::static_layout == layout_type::dynamic && layout == layout_type::dynamic)
1126 {
1127 layout = XTENSOR_DEFAULT_LAYOUT; // fall back to default layout
1128 }
1129 if (D::static_layout != layout_type::dynamic && layout != D::static_layout)
1130 {
1131 XTENSOR_THROW(std::runtime_error, "Cannot reshape with different layout if static layout != dynamic.");
1132 }
1133 m_layout = layout;
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());
1137 compute_strides<D::static_layout>(m_shape, m_layout, m_strides, m_backstrides);
1138 }
1139
1140 template <class D>
1141 template <class S>
1142 inline void
1143 xstrided_container<D>::reshape_impl(S&& _shape, std::true_type /* is signed */, layout_type layout)
1144 {
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)
1148 {
1149 XTENSOR_THROW(std::runtime_error, "Negative axis size cannot be inferred. Shape mismatch.");
1150 }
1151 std::decay_t<S> shape = _shape;
1152 tmp_value_type accumulator = 1;
1153 std::size_t neg_idx = 0;
1154 std::size_t i = 0;
1155 for (auto it = shape.begin(); it != shape.end(); ++it, i++)
1156 {
1157 auto&& dim = *it;
1158 if (dim < 0)
1159 {
1160 XTENSOR_ASSERT(dim == -1 && !neg_idx);
1161 neg_idx = i;
1162 }
1163 accumulator *= dim;
1164 }
1165 if (accumulator < 0)
1166 {
1167 shape[neg_idx] = static_cast<tmp_value_type>(this->size()) / std::abs(accumulator);
1168 }
1169 else if (this->size() != new_size)
1170 {
1171 XTENSOR_THROW(
1172 std::runtime_error,
1173 "Cannot reshape with incorrect number of elements. Do you mean to resize?"
1174 );
1175 }
1176 m_layout = layout;
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());
1180 compute_strides<D::static_layout>(m_shape, m_layout, m_strides, m_backstrides);
1181 }
1182
1183 template <class D>
1184 inline auto xstrided_container<D>::mutable_layout() noexcept -> layout_type&
1185 {
1186 return m_layout;
1187 }
1188}
1189
1190#endif
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
pointer data() 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.
Definition xstrides.hpp:570
standard mathematical functions for xexpressions
layout_type
Definition xlayout.hpp:24