xtensor
 
Loading...
Searching...
No Matches
xview.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_VIEW_HPP
11#define XTENSOR_VIEW_HPP
12
13#include <algorithm>
14#include <array>
15#include <cstddef>
16#include <tuple>
17#include <type_traits>
18#include <utility>
19
20#include <xtl/xclosure.hpp>
21#include <xtl/xmeta_utils.hpp>
22#include <xtl/xsequence.hpp>
23#include <xtl/xtype_traits.hpp>
24
25#include "../containers/xarray.hpp"
26#include "../containers/xcontainer.hpp"
27#include "../containers/xtensor.hpp"
28#include "../core/xaccessible.hpp"
29#include "../core/xiterable.hpp"
30#include "../core/xsemantic.hpp"
31#include "../core/xtensor_config.hpp"
32#include "../core/xtensor_forward.hpp"
33#include "../views/xbroadcast.hpp"
34#include "../views/xslice.hpp"
35#include "../views/xview_utils.hpp"
36
37namespace xt
38{
39
40 /*******************
41 * xview extension *
42 *******************/
43
44 namespace extension
45 {
46 template <class Tag, class CT, class... S>
48
49 template <class CT, class... S>
51 {
52 using type = xtensor_empty_base;
53 };
54
55 template <class CT, class... S>
56 struct xview_base : xview_base_impl<xexpression_tag_t<CT>, CT, S...>
57 {
58 };
59
60 template <class CT, class... S>
61 using xview_base_t = typename xview_base<CT, S...>::type;
62 }
63
64 /*********************
65 * xview declaration *
66 *********************/
67
68 template <bool is_const, class CT, class... S>
69 class xview_stepper;
70
71 template <class ST, class... S>
72 struct xview_shape_type;
73
74 namespace detail
75 {
76
77 template <class T>
78 struct is_xrange : std::false_type
79 {
80 };
81
82 template <class T>
83 struct is_xrange<xrange<T>> : std::true_type
84 {
85 };
86
87 template <class S>
88 struct is_xall_slice : std::false_type
89 {
90 };
91
92 template <class T>
93 struct is_xall_slice<xall<T>> : std::true_type
94 {
95 };
96
97 template <layout_type L, bool valid, bool all_seen, bool range_seen, class V>
98 struct is_contiguous_view_impl
99 {
100 static constexpr bool value = false;
101 };
102
103 template <class T>
104 struct static_dimension
105 {
106 static constexpr std::ptrdiff_t value = -1;
107 };
108
109 template <class T, std::size_t N>
110 struct static_dimension<std::array<T, N>>
111 {
112 static constexpr std::ptrdiff_t value = static_cast<std::ptrdiff_t>(N);
113 };
114
115 template <class T, std::size_t N>
116 struct static_dimension<xt::const_array<T, N>>
117 {
118 static constexpr std::ptrdiff_t value = static_cast<std::ptrdiff_t>(N);
119 };
120
121 template <std::size_t... I>
122 struct static_dimension<xt::fixed_shape<I...>>
123 {
124 static constexpr std::ptrdiff_t value = sizeof...(I);
125 };
126
127 // if we have the same number of integers as we have static dimensions
128 // this can be interpreted like a xscalar
129 template <class CT, class... S>
130 struct is_xscalar_impl<xview<CT, S...>>
131 {
132 static constexpr bool value = static_cast<std::ptrdiff_t>(integral_count<S...>()
133 ) == static_dimension<typename std::decay_t<CT>::shape_type>::value
134 ? true
135 : false;
136 };
137
138 template <class S>
139 struct is_strided_slice_impl : std::true_type
140 {
141 };
142
143 template <class T>
144 struct is_strided_slice_impl<xkeep_slice<T>> : std::false_type
145 {
146 };
147
148 template <class T>
149 struct is_strided_slice_impl<xdrop_slice<T>> : std::false_type
150 {
151 };
152
153 // If we have no discontiguous slices, we can calculate strides for this view.
154 template <class E, class... S>
155 struct is_strided_view
156 : std::integral_constant<
157 bool,
158 std::conjunction<has_data_interface<E>, is_strided_slice_impl<std::decay_t<S>>...>::value>
159 {
160 };
161
162 // if row major the view can only be (statically) computed as contiguous if:
163 // any number of integers is followed by either one or no range which
164 // are followed by explicit (or implicit) all's
165 //
166 // e.g.
167 // (i, j, all(), all()) == contiguous
168 // (i, range(0, 2), all()) == contiguous
169 // (i) == contiguous (implicit all slices)
170 // (i, all(), j) == *not* contiguous
171 // (i, range(0, 2), range(0, 2)) == *not* contiguous etc.
172 template <bool valid, bool all_seen, bool range_seen, class V>
173 struct is_contiguous_view_impl<layout_type::row_major, valid, all_seen, range_seen, V>
174 {
175 using slice = xtl::mpl::front_t<V>;
176 static constexpr bool is_range_slice = is_xrange<slice>::value;
177 static constexpr bool is_int_slice = xtl::is_integral<slice>::value;
178 static constexpr bool is_all_slice = is_xall_slice<slice>::value;
179 static constexpr bool have_all_seen = all_seen || is_all_slice;
180 static constexpr bool have_range_seen = is_range_slice;
181
182 static constexpr bool is_valid = valid
183 && (have_all_seen
184 ? is_all_slice
185 : (!range_seen && (is_int_slice || is_range_slice)));
186
187 static constexpr bool value = is_contiguous_view_impl < layout_type::row_major, is_valid,
188 have_all_seen, range_seen || is_range_slice,
189 xtl::mpl::pop_front_t < V >> ::value;
190 };
191
192 template <bool valid, bool all_seen, bool range_seen>
193 struct is_contiguous_view_impl<layout_type::row_major, valid, all_seen, range_seen, xtl::mpl::vector<>>
194 {
195 static constexpr bool value = valid;
196 };
197
198 // For column major the *same* but reverse is true -- with the additional
199 // constraint that we have to know the dimension at compile time otherwise
200 // we cannot make the decision as there might be implicit all's following.
201 template <bool valid, bool int_seen, bool range_seen, class V>
202 struct is_contiguous_view_impl<layout_type::column_major, valid, int_seen, range_seen, V>
203 {
204 using slice = xtl::mpl::front_t<V>;
205 static constexpr bool is_range_slice = is_xrange<slice>::value;
206 static constexpr bool is_int_slice = xtl::is_integral<slice>::value;
207 static constexpr bool is_all_slice = is_xall_slice<slice>::value;
208
209 static constexpr bool have_int_seen = int_seen || is_int_slice;
210
211 static constexpr bool is_valid = valid
212 && (have_int_seen
213 ? is_int_slice
214 : (!range_seen && (is_all_slice || is_range_slice)));
215 static constexpr bool value = is_contiguous_view_impl < layout_type::column_major, is_valid,
216 have_int_seen, is_range_slice || range_seen,
217 xtl::mpl::pop_front_t < V >> ::value;
218 };
219
220 template <bool valid, bool int_seen, bool range_seen>
221 struct is_contiguous_view_impl<layout_type::column_major, valid, int_seen, range_seen, xtl::mpl::vector<>>
222 {
223 static constexpr bool value = valid;
224 };
225
226 // TODO relax has_data_interface constraint here!
227 template <class E, class... S>
228 struct is_contiguous_view
229 : std::integral_constant<
230 bool,
231 has_data_interface<E>::value
232 && !(
233 E::static_layout == layout_type::column_major
234 && static_cast<std::size_t>(static_dimension<typename E::shape_type>::value) != sizeof...(S)
235 )
236 && is_contiguous_view_impl<E::static_layout, true, false, false, xtl::mpl::vector<S...>>::value>
237 {
238 };
239
240 template <layout_type L, class T, std::ptrdiff_t offset>
241 struct unwrap_offset_container
242 {
243 using type = void;
244 };
245
246 template <class T, std::ptrdiff_t offset>
247 struct unwrap_offset_container<layout_type::row_major, T, offset>
248 {
249 using type = sequence_view<T, offset, static_dimension<T>::value>;
250 };
251
252 template <class T, std::ptrdiff_t start, std::ptrdiff_t end, std::ptrdiff_t offset>
253 struct unwrap_offset_container<layout_type::row_major, sequence_view<T, start, end>, offset>
254 {
255 using type = sequence_view<T, start + offset, end>;
256 };
257
258 template <class T, std::ptrdiff_t offset>
259 struct unwrap_offset_container<layout_type::column_major, T, offset>
260 {
261 using type = sequence_view<T, 0, static_dimension<T>::value - offset>;
262 };
263
264 template <class T, std::ptrdiff_t start, std::ptrdiff_t end, std::ptrdiff_t offset>
265 struct unwrap_offset_container<layout_type::column_major, sequence_view<T, start, end>, offset>
266 {
267 using type = sequence_view<T, start, end - offset>;
268 };
269
270 template <class E, class... S>
271 struct get_contigous_shape_type
272 {
273 // if we have no `range` in the slices we can re-use the shape with an offset
274 using type = std::conditional_t<
275 std::disjunction<is_xrange<S>...>::value,
276 typename xview_shape_type<typename E::shape_type, S...>::type,
277 // In the false branch we know that we have only integers at the front OR end, and NO range
278 typename unwrap_offset_container<E::static_layout, typename E::inner_shape_type, integral_count<S...>()>::type>;
279 };
280
281 template <class T>
282 struct is_sequence_view : std::integral_constant<bool, false>
283 {
284 };
285
286 template <class T, std::ptrdiff_t S, std::ptrdiff_t E>
287 struct is_sequence_view<sequence_view<T, S, E>> : std::integral_constant<bool, true>
288 {
289 };
290 }
291
292 template <class E, class... S>
293 concept contiguous_view_concept = detail::is_contiguous_view<E, S...>::value;
294 template <class E, class... S>
295 concept strided_view_concept = detail::is_strided_view<std::decay_t<E>, S...>::value;
296
297 template <class CT, class... S>
299 {
300 using xexpression_type = std::decay_t<CT>;
301 using reference = inner_reference_t<CT>;
302 using const_reference = typename xexpression_type::const_reference;
303 using size_type = typename xexpression_type::size_type;
304 using temporary_type = view_temporary_type_t<xexpression_type, S...>;
305
306 static constexpr layout_type layout = detail::is_contiguous_view<xexpression_type, S...>::value
307 ? xexpression_type::static_layout
309
310 static constexpr bool is_const = std::is_const<std::remove_reference_t<CT>>::value;
311
312 using extract_storage_type = xtl::mpl::eval_if_t<
314 detail::expr_storage_type<xexpression_type>,
316 using storage_type = std::conditional_t<is_const, const extract_storage_type, extract_storage_type>;
317 };
318
319 template <class CT, class... S>
320 struct xiterable_inner_types<xview<CT, S...>>
321 {
322 using xexpression_type = std::decay_t<CT>;
323
324 static constexpr bool is_strided_view = detail::is_strided_view<xexpression_type, S...>::value;
325 static constexpr bool is_contiguous_view = detail::is_contiguous_view<xexpression_type, S...>::value;
326
327 using inner_shape_type = std::conditional_t<
328 is_contiguous_view,
329 typename detail::get_contigous_shape_type<xexpression_type, S...>::type,
330 typename xview_shape_type<typename xexpression_type::shape_type, S...>::type>;
331
332 using stepper = std::conditional_t<
333 is_strided_view,
334 xstepper<xview<CT, S...>>,
336
337 using const_stepper = std::conditional_t<
338 is_strided_view,
339 xstepper<const xview<CT, S...>>,
341 };
342
357 template <class CT, class... S>
358 class xview : public xview_semantic<xview<CT, S...>>,
359 public std::conditional_t<
360 detail::is_contiguous_view<std::decay_t<CT>, S...>::value,
361 xcontiguous_iterable<xview<CT, S...>>,
362 xiterable<xview<CT, S...>>>,
363 public xaccessible<xview<CT, S...>>,
364 public extension::xview_base_t<CT, S...>
365 {
366 public:
367
368 using self_type = xview<CT, S...>;
369 using inner_types = xcontainer_inner_types<self_type>;
370 using xexpression_type = std::decay_t<CT>;
371 using semantic_base = xview_semantic<self_type>;
372 using temporary_type = typename xcontainer_inner_types<self_type>::temporary_type;
373
374 using accessible_base = xaccessible<self_type>;
375 using extension_base = extension::xview_base_t<CT, S...>;
376 using expression_tag = typename extension_base::expression_tag;
377
378 static constexpr bool is_const = std::is_const<std::remove_reference_t<CT>>::value;
379 using value_type = typename xexpression_type::value_type;
380 using simd_value_type = xt_simd::simd_type<value_type>;
381 using bool_load_type = typename xexpression_type::bool_load_type;
382 using reference = typename inner_types::reference;
383 using const_reference = typename inner_types::const_reference;
384 using pointer = std::
385 conditional_t<is_const, typename xexpression_type::const_pointer, typename xexpression_type::pointer>;
386 using const_pointer = typename xexpression_type::const_pointer;
387 using size_type = typename inner_types::size_type;
388 using difference_type = typename xexpression_type::difference_type;
389
390 static constexpr layout_type static_layout = inner_types::layout;
391 static constexpr bool contiguous_layout = static_layout != layout_type::dynamic;
392
393 static constexpr bool is_strided_view = detail::is_strided_view<xexpression_type, S...>::value;
394 static constexpr bool is_contiguous_view = contiguous_layout;
395
396 using iterable_base = xiterable<self_type>;
397 using inner_shape_type = typename iterable_base::inner_shape_type;
398 using shape_type = typename xview_shape_type<typename xexpression_type::shape_type, S...>::type;
399
400 using xexpression_inner_strides_type = xtl::mpl::eval_if_t<
402 detail::expr_inner_strides_type<xexpression_type>,
404
405 using xexpression_inner_backstrides_type = xtl::mpl::eval_if_t<
407 detail::expr_inner_backstrides_type<xexpression_type>,
409
410 using storage_type = typename inner_types::storage_type;
411
412 static constexpr bool has_trivial_strides = is_contiguous_view
413 && !std::disjunction<detail::is_xrange<S>...>::value;
414 using inner_strides_type = std::conditional_t<
415 has_trivial_strides,
416 typename detail::unwrap_offset_container<
417 xexpression_type::static_layout,
418 xexpression_inner_strides_type,
419 integral_count<S...>()>::type,
420 get_strides_t<shape_type>>;
421
422 using inner_backstrides_type = std::conditional_t<
423 has_trivial_strides,
424 typename detail::unwrap_offset_container<
425 xexpression_type::static_layout,
426 xexpression_inner_backstrides_type,
427 integral_count<S...>()>::type,
428 get_strides_t<shape_type>>;
429
430 using strides_type = get_strides_t<shape_type>;
431 using backstrides_type = strides_type;
432
433
434 using slice_type = std::tuple<S...>;
435
436 using stepper = typename iterable_base::stepper;
437 using const_stepper = typename iterable_base::const_stepper;
438
439 using linear_iterator = std::conditional_t<
441 std::conditional_t<is_const, typename xexpression_type::const_linear_iterator, typename xexpression_type::linear_iterator>,
442 typename iterable_base::linear_iterator>;
443 using const_linear_iterator = std::conditional_t<
445 typename xexpression_type::const_linear_iterator,
446 typename iterable_base::const_linear_iterator>;
447
448 using reverse_linear_iterator = std::reverse_iterator<linear_iterator>;
449 using const_reverse_linear_iterator = std::reverse_iterator<const_linear_iterator>;
450
451 using container_iterator = pointer;
452 using const_container_iterator = const_pointer;
453 static constexpr std::size_t rank = SIZE_MAX;
454
455 // The FSL argument prevents the compiler from calling this constructor
456 // instead of the copy constructor when sizeof...(SL) == 0.
457 template <class CTA, class FSL, class... SL>
458 explicit xview(CTA&& e, FSL&& first_slice, SL&&... slices) noexcept;
459
460 xview(const xview&) = default;
461 self_type& operator=(const xview& rhs);
462
463 template <class E>
464 self_type& operator=(const xexpression<E>& e);
465
466 template <class E>
467 disable_xexpression<E, self_type>& operator=(const E& e);
468
469 const inner_shape_type& shape() const noexcept;
470 const slice_type& slices() const noexcept;
471 layout_type layout() const noexcept;
472 bool is_contiguous() const noexcept;
473 using accessible_base::shape;
474
475 template <class T>
476 void fill(const T& value);
477
478 template <class... Args>
479 reference operator()(Args... args);
480 template <class... Args>
481 reference unchecked(Args... args);
482 template <class It>
483 reference element(It first, It last);
484
485 template <class... Args>
486 const_reference operator()(Args... args) const;
487 template <class... Args>
488 const_reference unchecked(Args... args) const;
489 template <class It>
490 const_reference element(It first, It last) const;
491
492 xexpression_type& expression() noexcept;
493 const xexpression_type& expression() const noexcept;
494
495 template <class ST>
496 bool broadcast_shape(ST& shape, bool reuse_cache = false) const;
497
498 template <class ST>
499 bool has_linear_assign(const ST& strides) const;
500
501 template <class ST>
502 stepper stepper_begin(const ST& shape);
503 template <class ST>
504 stepper stepper_end(const ST& shape, layout_type l);
505
506 template <class ST>
507 const_stepper stepper_begin(const ST& shape) const;
508 template <class ST>
509 const_stepper stepper_end(const ST& shape, layout_type l) const;
510
511 template <class T = xexpression_type>
512 storage_type& storage()
513 requires(has_data_interface_concept<T>);
514
515 template <class T = xexpression_type>
516 const storage_type& storage() const
517 requires(has_data_interface_concept<T>);
518
519 template <class T = xexpression_type>
520 linear_iterator linear_begin()
521 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>);
522
523 template <class T = xexpression_type>
524 linear_iterator linear_end()
525 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>);
526
527 template <class T = xexpression_type>
528 const_linear_iterator linear_begin() const
529 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>);
530
531 template <class T = xexpression_type>
532 const_linear_iterator linear_end() const
533 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>);
534
535 template <class T = xexpression_type>
536 const_linear_iterator linear_cbegin() const
537 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>);
538
539 template <class T = xexpression_type>
540 const_linear_iterator linear_cend() const
541 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>);
542
543 template <class T = xexpression_type>
544 reverse_linear_iterator linear_rbegin()
545 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>);
546
547 template <class T = xexpression_type>
548 reverse_linear_iterator linear_rend()
549 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>);
550
551 template <class T = xexpression_type>
552 const_reverse_linear_iterator linear_rbegin() const
553 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>);
554
555 template <class T = xexpression_type>
556 const_reverse_linear_iterator linear_rend() const
557 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>);
558
559 template <class T = xexpression_type>
560 const_reverse_linear_iterator linear_crbegin() const
561 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>);
562
563 template <class T = xexpression_type>
564 const_reverse_linear_iterator linear_crend() const
565 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>);
566
567 template <class T = xexpression_type>
568 const inner_strides_type& strides() const
569 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>);
570
571 template <class T = xexpression_type>
572 const inner_strides_type& backstrides() const
573 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>);
574
575 template <class T = xexpression_type>
576 const_pointer data() const
577 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>);
578
579 template <class T = xexpression_type>
580 pointer data()
581 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>);
582
583 template <class T = xexpression_type>
584 std::size_t data_offset() const noexcept
585 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>);
586
587 template <class It>
588 inline It data_xbegin_impl(It begin) const noexcept;
589
590 template <class It>
591 inline It data_xend_impl(It begin, layout_type l, size_type offset) const noexcept;
592 inline container_iterator data_xbegin() noexcept;
593 inline const_container_iterator data_xbegin() const noexcept;
594 inline container_iterator data_xend(layout_type l, size_type offset) noexcept;
595
596 inline const_container_iterator data_xend(layout_type l, size_type offset) const noexcept;
597
598 // Conversion operator enabled for statically "scalar" views
599 template <xscalar_concept ST = self_type>
600 operator reference()
601 {
602 return (*this)();
603 }
604
605 template <xscalar_concept ST = self_type>
606 operator const_reference() const
607 {
608 return (*this)();
609 }
610
611 size_type underlying_size(size_type dim) const;
612
613 xtl::xclosure_pointer<self_type&> operator&() &;
614 xtl::xclosure_pointer<const self_type&> operator&() const&;
615 xtl::xclosure_pointer<self_type> operator&() &&;
616
617 template <class E, class T = xexpression_type>
618 void assign_to(xexpression<E>& e, bool force_resize) const
619 requires(has_data_interface_concept<T> and contiguous_view_concept<E, S...>);
620
621 template <class E>
622 using rebind_t = xview<E, S...>;
623
624 template <class E>
625 rebind_t<E> build_view(E&& e) const;
626
627 //
628 // SIMD interface
629 //
630
631 template <class requested_type>
632 using simd_return_type = xt_simd::simd_return_type<value_type, requested_type>;
633
634 template <class align, class simd, class T = xexpression_type>
635 void store_simd(size_type i, const simd& e)
636 requires(has_simd_interface_concept<T> and strided_view_concept<CT, S...>);
637
638 template <
639 class align,
640 class requested_type = value_type,
641 std::size_t N = xt_simd::simd_traits<requested_type>::size,
642 class T = xexpression_type>
643 simd_return_type<requested_type> load_simd(size_type i) const
644 requires(has_simd_interface_concept<T> and strided_view_concept<CT, S...>);
645
646 template <class T = xexpression_type>
647 reference data_element(size_type i)
648 requires(has_simd_interface_concept<T> and strided_view_concept<CT, S...>);
649
650 template <class T = xexpression_type>
651 const_reference data_element(size_type i) const
652 requires(has_simd_interface_concept<T> and strided_view_concept<CT, S...>);
653
654 template <class T = xexpression_type>
655 reference flat(size_type i)
656 requires(has_simd_interface_concept<T> and strided_view_concept<CT, S...>);
657
658 template <class T = xexpression_type>
659 const_reference flat(size_type i) const
660 requires(has_simd_interface_concept<T> and strided_view_concept<CT, S...>);
661
662 private:
663
664 // VS 2015 workaround (yes, really)
665 template <std::size_t I>
666 struct lesser_condition
667 {
668 static constexpr bool value = (I + newaxis_count_before<S...>(I + 1) < sizeof...(S));
669 };
670
671 CT m_e;
672 slice_type m_slices;
673 inner_shape_type m_shape;
674 mutable inner_strides_type m_strides;
675 mutable inner_backstrides_type m_backstrides;
676 mutable std::size_t m_data_offset;
677 mutable bool m_strides_computed;
678
679 template <class CTA, class FSL, class... SL>
680 explicit xview(std::true_type, CTA&& e, FSL&& first_slice, SL&&... slices) noexcept;
681
682 template <class CTA, class FSL, class... SL>
683 explicit xview(std::false_type, CTA&& e, FSL&& first_slice, SL&&... slices) noexcept;
684
685 template <class... Args>
686 auto make_index_sequence(Args... args) const noexcept;
687
688 void compute_strides(std::true_type) const;
689 void compute_strides(std::false_type) const;
690
691 reference access();
692
693 template <class Arg, class... Args>
694 reference access(Arg arg, Args... args);
695
696 const_reference access() const;
697
698 template <class Arg, class... Args>
699 const_reference access(Arg arg, Args... args) const;
700
701 template <typename std::decay_t<CT>::size_type... I, class... Args>
702 reference unchecked_impl(std::index_sequence<I...>, Args... args);
703
704 template <typename std::decay_t<CT>::size_type... I, class... Args>
705 const_reference unchecked_impl(std::index_sequence<I...>, Args... args) const;
706
707 template <typename std::decay_t<CT>::size_type... I, class... Args>
708 reference access_impl(std::index_sequence<I...>, Args... args);
709
710 template <typename std::decay_t<CT>::size_type... I, class... Args>
711 const_reference access_impl(std::index_sequence<I...>, Args... args) const;
712
713 template <typename std::decay_t<CT>::size_type I, class... Args>
714 size_type index(Args... args) const;
715
716 template <typename std::decay_t<CT>::size_type, class T>
717 size_type sliced_access(const xslice<T>& slice) const;
718
719 template <typename std::decay_t<CT>::size_type I, class T, class Arg, class... Args>
720 size_type sliced_access(const xslice<T>& slice, Arg arg, Args... args) const;
721
722 template <typename std::decay_t<CT>::size_type I, class T, class... Args>
723 disable_xslice<T, size_type> sliced_access(const T& squeeze, Args...) const;
724
725 using base_index_type = xindex_type_t<typename xexpression_type::shape_type>;
726
727 template <class It>
728 base_index_type make_index(It first, It last) const;
729
730 void assign_temporary_impl(temporary_type&& tmp);
731
732 template <std::size_t... I>
733 std::size_t data_offset_impl(std::index_sequence<I...>) const noexcept;
734
735 template <std::size_t... I>
736 auto compute_strides_impl(std::index_sequence<I...>) const noexcept;
737
738 inner_shape_type compute_shape(std::true_type) const;
739 inner_shape_type compute_shape(std::false_type) const;
740
741 template <class E, std::size_t... I>
742 rebind_t<E> build_view_impl(E&& e, std::index_sequence<I...>) const;
743
744 friend class xview_semantic<xview<CT, S...>>;
745 };
746
747 template <class E, class... S>
748 auto view(E&& e, S&&... slices);
749
750 template <class E>
751 auto row(E&& e, std::ptrdiff_t index);
752
753 template <class E>
754 auto col(E&& e, std::ptrdiff_t index);
755
756 /*****************************
757 * xview_stepper declaration *
758 *****************************/
759
760 namespace detail
761 {
762 template <class V>
763 struct get_stepper_impl
764 {
765 using xexpression_type = typename V::xexpression_type;
766 using type = typename xexpression_type::stepper;
767 };
768
769 template <class V>
770 struct get_stepper_impl<const V>
771 {
772 using xexpression_type = typename V::xexpression_type;
773 using type = typename xexpression_type::const_stepper;
774 };
775 }
776
777 template <class V>
778 using get_stepper = typename detail::get_stepper_impl<V>::type;
779
780 template <bool is_const, class CT, class... S>
781 class xview_stepper
782 {
783 public:
784
785 using view_type = std::conditional_t<is_const, const xview<CT, S...>, xview<CT, S...>>;
786 using substepper_type = get_stepper<view_type>;
787
788 using value_type = typename substepper_type::value_type;
789 using reference = typename substepper_type::reference;
790 using pointer = typename substepper_type::pointer;
791 using difference_type = typename substepper_type::difference_type;
792 using size_type = typename view_type::size_type;
793
794 using shape_type = typename substepper_type::shape_type;
795
796 xview_stepper() = default;
797 xview_stepper(
798 view_type* view,
799 substepper_type it,
800 size_type offset,
801 bool end = false,
802 layout_type l = XTENSOR_DEFAULT_TRAVERSAL
803 );
804
805 reference operator*() const;
806
807 void step(size_type dim);
808 void step_back(size_type dim);
809 void step(size_type dim, size_type n);
810 void step_back(size_type dim, size_type n);
811 void reset(size_type dim);
812 void reset_back(size_type dim);
813
814 void to_begin();
815 void to_end(layout_type l);
816
817 private:
818
819 bool is_newaxis_slice(size_type index) const noexcept;
820 void to_end_impl(layout_type l);
821
822 template <class F>
823 void common_step_forward(size_type dim, F f);
824 template <class F>
825 void common_step_backward(size_type dim, F f);
826
827 template <class F>
828 void common_step_forward(size_type dim, size_type n, F f);
829 template <class F>
830 void common_step_backward(size_type dim, size_type n, F f);
831
832 template <class F>
833 void common_reset(size_type dim, F f, bool backwards);
834
835 view_type* p_view;
836 substepper_type m_it;
837 size_type m_offset;
838 std::array<std::size_t, sizeof...(S)> m_index_keeper;
839 };
840
841 // meta-function returning the shape type for an xview
842 template <class ST, class... S>
844 {
845 using type = ST;
846 };
847
848 template <class I, std::size_t L, class... S>
849 struct xview_shape_type<std::array<I, L>, S...>
850 {
851 using type = std::array<I, L - integral_count<S...>() + newaxis_count<S...>()>;
852 };
853
854 template <std::size_t... I, class... S>
855 struct xview_shape_type<fixed_shape<I...>, S...>
856 {
857 using type = typename xview_shape_type<std::array<std::size_t, sizeof...(I)>, S...>::type;
858 };
859
860 /************************
861 * xview implementation *
862 ************************/
863
867
869
878 template <class CT, class... S>
879 template <class CTA, class FSL, class... SL>
880 xview<CT, S...>::xview(CTA&& e, FSL&& first_slice, SL&&... slices) noexcept
881 : xview(
882 std::integral_constant<bool, has_trivial_strides>{},
883 std::forward<CTA>(e),
884 std::forward<FSL>(first_slice),
885 std::forward<SL>(slices)...
886 )
887 {
888 }
889
890 // trivial strides initializer
891 template <class CT, class... S>
892 template <class CTA, class FSL, class... SL>
893 xview<CT, S...>::xview(std::true_type, CTA&& e, FSL&& first_slice, SL&&... slices) noexcept
894 : m_e(std::forward<CTA>(e))
895 , m_slices(std::forward<FSL>(first_slice), std::forward<SL>(slices)...)
896 , m_shape(compute_shape(detail::is_sequence_view<inner_shape_type>{}))
897 , m_strides(m_e.strides())
898 , m_backstrides(m_e.backstrides())
899 , m_data_offset(data_offset_impl(std::make_index_sequence<sizeof...(S)>()))
900 , m_strides_computed(true)
901 {
902 }
903
904 template <class CT, class... S>
905 template <class CTA, class FSL, class... SL>
906 xview<CT, S...>::xview(std::false_type, CTA&& e, FSL&& first_slice, SL&&... slices) noexcept
907 : m_e(std::forward<CTA>(e))
908 , m_slices(std::forward<FSL>(first_slice), std::forward<SL>(slices)...)
909 , m_shape(compute_shape(std::false_type{}))
910 , m_strides_computed(false)
911 {
912 }
913
915
916 template <class CT, class... S>
917 inline auto xview<CT, S...>::operator=(const xview& rhs) -> self_type&
918 {
919 temporary_type tmp(rhs);
920 return this->assign_temporary(std::move(tmp));
921 }
922
927
930 template <class CT, class... S>
931 template <class E>
932 inline auto xview<CT, S...>::operator=(const xexpression<E>& e) -> self_type&
933 {
934 return semantic_base::operator=(e);
935 }
936
938
939 template <class CT, class... S>
940 template <class E>
941 inline auto xview<CT, S...>::operator=(const E& e) -> disable_xexpression<E, self_type>&
942 {
943 this->fill(e);
944 return *this;
945 }
946
951
954 template <class CT, class... S>
955 inline auto xview<CT, S...>::shape() const noexcept -> const inner_shape_type&
956 {
957 return m_shape;
958 }
959
963 template <class CT, class... S>
964 inline auto xview<CT, S...>::slices() const noexcept -> const slice_type&
965 {
966 return m_slices;
967 }
968
972 template <class CT, class... S>
973 inline layout_type xview<CT, S...>::layout() const noexcept
974 {
975 if constexpr (is_strided_view)
976 {
977 if (static_layout != layout_type::dynamic)
978 {
979 return static_layout;
980 }
981 else
982 {
983 bool strides_match = do_strides_match(shape(), strides(), m_e.layout(), true);
984 return strides_match ? m_e.layout() : layout_type::dynamic;
985 }
986 }
987 else
988 {
990 }
991 }
992
993 template <class CT, class... S>
994 inline bool xview<CT, S...>::is_contiguous() const noexcept
995 {
996 return layout() != layout_type::dynamic;
997 }
998
1000
1005
1010 template <class CT, class... S>
1011 template <class T>
1012 inline void xview<CT, S...>::fill(const T& value)
1013 {
1014 if constexpr (static_layout != layout_type::dynamic)
1015 {
1016 std::fill(linear_begin(), linear_end(), value);
1017 }
1018 else
1019 {
1020 std::fill(this->begin(), this->end(), value);
1021 }
1022 }
1023
1030 template <class CT, class... S>
1031 template <class... Args>
1032 inline auto xview<CT, S...>::operator()(Args... args) -> reference
1033 {
1034 XTENSOR_TRY(check_index(shape(), args...));
1035 XTENSOR_CHECK_DIMENSION(shape(), args...);
1036 // The static cast prevents the compiler from instantiating the template methods with signed integers,
1037 // leading to warning about signed/unsigned conversions in the deeper layers of the access methods
1038 return access(static_cast<size_type>(args)...);
1039 }
1040
1060 template <class CT, class... S>
1061 template <class... Args>
1062 inline auto xview<CT, S...>::unchecked(Args... args) -> reference
1063 {
1064 return unchecked_impl(make_index_sequence(args...), static_cast<size_type>(args)...);
1065 }
1066
1067 template <class CT, class... S>
1068 template <class It>
1069 inline auto xview<CT, S...>::element(It first, It last) -> reference
1070 {
1071 XTENSOR_TRY(check_element_index(shape(), first, last));
1072 // TODO: avoid memory allocation
1073 auto index = make_index(first, last);
1074 return m_e.element(index.cbegin(), index.cend());
1075 }
1076
1083 template <class CT, class... S>
1084 template <class... Args>
1085 inline auto xview<CT, S...>::operator()(Args... args) const -> const_reference
1086 {
1087 XTENSOR_TRY(check_index(shape(), args...));
1088 XTENSOR_CHECK_DIMENSION(shape(), args...);
1089 // The static cast prevents the compiler from instantiating the template methods with signed integers,
1090 // leading to warning about signed/unsigned conversions in the deeper layers of the access methods
1091 return access(static_cast<size_type>(args)...);
1092 }
1093
1113 template <class CT, class... S>
1114 template <class... Args>
1115 inline auto xview<CT, S...>::unchecked(Args... args) const -> const_reference
1116 {
1117 return unchecked_impl(make_index_sequence(args...), static_cast<size_type>(args)...);
1118 }
1119
1120 template <class CT, class... S>
1121 template <class It>
1122 inline auto xview<CT, S...>::element(It first, It last) const -> const_reference
1123 {
1124 // TODO: avoid memory allocation
1125 auto index = make_index(first, last);
1126 return m_e.element(index.cbegin(), index.cend());
1127 }
1128
1132 template <class CT, class... S>
1133 inline auto xview<CT, S...>::expression() noexcept -> xexpression_type&
1134 {
1135 return m_e;
1136 }
1137
1141 template <class CT, class... S>
1142 inline auto xview<CT, S...>::expression() const noexcept -> const xexpression_type&
1143 {
1144 return m_e;
1145 }
1146
1152 template <class CT, class... S>
1153 template <class T>
1154 inline auto xview<CT, S...>::storage() -> storage_type&
1156 {
1157 return m_e.storage();
1158 }
1159
1160 template <class CT, class... S>
1161 template <class T>
1162 inline auto xview<CT, S...>::storage() const -> const storage_type&
1163 requires(has_data_interface_concept<T>)
1164 {
1165 return m_e.storage();
1166 }
1167
1168 template <class CT, class... S>
1169 template <class T>
1170 auto xview<CT, S...>::linear_begin() -> linear_iterator
1171 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>)
1172 {
1173 return m_e.storage().begin() + data_offset();
1174 }
1175
1176 template <class CT, class... S>
1177 template <class T>
1178 auto xview<CT, S...>::linear_end() -> linear_iterator
1180 {
1181 return m_e.storage().begin() + data_offset() + this->size();
1182 }
1183
1184 template <class CT, class... S>
1185 template <class T>
1186 auto xview<CT, S...>::linear_begin() const -> const_linear_iterator
1187 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>)
1188 {
1189 return linear_cbegin();
1190 }
1191
1192 template <class CT, class... S>
1193 template <class T>
1194 auto xview<CT, S...>::linear_end() const -> const_linear_iterator
1195 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>)
1196 {
1197 return linear_cend();
1198 }
1199
1200 template <class CT, class... S>
1201 template <class T>
1202 auto xview<CT, S...>::linear_cbegin() const -> const_linear_iterator
1203 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>)
1204 {
1205 return m_e.storage().cbegin() + data_offset();
1206 }
1207
1208 template <class CT, class... S>
1209 template <class T>
1210 auto xview<CT, S...>::linear_cend() const -> const_linear_iterator
1211 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>)
1212 {
1213 return m_e.storage().cbegin() + data_offset() + this->size();
1214 }
1215
1216 template <class CT, class... S>
1217 template <class T>
1218 auto xview<CT, S...>::linear_rbegin() -> reverse_linear_iterator
1220 {
1221 return reverse_linear_iterator(linear_end());
1222 }
1223
1224 template <class CT, class... S>
1225 template <class T>
1226 auto xview<CT, S...>::linear_rend() -> reverse_linear_iterator
1228 {
1229 return reverse_linear_iterator(linear_begin());
1230 }
1231
1232 template <class CT, class... S>
1233 template <class T>
1234 auto xview<CT, S...>::linear_rbegin() const -> const_reverse_linear_iterator
1235 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>)
1236 {
1237 return linear_crbegin();
1238 }
1239
1240 template <class CT, class... S>
1241 template <class T>
1242 auto xview<CT, S...>::linear_rend() const -> const_reverse_linear_iterator
1243 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>)
1244 {
1245 return linear_crend();
1246 }
1247
1248 template <class CT, class... S>
1249 template <class T>
1250 auto xview<CT, S...>::linear_crbegin() const -> const_reverse_linear_iterator
1251 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>)
1252 {
1253 return const_reverse_linear_iterator(linear_end());
1254 }
1255
1256 template <class CT, class... S>
1257 template <class T>
1258 auto xview<CT, S...>::linear_crend() const -> const_reverse_linear_iterator
1259 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>)
1260 {
1261 return const_reverse_linear_iterator(linear_begin());
1262 }
1263
1267 template <class CT, class... S>
1268 template <class T>
1269 inline auto xview<CT, S...>::strides() const
1270 -> const inner_strides_type& requires(has_data_interface_concept<T>and strided_view_concept<CT, S...>) {
1271 if (!m_strides_computed)
1272 {
1273 compute_strides(std::integral_constant<bool, has_trivial_strides>{});
1274 m_strides_computed = true;
1275 }
1276 return m_strides;
1277 }
1278
1279 template <class CT, class... S>
1280 template <class T>
1281 inline auto xview<CT, S...>::backstrides() const
1282 -> const inner_strides_type& requires(has_data_interface_concept<T>and strided_view_concept<CT, S...>) {
1283 if (!m_strides_computed)
1284 {
1285 compute_strides(std::integral_constant<bool, has_trivial_strides>{});
1286 m_strides_computed = true;
1287 }
1288 return m_backstrides;
1289 }
1290
1294 template <class CT, class... S>
1295 template <class T>
1296 inline auto xview<CT, S...>::data() const -> const_pointer
1297 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>)
1298 {
1299 return m_e.data();
1300 }
1301
1302 template <class CT, class... S>
1303 template <class T>
1304 inline auto xview<CT, S...>::data() -> pointer
1306 {
1307 return m_e.data();
1308 }
1309
1310 template <class CT, class... S>
1311 template <std::size_t... I>
1312 inline std::size_t xview<CT, S...>::data_offset_impl(std::index_sequence<I...>) const noexcept
1313 {
1314 auto temp = std::array<std::ptrdiff_t, sizeof...(S)>(
1315 {(static_cast<ptrdiff_t>(xt::value(std::get<I>(m_slices), 0)))...}
1316 );
1317
1318 std::ptrdiff_t result = 0;
1319 std::size_t i = 0;
1320 for (; i < std::min(sizeof...(S), m_e.strides().size()); ++i)
1321 {
1322 result += temp[i] * m_e.strides()[i - newaxis_count_before<S...>(i)];
1323 }
1324 for (; i < sizeof...(S); ++i)
1325 {
1326 result += temp[i];
1327 }
1328 return static_cast<std::size_t>(result) + m_e.data_offset();
1329 }
1330
1334 template <class CT, class... S>
1335 template <class T>
1336 inline std::size_t xview<CT, S...>::data_offset() const noexcept
1337 requires(has_data_interface_concept<T> and strided_view_concept<CT, S...>)
1338 {
1339 if (!m_strides_computed)
1340 {
1341 compute_strides(std::integral_constant<bool, has_trivial_strides>{});
1342 m_strides_computed = true;
1343 }
1344 return m_data_offset;
1345 }
1346
1348
1349 template <class CT, class... S>
1350 inline auto xview<CT, S...>::underlying_size(size_type dim) const -> size_type
1351 {
1352 return m_e.shape()[dim];
1353 }
1354
1355 template <class CT, class... S>
1356 inline auto xview<CT, S...>::operator&() & -> xtl::xclosure_pointer<self_type&>
1357 {
1358 return xtl::closure_pointer(*this);
1359 }
1360
1361 template <class CT, class... S>
1362 inline auto xview<CT, S...>::operator&() const& -> xtl::xclosure_pointer<const self_type&>
1363 {
1364 return xtl::closure_pointer(*this);
1365 }
1366
1367 template <class CT, class... S>
1368 inline auto xview<CT, S...>::operator&() && -> xtl::xclosure_pointer<self_type>
1369 {
1370 return xtl::closure_pointer(std::move(*this));
1371 }
1372
1377
1383 template <class CT, class... S>
1384 template <class ST>
1385 inline bool xview<CT, S...>::broadcast_shape(ST& shape, bool) const
1386 {
1387 return xt::broadcast_shape(m_shape, shape);
1388 }
1389
1395 template <class CT, class... S>
1396 template <class ST>
1397 inline bool xview<CT, S...>::has_linear_assign(const ST& str) const
1398 {
1399 if constexpr (is_strided_view)
1400 {
1401 return str.size() == strides().size() && std::equal(str.cbegin(), str.cend(), strides().begin());
1402 }
1403 else
1404 {
1405 return false;
1406 }
1407 }
1408
1410
1411 template <class CT, class... S>
1412 template <class It>
1413 inline It xview<CT, S...>::data_xbegin_impl(It begin) const noexcept
1414 {
1415 return begin + data_offset();
1416 }
1417
1418 template <class CT, class... S>
1419 template <class It>
1420 inline It xview<CT, S...>::data_xend_impl(It begin, layout_type l, size_type offset) const noexcept
1421 {
1422 return strided_data_end(*this, begin, l, offset);
1423 }
1424
1425 template <class CT, class... S>
1426 inline auto xview<CT, S...>::data_xbegin() noexcept -> container_iterator
1427 {
1428 return data_xbegin_impl(data());
1429 }
1430
1431 template <class CT, class... S>
1432 inline auto xview<CT, S...>::data_xbegin() const noexcept -> const_container_iterator
1433 {
1434 return data_xbegin_impl(data());
1435 }
1436
1437 template <class CT, class... S>
1438 inline auto xview<CT, S...>::data_xend(layout_type l, size_type offset) noexcept -> container_iterator
1439 {
1440 return data_xend_impl(data() + data_offset(), l, offset);
1441 }
1442
1443 template <class CT, class... S>
1444 inline auto xview<CT, S...>::data_xend(layout_type l, size_type offset) const noexcept
1445 -> const_container_iterator
1446 {
1447 return data_xend_impl(data() + data_offset(), l, offset);
1448 }
1449
1450 // Assign to operator enabled for contigous views
1451 template <class CT, class... S>
1452 template <class E, class T>
1453 void xview<CT, S...>::assign_to(xexpression<E>& e, bool force_resize) const
1455 {
1456 auto& de = e.derived_cast();
1457 de.resize(shape(), force_resize);
1458 std::copy(data() + data_offset(), data() + data_offset() + de.size(), de.template begin<static_layout>());
1459 }
1460
1461 template <class CT, class... S>
1462 template <class E, std::size_t... I>
1463 inline auto xview<CT, S...>::build_view_impl(E&& e, std::index_sequence<I...>) const -> rebind_t<E>
1464 {
1465 return rebind_t<E>(std::forward<E>(e), std::get<I>(m_slices)...);
1466 }
1467
1468 template <class CT, class... S>
1469 template <class E>
1470 inline auto xview<CT, S...>::build_view(E&& e) const -> rebind_t<E>
1471 {
1472 return build_view_impl(std::forward<E>(e), std::make_index_sequence<sizeof...(S)>());
1473 }
1474
1475 template <class CT, class... S>
1476 template <class align, class simd, class T>
1477 inline auto xview<CT, S...>::store_simd(size_type i, const simd& e) -> void
1479 {
1480 return m_e.template store_simd<xt_simd::unaligned_mode>(data_offset() + i, e);
1481 }
1482
1483 template <class CT, class... S>
1484 template <class align, class requested_type, std::size_t N, class T>
1485 inline auto xview<CT, S...>::load_simd(size_type i) const -> simd_return_type<requested_type>
1487 {
1488 return m_e.template load_simd<xt_simd::unaligned_mode, requested_type>(data_offset() + i);
1489 }
1490
1491 template <class CT, class... S>
1492 template <class T>
1493 inline auto xview<CT, S...>::data_element(size_type i) -> reference
1495 {
1496 return m_e.data_element(data_offset() + i);
1497 }
1498
1499 template <class CT, class... S>
1500 template <class T>
1501 inline auto xview<CT, S...>::data_element(size_type i) const -> const_reference
1503 {
1504 return m_e.data_element(data_offset() + i);
1505 }
1506
1507 template <class CT, class... S>
1508 template <class T>
1509 inline auto xview<CT, S...>::flat(size_type i) -> reference
1511 {
1512 XTENSOR_ASSERT(is_contiguous());
1513 return m_e.flat(data_offset() + i);
1514 }
1515
1516 template <class CT, class... S>
1517 template <class T>
1518 inline auto xview<CT, S...>::flat(size_type i) const -> const_reference
1520 {
1521 XTENSOR_ASSERT(is_contiguous());
1522 return m_e.flat(data_offset() + i);
1523 }
1524
1525 template <class CT, class... S>
1526 template <class... Args>
1527 inline auto xview<CT, S...>::make_index_sequence(Args...) const noexcept
1528 {
1529 return std::make_index_sequence<
1530 (sizeof...(Args) + integral_count<S...>() > newaxis_count<S...>()
1531 ? sizeof...(Args) + integral_count<S...>() - newaxis_count<S...>()
1532 : 0)>();
1533 }
1534
1535 template <class CT, class... S>
1536 template <std::size_t... I>
1537 inline auto xview<CT, S...>::compute_strides_impl(std::index_sequence<I...>) const noexcept
1538 {
1539 std::size_t original_dim = m_e.dimension();
1540 return std::array<std::ptrdiff_t, sizeof...(I)>(
1541 {(static_cast<std::ptrdiff_t>(xt::step_size(std::get<integral_skip<S...>(I)>(m_slices), 1))
1542 * ((integral_skip<S...>(I) - newaxis_count_before<S...>(integral_skip<S...>(I))) < original_dim
1543 ? m_e.strides()[integral_skip<S...>(I) - newaxis_count_before<S...>(integral_skip<S...>(I))]
1544 : 1))...}
1545 );
1546 }
1547
1548 template <class CT, class... S>
1549 inline void xview<CT, S...>::compute_strides(std::false_type) const
1550 {
1551 m_strides = xtl::make_sequence<inner_strides_type>(this->dimension(), 0);
1552 m_backstrides = xtl::make_sequence<inner_strides_type>(this->dimension(), 0);
1553
1554 constexpr std::size_t n_strides = sizeof...(S) - integral_count<S...>();
1555
1556 auto slice_strides = compute_strides_impl(std::make_index_sequence<n_strides>());
1557
1558 for (std::size_t i = 0; i < n_strides; ++i)
1559 {
1560 m_strides[i] = slice_strides[i];
1561 // adapt strides for shape[i] == 1 to make consistent with rest of xtensor
1562 detail::adapt_strides(shape(), m_strides, &m_backstrides, i);
1563 }
1564 for (std::size_t i = n_strides; i < this->dimension(); ++i)
1565 {
1566 m_strides[i] = m_e.strides()[i + integral_count<S...>() - newaxis_count<S...>()];
1567 detail::adapt_strides(shape(), m_strides, &m_backstrides, i);
1568 }
1569
1570 m_data_offset = data_offset_impl(std::make_index_sequence<sizeof...(S)>());
1571 }
1572
1573 template <class CT, class... S>
1574 inline void xview<CT, S...>::compute_strides(std::true_type) const
1575 {
1576 }
1577
1578 template <class CT, class... S>
1579 inline auto xview<CT, S...>::access() -> reference
1580 {
1581 return access_impl(make_index_sequence());
1582 }
1583
1584 template <class CT, class... S>
1585 template <class Arg, class... Args>
1586 inline auto xview<CT, S...>::access(Arg arg, Args... args) -> reference
1587 {
1588 if (sizeof...(Args) >= this->dimension())
1589 {
1590 return access(args...);
1591 }
1592 return access_impl(make_index_sequence(arg, args...), arg, args...);
1593 }
1594
1595 template <class CT, class... S>
1596 inline auto xview<CT, S...>::access() const -> const_reference
1597 {
1598 return access_impl(make_index_sequence());
1599 }
1600
1601 template <class CT, class... S>
1602 template <class Arg, class... Args>
1603 inline auto xview<CT, S...>::access(Arg arg, Args... args) const -> const_reference
1604 {
1605 if (sizeof...(Args) >= this->dimension())
1606 {
1607 return access(args...);
1608 }
1609 return access_impl(make_index_sequence(arg, args...), arg, args...);
1610 }
1611
1612 template <class CT, class... S>
1613 template <typename std::decay_t<CT>::size_type... I, class... Args>
1614 inline auto xview<CT, S...>::unchecked_impl(std::index_sequence<I...>, Args... args) -> reference
1615 {
1616 return m_e.unchecked(index<I>(args...)...);
1617 }
1618
1619 template <class CT, class... S>
1620 template <typename std::decay_t<CT>::size_type... I, class... Args>
1621 inline auto xview<CT, S...>::unchecked_impl(std::index_sequence<I...>, Args... args) const
1622 -> const_reference
1623 {
1624 return m_e.unchecked(index<I>(args...)...);
1625 }
1626
1627 template <class CT, class... S>
1628 template <typename std::decay_t<CT>::size_type... I, class... Args>
1629 inline auto xview<CT, S...>::access_impl(std::index_sequence<I...>, Args... args) -> reference
1630 {
1631 return m_e(index<I>(args...)...);
1632 }
1633
1634 template <class CT, class... S>
1635 template <typename std::decay_t<CT>::size_type... I, class... Args>
1636 inline auto xview<CT, S...>::access_impl(std::index_sequence<I...>, Args... args) const -> const_reference
1637 {
1638 return m_e(index<I>(args...)...);
1639 }
1640
1641 template <class CT, class... S>
1642 template <typename std::decay_t<CT>::size_type I, class... Args>
1643 inline auto xview<CT, S...>::index(Args... args) const -> size_type
1644 {
1645 if constexpr (lesser_condition<I>::value)
1646 {
1647 return sliced_access<I - integral_count_before<S...>(I) + newaxis_count_before<S...>(I + 1)>(
1648 std::get<I + newaxis_count_before<S...>(I + 1)>(m_slices),
1649 args...
1650 );
1651 }
1652 else
1653 {
1654 return argument<I - integral_count<S...>() + newaxis_count<S...>()>(args...);
1655 }
1656 }
1657
1658 template <class CT, class... S>
1659 template <typename std::decay_t<CT>::size_type I, class T>
1660 inline auto xview<CT, S...>::sliced_access(const xslice<T>& slice) const -> size_type
1661 {
1662 return static_cast<size_type>(slice.derived_cast()(0));
1663 }
1664
1665 template <class CT, class... S>
1666 template <typename std::decay_t<CT>::size_type I, class T, class Arg, class... Args>
1667 inline auto xview<CT, S...>::sliced_access(const xslice<T>& slice, Arg arg, Args... args) const -> size_type
1668 {
1669 using ST = typename T::size_type;
1670 return static_cast<size_type>(
1671 slice.derived_cast()(argument<I>(static_cast<ST>(arg), static_cast<ST>(args)...))
1672 );
1673 }
1674
1675 template <class CT, class... S>
1676 template <typename std::decay_t<CT>::size_type I, class T, class... Args>
1677 inline auto xview<CT, S...>::sliced_access(const T& squeeze, Args...) const -> disable_xslice<T, size_type>
1678 {
1679 return static_cast<size_type>(squeeze);
1680 }
1681
1682 template <class CT, class... S>
1683 template <class It>
1684 inline auto xview<CT, S...>::make_index(It first, It last) const -> base_index_type
1685 {
1686 auto index = xtl::make_sequence<base_index_type>(m_e.dimension(), 0);
1687 using diff_type = typename std::iterator_traits<It>::difference_type;
1688 using ivalue_type = typename base_index_type::value_type;
1689 auto func1 = [&first](const auto& s) noexcept
1690 {
1691 return get_slice_value(s, first);
1692 };
1693 auto func2 = [](const auto& s) noexcept
1694 {
1695 return xt::value(s, 0);
1696 };
1697
1698 auto s = static_cast<diff_type>(
1699 (std::min)(static_cast<size_type>(std::distance(first, last)), this->dimension())
1700 );
1701 auto first_copy = last - s;
1702 for (size_type i = 0; i != m_e.dimension(); ++i)
1703 {
1704 size_type k = newaxis_skip<S...>(i);
1705
1706 // need to advance captured `first`
1707 first = first_copy;
1708 std::advance(first, static_cast<diff_type>(k - xt::integral_count_before<S...>(i)));
1709
1710 if (first < last)
1711 {
1712 index[i] = k < sizeof...(S) ? apply<size_type>(k, func1, m_slices)
1713 : static_cast<ivalue_type>(*first);
1714 }
1715 else
1716 {
1717 index[i] = k < sizeof...(S) ? apply<size_type>(k, func2, m_slices) : ivalue_type(0);
1718 }
1719 }
1720 return index;
1721 }
1722
1723 template <class CT, class... S>
1724 inline auto xview<CT, S...>::compute_shape(std::true_type) const -> inner_shape_type
1725 {
1726 return inner_shape_type(m_e.shape());
1727 }
1728
1729 template <class CT, class... S>
1730 inline auto xview<CT, S...>::compute_shape(std::false_type) const -> inner_shape_type
1731 {
1732 std::size_t dim = m_e.dimension() - integral_count<S...>() + newaxis_count<S...>();
1733 auto shape = xtl::make_sequence<inner_shape_type>(dim, 0);
1734 auto func = [](const auto& s) noexcept
1735 {
1736 return get_size(s);
1737 };
1738 for (size_type i = 0; i != dim; ++i)
1739 {
1740 size_type index = integral_skip<S...>(i);
1741 shape[i] = index < sizeof...(S) ? apply<size_type>(index, func, m_slices)
1742 : m_e.shape()[index - newaxis_count_before<S...>(index)];
1743 }
1744 return shape;
1745 }
1746
1747 namespace xview_detail
1748 {
1749 template <class V, class T>
1750 inline void run_assign_temporary_impl(V& v, const T& t, std::true_type /* enable strided assign */)
1751 {
1752 strided_loop_assigner<true>::run(v, t);
1753 }
1754
1755 template <class V, class T>
1756 inline void
1757 run_assign_temporary_impl(V& v, const T& t, std::false_type /* fallback to iterator assign */)
1758 {
1759 std::copy(t.cbegin(), t.cend(), v.begin());
1760 }
1761 }
1762
1763 template <class CT, class... S>
1764 inline void xview<CT, S...>::assign_temporary_impl(temporary_type&& tmp)
1765 {
1766 constexpr bool fast_assign = detail::is_strided_view<xexpression_type, S...>::value
1767 && xassign_traits<xview<CT, S...>, temporary_type>::simd_strided_assign();
1768 xview_detail::run_assign_temporary_impl(*this, tmp, std::integral_constant<bool, fast_assign>{});
1769 }
1770
1771 namespace detail
1772 {
1773 template <class E, class... S>
1774 inline std::size_t get_underlying_shape_index(std::size_t I)
1775 {
1776 return I - newaxis_count_before<get_slice_type<E, S>...>(I);
1777 }
1778
1779 template <class... S>
1780 struct check_slice;
1781
1782 template <>
1783 struct check_slice<>
1784 {
1785 using type = void_t<>;
1786 };
1787
1788 template <class S, class... SL>
1789 struct check_slice<S, SL...>
1790 {
1791 static_assert(!std::is_same<S, xellipsis_tag>::value, "ellipsis not supported vith xview");
1792 using type = typename check_slice<SL...>::type;
1793 };
1794
1795 template <class E, std::size_t... I, class... S>
1796 inline auto make_view_impl(E&& e, std::index_sequence<I...>, S&&... slices)
1797 {
1798 // Checks that no ellipsis slice is used
1799 using view_type = xview<xtl::closure_type_t<E>, get_slice_type<std::decay_t<E>, S>...>;
1800 return view_type(
1801 std::forward<E>(e),
1802 get_slice_implementation(
1803 e,
1804 std::forward<S>(slices),
1805 get_underlying_shape_index<std::decay_t<E>, S...>(I)
1806 )...
1807 );
1808 }
1809 }
1810
1820 template <class E, class... S>
1821 inline auto view(E&& e, S&&... slices)
1822 {
1823 return detail::make_view_impl(
1824 std::forward<E>(e),
1825 std::make_index_sequence<sizeof...(S)>(),
1826 std::forward<S>(slices)...
1827 );
1828 }
1829
1830 namespace detail
1831 {
1832 class row_impl
1833 {
1834 public:
1835
1836 template <class E>
1837 inline static auto make(E&& e, const std::ptrdiff_t index)
1838 {
1839 const auto shape = e.shape();
1840 check_dimension(shape);
1841 return view(e, index, xt::all());
1842 }
1843
1844 private:
1845
1846 template <class S>
1847 inline static void check_dimension(const S& shape)
1848 {
1849 if (shape.size() != 2)
1850 {
1851 XTENSOR_THROW(
1852 std::invalid_argument,
1853 "A row can only be accessed on an expression with exact two dimensions"
1854 );
1855 }
1856 }
1857
1858 template <class T, std::size_t N>
1859 inline static void check_dimension(const std::array<T, N>&)
1860 {
1861 static_assert(N == 2, "A row can only be accessed on an expression with exact two dimensions");
1862 }
1863 };
1864
1865 class column_impl
1866 {
1867 public:
1868
1869 template <class E>
1870 inline static auto make(E&& e, const std::ptrdiff_t index)
1871 {
1872 const auto shape = e.shape();
1873 check_dimension(shape);
1874 return view(e, xt::all(), index);
1875 }
1876
1877 private:
1878
1879 template <class S>
1880 inline static void check_dimension(const S& shape)
1881 {
1882 if (shape.size() != 2)
1883 {
1884 XTENSOR_THROW(
1885 std::invalid_argument,
1886 "A column can only be accessed on an expression with exact two dimensions"
1887 );
1888 }
1889 }
1890
1891 template <class T, std::size_t N>
1892 inline static void check_dimension(const std::array<T, N>&)
1893 {
1894 static_assert(N == 2, "A column can only be accessed on an expression with exact two dimensions");
1895 }
1896 };
1897 }
1898
1908 template <class E>
1909 inline auto row(E&& e, std::ptrdiff_t index)
1910 {
1911 return detail::row_impl::make(e, index);
1912 }
1913
1923 template <class E>
1924 inline auto col(E&& e, std::ptrdiff_t index)
1925 {
1926 return detail::column_impl::make(e, index);
1927 }
1928
1929 /***************
1930 * stepper api *
1931 ***************/
1932
1933 template <class CT, class... S>
1934 template <class ST>
1935 inline auto xview<CT, S...>::stepper_begin(const ST& shape) -> stepper
1936 {
1937 const size_type offset = shape.size() - this->dimension();
1938 if constexpr (is_strided_view)
1939 {
1940 return stepper(this, data_xbegin(), offset);
1941 }
1942 else
1943 {
1944 return stepper(this, m_e.stepper_begin(m_e.shape()), offset);
1945 }
1946 }
1947
1948 template <class CT, class... S>
1949 template <class ST>
1950 inline auto xview<CT, S...>::stepper_end(const ST& shape, layout_type l) -> stepper
1951 {
1952 const size_type offset = shape.size() - this->dimension();
1953 if constexpr (is_strided_view)
1954 {
1955 return stepper(this, data_xend(l, offset), offset);
1956 }
1957 else
1958 {
1959 return stepper(this, m_e.stepper_end(m_e.shape(), l), offset, true, l);
1960 }
1961 }
1962
1963 template <class CT, class... S>
1964 template <class ST>
1965 inline auto xview<CT, S...>::stepper_begin(const ST& shape) const -> const_stepper
1966 {
1967 const size_type offset = shape.size() - this->dimension();
1968 if constexpr (is_strided_view)
1969 {
1970 return const_stepper(this, data_xbegin(), offset);
1971 }
1972 else
1973 {
1974 const xexpression_type& e = m_e;
1975 return const_stepper(this, e.stepper_begin(m_e.shape()), offset);
1976 }
1977 }
1978
1979 template <class CT, class... S>
1980 template <class ST>
1981 inline auto xview<CT, S...>::stepper_end(const ST& shape, layout_type l) const -> const_stepper
1982 {
1983 const size_type offset = shape.size() - this->dimension();
1984 if constexpr (is_strided_view)
1985 {
1986 return const_stepper(this, data_xend(l, offset), offset);
1987 }
1988 else
1989 {
1990 const xexpression_type& e = m_e;
1991 return const_stepper(this, e.stepper_end(m_e.shape(), l), offset, true, l);
1992 }
1993 }
1994
1995 /********************************
1996 * xview_stepper implementation *
1997 ********************************/
1998
1999 template <bool is_const, class CT, class... S>
2000 inline xview_stepper<is_const, CT, S...>::xview_stepper(
2001 view_type* view,
2002 substepper_type it,
2003 size_type offset,
2004 bool end,
2005 layout_type l
2006 )
2007 : p_view(view)
2008 , m_it(it)
2009 , m_offset(offset)
2010 {
2011 if (!end)
2012 {
2013 std::fill(m_index_keeper.begin(), m_index_keeper.end(), 0);
2014 auto func = [](const auto& s) noexcept
2015 {
2016 return xt::value(s, 0);
2017 };
2018 for (size_type i = 0; i < sizeof...(S); ++i)
2019 {
2020 if (!is_newaxis_slice(i))
2021 {
2022 size_type s = apply<size_type>(i, func, p_view->slices());
2023 size_type index = i - newaxis_count_before<S...>(i);
2024 m_it.step(index, s);
2025 }
2026 }
2027 }
2028 else
2029 {
2030 to_end_impl(l);
2031 }
2032 }
2033
2034 template <bool is_const, class CT, class... S>
2035 inline auto xview_stepper<is_const, CT, S...>::operator*() const -> reference
2036 {
2037 return *m_it;
2038 }
2039
2040 template <bool is_const, class CT, class... S>
2041 inline void xview_stepper<is_const, CT, S...>::step(size_type dim)
2042 {
2043 auto func = [this](size_type index, size_type offset)
2044 {
2045 m_it.step(index, offset);
2046 };
2047 common_step_forward(dim, func);
2048 }
2049
2050 template <bool is_const, class CT, class... S>
2051 inline void xview_stepper<is_const, CT, S...>::step_back(size_type dim)
2052 {
2053 auto func = [this](size_type index, size_type offset)
2054 {
2055 m_it.step_back(index, offset);
2056 };
2057 common_step_backward(dim, func);
2058 }
2059
2060 template <bool is_const, class CT, class... S>
2061 inline void xview_stepper<is_const, CT, S...>::step(size_type dim, size_type n)
2062 {
2063 auto func = [this](size_type index, size_type offset)
2064 {
2065 m_it.step(index, offset);
2066 };
2067 common_step_forward(dim, n, func);
2068 }
2069
2070 template <bool is_const, class CT, class... S>
2071 inline void xview_stepper<is_const, CT, S...>::step_back(size_type dim, size_type n)
2072 {
2073 auto func = [this](size_type index, size_type offset)
2074 {
2075 m_it.step_back(index, offset);
2076 };
2077 common_step_backward(dim, n, func);
2078 }
2079
2080 template <bool is_const, class CT, class... S>
2081 inline void xview_stepper<is_const, CT, S...>::reset(size_type dim)
2082 {
2083 auto func = [this](size_type index, size_type offset)
2084 {
2085 m_it.step_back(index, offset);
2086 };
2087 common_reset(dim, func, false);
2088 }
2089
2090 template <bool is_const, class CT, class... S>
2091 inline void xview_stepper<is_const, CT, S...>::reset_back(size_type dim)
2092 {
2093 auto func = [this](size_type index, size_type offset)
2094 {
2095 m_it.step(index, offset);
2096 };
2097 common_reset(dim, func, true);
2098 }
2099
2100 template <bool is_const, class CT, class... S>
2101 inline void xview_stepper<is_const, CT, S...>::to_begin()
2102 {
2103 std::fill(m_index_keeper.begin(), m_index_keeper.end(), 0);
2104 m_it.to_begin();
2105 }
2106
2107 template <bool is_const, class CT, class... S>
2108 inline void xview_stepper<is_const, CT, S...>::to_end(layout_type l)
2109 {
2110 m_it.to_end(l);
2111 to_end_impl(l);
2112 }
2113
2114 template <bool is_const, class CT, class... S>
2115 inline bool xview_stepper<is_const, CT, S...>::is_newaxis_slice(size_type index) const noexcept
2116 {
2117 // A bit tricky but avoids a lot of template instantiations
2118 return newaxis_count_before<S...>(index + 1) != newaxis_count_before<S...>(index);
2119 }
2120
2121 template <bool is_const, class CT, class... S>
2122 inline void xview_stepper<is_const, CT, S...>::to_end_impl(layout_type l)
2123 {
2124 auto func = [](const auto& s) noexcept
2125 {
2126 return xt::value(s, get_size(s) - 1);
2127 };
2128 auto size_func = [](const auto& s) noexcept
2129 {
2130 return get_size(s);
2131 };
2132
2133 for (size_type i = 0; i < sizeof...(S); ++i)
2134 {
2135 if (!is_newaxis_slice(i))
2136 {
2137 size_type s = apply<size_type>(i, func, p_view->slices());
2138 size_type ix = apply<size_type>(i, size_func, p_view->slices());
2139 m_index_keeper[i] = ix - size_type(1);
2140 size_type index = i - newaxis_count_before<S...>(i);
2141 s = p_view->underlying_size(index) - 1 - s;
2142 m_it.step_back(index, s);
2143 }
2144 }
2145 if (l == layout_type::row_major)
2146 {
2147 for (size_type i = sizeof...(S); i > 0; --i)
2148 {
2149 if (!is_newaxis_slice(i - 1))
2150 {
2151 m_index_keeper[i - 1]++;
2152 break;
2153 }
2154 }
2155 }
2156 else if (l == layout_type::column_major)
2157 {
2158 for (size_type i = 0; i < sizeof...(S); ++i)
2159 {
2160 if (!is_newaxis_slice(i))
2161 {
2162 m_index_keeper[i]++;
2163 break;
2164 }
2165 }
2166 }
2167 else
2168 {
2169 XTENSOR_THROW(std::runtime_error, "Iteration only allowed in row or column major.");
2170 }
2171 }
2172
2173 template <bool is_const, class CT, class... S>
2174 template <class F>
2175 void xview_stepper<is_const, CT, S...>::common_step_forward(size_type dim, F f)
2176 {
2177 if (dim >= m_offset)
2178 {
2179 auto func = [&dim, this](const auto& s) noexcept
2180 {
2181 return step_size(s, this->m_index_keeper[dim]++, 1);
2182 };
2183 size_type index = integral_skip<S...>(dim);
2184 if (!is_newaxis_slice(index))
2185 {
2186 size_type step_size = index < sizeof...(S) ? apply<size_type>(index, func, p_view->slices())
2187 : 1;
2188 index -= newaxis_count_before<S...>(index);
2189 f(index, step_size);
2190 }
2191 }
2192 }
2193
2194 template <bool is_const, class CT, class... S>
2195 template <class F>
2196 void xview_stepper<is_const, CT, S...>::common_step_forward(size_type dim, size_type n, F f)
2197 {
2198 if (dim >= m_offset)
2199 {
2200 auto func = [&dim, &n, this](const auto& s) noexcept
2201 {
2202 auto st_size = step_size(s, this->m_index_keeper[dim], n);
2203 this->m_index_keeper[dim] += n;
2204 return size_type(st_size);
2205 };
2206
2207 size_type index = integral_skip<S...>(dim);
2208 if (!is_newaxis_slice(index))
2209 {
2210 size_type step_size = index < sizeof...(S) ? apply<size_type>(index, func, p_view->slices())
2211 : n;
2212 index -= newaxis_count_before<S...>(index);
2213 f(index, step_size);
2214 }
2215 }
2216 }
2217
2218 template <bool is_const, class CT, class... S>
2219 template <class F>
2220 void xview_stepper<is_const, CT, S...>::common_step_backward(size_type dim, F f)
2221 {
2222 if (dim >= m_offset)
2223 {
2224 auto func = [&dim, this](const auto& s) noexcept
2225 {
2226 this->m_index_keeper[dim]--;
2227 return step_size(s, this->m_index_keeper[dim], 1);
2228 };
2229 size_type index = integral_skip<S...>(dim);
2230 if (!is_newaxis_slice(index))
2231 {
2232 size_type step_size = index < sizeof...(S) ? apply<size_type>(index, func, p_view->slices())
2233 : 1;
2234 index -= newaxis_count_before<S...>(index);
2235 f(index, step_size);
2236 }
2237 }
2238 }
2239
2240 template <bool is_const, class CT, class... S>
2241 template <class F>
2242 void xview_stepper<is_const, CT, S...>::common_step_backward(size_type dim, size_type n, F f)
2243 {
2244 if (dim >= m_offset)
2245 {
2246 auto func = [&dim, &n, this](const auto& s) noexcept
2247 {
2248 this->m_index_keeper[dim] -= n;
2249 return step_size(s, this->m_index_keeper[dim], n);
2250 };
2251
2252 size_type index = integral_skip<S...>(dim);
2253 if (!is_newaxis_slice(index))
2254 {
2255 size_type step_size = index < sizeof...(S) ? apply<size_type>(index, func, p_view->slices())
2256 : n;
2257 index -= newaxis_count_before<S...>(index);
2258 f(index, step_size);
2259 }
2260 }
2261 }
2262
2263 template <bool is_const, class CT, class... S>
2264 template <class F>
2265 void xview_stepper<is_const, CT, S...>::common_reset(size_type dim, F f, bool backwards)
2266 {
2267 auto size_func = [](const auto& s) noexcept
2268 {
2269 return get_size(s);
2270 };
2271 auto end_func = [](const auto& s) noexcept
2272 {
2273 return xt::value(s, get_size(s) - 1) - xt::value(s, 0);
2274 };
2275
2276 size_type index = integral_skip<S...>(dim);
2277 if (!is_newaxis_slice(index))
2278 {
2279 if (dim < m_index_keeper.size())
2280 {
2281 size_type size = index < sizeof...(S) ? apply<size_type>(index, size_func, p_view->slices())
2282 : p_view->shape()[dim];
2283 m_index_keeper[dim] = backwards ? size - 1 : 0;
2284 }
2285
2286 size_type reset_n = index < sizeof...(S) ? apply<size_type>(index, end_func, p_view->slices())
2287 : p_view->shape()[dim] - 1;
2288 index -= newaxis_count_before<S...>(index);
2289 f(index, reset_n);
2290 }
2291 }
2292}
2293
2294#endif
Fixed shape implementation for compile time defined arrays.
Base class for xexpressions.
Base class for multidimensional iterable expressions.
Multidimensional view with tensor semantic.
Definition xview.hpp:365
xview(CTA &&e, FSL &&first_slice, SL &&... slices) noexcept
Constructs a view on the specified xexpression.
Definition xview.hpp:880
const slice_type & slices() const noexcept
bool has_linear_assign(const ST &strides) const
const inner_shape_type & shape() const noexcept
Returns the shape of the view.
Definition xview.hpp:955
bool broadcast_shape(ST &shape, bool reuse_cache=false) const
xexpression_type & expression() noexcept
void fill(const T &value)
std::size_t data_offset() const noexcept
layout_type layout() const noexcept
bool all(E &&e)
Any.
auto arg(E &&e) noexcept
Calculates the phase angle (in radians) elementwise for the complex numbers in e.
Definition xcomplex.hpp:221
auto squeeze(E &&e)
Returns a squeeze view of the given expression.
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:566
auto strides(const E &e, stride_type type=stride_type::normal) noexcept
Get strides of an object.
Definition xstrides.hpp:248
standard mathematical functions for xexpressions
auto all() noexcept
Returns a slice representing a full dimension, to be used as an argument of view function.
Definition xslice.hpp:234
auto row(E &&e, std::ptrdiff_t index)
Constructs and returns a row (sliced view) on the specified expression.
Definition xview.hpp:1909
layout_type
Definition xlayout.hpp:24
auto col(E &&e, std::ptrdiff_t index)
Constructs and returns a column (sliced view) on the specified expression.
Definition xview.hpp:1924
auto view(E &&e, S &&... slices)
Constructs and returns a view on the specified xexpression.
Definition xview.hpp:1821