xtensor
Loading...
Searching...
No Matches
xsort.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_SORT_HPP
11#define XTENSOR_SORT_HPP
12
13#include <algorithm>
14#include <cmath>
15#include <iterator>
16#include <utility>
17
18#include <xtl/xcompare.hpp>
19
20#include "xadapt.hpp"
21#include "xarray.hpp"
22#include "xeval.hpp"
23#include "xindex_view.hpp"
24#include "xmanipulation.hpp"
25#include "xmath.hpp"
26#include "xslice.hpp" // for xnone
27#include "xtensor.hpp"
28#include "xtensor_config.hpp"
29#include "xtensor_forward.hpp"
30#include "xview.hpp"
31
32namespace xt
33{
41 namespace detail
42 {
43 template <class T>
44 std::ptrdiff_t adjust_secondary_stride(std::ptrdiff_t stride, T shape)
45 {
46 return stride != 0 ? stride : static_cast<std::ptrdiff_t>(shape);
47 }
48
49 template <class E>
50 inline std::ptrdiff_t get_secondary_stride(const E& ev)
51 {
52 if (ev.layout() == layout_type::row_major)
53 {
54 return adjust_secondary_stride(ev.strides()[ev.dimension() - 2], *(ev.shape().end() - 1));
55 }
56
57 return adjust_secondary_stride(ev.strides()[1], *(ev.shape().begin()));
58 }
59
60 template <class E>
61 inline std::size_t leading_axis_n_iters(const E& ev)
62 {
63 if (ev.layout() == layout_type::row_major)
64 {
65 return std::accumulate(
66 ev.shape().begin(),
67 ev.shape().end() - 1,
68 std::size_t(1),
69 std::multiplies<>()
70 );
71 }
72 return std::accumulate(ev.shape().begin() + 1, ev.shape().end(), std::size_t(1), std::multiplies<>());
73 }
74
75 template <class E, class F>
76 inline void call_over_leading_axis(E& ev, F&& fct)
77 {
78 XTENSOR_ASSERT(ev.dimension() >= 2);
79
80 const std::size_t n_iters = leading_axis_n_iters(ev);
81 const std::ptrdiff_t secondary_stride = get_secondary_stride(ev);
82
83 const auto begin = ev.data();
84 const auto end = begin + n_iters * secondary_stride;
85 for (auto iter = begin; iter != end; iter += secondary_stride)
86 {
87 fct(iter, iter + secondary_stride);
88 }
89 }
90
91 template <class E1, class E2, class F>
92 inline void call_over_leading_axis(E1& e1, E2& e2, F&& fct)
93 {
94 XTENSOR_ASSERT(e1.dimension() >= 2);
95 XTENSOR_ASSERT(e1.dimension() == e2.dimension());
96
97 const std::size_t n_iters = leading_axis_n_iters(e1);
98 const std::ptrdiff_t secondary_stride1 = get_secondary_stride(e1);
99 const std::ptrdiff_t secondary_stride2 = get_secondary_stride(e2);
100 XTENSOR_ASSERT(secondary_stride1 == secondary_stride2);
101
102 const auto begin1 = e1.data();
103 const auto end1 = begin1 + n_iters * secondary_stride1;
104 const auto begin2 = e2.data();
105 const auto end2 = begin2 + n_iters * secondary_stride2;
106 auto iter1 = begin1;
107 auto iter2 = begin2;
108 for (; (iter1 != end1) && (iter2 != end2); iter1 += secondary_stride1, iter2 += secondary_stride2)
109 {
110 fct(iter1, iter1 + secondary_stride1, iter2, iter2 + secondary_stride2);
111 }
112 }
113
114 template <class E>
115 inline std::size_t leading_axis(const E& e)
116 {
117 if (e.layout() == layout_type::row_major)
118 {
119 return e.dimension() - 1;
120 }
121 else if (e.layout() == layout_type::column_major)
122 {
123 return 0;
124 }
125 XTENSOR_THROW(std::runtime_error, "Layout not supported.");
126 }
127
128 // get permutations to transpose and reverse-transpose array
129 inline std::pair<dynamic_shape<std::size_t>, dynamic_shape<std::size_t>>
130 get_permutations(std::size_t dim, std::size_t ax, layout_type layout)
131 {
132 dynamic_shape<std::size_t> permutation(dim);
133 std::iota(permutation.begin(), permutation.end(), std::size_t(0));
134 permutation.erase(permutation.begin() + std::ptrdiff_t(ax));
135
136 if (layout == layout_type::row_major)
137 {
138 permutation.push_back(ax);
139 }
140 else
141 {
142 permutation.insert(permutation.begin(), ax);
143 }
144
145 // TODO find a more clever way to get reverse permutation?
146 dynamic_shape<std::size_t> reverse_permutation;
147 for (std::size_t i = 0; i < dim; ++i)
148 {
149 auto it = std::find(permutation.begin(), permutation.end(), i);
150 reverse_permutation.push_back(std::size_t(std::distance(permutation.begin(), it)));
151 }
152
153 return std::make_pair(std::move(permutation), std::move(reverse_permutation));
154 }
155
156 template <class R, class E, class F>
157 inline R map_axis(const E& e, std::ptrdiff_t axis, F&& lambda)
158 {
159 if (e.dimension() == 1)
160 {
161 R res = e;
162 lambda(res.begin(), res.end());
163 return res;
164 }
165
166 const std::size_t ax = normalize_axis(e.dimension(), axis);
167 if (ax == detail::leading_axis(e))
168 {
169 R res = e;
170 detail::call_over_leading_axis(res, std::forward<F>(lambda));
171 return res;
172 }
173
174 dynamic_shape<std::size_t> permutation, reverse_permutation;
175 std::tie(permutation, reverse_permutation) = get_permutations(e.dimension(), ax, e.layout());
176 R res = transpose(e, permutation);
177 detail::call_over_leading_axis(res, std::forward<F>(lambda));
178 res = transpose(res, reverse_permutation);
179 return res;
180 }
181
182 template <class VT>
183 struct flatten_sort_result_type_impl
184 {
185 using type = VT;
186 };
187
188 template <class VT, std::size_t N, layout_type L>
189 struct flatten_sort_result_type_impl<xtensor<VT, N, L>>
190 {
191 using type = xtensor<VT, 1, L>;
192 };
193
194 template <class VT, class S, layout_type L>
195 struct flatten_sort_result_type_impl<xtensor_fixed<VT, S, L>>
196 {
197 using type = xtensor_fixed<VT, xshape<fixed_compute_size<S>::value>, L>;
198 };
199
200 template <class VT>
201 struct flatten_sort_result_type : flatten_sort_result_type_impl<common_tensor_type_t<VT>>
202 {
203 };
204
205 template <class VT>
206 using flatten_sort_result_type_t = typename flatten_sort_result_type<VT>::type;
207
208 template <class E, class R = flatten_sort_result_type_t<E>>
209 inline auto flat_sort_impl(const xexpression<E>& e)
210 {
211 const auto& de = e.derived_cast();
212 R ev;
213 ev.resize({static_cast<typename R::shape_type::value_type>(de.size())});
214
215 std::copy(de.cbegin(), de.cend(), ev.begin());
216 std::sort(ev.begin(), ev.end());
217
218 return ev;
219 }
220 }
221
222 template <class E>
223 inline auto sort(const xexpression<E>& e, placeholders::xtuph /*t*/)
224 {
225 return detail::flat_sort_impl(e);
226 }
227
228 namespace detail
229 {
230 template <class T>
231 struct sort_eval_type
232 {
233 using type = typename T::temporary_type;
234 };
235
236 template <class T, std::size_t... I, layout_type L>
237 struct sort_eval_type<xtensor_fixed<T, fixed_shape<I...>, L>>
238 {
239 using type = xtensor<T, sizeof...(I), L>;
240 };
241 }
242
254 template <class E>
255 inline auto sort(const xexpression<E>& e, std::ptrdiff_t axis = -1)
256 {
257 using eval_type = typename detail::sort_eval_type<E>::type;
258
259 return detail::map_axis<eval_type>(
260 e.derived_cast(),
261 axis,
262 [](auto begin, auto end)
263 {
264 std::sort(begin, end);
265 }
266 );
267 }
268
269 /*****************************
270 * Implementation of argsort *
271 *****************************/
272
278 enum class sorting_method
279 {
284 quick,
289 stable,
290 };
291
292 namespace detail
293 {
294 template <class ConstRandomIt, class RandomIt, class Compare, class Method>
295 inline void argsort_iter(
296 ConstRandomIt data_begin,
297 ConstRandomIt data_end,
298 RandomIt idx_begin,
299 RandomIt idx_end,
300 Compare comp,
301 Method method
302 )
303 {
304 XTENSOR_ASSERT(std::distance(data_begin, data_end) >= 0);
305 XTENSOR_ASSERT(std::distance(idx_begin, idx_end) == std::distance(data_begin, data_end));
306 (void) idx_end; // TODO(C++17) [[maybe_unused]] only used in assertion.
307
308 std::iota(idx_begin, idx_end, 0);
309 switch (method)
310 {
312 {
313 std::sort(
314 idx_begin,
315 idx_end,
316 [&](const auto i, const auto j)
317 {
318 return comp(*(data_begin + i), *(data_begin + j));
319 }
320 );
321 }
323 {
324 std::stable_sort(
325 idx_begin,
326 idx_end,
327 [&](const auto i, const auto j)
328 {
329 return comp(*(data_begin + i), *(data_begin + j));
330 }
331 );
332 }
333 }
334 }
335
336 template <class ConstRandomIt, class RandomIt, class Method>
337 inline void
338 argsort_iter(ConstRandomIt data_begin, ConstRandomIt data_end, RandomIt idx_begin, RandomIt idx_end, Method method)
339 {
340 return argsort_iter(
341 std::move(data_begin),
342 std::move(data_end),
343 std::move(idx_begin),
344 std::move(idx_end),
345 [](const auto& x, const auto& y) -> bool
346 {
347 return x < y;
348 },
349 method
350 );
351 }
352
353 template <class VT, class T>
354 struct rebind_value_type
355 {
356 using type = xarray<VT, xt::layout_type::dynamic>;
357 };
358
359 template <class VT, class EC, layout_type L>
360 struct rebind_value_type<VT, xarray<EC, L>>
361 {
362 using type = xarray<VT, L>;
363 };
364
365 template <class VT, class EC, std::size_t N, layout_type L>
366 struct rebind_value_type<VT, xtensor<EC, N, L>>
367 {
368 using type = xtensor<VT, N, L>;
369 };
370
371 template <class VT, class ET, class S, layout_type L>
372 struct rebind_value_type<VT, xtensor_fixed<ET, S, L>>
373 {
374 using type = xtensor_fixed<VT, S, L>;
375 };
376
377 template <class VT, class T>
378 struct flatten_rebind_value_type
379 {
380 using type = typename rebind_value_type<VT, T>::type;
381 };
382
383 template <class VT, class EC, std::size_t N, layout_type L>
384 struct flatten_rebind_value_type<VT, xtensor<EC, N, L>>
385 {
386 using type = xtensor<VT, 1, L>;
387 };
388
389 template <class VT, class ET, class S, layout_type L>
390 struct flatten_rebind_value_type<VT, xtensor_fixed<ET, S, L>>
391 {
392 using type = xtensor_fixed<VT, xshape<fixed_compute_size<S>::value>, L>;
393 };
394
395 template <class T>
396 struct argsort_result_type
397 {
398 using type = typename rebind_value_type<typename T::temporary_type::size_type, typename T::temporary_type>::type;
399 };
400
401 template <class T>
402 struct linear_argsort_result_type
403 {
404 using type = typename flatten_rebind_value_type<
405 typename T::temporary_type::size_type,
406 typename T::temporary_type>::type;
407 };
408
409 template <class E, class R = typename detail::linear_argsort_result_type<E>::type, class Method>
410 inline auto flatten_argsort_impl(const xexpression<E>& e, Method method)
411 {
412 const auto& de = e.derived_cast();
413
414 auto cit = de.template begin<layout_type::row_major>();
415 using const_iterator = decltype(cit);
416 auto ad = xiterator_adaptor<const_iterator, const_iterator>(cit, cit, de.size());
417
418 using result_type = R;
419 result_type result;
420 result.resize({de.size()});
421
422 detail::argsort_iter(de.cbegin(), de.cend(), result.begin(), result.end(), method);
423
424 return result;
425 }
426 }
427
428 template <class E>
429 inline auto
430 argsort(const xexpression<E>& e, placeholders::xtuph /*t*/, sorting_method method = sorting_method::quick)
431 {
432 return detail::flatten_argsort_impl(e, method);
433 }
434
450 template <class E>
451 inline auto
452 argsort(const xexpression<E>& e, std::ptrdiff_t axis = -1, sorting_method method = sorting_method::quick)
453 {
454 using eval_type = typename detail::sort_eval_type<E>::type;
455 using result_type = typename detail::argsort_result_type<eval_type>::type;
456
457 const auto& de = e.derived_cast();
458
459 std::size_t ax = normalize_axis(de.dimension(), axis);
460
461 if (de.dimension() == 1)
462 {
463 return detail::flatten_argsort_impl<E, result_type>(e, method);
464 }
465
466 const auto argsort = [&method](auto res_begin, auto res_end, auto ev_begin, auto ev_end)
467 {
468 detail::argsort_iter(ev_begin, ev_end, res_begin, res_end, method);
469 };
470
471 if (ax == detail::leading_axis(de))
472 {
473 result_type res = result_type::from_shape(de.shape());
474 detail::call_over_leading_axis(res, de, argsort);
475 return res;
476 }
477
479 std::tie(permutation, reverse_permutation) = detail::get_permutations(de.dimension(), ax, de.layout());
480 eval_type ev = transpose(de, permutation);
481 result_type res = result_type::from_shape(ev.shape());
482 detail::call_over_leading_axis(res, ev, argsort);
484 return res;
485 }
486
487 /************************************************
488 * Implementation of partition and argpartition *
489 ************************************************/
490
491 namespace detail
492 {
504 template <class RandomIt, class Iter, class Compare>
505 inline void
506 partition_iter(RandomIt data_begin, RandomIt data_end, Iter kth_begin, Iter kth_end, Compare comp)
507 {
508 XTENSOR_ASSERT(std::distance(data_begin, data_end) >= 0);
509 XTENSOR_ASSERT(std::distance(kth_begin, kth_end) >= 0);
510
511 using idx_type = typename std::iterator_traits<Iter>::value_type;
512
513 idx_type k_last = static_cast<idx_type>(std::distance(data_begin, data_end));
514 for (; kth_begin != kth_end; ++kth_begin)
515 {
516 std::nth_element(data_begin, data_begin + *kth_begin, data_begin + k_last, std::move(comp));
517 k_last = *kth_begin;
518 }
519 }
520
521 template <class RandomIt, class Iter>
522 inline void partition_iter(RandomIt data_begin, RandomIt data_end, Iter kth_begin, Iter kth_end)
523 {
524 return partition_iter(
525 std::move(data_begin),
526 std::move(data_end),
527 std::move(kth_begin),
528 std::move(kth_end),
529 [](const auto& x, const auto& y) -> bool
530 {
531 return x < y;
532 }
533 );
534 }
535 }
536
564 template <
565 class E,
566 class C,
567 class R = detail::flatten_sort_result_type_t<E>,
568 class = std::enable_if_t<!xtl::is_integral<C>::value, int>>
570 {
571 const auto& de = e.derived_cast();
572
573 R ev = R::from_shape({de.size()});
574 std::sort(kth_container.begin(), kth_container.end());
575
576 std::copy(de.linear_cbegin(), de.linear_cend(), ev.linear_begin()); // flatten
577
578 detail::partition_iter(ev.linear_begin(), ev.linear_end(), kth_container.rbegin(), kth_container.rend());
579
580 return ev;
581 }
582
583 template <class E, class I, std::size_t N, class R = detail::flatten_sort_result_type_t<E>>
584 inline R partition(const xexpression<E>& e, const I (&kth_container)[N], placeholders::xtuph tag)
585 {
586 return partition(
587 e,
588 xtl::forward_sequence<std::array<std::size_t, N>, decltype(kth_container)>(kth_container),
589 tag
590 );
591 }
592
593 template <class E, class R = detail::flatten_sort_result_type_t<E>>
594 inline R partition(const xexpression<E>& e, std::size_t kth, placeholders::xtuph tag)
595 {
596 return partition(e, std::array<std::size_t, 1>({kth}), tag);
597 }
598
599 template <class E, class C, class = std::enable_if_t<!xtl::is_integral<C>::value, int>>
600 inline auto partition(const xexpression<E>& e, C kth_container, std::ptrdiff_t axis = -1)
601 {
602 using eval_type = typename detail::sort_eval_type<E>::type;
603
604 std::sort(kth_container.begin(), kth_container.end());
605
606 return detail::map_axis<eval_type>(
607 e.derived_cast(),
608 axis,
609 [&kth_container](auto begin, auto end)
610 {
611 detail::partition_iter(begin, end, kth_container.rbegin(), kth_container.rend());
612 }
613 );
614 }
615
616 template <class E, class T, std::size_t N>
617 inline auto partition(const xexpression<E>& e, const T (&kth_container)[N], std::ptrdiff_t axis = -1)
618 {
619 return partition(
620 e,
621 xtl::forward_sequence<std::array<std::size_t, N>, decltype(kth_container)>(kth_container),
622 axis
623 );
624 }
625
626 template <class E>
627 inline auto partition(const xexpression<E>& e, std::size_t kth, std::ptrdiff_t axis = -1)
628 {
629 return partition(e, std::array<std::size_t, 1>({kth}), axis);
630 }
631
659 template <
660 class E,
661 class C,
662 class R = typename detail::linear_argsort_result_type<typename detail::sort_eval_type<E>::type>::type,
663 class = std::enable_if_t<!xtl::is_integral<C>::value, int>>
664 inline R argpartition(const xexpression<E>& e, C kth_container, placeholders::xtuph)
665 {
666 using eval_type = typename detail::sort_eval_type<E>::type;
667 using result_type = typename detail::linear_argsort_result_type<eval_type>::type;
668
669 const auto& de = e.derived_cast();
670
671 result_type res = result_type::from_shape({de.size()});
672
673 std::sort(kth_container.begin(), kth_container.end());
674
675 std::iota(res.linear_begin(), res.linear_end(), 0);
676
677 detail::partition_iter(
678 res.linear_begin(),
679 res.linear_end(),
680 kth_container.rbegin(),
681 kth_container.rend(),
682 [&de](std::size_t a, std::size_t b)
683 {
684 return de[a] < de[b];
685 }
686 );
687
688 return res;
689 }
690
691 template <class E, class I, std::size_t N>
692 inline auto argpartition(const xexpression<E>& e, const I (&kth_container)[N], placeholders::xtuph tag)
693 {
694 return argpartition(
695 e,
696 xtl::forward_sequence<std::array<std::size_t, N>, decltype(kth_container)>(kth_container),
697 tag
698 );
699 }
700
701 template <class E>
702 inline auto argpartition(const xexpression<E>& e, std::size_t kth, placeholders::xtuph tag)
703 {
704 return argpartition(e, std::array<std::size_t, 1>({kth}), tag);
705 }
706
707 template <class E, class C, class = std::enable_if_t<!xtl::is_integral<C>::value, int>>
708 inline auto argpartition(const xexpression<E>& e, C kth_container, std::ptrdiff_t axis = -1)
709 {
710 using eval_type = typename detail::sort_eval_type<E>::type;
711 using result_type = typename detail::argsort_result_type<eval_type>::type;
712
713 const auto& de = e.derived_cast();
714
715 if (de.dimension() == 1)
716 {
717 return argpartition<E, C, result_type>(e, std::forward<C>(kth_container), xnone());
718 }
719
720 std::sort(kth_container.begin(), kth_container.end());
721 const auto argpartition_w_kth =
722 [&kth_container](auto res_begin, auto res_end, auto ev_begin, auto /*ev_end*/)
723 {
724 std::iota(res_begin, res_end, 0);
725 detail::partition_iter(
726 res_begin,
727 res_end,
728 kth_container.rbegin(),
729 kth_container.rend(),
730 [&ev_begin](auto const& i, auto const& j)
731 {
732 return *(ev_begin + i) < *(ev_begin + j);
733 }
734 );
735 };
736
737 const std::size_t ax = normalize_axis(de.dimension(), axis);
738 if (ax == detail::leading_axis(de))
739 {
740 result_type res = result_type::from_shape(de.shape());
741 detail::call_over_leading_axis(res, de, argpartition_w_kth);
742 return res;
743 }
744
745 dynamic_shape<std::size_t> permutation, reverse_permutation;
746 std::tie(permutation, reverse_permutation) = detail::get_permutations(de.dimension(), ax, de.layout());
747 eval_type ev = transpose(de, permutation);
748 result_type res = result_type::from_shape(ev.shape());
749 detail::call_over_leading_axis(res, ev, argpartition_w_kth);
750 res = transpose(res, reverse_permutation);
751 return res;
752 }
753
754 template <class E, class I, std::size_t N>
755 inline auto argpartition(const xexpression<E>& e, const I (&kth_container)[N], std::ptrdiff_t axis = -1)
756 {
757 return argpartition(
758 e,
759 xtl::forward_sequence<std::array<std::size_t, N>, decltype(kth_container)>(kth_container),
760 axis
761 );
762 }
763
764 template <class E>
765 inline auto argpartition(const xexpression<E>& e, std::size_t kth, std::ptrdiff_t axis = -1)
766 {
767 return argpartition(e, std::array<std::size_t, 1>({kth}), axis);
768 }
769
770 /******************
771 * xt::quantile *
772 ******************/
773
774 namespace detail
775 {
776 template <class S, class I, class K, class O>
777 inline void select_indices_impl(
778 const S& shape,
779 const I& indices,
780 std::size_t axis,
781 std::size_t current_dim,
782 const K& current_index,
783 O& out
784 )
785 {
786 using id_t = typename K::value_type;
787 if ((current_dim < shape.size() - 1) && (current_dim == axis))
788 {
789 for (auto i : indices)
790 {
791 auto idx = current_index;
792 idx[current_dim] = i;
793 select_indices_impl(shape, indices, axis, current_dim + 1, idx, out);
794 }
795 }
796 else if ((current_dim < shape.size() - 1) && (current_dim != axis))
797 {
798 for (id_t i = 0; xtl::cmp_less(i, shape[current_dim]); ++i)
799 {
800 auto idx = current_index;
801 idx[current_dim] = i;
802 select_indices_impl(shape, indices, axis, current_dim + 1, idx, out);
803 }
804 }
805 else if ((current_dim == shape.size() - 1) && (current_dim == axis))
806 {
807 for (auto i : indices)
808 {
809 auto idx = current_index;
810 idx[current_dim] = i;
811 out.push_back(std::move(idx));
812 }
813 }
814 else if ((current_dim == shape.size() - 1) && (current_dim != axis))
815 {
816 for (id_t i = 0; xtl::cmp_less(i, shape[current_dim]); ++i)
817 {
818 auto idx = current_index;
819 idx[current_dim] = i;
820 out.push_back(std::move(idx));
821 }
822 }
823 }
824
825 template <class S, class I>
826 inline auto select_indices(const S& shape, const I& indices, std::size_t axis)
827 {
828 using index_type = get_strides_t<S>;
829 auto out = std::vector<index_type>();
830 select_indices_impl(shape, indices, axis, 0, xtl::make_sequence<index_type>(shape.size()), out);
831 return out;
832 }
833
834 // TODO remove when fancy index views are implemented
835 // Poor man's indexing along a single axis as in NumPy a[:, [1, 3, 4]]
836 template <class E, class I>
837 inline auto fancy_indexing(E&& e, const I& indices, std::ptrdiff_t axis)
838 {
839 const std::size_t ax = normalize_axis(e.dimension(), axis);
840 using shape_t = get_strides_t<typename std::decay_t<E>::shape_type>;
841 auto shape = xtl::forward_sequence<shape_t, decltype(e.shape())>(e.shape());
842 shape[ax] = indices.size();
843 return reshape_view(
844 index_view(std::forward<E>(e), select_indices(e.shape(), indices, ax)),
845 std::move(shape)
846 );
847 }
848
849 template <class T, class I, class P>
850 inline auto quantile_kth_gamma(std::size_t n, const P& probas, T alpha, T beta)
851 {
852 const auto m = alpha + probas * (T(1) - alpha - beta);
853 // Evaluting since reused a lot
854 const auto p_n_m = eval(probas * static_cast<T>(n) + m - 1);
855 // Previous (virtual) index, may be out of bounds
856 const auto j = floor(p_n_m);
857 const auto j_jp1 = concatenate(xtuple(j, j + 1));
858 // Both interpolation indices, k and k+1
859 const auto k_kp1 = xt::cast<std::size_t>(clip(j_jp1, T(0), T(n - 1)));
860 // Both interpolation coefficients, 1-gamma and gamma
861 const auto omg_g = concatenate(xtuple(T(1) - (p_n_m - j), p_n_m - j));
862 return std::make_pair(eval(k_kp1), eval(omg_g));
863 }
864
865 // TODO should implement unsqueeze rather
866 template <class S>
867 inline auto unsqueeze_shape(const S& shape, std::size_t axis)
868 {
869 XTENSOR_ASSERT(axis <= shape.size());
870 auto new_shape = xtl::forward_sequence<xt::svector<std::size_t>, decltype(shape)>(shape);
871 new_shape.insert(new_shape.begin() + axis, 1);
872 return new_shape;
873 }
874 }
875
905 template <class T = double, class E, class P>
906 inline auto quantile(E&& e, const P& probas, std::ptrdiff_t axis, T alpha, T beta)
907 {
908 XTENSOR_ASSERT(all(0. <= probas));
909 XTENSOR_ASSERT(all(probas <= 1.));
910 XTENSOR_ASSERT(0. <= alpha);
911 XTENSOR_ASSERT(alpha <= 1.);
912 XTENSOR_ASSERT(0. <= beta);
913 XTENSOR_ASSERT(beta <= 1.);
914
916 using id_t = typename tmp_shape_t::value_type;
917
918 const std::size_t ax = normalize_axis(e.dimension(), axis);
919 const std::size_t n = e.shape()[ax];
920 auto kth_gamma = detail::quantile_kth_gamma<T, id_t, P>(n, probas, alpha, beta);
921
922 // Select relevant values for computing interpolating quantiles
923 auto e_partition = xt::partition(std::forward<E>(e), kth_gamma.first, ax);
924 auto e_kth = detail::fancy_indexing(std::move(e_partition), std::move(kth_gamma.first), ax);
925
926 // Reshape interpolation coefficients
927 auto gm1_g_shape = xtl::make_sequence<tmp_shape_t>(e.dimension(), 1);
928 gm1_g_shape[ax] = kth_gamma.second.size();
929 auto gm1_g_reshaped = reshape_view(std::move(kth_gamma.second), std::move(gm1_g_shape));
930
931 // Compute interpolation
932 // TODO(C++20) use (and create) xt::lerp in C++
933 auto e_kth_g = std::move(e_kth) * std::move(gm1_g_reshaped);
934 // Reshape pairwise interpolate for suming along new axis
935 auto e_kth_g_shape = detail::unsqueeze_shape(e_kth_g.shape(), ax);
936 e_kth_g_shape[ax] = 2;
937 e_kth_g_shape[ax + 1] /= 2;
938 auto quantiles = xt::sum(reshape_view(std::move(e_kth_g), std::move(e_kth_g_shape)), ax);
939 // Cannot do a transpose on a non-strided expression so we have to eval
940 return moveaxis(eval(std::move(quantiles)), ax, 0);
941 }
942
943 // Static proba array overload
944 template <class T = double, class E, std::size_t N>
945 inline auto quantile(E&& e, const T (&probas)[N], std::ptrdiff_t axis, T alpha, T beta)
946 {
947 return quantile(std::forward<E>(e), adapt(probas, {N}), axis, alpha, beta);
948 }
949
959 template <class T = double, class E, class P>
960 inline auto quantile(E&& e, const P& probas, T alpha, T beta)
961 {
962 return quantile(xt::ravel(std::forward<E>(e)), probas, 0, alpha, beta);
963 }
964
965 // Static proba array overload
966 template <class T = double, class E, std::size_t N>
967 inline auto quantile(E&& e, const T (&probas)[N], T alpha, T beta)
968 {
969 return quantile(std::forward<E>(e), adapt(probas, {N}), alpha, beta);
970 }
971
984 {
988 hazen,
990 weibull,
992 linear,
997 };
998
1008 template <class T = double, class E, class P>
1009 inline auto
1010 quantile(E&& e, const P& probas, std::ptrdiff_t axis, quantile_method method = quantile_method::linear)
1011 {
1012 T alpha = 0.;
1013 T beta = 0.;
1014 switch (method)
1015 {
1017 {
1018 alpha = 0.;
1019 beta = 1.;
1020 break;
1021 }
1023 {
1024 alpha = 0.5;
1025 beta = 0.5;
1026 break;
1027 }
1029 {
1030 alpha = 0.;
1031 beta = 0.;
1032 break;
1033 }
1035 {
1036 alpha = 1.;
1037 beta = 1.;
1038 break;
1039 }
1041 {
1042 alpha = 1. / 3.;
1043 beta = 1. / 3.;
1044 break;
1045 }
1047 {
1048 alpha = 3. / 8.;
1049 beta = 3. / 8.;
1050 break;
1051 }
1052 }
1053 return quantile(std::forward<E>(e), probas, axis, alpha, beta);
1054 }
1055
1056 // Static proba array overload
1057 template <class T = double, class E, std::size_t N>
1058 inline auto
1059 quantile(E&& e, const T (&probas)[N], std::ptrdiff_t axis, quantile_method method = quantile_method::linear)
1060 {
1061 return quantile(std::forward<E>(e), adapt(probas, {N}), axis, method);
1062 }
1063
1075 template <class T = double, class E, class P>
1077 {
1078 return quantile(xt::ravel(std::forward<E>(e)), probas, 0, method);
1079 }
1080
1081 // Static proba array overload
1082 template <class T = double, class E, std::size_t N>
1083 inline auto quantile(E&& e, const T (&probas)[N], quantile_method method = quantile_method::linear)
1084 {
1085 return quantile(std::forward<E>(e), adapt(probas, {N}), method);
1086 }
1087
1088 /****************
1089 * xt::median *
1090 ****************/
1091
1092 template <class E>
1093 inline typename std::decay_t<E>::value_type median(E&& e)
1094 {
1095 using value_type = typename std::decay_t<E>::value_type;
1096 auto sz = e.size();
1097 if (sz % 2 == 0)
1098 {
1099 std::size_t szh = sz / 2; // integer floor div
1100 std::array<std::size_t, 2> kth = {szh - 1, szh};
1101 auto values = xt::partition(xt::flatten(e), kth);
1102 return (values[kth[0]] + values[kth[1]]) / value_type(2);
1103 }
1104 else
1105 {
1106 std::array<std::size_t, 1> kth = {(sz - 1) / 2};
1107 auto values = xt::partition(xt::flatten(e), kth);
1108 return values[kth[0]];
1109 }
1110 }
1111
1125 template <class E>
1126 inline auto median(E&& e, std::ptrdiff_t axis)
1127 {
1128 std::size_t ax = normalize_axis(e.dimension(), axis);
1129 std::size_t sz = e.shape()[ax];
1130 xstrided_slice_vector sv(e.dimension(), xt::all());
1131
1132 if (sz % 2 == 0)
1133 {
1134 std::size_t szh = sz / 2; // integer floor div
1135 std::array<std::size_t, 2> kth = {szh - 1, szh};
1136 auto values = xt::partition(std::forward<E>(e), kth, static_cast<ptrdiff_t>(ax));
1137 sv[ax] = xt::range(szh - 1, szh + 1);
1138 return xt::mean(xt::strided_view(std::move(values), std::move(sv)), {ax});
1139 }
1140 else
1141 {
1142 std::size_t szh = (sz - 1) / 2;
1143 std::array<std::size_t, 1> kth = {(sz - 1) / 2};
1144 auto values = xt::partition(std::forward<E>(e), kth, static_cast<ptrdiff_t>(ax));
1145 sv[ax] = xt::range(szh, szh + 1);
1146 return xt::mean(xt::strided_view(std::move(values), std::move(sv)), {ax});
1147 }
1148 }
1149
1150 namespace detail
1151 {
1152 template <class T>
1153 struct argfunc_result_type
1154 {
1155 using type = xarray<std::size_t>;
1156 };
1157
1158 template <class T, std::size_t N>
1159 struct argfunc_result_type<xtensor<T, N>>
1160 {
1161 using type = xtensor<std::size_t, N - 1>;
1162 };
1163
1164 template <layout_type L, class E, class F>
1165 inline typename argfunc_result_type<E>::type arg_func_impl(const E& e, std::size_t axis, F&& cmp)
1166 {
1167 using eval_type = typename detail::sort_eval_type<E>::type;
1168 using value_type = typename E::value_type;
1169 using result_type = typename argfunc_result_type<E>::type;
1170 using result_shape_type = typename result_type::shape_type;
1171
1172 if (e.dimension() == 1)
1173 {
1174 auto begin = e.template begin<L>();
1175 auto end = e.template end<L>();
1176 // todo C++17 : constexpr
1177 if (std::is_same<F, std::less<value_type>>::value)
1178 {
1179 std::size_t i = static_cast<std::size_t>(std::distance(begin, std::min_element(begin, end)));
1180 return xtensor<size_t, 0>{i};
1181 }
1182 else
1183 {
1184 std::size_t i = static_cast<std::size_t>(std::distance(begin, std::max_element(begin, end)));
1185 return xtensor<size_t, 0>{i};
1186 }
1187 }
1188
1189 result_shape_type alt_shape;
1190 xt::resize_container(alt_shape, e.dimension() - 1);
1191
1192 // Excluding copy, copy all of shape except for axis
1193 std::copy(e.shape().cbegin(), e.shape().cbegin() + std::ptrdiff_t(axis), alt_shape.begin());
1194 std::copy(
1195 e.shape().cbegin() + std::ptrdiff_t(axis) + 1,
1196 e.shape().cend(),
1197 alt_shape.begin() + std::ptrdiff_t(axis)
1198 );
1199
1200 result_type result = result_type::from_shape(std::move(alt_shape));
1201 auto result_iter = result.template begin<L>();
1202
1203 auto arg_func_lambda = [&result_iter, &cmp](auto begin, auto end)
1204 {
1205 std::size_t idx = 0;
1206 value_type val = *begin;
1207 ++begin;
1208 for (std::size_t i = 1; begin != end; ++begin, ++i)
1209 {
1210 if (cmp(*begin, val))
1211 {
1212 val = *begin;
1213 idx = i;
1214 }
1215 }
1216 *result_iter = idx;
1217 ++result_iter;
1218 };
1219
1220 if (axis != detail::leading_axis(e))
1221 {
1222 dynamic_shape<std::size_t> permutation, reverse_permutation;
1223 std::tie(
1224 permutation,
1225 reverse_permutation
1226 ) = detail::get_permutations(e.dimension(), axis, e.layout());
1227
1228 // note: creating copy
1229 eval_type input = transpose(e, permutation);
1230 detail::call_over_leading_axis(input, arg_func_lambda);
1231 return result;
1232 }
1233 else
1234 {
1235 auto&& input = eval(e);
1236 detail::call_over_leading_axis(input, arg_func_lambda);
1237 return result;
1238 }
1239 }
1240 }
1241
1242 template <layout_type L = XTENSOR_DEFAULT_TRAVERSAL, class E>
1243 inline auto argmin(const xexpression<E>& e)
1244 {
1245 using value_type = typename E::value_type;
1246 auto&& ed = eval(e.derived_cast());
1247 auto begin = ed.template begin<L>();
1248 auto end = ed.template end<L>();
1249 std::size_t i = static_cast<std::size_t>(std::distance(begin, std::min_element(begin, end)));
1250 return xtensor<size_t, 0>{i};
1251 }
1252
1263 template <layout_type L = XTENSOR_DEFAULT_TRAVERSAL, class E>
1264 inline auto argmin(const xexpression<E>& e, std::ptrdiff_t axis)
1265 {
1266 using value_type = typename E::value_type;
1267 auto&& ed = eval(e.derived_cast());
1268 std::size_t ax = normalize_axis(ed.dimension(), axis);
1269 return detail::arg_func_impl<L>(ed, ax, std::less<value_type>());
1270 }
1271
1272 template <layout_type L = XTENSOR_DEFAULT_TRAVERSAL, class E>
1273 inline auto argmax(const xexpression<E>& e)
1274 {
1275 using value_type = typename E::value_type;
1276 auto&& ed = eval(e.derived_cast());
1277 auto begin = ed.template begin<L>();
1278 auto end = ed.template end<L>();
1279 std::size_t i = static_cast<std::size_t>(std::distance(begin, std::max_element(begin, end)));
1280 return xtensor<size_t, 0>{i};
1281 }
1282
1294 template <layout_type L = XTENSOR_DEFAULT_TRAVERSAL, class E>
1295 inline auto argmax(const xexpression<E>& e, std::ptrdiff_t axis)
1296 {
1297 using value_type = typename E::value_type;
1298 auto&& ed = eval(e.derived_cast());
1299 std::size_t ax = normalize_axis(ed.dimension(), axis);
1300 return detail::arg_func_impl<L>(ed, ax, std::greater<value_type>());
1301 }
1302
1310 template <class E>
1311 inline auto unique(const xexpression<E>& e)
1312 {
1313 auto sorted = sort(e, xnone());
1314 auto end = std::unique(sorted.begin(), sorted.end());
1315 std::size_t sz = static_cast<std::size_t>(std::distance(sorted.begin(), end));
1316 // TODO check if we can shrink the vector without reallocation
1317 using value_type = typename E::value_type;
1318 auto result = xtensor<value_type, 1>::from_shape({sz});
1319 std::copy(sorted.begin(), end, result.begin());
1320 return result;
1321 }
1322
1331 template <class E1, class E2>
1333 {
1334 using value_type = typename E1::value_type;
1335
1336 auto unique1 = unique(ar1);
1337 auto unique2 = unique(ar2);
1338
1340
1341 auto end = std::set_difference(unique1.begin(), unique1.end(), unique2.begin(), unique2.end(), tmp.begin());
1342
1343 std::size_t sz = static_cast<std::size_t>(std::distance(tmp.begin(), end));
1344
1345 auto result = xtensor<value_type, 1>::from_shape({sz});
1346
1347 std::copy(tmp.begin(), end, result.begin());
1348
1349 return result;
1350 }
1351}
1352
1353#endif
auto clip(E1 &&e1, E2 &&lo, E3 &&hi) noexcept -> detail::xfunction_type_t< math::clamp_fun, E1, E2, E3 >
Clip values between hi and lo.
Definition xmath.hpp:815
auto floor(E &&e) noexcept -> detail::xfunction_type_t< math::floor_fun, E >
floor function.
Definition xmath.hpp:1592
auto sum(E &&e, X &&axes, EVS es=EVS())
Sum of elements over given axes.
Definition xmath.hpp:1841
auto mean(E &&e, X &&axes, EVS es=EVS())
Mean of elements over given axes.
Definition xmath.hpp:1935
auto adapt(C &&container, const SC &shape, layout_type l=L)
Constructs:
auto eval(T &&t) -> std::enable_if_t< detail::is_container< std::decay_t< T > >::value, T && >
Force evaluation of xexpression.
Definition xeval.hpp:46
auto flatten(E &&e)
Return a flatten view of the given expression.
auto moveaxis(E &&e, std::ptrdiff_t src, std::ptrdiff_t dest)
Return a new expression with an axis move to a new position.
auto ravel(E &&e)
Return a flatten view of the given expression.
auto transpose(E &&e) noexcept
Returns a transpose view by reversing the dimensions of xexpression e.
quantile_method
Quantile interpolation method.
Definition xsort.hpp:984
auto unique(const xexpression< E > &e)
Find unique elements of a xexpression.
Definition xsort.hpp:1311
R partition(const xexpression< E > &e, C kth_container, placeholders::xtuph)
Partially sort xexpression.
Definition xsort.hpp:569
auto quantile(E &&e, const P &probas, std::ptrdiff_t axis, T alpha, T beta)
Compute quantiles over the given axis.
Definition xsort.hpp:906
auto setdiff1d(const xexpression< E1 > &ar1, const xexpression< E2 > &ar2)
Find the set difference of two xexpressions.
Definition xsort.hpp:1332
R argpartition(const xexpression< E > &e, C kth_container, placeholders::xtuph)
Partially sort arguments.
Definition xsort.hpp:664
@ weibull
Method 6 of (Hyndman and Fan, 1996) with alpha=0 and beta=0.
@ interpolated_inverted_cdf
Method 4 of (Hyndman and Fan, 1996) with alpha=0 and beta=1.
@ linear
Method 7 of (Hyndman and Fan, 1996) with alpha=1 and beta=1.
@ normal_unbiased
Method 9 of (Hyndman and Fan, 1996) with alpha=3/8 and beta=3/8.
@ median_unbiased
Method 8 of (Hyndman and Fan, 1996) with alpha=1/3 and beta=1/3.
@ hazen
Method 5 of (Hyndman and Fan, 1996) with alpha=1/2 and beta=1/2.
standard mathematical functions for xexpressions
auto range(A start_val, B stop_val)
Select a range from start_val to stop_val (excluded).
Definition xslice.hpp:818
auto all() noexcept
Returns a slice representing a full dimension, to be used as an argument of view function.
Definition xslice.hpp:234
xtensor_container< uvector< T, A >, N, L > xtensor
Alias template on xtensor_container with default parameters for data container type.
auto concatenate(std::tuple< CT... > &&t, std::size_t axis=0)
Concatenates xexpressions along axis.
Definition xbuilder.hpp:830
xarray_container< uvector< T, A >, L, xt::svector< typename uvector< T, A >::size_type, 4, SA, true > > xarray
Alias template on xarray_container with default parameters for data container type and shape / stride...
layout_type
Definition xlayout.hpp:24
std::vector< xstrided_slice< std::ptrdiff_t > > xstrided_slice_vector
vector of slices used to build a xstrided_view
auto index_view(E &&e, I &&indices) noexcept
creates an indexview from a container of indices.
sorting_method
Sorting method.
Definition xsort.hpp:279
@ quick
Faster method but with no guarantee on preservation of order of equal elements https://en....
@ stable
Slower method but with guarantee on preservation of order of equal elements https://en....
auto strided_view(E &&e, S &&shape, X &&stride, std::size_t offset=0, layout_type layout=L) noexcept
Construct a strided view from an xexpression, shape, strides and offset.
auto xtuple(Types &&... args)
Creates tuples from arguments for concatenate and stack.
Definition xbuilder.hpp:707
xfixed_container< T, FSH, L, Sharable > xtensor_fixed
Alias template on xfixed_container with default parameters for layout type.