xtensor
Loading...
Searching...
No Matches
xdynamic_view.hpp
1/***************************************************************************
2 * Copyright (c) 2016, Johan Mabille, Sylvain Corlay and Wolf Vollprecht *
3 * *
4 * Distributed under the terms of the BSD 3-Clause License. *
5 * *
6 * The full license is in the file LICENSE, distributed with this software. *
7 ****************************************************************************/
8
9#ifndef XTENSOR_DYNAMIC_VIEW_HPP
10#define XTENSOR_DYNAMIC_VIEW_HPP
11
12#include <xtl/xsequence.hpp>
13#include <xtl/xvariant.hpp>
14
15#include "xexpression.hpp"
16#include "xiterable.hpp"
17#include "xlayout.hpp"
18#include "xsemantic.hpp"
19#include "xstrided_view_base.hpp"
20
21namespace xt
22{
23
24 template <class CT, class S, layout_type L, class FST>
25 class xdynamic_view;
26
27 template <class CT, class S, layout_type L, class FST>
29 {
30 using xexpression_type = std::decay_t<CT>;
31 using undecay_expression = CT;
33 using const_reference = typename xexpression_type::const_reference;
34 using size_type = typename xexpression_type::size_type;
35 using shape_type = std::decay_t<S>;
36 using undecay_shape = S;
37 using storage_getter = FST;
38 using inner_storage_type = typename storage_getter::type;
39 using temporary_type = xarray<std::decay_t<typename xexpression_type::value_type>, xexpression_type::static_layout>;
40 static constexpr layout_type layout = L;
41 };
42
43 template <class CT, class S, layout_type L, class FST>
45 {
46 using inner_shape_type = S;
49
50#if defined(__GNUC__) && !defined(__clang__) && __GNUC__ == 8
51 static constexpr auto
53 static constexpr auto
55#endif
56
57 // TODO: implement efficient stepper specific to the dynamic_view
60 };
61
62 /****************************
63 * xdynamic_view extensions *
64 ****************************/
65
66 namespace extension
67 {
68 template <class Tag, class CT, class S, layout_type L, class FST>
70
71 template <class CT, class S, layout_type L, class FST>
76
77 template <class CT, class S, layout_type L, class FST>
78 struct xdynamic_view_base : xdynamic_view_base_impl<xexpression_tag_t<CT>, CT, S, L, FST>
79 {
80 };
81
82 template <class CT, class S, layout_type L, class FST>
83 using xdynamic_view_base_t = typename xdynamic_view_base<CT, S, L, FST>::type;
84 }
85
86 /*****************
87 * xdynamic_view *
88 *****************/
89
90 namespace detail
91 {
92 template <class T>
93 class xfake_slice;
94 }
95
96 template <class CT, class S, layout_type L = layout_type::dynamic, class FST = detail::flat_storage_getter<CT, XTENSOR_DEFAULT_TRAVERSAL>>
97 class xdynamic_view : public xview_semantic<xdynamic_view<CT, S, L, FST>>,
98 public xiterable<xdynamic_view<CT, S, L, FST>>,
99 public extension::xdynamic_view_base_t<CT, S, L, FST>,
100 private xstrided_view_base<xdynamic_view<CT, S, L, FST>>
101 {
102 public:
103
107 using extension_base = extension::xdynamic_view_base_t<CT, S, L, FST>;
108 using expression_tag = typename extension_base::expression_tag;
109
110 using xexpression_type = typename base_type::xexpression_type;
111 using base_type::is_const;
112
113 using value_type = typename base_type::value_type;
114 using reference = typename base_type::reference;
115 using const_reference = typename base_type::const_reference;
116 using pointer = typename base_type::pointer;
117 using const_pointer = typename base_type::const_pointer;
118 using size_type = typename base_type::size_type;
119 using difference_type = typename base_type::difference_type;
120
121 using inner_storage_type = typename base_type::inner_storage_type;
122 using storage_type = typename base_type::storage_type;
123
125 using inner_shape_type = typename iterable_base::inner_shape_type;
126 using inner_strides_type = typename base_type::inner_strides_type;
127 using inner_backstrides_type = typename base_type::inner_backstrides_type;
128
129 using shape_type = typename base_type::shape_type;
130 using strides_type = typename base_type::strides_type;
131 using backstrides_type = typename base_type::backstrides_type;
132
133 using stepper = typename iterable_base::stepper;
134 using const_stepper = typename iterable_base::const_stepper;
135
136 using base_type::contiguous_layout;
137 using base_type::static_layout;
138
139 using temporary_type = typename xcontainer_inner_types<self_type>::temporary_type;
141
142 using simd_value_type = typename base_type::simd_value_type;
143 using bool_load_type = typename base_type::bool_load_type;
144
145 using strides_vt = typename strides_type::value_type;
146 using slice_type = xtl::variant<detail::xfake_slice<strides_vt>, xkeep_slice<strides_vt>, xdrop_slice<strides_vt>>;
147 using slice_vector_type = std::vector<slice_type>;
148
149 template <class CTA, class SA>
151 CTA&& e,
152 SA&& shape,
153 get_strides_t<S>&& strides,
154 std::size_t offset,
156 slice_vector_type&& slices,
158 ) noexcept;
159
160 template <class E>
161 self_type& operator=(const xexpression<E>& e);
162
163 template <class E>
164 disable_xexpression<E, self_type>& operator=(const E& e);
165
167 using base_type::is_contiguous;
168 using base_type::layout;
169 using base_type::shape;
170 using base_type::size;
171
172 // Explicitly deleting strides method to avoid compilers complaining
173 // about not being able to call the strides method from xstrided_view_base
174 // private base
175 const inner_strides_type& strides() const noexcept = delete;
176
177 reference operator()();
178 const_reference operator()() const;
179
180 template <class... Args>
181 reference operator()(Args... args);
182
183 template <class... Args>
184 const_reference operator()(Args... args) const;
185
186 template <class... Args>
187 reference unchecked(Args... args);
188
189 template <class... Args>
190 const_reference unchecked(Args... args) const;
191
192 reference flat(size_type index);
193 const_reference flat(size_type index) const;
194
195 using base_type::operator[];
196 using base_type::at;
197 using base_type::back;
198 using base_type::front;
200 using base_type::periodic;
201
202 template <class It>
203 reference element(It first, It last);
204
205 template <class It>
206 const_reference element(It first, It last) const;
207
208 size_type data_offset() const noexcept;
209
210 // Explicitly deleting data methods so has_data_interface results
211 // to false instead of having compilers complaining about not being
212 // able to call the methods from the private base
213 value_type* data() noexcept = delete;
214 const value_type* data() const noexcept = delete;
215
218 using base_type::storage;
219
220 template <class O>
221 bool has_linear_assign(const O& str) const noexcept;
222
223 template <class T>
224 void fill(const T& value);
225
226 template <class ST>
227 stepper stepper_begin(const ST& shape);
228 template <class ST>
229 stepper stepper_end(const ST& shape, layout_type l);
230
231 template <class ST>
232 const_stepper stepper_begin(const ST& shape) const;
233 template <class ST>
234 const_stepper stepper_end(const ST& shape, layout_type l) const;
235
236 using container_iterator = std::
237 conditional_t<is_const, typename storage_type::const_iterator, typename storage_type::iterator>;
238 using const_container_iterator = typename storage_type::const_iterator;
239
240 template <class E>
242
243 template <class E>
244 rebind_t<E> build_view(E&& e) const;
245
246 private:
247
248 using offset_type = typename base_type::offset_type;
249
250 slice_vector_type m_slices;
251 inner_strides_type m_adj_strides;
252
253 container_iterator data_xbegin() noexcept;
254 const_container_iterator data_xbegin() const noexcept;
255 container_iterator data_xend(layout_type l, size_type offset) noexcept;
256 const_container_iterator data_xend(layout_type l, size_type offset) const noexcept;
257
258 template <class It>
259 It data_xbegin_impl(It begin) const noexcept;
260
261 template <class It>
262 It data_xend_impl(It end, layout_type l, size_type offset) const noexcept;
263
264 void assign_temporary_impl(temporary_type&& tmp);
265
266 template <class T, class... Args>
267 offset_type adjust_offset(offset_type offset, T idx, Args... args) const noexcept;
268 offset_type adjust_offset(offset_type offset) const noexcept;
269
270 template <class T, class... Args>
271 offset_type
272 adjust_offset_impl(offset_type offset, size_type idx_offset, T idx, Args... args) const noexcept;
273 offset_type adjust_offset_impl(offset_type offset, size_type idx_offset) const noexcept;
274
275 template <class It>
276 offset_type adjust_element_offset(offset_type offset, It first, It last) const noexcept;
277
278 template <class C>
279 friend class xstepper;
280 friend class xview_semantic<self_type>;
281 friend class xaccessible<self_type>;
282 friend class xconst_accessible<self_type>;
283 };
284
285 /**************************
286 * xdynamic_view builders *
287 **************************/
288
289 template <class T>
290 using xdynamic_slice = xtl::variant<
291 T,
292
296
300
303
304 xrange<T>,
306
309
310 xall_tag,
313
314 using xdynamic_slice_vector = std::vector<xdynamic_slice<std::ptrdiff_t>>;
315
316 template <class E>
317 auto dynamic_view(E&& e, const xdynamic_slice_vector& slices);
318
319 /******************************
320 * xfake_slice implementation *
321 ******************************/
322
323 namespace detail
324 {
325 template <class T>
326 class xfake_slice : public xslice<xfake_slice<T>>
327 {
328 public:
329
330 using size_type = T;
331 using self_type = xfake_slice<T>;
332
333 xfake_slice() = default;
334
335 size_type operator()(size_type /*i*/) const noexcept
336 {
337 return size_type(0);
338 }
339
340 size_type size() const noexcept
341 {
342 return size_type(1);
343 }
344
345 size_type step_size() const noexcept
346 {
347 return size_type(0);
348 }
349
350 size_type step_size(std::size_t /*i*/, std::size_t /*n*/ = 1) const noexcept
351 {
352 return size_type(0);
353 }
354
355 size_type revert_index(std::size_t i) const noexcept
356 {
357 return i;
358 }
359
360 bool contains(size_type /*i*/) const noexcept
361 {
362 return true;
363 }
364
365 bool operator==(const self_type& /*rhs*/) const noexcept
366 {
367 return true;
368 }
369
370 bool operator!=(const self_type& /*rhs*/) const noexcept
371 {
372 return false;
373 }
374 };
375 }
376
377 /********************************
378 * xdynamic_view implementation *
379 ********************************/
380
381 template <class CT, class S, layout_type L, class FST>
382 template <class CTA, class SA>
383 inline xdynamic_view<CT, S, L, FST>::xdynamic_view(
384 CTA&& e,
385 SA&& shape,
386 get_strides_t<S>&& strides,
387 std::size_t offset,
388 layout_type layout,
389 slice_vector_type&& slices,
390 get_strides_t<S>&& adj_strides
391 ) noexcept
392 : base_type(std::forward<CTA>(e), std::forward<SA>(shape), std::move(strides), offset, layout)
393 , m_slices(std::move(slices))
394 , m_adj_strides(std::move(adj_strides))
395 {
396 }
397
398 template <class CT, class S, layout_type L, class FST>
399 template <class E>
400 inline auto xdynamic_view<CT, S, L, FST>::operator=(const xexpression<E>& e) -> self_type&
401 {
402 return semantic_base::operator=(e);
403 }
404
405 template <class CT, class S, layout_type L, class FST>
406 template <class E>
407 inline auto xdynamic_view<CT, S, L, FST>::operator=(const E& e) -> disable_xexpression<E, self_type>&
408 {
409 std::fill(this->begin(), this->end(), e);
410 return *this;
411 }
412
413 template <class CT, class S, layout_type L, class FST>
414 inline auto xdynamic_view<CT, S, L, FST>::operator()() -> reference
415 {
416 return base_type::storage()[data_offset()];
417 }
418
419 template <class CT, class S, layout_type L, class FST>
420 inline auto xdynamic_view<CT, S, L, FST>::operator()() const -> const_reference
421 {
422 return base_type::storage()[data_offset()];
423 }
424
425 template <class CT, class S, layout_type L, class FST>
426 template <class... Args>
427 inline auto xdynamic_view<CT, S, L, FST>::operator()(Args... args) -> reference
428 {
429 XTENSOR_TRY(check_index(base_type::shape(), args...));
430 XTENSOR_CHECK_DIMENSION(base_type::shape(), args...);
431 offset_type offset = base_type::compute_index(args...);
432 offset = adjust_offset(offset, args...);
433 return base_type::storage()[static_cast<size_type>(offset)];
434 }
435
436 template <class CT, class S, layout_type L, class FST>
437 template <class... Args>
438 inline auto xdynamic_view<CT, S, L, FST>::operator()(Args... args) const -> const_reference
439 {
440 XTENSOR_TRY(check_index(base_type::shape(), args...));
441 XTENSOR_CHECK_DIMENSION(base_type::shape(), args...);
442 offset_type offset = base_type::compute_index(args...);
443 offset = adjust_offset(offset, args...);
444 return base_type::storage()[static_cast<size_type>(offset)];
445 }
446
447 template <class CT, class S, layout_type L, class FST>
448 template <class O>
449 inline bool xdynamic_view<CT, S, L, FST>::has_linear_assign(const O&) const noexcept
450 {
451 return false;
452 }
453
454 template <class CT, class S, layout_type L, class FST>
455 template <class... Args>
456 inline auto xdynamic_view<CT, S, L, FST>::unchecked(Args... args) -> reference
457 {
458 offset_type offset = base_type::compute_unchecked_index(args...);
459 offset = adjust_offset(args...);
460 return base_type::storage()[static_cast<size_type>(offset)];
461 }
462
463 template <class CT, class S, layout_type L, class FST>
464 template <class... Args>
465 inline auto xdynamic_view<CT, S, L, FST>::unchecked(Args... args) const -> const_reference
466 {
467 offset_type offset = base_type::compute_unchecked_index(args...);
468 offset = adjust_offset(args...);
469 return base_type::storage()[static_cast<size_type>(offset)];
470 }
471
472 template <class CT, class S, layout_type L, class FST>
473 inline auto xdynamic_view<CT, S, L, FST>::flat(size_type i) -> reference
474 {
475 return base_type::storage()[data_offset() + i];
476 }
477
478 template <class CT, class S, layout_type L, class FST>
479 inline auto xdynamic_view<CT, S, L, FST>::flat(size_type i) const -> const_reference
480 {
481 return base_type::storage()[data_offset() + i];
482 }
483
484 template <class CT, class S, layout_type L, class FST>
485 template <class It>
486 inline auto xdynamic_view<CT, S, L, FST>::element(It first, It last) -> reference
487 {
488 XTENSOR_TRY(check_element_index(base_type::shape(), first, last));
489 offset_type offset = base_type::compute_element_index(first, last);
490 offset = adjust_element_offset(offset, first, last);
491 return base_type::storage()[static_cast<size_type>(offset)];
492 }
493
494 template <class CT, class S, layout_type L, class FST>
495 template <class It>
496 inline auto xdynamic_view<CT, S, L, FST>::element(It first, It last) const -> const_reference
497 {
498 XTENSOR_TRY(check_element_index(base_type::shape(), first, last));
499 offset_type offset = base_type::compute_element_index(first, last);
500 offset = adjust_element_offset(offset, first, last);
501 return base_type::storage()[static_cast<size_type>(offset)];
502 }
503
504 template <class CT, class S, layout_type L, class FST>
505 inline auto xdynamic_view<CT, S, L, FST>::data_offset() const noexcept -> size_type
506 {
507 size_type offset = base_type::data_offset();
508 size_type sl_offset = xtl::visit(
509 [](const auto& sl)
510 {
511 return sl(size_type(0));
512 },
513 m_slices[0]
514 );
515 return offset + sl_offset * m_adj_strides[0];
516 }
517
518 template <class CT, class S, layout_type L, class FST>
519 template <class T>
520 inline void xdynamic_view<CT, S, L, FST>::fill(const T& value)
521 {
522 return std::fill(this->linear_begin(), this->linear_end(), value);
523 }
524
525 template <class CT, class S, layout_type L, class FST>
526 template <class ST>
527 inline auto xdynamic_view<CT, S, L, FST>::stepper_begin(const ST& shape) -> stepper
528 {
529 size_type offset = shape.size() - dimension();
530 return stepper(this, offset);
531 }
532
533 template <class CT, class S, layout_type L, class FST>
534 template <class ST>
535 inline auto xdynamic_view<CT, S, L, FST>::stepper_end(const ST& shape, layout_type /*l*/) -> stepper
536 {
537 size_type offset = shape.size() - dimension();
538 return stepper(this, offset, true);
539 }
540
541 template <class CT, class S, layout_type L, class FST>
542 template <class ST>
543 inline auto xdynamic_view<CT, S, L, FST>::stepper_begin(const ST& shape) const -> const_stepper
544 {
545 size_type offset = shape.size() - dimension();
546 return const_stepper(this, offset);
547 }
548
549 template <class CT, class S, layout_type L, class FST>
550 template <class ST>
551 inline auto xdynamic_view<CT, S, L, FST>::stepper_end(const ST& shape, layout_type /*l*/) const
552 -> const_stepper
553 {
554 size_type offset = shape.size() - dimension();
555 return const_stepper(this, offset, true);
556 }
557
558 template <class CT, class S, layout_type L, class FST>
559 template <class E>
560 inline auto xdynamic_view<CT, S, L, FST>::build_view(E&& e) const -> rebind_t<E>
561 {
562 inner_shape_type sh(this->shape());
563 inner_strides_type str(base_type::strides());
564 slice_vector_type svt(m_slices);
565 inner_strides_type adj_str(m_adj_strides);
566 return rebind_t<E>(
567 std::forward<E>(e),
568 std::move(sh),
569 std::move(str),
570 base_type::data_offset(),
571 this->layout(),
572 std::move(svt),
573 std::move(adj_str)
574 );
575 }
576
577 template <class CT, class S, layout_type L, class FST>
578 inline auto xdynamic_view<CT, S, L, FST>::data_xbegin() noexcept -> container_iterator
579 {
580 return data_xbegin_impl(this->storage().begin());
581 }
582
583 template <class CT, class S, layout_type L, class FST>
584 inline auto xdynamic_view<CT, S, L, FST>::data_xbegin() const noexcept -> const_container_iterator
585 {
586 return data_xbegin_impl(this->storage().cbegin());
587 }
588
589 template <class CT, class S, layout_type L, class FST>
590 inline auto xdynamic_view<CT, S, L, FST>::data_xend(layout_type l, size_type offset) noexcept
591 -> container_iterator
592 {
593 return data_xend_impl(this->storage().begin(), l, offset);
594 }
595
596 template <class CT, class S, layout_type L, class FST>
597 inline auto xdynamic_view<CT, S, L, FST>::data_xend(layout_type l, size_type offset) const noexcept
598 -> const_container_iterator
599 {
600 return data_xend_impl(this->storage().cbegin(), l, offset);
601 }
602
603 template <class CT, class S, layout_type L, class FST>
604 template <class It>
605 inline It xdynamic_view<CT, S, L, FST>::data_xbegin_impl(It begin) const noexcept
606 {
607 return begin + static_cast<std::ptrdiff_t>(data_offset());
608 }
609
610 // TODO: fix the data_xend implementation and assign_temporary_impl
611
612 template <class CT, class S, layout_type L, class FST>
613 template <class It>
614 inline It
615 xdynamic_view<CT, S, L, FST>::data_xend_impl(It begin, layout_type l, size_type offset) const noexcept
616 {
617 return strided_data_end(*this, begin + std::ptrdiff_t(data_offset()), l, offset);
618 }
619
620 template <class CT, class S, layout_type L, class FST>
621 inline void xdynamic_view<CT, S, L, FST>::assign_temporary_impl(temporary_type&& tmp)
622 {
623 std::copy(tmp.cbegin(), tmp.cend(), this->begin());
624 }
625
626 template <class CT, class S, layout_type L, class FST>
627 template <class T, class... Args>
628 inline auto
629 xdynamic_view<CT, S, L, FST>::adjust_offset(offset_type offset, T idx, Args... args) const noexcept
630 -> offset_type
631 {
632 constexpr size_type nb_args = sizeof...(Args) + 1;
633 size_type dim = base_type::dimension();
634 offset_type res = nb_args > dim ? adjust_offset(offset, args...)
635 : adjust_offset_impl(offset, dim - nb_args, idx, args...);
636 return res;
637 }
638
639 template <class CT, class S, layout_type L, class FST>
640 inline auto xdynamic_view<CT, S, L, FST>::adjust_offset(offset_type offset) const noexcept -> offset_type
641 {
642 return offset;
643 }
644
645 template <class CT, class S, layout_type L, class FST>
646 template <class T, class... Args>
647 inline auto
648 xdynamic_view<CT, S, L, FST>::adjust_offset_impl(offset_type offset, size_type idx_offset, T idx, Args... args)
649 const noexcept -> offset_type
650 {
651 offset_type sl_offset = xtl::visit(
652 [idx](const auto& sl)
653 {
654 using type = typename std::decay_t<decltype(sl)>::size_type;
655 return sl(type(idx));
656 },
657 m_slices[idx_offset]
658 );
659 offset_type res = offset + sl_offset * m_adj_strides[idx_offset];
660 return adjust_offset_impl(res, idx_offset + 1, args...);
661 }
662
663 template <class CT, class S, layout_type L, class FST>
664 inline auto xdynamic_view<CT, S, L, FST>::adjust_offset_impl(offset_type offset, size_type) const noexcept
665 -> offset_type
666 {
667 return offset;
668 }
669
670 template <class CT, class S, layout_type L, class FST>
671 template <class It>
672 inline auto
673 xdynamic_view<CT, S, L, FST>::adjust_element_offset(offset_type offset, It first, It last) const noexcept
674 -> offset_type
675 {
676 auto dst = std::distance(first, last);
677 offset_type dim = static_cast<offset_type>(dimension());
678 offset_type loop_offset = dst < dim ? dim - dst : offset_type(0);
679 offset_type idx_offset = dim < dst ? dst - dim : offset_type(0);
680 offset_type res = offset;
681 for (offset_type i = loop_offset; i < dim; ++i, ++first)
682 {
683 offset_type j = static_cast<offset_type>(first[idx_offset]);
684 offset_type sl_offset = xtl::visit(
685 [j](const auto& sl)
686 {
687 return static_cast<offset_type>(sl(j));
688 },
689 m_slices[static_cast<std::size_t>(i)]
690 );
691 res += sl_offset * m_adj_strides[static_cast<std::size_t>(i)];
692 }
693 return res;
694 }
695
696 /*****************************************
697 * xdynamic_view builders implementation *
698 *****************************************/
699
700 namespace detail
701 {
702 template <class V>
703 struct adj_strides_policy
704 {
705 using slice_vector = V;
706 using strides_type = dynamic_shape<std::ptrdiff_t>;
707
708 slice_vector new_slices;
709 strides_type new_adj_strides;
710
711 protected:
712
713 inline void resize(std::size_t size)
714 {
715 new_slices.resize(size);
716 new_adj_strides.resize(size);
717 }
718
719 inline void set_fake_slice(std::size_t idx)
720 {
721 new_slices[idx] = xfake_slice<std::ptrdiff_t>();
722 new_adj_strides[idx] = std::ptrdiff_t(0);
723 }
724
725 template <class ST, class S>
726 bool fill_args(
727 const xdynamic_slice_vector& slices,
728 std::size_t sl_idx,
729 std::size_t i,
730 std::size_t old_shape,
731 const ST& old_stride,
732 S& shape,
733 get_strides_t<S>& strides
734 )
735 {
736 return fill_args_impl<xkeep_slice<std::ptrdiff_t>>(
737 slices,
738 sl_idx,
739 i,
740 old_shape,
741 old_stride,
742 shape,
743 strides
744 )
745 || fill_args_impl<xdrop_slice<std::ptrdiff_t>>(
746 slices,
747 sl_idx,
748 i,
749 old_shape,
750 old_stride,
751 shape,
752 strides
753 );
754 }
755
756 template <class SL, class ST, class S>
757 bool fill_args_impl(
758 const xdynamic_slice_vector& slices,
759 std::size_t sl_idx,
760 std::size_t i,
761 std::size_t old_shape,
762 const ST& old_stride,
763 S& shape,
764 get_strides_t<S>& strides
765 )
766 {
767 auto* sl = xtl::get_if<SL>(&slices[sl_idx]);
768 if (sl != nullptr)
769 {
770 new_slices[i] = *sl;
771 auto& ns = xtl::get<SL>(new_slices[i]);
772 ns.normalize(old_shape);
773 shape[i] = static_cast<std::size_t>(ns.size());
774 strides[i] = std::ptrdiff_t(0);
775 new_adj_strides[i] = static_cast<std::ptrdiff_t>(old_stride);
776 }
777 return sl != nullptr;
778 }
779 };
780 }
781
782 template <class E>
783 inline auto dynamic_view(E&& e, const xdynamic_slice_vector& slices)
784 {
785 using view_type = xdynamic_view<xclosure_t<E>, dynamic_shape<std::size_t>>;
786 using slice_vector = typename view_type::slice_vector_type;
787 using policy = detail::adj_strides_policy<slice_vector>;
788 detail::strided_view_args<policy> args;
789 args.fill_args(
790 e.shape(),
791 detail::get_strides<XTENSOR_DEFAULT_TRAVERSAL>(e),
792 detail::get_offset<XTENSOR_DEFAULT_TRAVERSAL>(e),
793 e.layout(),
794 slices
795 );
796 return view_type(
797 std::forward<E>(e),
798 std::move(args.new_shape),
799 std::move(args.new_strides),
800 args.new_offset,
801 args.new_layout,
802 std::move(args.new_slices),
803 std::move(args.new_adj_strides)
804 );
805 }
806}
807
808#endif
Base class for implementation of common expression access methods.
Base class for implementation of common expression constant access methods.
const_reference front() const
Returns a constant reference to first the element of the expression.
size_type size() const noexcept
Returns the size of the expression.
size_type dimension() const noexcept
Returns the number of dimensions of the expression.
bool in_bounds(Args... args) const
Returns true only if the the specified position is a valid entry in the expression.
const_reference back() const
Returns a constant reference to last the element of the expression.
layout_type layout() const noexcept
Returns the layout of the xtrided_view_base.
const inner_shape_type & shape() const noexcept
Returns the shape of the xtrided_view_base.
Base class for multidimensional iterable expressions.
layout_type layout() const noexcept
Returns the layout of the xtrided_view_base.
bool broadcast_shape(O &shape, bool reuse_cache=false) const
Broadcast the shape of the view to the specified parameter.
const inner_shape_type & shape() const noexcept
Returns the shape of the xtrided_view_base.
storage_type & storage() noexcept
Returns a reference to the buffer containing the elements of the view.
xexpression_type & expression() noexcept
Returns a reference to the underlying expression of the view.
Implementation of the xsemantic_base interface for multidimensional views.
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
bool operator==(const xaxis_iterator< CT > &lhs, const xaxis_iterator< CT > &rhs)
Checks equality of the iterators.
layout_type
Definition xlayout.hpp:24
bool operator!=(const xaxis_iterator< CT > &lhs, const xaxis_iterator< CT > &rhs)
Checks inequality of the iterators.