Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RVec.hxx
Go to the documentation of this file.
1// Author: Enrico Guiraud, Enric Tejedor, Danilo Piparo CERN 04/2021
2// Implementation adapted from from llvm::SmallVector.
3// See /math/vecops/ARCHITECTURE.md for more information.
4
5/*************************************************************************
6 * Copyright (C) 1995-2021, Rene Brun and Fons Rademakers. *
7 * All rights reserved. *
8 * *
9 * For the licensing terms see $ROOTSYS/LICENSE. *
10 * For the list of contributors see $ROOTSYS/README/CREDITS. *
11 *************************************************************************/
12
13#ifndef ROOT_RVEC
14#define ROOT_RVEC
15
16#if __cplusplus > 201402L
17#define R__RVEC_NODISCARD [[nodiscard]]
18#else
19#define R__RVEC_NODISCARD
20#endif
21
22#ifdef _WIN32
23 #ifndef M_PI
24 #ifndef _USE_MATH_DEFINES
25 #define _USE_MATH_DEFINES
26 #endif
27 #include <math.h>
28 #undef _USE_MATH_DEFINES
29 #endif
30 #define _VECOPS_USE_EXTERN_TEMPLATES false
31#else
32 #define _VECOPS_USE_EXTERN_TEMPLATES true
33#endif
34
35#include <Rtypes.h> // R__CLING_PTRCHECK
36#include <TError.h> // R__ASSERT
37
38#include <algorithm>
39#include <cmath>
40#include <cstring>
41#include <limits> // for numeric_limits
42#include <memory> // uninitialized_value_construct
43#include <new>
44#include <numeric> // for inner_product
45#include <sstream>
46#include <stdexcept>
47#include <string>
48#include <tuple>
49#include <type_traits>
50#include <utility>
51#include <vector>
52
53#ifdef R__HAS_VDT
54#include <vdt/vdtMath.h>
55#endif
56
57
58namespace ROOT {
59
60namespace VecOps {
61template<typename T>
62class RVec;
63}
64
65namespace Internal {
66namespace VecOps {
67
68template<typename T>
70
71// clang-format off
72template <typename>
73struct IsRVec : std::false_type {};
74
75template <typename T>
76struct IsRVec<ROOT::VecOps::RVec<T>> : std::true_type {};
77// clang-format on
78
79constexpr bool All(const bool *vals, std::size_t size)
80{
81 for (auto i = 0u; i < size; ++i)
82 if (!vals[i])
83 return false;
84 return true;
85}
86
87template <typename... T>
88std::size_t GetVectorsSize(const std::string &id, const RVec<T> &... vs)
89{
90 constexpr const auto nArgs = sizeof...(T);
91 const std::size_t sizes[] = {vs.size()...};
92 if (nArgs > 1) {
93 for (auto i = 1UL; i < nArgs; i++) {
94 if (sizes[0] == sizes[i])
95 continue;
96 std::string msg(id);
97 msg += ": input RVec instances have different lengths!";
98 throw std::runtime_error(msg);
99 }
100 }
101 return sizes[0];
102}
103
104template <typename F, typename... RVecs>
105auto MapImpl(F &&f, RVecs &&... vs) -> RVec<decltype(f(vs[0]...))>
106{
107 const auto size = GetVectorsSize("Map", vs...);
108 RVec<decltype(f(vs[0]...))> ret(size);
109
110 for (auto i = 0UL; i < size; i++)
111 ret[i] = f(vs[i]...);
112
113 return ret;
114}
115
116template <typename Tuple_t, std::size_t... Is>
117auto MapFromTuple(Tuple_t &&t, std::index_sequence<Is...>)
118 -> decltype(MapImpl(std::get<std::tuple_size<Tuple_t>::value - 1>(t), std::get<Is>(t)...))
119{
120 constexpr const auto tupleSizeM1 = std::tuple_size<Tuple_t>::value - 1;
121 return MapImpl(std::get<tupleSizeM1>(t), std::get<Is>(t)...);
122}
123
124/// Return the next power of two (in 64-bits) that is strictly greater than A.
125/// Return zero on overflow.
126inline uint64_t NextPowerOf2(uint64_t A)
127{
128 A |= (A >> 1);
129 A |= (A >> 2);
130 A |= (A >> 4);
131 A |= (A >> 8);
132 A |= (A >> 16);
133 A |= (A >> 32);
134 return A + 1;
135}
136
137/// This is all the stuff common to all SmallVectors.
138class R__CLING_PTRCHECK(off) SmallVectorBase {
139public:
140 // This limits the maximum size of an RVec<char> to ~4GB but we don't expect this to ever be a problem,
141 // and we prefer the smaller Size_T to reduce the size of each RVec object.
142 using Size_T = int32_t;
143
144protected:
145 void *fBeginX;
146 /// Always >= 0.
147 // Type is signed only for consistency with fCapacity.
149 /// Always >= -1. fCapacity == -1 indicates the RVec is in "memory adoption" mode.
151
152 /// The maximum value of the Size_T used.
153 static constexpr size_t SizeTypeMax() { return std::numeric_limits<Size_T>::max(); }
154
155 SmallVectorBase() = delete;
156 SmallVectorBase(void *FirstEl, size_t TotalCapacity) : fBeginX(FirstEl), fCapacity(TotalCapacity) {}
157
158 /// This is an implementation of the grow() method which only works
159 /// on POD-like data types and is out of line to reduce code duplication.
160 /// This function will report a fatal error if it cannot increase capacity.
161 void grow_pod(void *FirstEl, size_t MinSize, size_t TSize);
162
163 /// Report that MinSize doesn't fit into this vector's size type. Throws
164 /// std::length_error or calls report_fatal_error.
165 static void report_size_overflow(size_t MinSize);
166 /// Report that this vector is already at maximum capacity. Throws
167 /// std::length_error or calls report_fatal_error.
168 static void report_at_maximum_capacity();
169
170 /// If false, the RVec is in "memory adoption" mode, i.e. it is acting as a view on a memory buffer it does not own.
171 bool Owns() const { return fCapacity != -1; }
172
173public:
174 size_t size() const { return fSize; }
175 size_t capacity() const noexcept { return Owns() ? fCapacity : fSize; }
176
177 R__RVEC_NODISCARD bool empty() const { return !fSize; }
178
179 /// Set the array size to \p N, which the current array must have enough
180 /// capacity for.
181 ///
182 /// This does not construct or destroy any elements in the vector.
183 ///
184 /// Clients can use this in conjunction with capacity() to write past the end
185 /// of the buffer when they know that more elements are available, and only
186 /// update the size later. This avoids the cost of value initializing elements
187 /// which will only be overwritten.
188 void set_size(size_t N)
189 {
190 if (N > capacity()) {
191 throw std::runtime_error("Setting size to a value greater than capacity.");
192 }
193 fSize = N;
194 }
195};
196
197/// Used to figure out the offset of the first element of an RVec
198template <class T>
200 alignas(SmallVectorBase) char Base[sizeof(SmallVectorBase)];
201 alignas(T) char FirstEl[sizeof(T)];
202};
203
204/// This is the part of SmallVectorTemplateBase which does not depend on whether the type T is a POD.
205template <typename T>
206class R__CLING_PTRCHECK(off) SmallVectorTemplateCommon : public SmallVectorBase {
208
209 /// Find the address of the first element. For this pointer math to be valid
210 /// with small-size of 0 for T with lots of alignment, it's important that
211 /// SmallVectorStorage is properly-aligned even for small-size of 0.
212 void *getFirstEl() const
213 {
214 return const_cast<void *>(reinterpret_cast<const void *>(reinterpret_cast<const char *>(this) +
215 offsetof(SmallVectorAlignmentAndSize<T>, FirstEl)));
216 }
217 // Space after 'FirstEl' is clobbered, do not add any instance vars after it.
218
219protected:
220 SmallVectorTemplateCommon(size_t Size) : Base(getFirstEl(), Size) {}
221
222 void grow_pod(size_t MinSize, size_t TSize) { Base::grow_pod(getFirstEl(), MinSize, TSize); }
223
224 /// Return true if this is a smallvector which has not had dynamic
225 /// memory allocated for it.
226 bool isSmall() const { return this->fBeginX == getFirstEl(); }
227
228 /// Put this vector in a state of being small.
230 {
231 this->fBeginX = getFirstEl();
232 // from the original LLVM implementation:
233 // FIXME: Setting fCapacity to 0 is suspect.
234 this->fSize = this->fCapacity = 0;
235 }
236
237public:
238 // note that fSize is a _signed_ integer, but we expose it as an unsigned integer for consistency with STL containers
239 // as well as backward-compatibility
240 using size_type = size_t;
241 using difference_type = ptrdiff_t;
242 using value_type = T;
243 using iterator = T *;
244 using const_iterator = const T *;
245
246 using const_reverse_iterator = std::reverse_iterator<const_iterator>;
247 using reverse_iterator = std::reverse_iterator<iterator>;
248
249 using reference = T &;
250 using const_reference = const T &;
251 using pointer = T *;
252 using const_pointer = const T *;
253
254 using Base::capacity;
255 using Base::empty;
256 using Base::size;
257
258 // forward iterator creation methods.
259 iterator begin() noexcept { return (iterator)this->fBeginX; }
260 const_iterator begin() const noexcept { return (const_iterator)this->fBeginX; }
261 const_iterator cbegin() const noexcept { return (const_iterator)this->fBeginX; }
262 iterator end() noexcept { return begin() + size(); }
263 const_iterator end() const noexcept { return begin() + size(); }
264 const_iterator cend() const noexcept { return begin() + size(); }
265
266 // reverse iterator creation methods.
267 reverse_iterator rbegin() noexcept { return reverse_iterator(end()); }
268 const_reverse_iterator rbegin() const noexcept { return const_reverse_iterator(end()); }
269 const_reverse_iterator crbegin() const noexcept { return const_reverse_iterator(end()); }
270 reverse_iterator rend() noexcept { return reverse_iterator(begin()); }
271 const_reverse_iterator rend() const noexcept { return const_reverse_iterator(begin()); }
272 const_reverse_iterator crend() const noexcept { return const_reverse_iterator(begin()); }
273
274 size_type size_in_bytes() const { return size() * sizeof(T); }
275 size_type max_size() const noexcept { return std::min(this->SizeTypeMax(), size_type(-1) / sizeof(T)); }
276
277 size_t capacity_in_bytes() const { return capacity() * sizeof(T); }
278
279 /// Return a pointer to the vector's buffer, even if empty().
280 pointer data() noexcept { return pointer(begin()); }
281 /// Return a pointer to the vector's buffer, even if empty().
282 const_pointer data() const noexcept { return const_pointer(begin()); }
283
285 {
286 if (empty()) {
287 throw std::runtime_error("`front` called on an empty RVec");
288 }
289 return begin()[0];
290 }
291
293 {
294 if (empty()) {
295 throw std::runtime_error("`front` called on an empty RVec");
296 }
297 return begin()[0];
298 }
299
301 {
302 if (empty()) {
303 throw std::runtime_error("`back` called on an empty RVec");
304 }
305 return end()[-1];
306 }
307
309 {
310 if (empty()) {
311 throw std::runtime_error("`back` called on an empty RVec");
312 }
313 return end()[-1];
314 }
315};
316
317/// SmallVectorTemplateBase<TriviallyCopyable = false> - This is where we put
318/// method implementations that are designed to work with non-trivial T's.
319///
320/// We approximate is_trivially_copyable with trivial move/copy construction and
321/// trivial destruction. While the standard doesn't specify that you're allowed
322/// copy these types with memcpy, there is no way for the type to observe this.
323/// This catches the important case of std::pair<POD, POD>, which is not
324/// trivially assignable.
325template <typename T, bool = (std::is_trivially_copy_constructible<T>::value) &&
326 (std::is_trivially_move_constructible<T>::value) &&
327 std::is_trivially_destructible<T>::value>
328class R__CLING_PTRCHECK(off) SmallVectorTemplateBase : public SmallVectorTemplateCommon<T> {
329protected:
331
332 static void destroy_range(T *S, T *E)
333 {
334 while (S != E) {
335 --E;
336 E->~T();
337 }
338 }
339
340 /// Move the range [I, E) into the uninitialized memory starting with "Dest",
341 /// constructing elements as needed.
342 template <typename It1, typename It2>
343 static void uninitialized_move(It1 I, It1 E, It2 Dest)
344 {
345 std::uninitialized_copy(std::make_move_iterator(I), std::make_move_iterator(E), Dest);
346 }
347
348 /// Copy the range [I, E) onto the uninitialized memory starting with "Dest",
349 /// constructing elements as needed.
350 template <typename It1, typename It2>
351 static void uninitialized_copy(It1 I, It1 E, It2 Dest)
352 {
353 std::uninitialized_copy(I, E, Dest);
354 }
355
356 /// Grow the allocated memory (without initializing new elements), doubling
357 /// the size of the allocated memory. Guarantees space for at least one more
358 /// element, or MinSize more elements if specified.
359 void grow(size_t MinSize = 0);
360
361public:
362 void push_back(const T &Elt)
363 {
364 if (R__unlikely(this->size() >= this->capacity()))
365 this->grow();
366 ::new ((void *)this->end()) T(Elt);
367 this->set_size(this->size() + 1);
368 }
369
370 void push_back(T &&Elt)
371 {
372 if (R__unlikely(this->size() >= this->capacity()))
373 this->grow();
374 ::new ((void *)this->end()) T(::std::move(Elt));
375 this->set_size(this->size() + 1);
376 }
377
378 void pop_back()
379 {
380 this->set_size(this->size() - 1);
381 this->end()->~T();
382 }
383};
384
385// Define this out-of-line to dissuade the C++ compiler from inlining it.
386template <typename T, bool TriviallyCopyable>
388{
389 // Ensure we can fit the new capacity.
390 // This is only going to be applicable when the capacity is 32 bit.
391 if (MinSize > this->SizeTypeMax())
392 this->report_size_overflow(MinSize);
393
394 // Ensure we can meet the guarantee of space for at least one more element.
395 // The above check alone will not catch the case where grow is called with a
396 // default MinSize of 0, but the current capacity cannot be increased.
397 // This is only going to be applicable when the capacity is 32 bit.
398 if (this->capacity() == this->SizeTypeMax())
399 this->report_at_maximum_capacity();
400
401 // Always grow, even from zero.
402 size_t NewCapacity = size_t(NextPowerOf2(this->capacity() + 2));
403 NewCapacity = std::min(std::max(NewCapacity, MinSize), this->SizeTypeMax());
404 T *NewElts = static_cast<T *>(malloc(NewCapacity * sizeof(T)));
405 R__ASSERT(NewElts != nullptr);
406
407 // Move the elements over.
408 this->uninitialized_move(this->begin(), this->end(), NewElts);
409
410 if (this->Owns()) {
411 // Destroy the original elements.
412 destroy_range(this->begin(), this->end());
413
414 // If this wasn't grown from the inline copy, deallocate the old space.
415 if (!this->isSmall())
416 free(this->begin());
417 }
418
419 this->fBeginX = NewElts;
420 this->fCapacity = NewCapacity;
421}
422
423/// SmallVectorTemplateBase<TriviallyCopyable = true> - This is where we put
424/// method implementations that are designed to work with trivially copyable
425/// T's. This allows using memcpy in place of copy/move construction and
426/// skipping destruction.
427template <typename T>
428class R__CLING_PTRCHECK(off) SmallVectorTemplateBase<T, true> : public SmallVectorTemplateCommon<T> {
430
431protected:
433
434 // No need to do a destroy loop for POD's.
435 static void destroy_range(T *, T *) {}
436
437 /// Move the range [I, E) onto the uninitialized memory
438 /// starting with "Dest", constructing elements into it as needed.
439 template <typename It1, typename It2>
440 static void uninitialized_move(It1 I, It1 E, It2 Dest)
441 {
442 // Just do a copy.
443 uninitialized_copy(I, E, Dest);
444 }
445
446 /// Copy the range [I, E) onto the uninitialized memory
447 /// starting with "Dest", constructing elements into it as needed.
448 template <typename It1, typename It2>
449 static void uninitialized_copy(It1 I, It1 E, It2 Dest)
450 {
451 // Arbitrary iterator types; just use the basic implementation.
452 std::uninitialized_copy(I, E, Dest);
453 }
454
455 /// Copy the range [I, E) onto the uninitialized memory
456 /// starting with "Dest", constructing elements into it as needed.
457 template <typename T1, typename T2>
459 T1 *I, T1 *E, T2 *Dest,
460 typename std::enable_if<std::is_same<typename std::remove_const<T1>::type, T2>::value>::type * = nullptr)
461 {
462 // Use memcpy for PODs iterated by pointers (which includes SmallVector
463 // iterators): std::uninitialized_copy optimizes to memmove, but we can
464 // use memcpy here. Note that I and E are iterators and thus might be
465 // invalid for memcpy if they are equal.
466 if (I != E)
467 memcpy(reinterpret_cast<void *>(Dest), I, (E - I) * sizeof(T));
468 }
469
470 /// Double the size of the allocated memory, guaranteeing space for at
471 /// least one more element or MinSize if specified.
472 void grow(size_t MinSize = 0)
473 {
474 this->grow_pod(MinSize, sizeof(T));
475 }
476
477public:
482
483 void push_back(const T &Elt)
484 {
485 if (R__unlikely(this->size() >= this->capacity()))
486 this->grow();
487 memcpy(reinterpret_cast<void *>(this->end()), &Elt, sizeof(T));
488 this->set_size(this->size() + 1);
489 }
490
491 void pop_back() { this->set_size(this->size() - 1); }
492};
493
494/// Storage for the SmallVector elements. This is specialized for the N=0 case
495/// to avoid allocating unnecessary storage.
496template <typename T, unsigned N>
498 alignas(T) char InlineElts[N * sizeof(T)]{};
499};
500
501/// We need the storage to be properly aligned even for small-size of 0 so that
502/// the pointer math in \a SmallVectorTemplateCommon::getFirstEl() is
503/// well-defined.
504template <typename T>
505struct alignas(T) SmallVectorStorage<T, 0> {
506};
507
508/// The size of the inline storage of an RVec.
509/// Our policy is to allocate at least 8 elements (or more if they all fit into one cacheline)
510/// unless the size of the buffer with 8 elements would be over a certain maximum size.
511template <typename T>
513private:
514#ifdef R__HAS_HARDWARE_INTERFERENCE_SIZE
515 constexpr std::size_t cacheLineSize = std::hardware_destructive_interference_size;
516#else
517 // safe bet: assume the typical 64 bytes
518 static constexpr std::size_t cacheLineSize = 64;
519#endif
520 static constexpr unsigned elementsPerCacheLine = (cacheLineSize - sizeof(SmallVectorBase)) / sizeof(T);
521 static constexpr unsigned maxInlineByteSize = 1024;
522
523public:
524 static constexpr unsigned value =
525 elementsPerCacheLine >= 8 ? elementsPerCacheLine : (sizeof(T) * 8 > maxInlineByteSize ? 0 : 8);
526};
527
528// A C++14-compatible implementation of std::uninitialized_value_construct
529template <typename ForwardIt>
530void UninitializedValueConstruct(ForwardIt first, ForwardIt last)
531{
532#if __cplusplus < 201703L
533 for (; first != last; ++first)
534 new (static_cast<void *>(std::addressof(*first))) typename std::iterator_traits<ForwardIt>::value_type();
535#else
536 std::uninitialized_value_construct(first, last);
537#endif
538}
539
540/// An unsafe function to reset the buffer for which this RVec is acting as a view.
541///
542/// \note This is a low-level method that _must_ be called on RVecs that are already non-owning:
543/// - it does not put the RVec in "non-owning mode" (fCapacity == -1)
544/// - it does not free any owned buffer
545template <typename T>
546void ResetView(RVec<T> &v, T* addr, std::size_t sz)
547{
548 v.fBeginX = addr;
549 v.fSize = sz;
550}
551
552} // namespace VecOps
553} // namespace Internal
554
555namespace Detail {
556namespace VecOps {
557
558/// This class consists of common code factored out of the SmallVector class to
559/// reduce code duplication based on the SmallVector 'N' template parameter.
560template <typename T>
561class R__CLING_PTRCHECK(off) RVecImpl : public Internal::VecOps::SmallVectorTemplateBase<T> {
563
564public:
569
570protected:
571 // Default ctor - Initialize to empty.
572 explicit RVecImpl(unsigned N) : ROOT::Internal::VecOps::SmallVectorTemplateBase<T>(N) {}
573
574public:
575 RVecImpl(const RVecImpl &) = delete;
576
578 {
579 // Subclass has already destructed this vector's elements.
580 // If this wasn't grown from the inline copy, deallocate the old space.
581 if (!this->isSmall() && this->Owns())
582 free(this->begin());
583 }
584
585 // also give up adopted memory if applicable
586 void clear()
587 {
588 if (this->Owns()) {
589 this->destroy_range(this->begin(), this->end());
590 this->fSize = 0;
591 } else {
592 this->resetToSmall();
593 }
594 }
595
597 {
598 if (N < this->size()) {
599 if (this->Owns())
600 this->destroy_range(this->begin() + N, this->end());
601 this->set_size(N);
602 } else if (N > this->size()) {
603 if (this->capacity() < N)
604 this->grow(N);
605 for (auto I = this->end(), E = this->begin() + N; I != E; ++I)
606 new (&*I) T();
607 this->set_size(N);
608 }
609 }
610
611 void resize(size_type N, const T &NV)
612 {
613 if (N < this->size()) {
614 if (this->Owns())
615 this->destroy_range(this->begin() + N, this->end());
616 this->set_size(N);
617 } else if (N > this->size()) {
618 if (this->capacity() < N)
619 this->grow(N);
620 std::uninitialized_fill(this->end(), this->begin() + N, NV);
621 this->set_size(N);
622 }
623 }
624
626 {
627 if (this->capacity() < N)
628 this->grow(N);
629 }
630
631 void pop_back_n(size_type NumItems)
632 {
633 if (this->size() < NumItems) {
634 throw std::runtime_error("Popping back more elements than those available.");
635 }
636 if (this->Owns())
637 this->destroy_range(this->end() - NumItems, this->end());
638 this->set_size(this->size() - NumItems);
639 }
640
642 {
643 T Result = ::std::move(this->back());
644 this->pop_back();
645 return Result;
646 }
647
648 void swap(RVecImpl &RHS);
649
650 /// Add the specified range to the end of the SmallVector.
651 template <typename in_iter,
652 typename = typename std::enable_if<std::is_convertible<
653 typename std::iterator_traits<in_iter>::iterator_category, std::input_iterator_tag>::value>::type>
654 void append(in_iter in_start, in_iter in_end)
655 {
656 size_type NumInputs = std::distance(in_start, in_end);
657 if (NumInputs > this->capacity() - this->size())
658 this->grow(this->size() + NumInputs);
659
660 this->uninitialized_copy(in_start, in_end, this->end());
661 this->set_size(this->size() + NumInputs);
662 }
663
664 /// Append \p NumInputs copies of \p Elt to the end.
665 void append(size_type NumInputs, const T &Elt)
666 {
667 if (NumInputs > this->capacity() - this->size())
668 this->grow(this->size() + NumInputs);
669
670 std::uninitialized_fill_n(this->end(), NumInputs, Elt);
671 this->set_size(this->size() + NumInputs);
672 }
673
674 void append(std::initializer_list<T> IL) { append(IL.begin(), IL.end()); }
675
676 // from the original LLVM implementation:
677 // FIXME: Consider assigning over existing elements, rather than clearing &
678 // re-initializing them - for all assign(...) variants.
679
680 void assign(size_type NumElts, const T &Elt)
681 {
682 clear();
683 if (this->capacity() < NumElts)
684 this->grow(NumElts);
685 this->set_size(NumElts);
686 std::uninitialized_fill(this->begin(), this->end(), Elt);
687 }
688
689 template <typename in_iter,
690 typename = typename std::enable_if<std::is_convertible<
691 typename std::iterator_traits<in_iter>::iterator_category, std::input_iterator_tag>::value>::type>
692 void assign(in_iter in_start, in_iter in_end)
693 {
694 clear();
695 append(in_start, in_end);
696 }
697
698 void assign(std::initializer_list<T> IL)
699 {
700 clear();
701 append(IL);
702 }
703
705 {
706 // Just cast away constness because this is a non-const member function.
707 iterator I = const_cast<iterator>(CI);
708
709 if (I < this->begin() || I >= this->end()) {
710 throw std::runtime_error("The iterator passed to `erase` is out of bounds.");
711 }
712
713 iterator N = I;
714 // Shift all elts down one.
715 std::move(I + 1, this->end(), I);
716 // Drop the last elt.
717 this->pop_back();
718 return (N);
719 }
720
722 {
723 // Just cast away constness because this is a non-const member function.
724 iterator S = const_cast<iterator>(CS);
725 iterator E = const_cast<iterator>(CE);
726
727 if (S < this->begin() || E > this->end() || S > E) {
728 throw std::runtime_error("Invalid start/end pair passed to `erase` (out of bounds or start > end).");
729 }
730
731 iterator N = S;
732 // Shift all elts down.
733 iterator I = std::move(E, this->end(), S);
734 // Drop the last elts.
735 if (this->Owns())
736 this->destroy_range(I, this->end());
737 this->set_size(I - this->begin());
738 return (N);
739 }
740
742 {
743 if (I == this->end()) { // Important special case for empty vector.
744 this->push_back(::std::move(Elt));
745 return this->end() - 1;
746 }
747
748 if (I < this->begin() || I > this->end()) {
749 throw std::runtime_error("The iterator passed to `insert` is out of bounds.");
750 }
751
752 if (this->size() >= this->capacity()) {
753 size_t EltNo = I - this->begin();
754 this->grow();
755 I = this->begin() + EltNo;
756 }
757
758 ::new ((void *)this->end()) T(::std::move(this->back()));
759 // Push everything else over.
760 std::move_backward(I, this->end() - 1, this->end());
761 this->set_size(this->size() + 1);
762
763 // If we just moved the element we're inserting, be sure to update
764 // the reference.
765 T *EltPtr = &Elt;
766 if (I <= EltPtr && EltPtr < this->end())
767 ++EltPtr;
768
769 *I = ::std::move(*EltPtr);
770 return I;
771 }
772
773 iterator insert(iterator I, const T &Elt)
774 {
775 if (I == this->end()) { // Important special case for empty vector.
776 this->push_back(Elt);
777 return this->end() - 1;
778 }
779
780 if (I < this->begin() || I > this->end()) {
781 throw std::runtime_error("The iterator passed to `insert` is out of bounds.");
782 }
783
784 if (this->size() >= this->capacity()) {
785 size_t EltNo = I - this->begin();
786 this->grow();
787 I = this->begin() + EltNo;
788 }
789 ::new ((void *)this->end()) T(std::move(this->back()));
790 // Push everything else over.
791 std::move_backward(I, this->end() - 1, this->end());
792 this->set_size(this->size() + 1);
793
794 // If we just moved the element we're inserting, be sure to update
795 // the reference.
796 const T *EltPtr = &Elt;
797 if (I <= EltPtr && EltPtr < this->end())
798 ++EltPtr;
799
800 *I = *EltPtr;
801 return I;
802 }
803
804 iterator insert(iterator I, size_type NumToInsert, const T &Elt)
805 {
806 // Convert iterator to elt# to avoid invalidating iterator when we reserve()
807 size_t InsertElt = I - this->begin();
808
809 if (I == this->end()) { // Important special case for empty vector.
810 append(NumToInsert, Elt);
811 return this->begin() + InsertElt;
812 }
813
814 if (I < this->begin() || I > this->end()) {
815 throw std::runtime_error("The iterator passed to `insert` is out of bounds.");
816 }
817
818 // Ensure there is enough space.
819 reserve(this->size() + NumToInsert);
820
821 // Uninvalidate the iterator.
822 I = this->begin() + InsertElt;
823
824 // If there are more elements between the insertion point and the end of the
825 // range than there are being inserted, we can use a simple approach to
826 // insertion. Since we already reserved space, we know that this won't
827 // reallocate the vector.
828 if (size_t(this->end() - I) >= NumToInsert) {
829 T *OldEnd = this->end();
830 append(std::move_iterator<iterator>(this->end() - NumToInsert), std::move_iterator<iterator>(this->end()));
831
832 // Copy the existing elements that get replaced.
833 std::move_backward(I, OldEnd - NumToInsert, OldEnd);
834
835 std::fill_n(I, NumToInsert, Elt);
836 return I;
837 }
838
839 // Otherwise, we're inserting more elements than exist already, and we're
840 // not inserting at the end.
841
842 // Move over the elements that we're about to overwrite.
843 T *OldEnd = this->end();
844 this->set_size(this->size() + NumToInsert);
845 size_t NumOverwritten = OldEnd - I;
846 this->uninitialized_move(I, OldEnd, this->end() - NumOverwritten);
847
848 // Replace the overwritten part.
849 std::fill_n(I, NumOverwritten, Elt);
850
851 // Insert the non-overwritten middle part.
852 std::uninitialized_fill_n(OldEnd, NumToInsert - NumOverwritten, Elt);
853 return I;
854 }
855
856 template <typename ItTy,
857 typename = typename std::enable_if<std::is_convertible<
858 typename std::iterator_traits<ItTy>::iterator_category, std::input_iterator_tag>::value>::type>
859 iterator insert(iterator I, ItTy From, ItTy To)
860 {
861 // Convert iterator to elt# to avoid invalidating iterator when we reserve()
862 size_t InsertElt = I - this->begin();
863
864 if (I == this->end()) { // Important special case for empty vector.
865 append(From, To);
866 return this->begin() + InsertElt;
867 }
868
869 if (I < this->begin() || I > this->end()) {
870 throw std::runtime_error("The iterator passed to `insert` is out of bounds.");
871 }
872
873 size_t NumToInsert = std::distance(From, To);
874
875 // Ensure there is enough space.
876 reserve(this->size() + NumToInsert);
877
878 // Uninvalidate the iterator.
879 I = this->begin() + InsertElt;
880
881 // If there are more elements between the insertion point and the end of the
882 // range than there are being inserted, we can use a simple approach to
883 // insertion. Since we already reserved space, we know that this won't
884 // reallocate the vector.
885 if (size_t(this->end() - I) >= NumToInsert) {
886 T *OldEnd = this->end();
887 append(std::move_iterator<iterator>(this->end() - NumToInsert), std::move_iterator<iterator>(this->end()));
888
889 // Copy the existing elements that get replaced.
890 std::move_backward(I, OldEnd - NumToInsert, OldEnd);
891
892 std::copy(From, To, I);
893 return I;
894 }
895
896 // Otherwise, we're inserting more elements than exist already, and we're
897 // not inserting at the end.
898
899 // Move over the elements that we're about to overwrite.
900 T *OldEnd = this->end();
901 this->set_size(this->size() + NumToInsert);
902 size_t NumOverwritten = OldEnd - I;
903 this->uninitialized_move(I, OldEnd, this->end() - NumOverwritten);
904
905 // Replace the overwritten part.
906 for (T *J = I; NumOverwritten > 0; --NumOverwritten) {
907 *J = *From;
908 ++J;
909 ++From;
910 }
911
912 // Insert the non-overwritten middle part.
913 this->uninitialized_copy(From, To, OldEnd);
914 return I;
915 }
916
917 void insert(iterator I, std::initializer_list<T> IL) { insert(I, IL.begin(), IL.end()); }
918
919 template <typename... ArgTypes>
920 reference emplace_back(ArgTypes &&...Args)
921 {
922 if (R__unlikely(this->size() >= this->capacity()))
923 this->grow();
924 ::new ((void *)this->end()) T(std::forward<ArgTypes>(Args)...);
925 this->set_size(this->size() + 1);
926 return this->back();
927 }
928
930
932};
933
934template <typename T>
936{
937 if (this == &RHS)
938 return;
939
940 // We can only avoid copying elements if neither vector is small.
941 if (!this->isSmall() && !RHS.isSmall()) {
942 std::swap(this->fBeginX, RHS.fBeginX);
943 std::swap(this->fSize, RHS.fSize);
944 std::swap(this->fCapacity, RHS.fCapacity);
945 return;
946 }
947
948 // This block handles the swap of a small and a non-owning vector
949 // It is more efficient to first move the non-owning vector, hence the 2 cases
950 if (this->isSmall() && !RHS.Owns()) { // the right vector is non-owning
951 RVecImpl<T> temp(0);
952 temp = std::move(RHS);
953 RHS = std::move(*this);
954 *this = std::move(temp);
955 return;
956 } else if (RHS.isSmall() && !this->Owns()) { // the left vector is non-owning
957 RVecImpl<T> temp(0);
958 temp = std::move(*this);
959 *this = std::move(RHS);
960 RHS = std::move(temp);
961 return;
962 }
963
964 if (RHS.size() > this->capacity())
965 this->grow(RHS.size());
966 if (this->size() > RHS.capacity())
967 RHS.grow(this->size());
968
969 // Swap the shared elements.
970 size_t NumShared = this->size();
971 if (NumShared > RHS.size())
972 NumShared = RHS.size();
973 for (size_type i = 0; i != NumShared; ++i)
974 std::iter_swap(this->begin() + i, RHS.begin() + i);
975
976 // Copy over the extra elts.
977 if (this->size() > RHS.size()) {
978 size_t EltDiff = this->size() - RHS.size();
979 this->uninitialized_copy(this->begin() + NumShared, this->end(), RHS.end());
980 RHS.set_size(RHS.size() + EltDiff);
981 if (this->Owns())
982 this->destroy_range(this->begin() + NumShared, this->end());
983 this->set_size(NumShared);
984 } else if (RHS.size() > this->size()) {
985 size_t EltDiff = RHS.size() - this->size();
986 this->uninitialized_copy(RHS.begin() + NumShared, RHS.end(), this->end());
987 this->set_size(this->size() + EltDiff);
988 if (RHS.Owns())
989 this->destroy_range(RHS.begin() + NumShared, RHS.end());
990 RHS.set_size(NumShared);
991 }
992}
993
994template <typename T>
996{
997 // Avoid self-assignment.
998 if (this == &RHS)
999 return *this;
1000
1001 // If we already have sufficient space, assign the common elements, then
1002 // destroy any excess.
1003 size_t RHSSize = RHS.size();
1004 size_t CurSize = this->size();
1005 if (CurSize >= RHSSize) {
1006 // Assign common elements.
1007 iterator NewEnd;
1008 if (RHSSize)
1009 NewEnd = std::copy(RHS.begin(), RHS.begin() + RHSSize, this->begin());
1010 else
1011 NewEnd = this->begin();
1012
1013 // Destroy excess elements.
1014 if (this->Owns())
1015 this->destroy_range(NewEnd, this->end());
1016
1017 // Trim.
1018 this->set_size(RHSSize);
1019 return *this;
1020 }
1021
1022 // If we have to grow to have enough elements, destroy the current elements.
1023 // This allows us to avoid copying them during the grow.
1024 // From the original LLVM implementation:
1025 // FIXME: don't do this if they're efficiently moveable.
1026 if (this->capacity() < RHSSize) {
1027 if (this->Owns()) {
1028 // Destroy current elements.
1029 this->destroy_range(this->begin(), this->end());
1030 }
1031 this->set_size(0);
1032 CurSize = 0;
1033 this->grow(RHSSize);
1034 } else if (CurSize) {
1035 // Otherwise, use assignment for the already-constructed elements.
1036 std::copy(RHS.begin(), RHS.begin() + CurSize, this->begin());
1037 }
1038
1039 // Copy construct the new elements in place.
1040 this->uninitialized_copy(RHS.begin() + CurSize, RHS.end(), this->begin() + CurSize);
1041
1042 // Set end.
1043 this->set_size(RHSSize);
1044 return *this;
1045}
1046
1047template <typename T>
1049{
1050 // Avoid self-assignment.
1051 if (this == &RHS)
1052 return *this;
1053
1054 // If the RHS isn't small, clear this vector and then steal its buffer.
1055 if (!RHS.isSmall()) {
1056 if (this->Owns()) {
1057 this->destroy_range(this->begin(), this->end());
1058 if (!this->isSmall())
1059 free(this->begin());
1060 }
1061 this->fBeginX = RHS.fBeginX;
1062 this->fSize = RHS.fSize;
1063 this->fCapacity = RHS.fCapacity;
1064 RHS.resetToSmall();
1065 return *this;
1066 }
1067
1068 // If we already have sufficient space, assign the common elements, then
1069 // destroy any excess.
1070 size_t RHSSize = RHS.size();
1071 size_t CurSize = this->size();
1072 if (CurSize >= RHSSize) {
1073 // Assign common elements.
1074 iterator NewEnd = this->begin();
1075 if (RHSSize)
1076 NewEnd = std::move(RHS.begin(), RHS.end(), NewEnd);
1077
1078 // Destroy excess elements and trim the bounds.
1079 if (this->Owns())
1080 this->destroy_range(NewEnd, this->end());
1081 this->set_size(RHSSize);
1082
1083 // Clear the RHS.
1084 RHS.clear();
1085
1086 return *this;
1087 }
1088
1089 // If we have to grow to have enough elements, destroy the current elements.
1090 // This allows us to avoid copying them during the grow.
1091 // From the original LLVM implementation:
1092 // FIXME: this may not actually make any sense if we can efficiently move
1093 // elements.
1094 if (this->capacity() < RHSSize) {
1095 if (this->Owns()) {
1096 // Destroy current elements.
1097 this->destroy_range(this->begin(), this->end());
1098 }
1099 this->set_size(0);
1100 CurSize = 0;
1101 this->grow(RHSSize);
1102 } else if (CurSize) {
1103 // Otherwise, use assignment for the already-constructed elements.
1104 std::move(RHS.begin(), RHS.begin() + CurSize, this->begin());
1105 }
1106
1107 // Move-construct the new elements in place.
1108 this->uninitialized_move(RHS.begin() + CurSize, RHS.end(), this->begin() + CurSize);
1109
1110 // Set end.
1111 this->set_size(RHSSize);
1112
1113 RHS.clear();
1114 return *this;
1115}
1116
1117template <typename T>
1119{
1120 return v.isSmall();
1121}
1122
1123template <typename T>
1125{
1126 return !v.Owns();
1127}
1128
1129} // namespace VecOps
1130} // namespace Detail
1131
1132namespace VecOps {
1133// Note that we open here with @{ the Doxygen group vecops and it is
1134// closed again at the end of the C++ namespace VecOps
1135/**
1136 * \defgroup vecops VecOps
1137 * A "std::vector"-like collection of values implementing handy operation to analyse them
1138 * @{
1139*/
1140
1141// From the original SmallVector code:
1142// This is a 'vector' (really, a variable-sized array), optimized
1143// for the case when the array is small. It contains some number of elements
1144// in-place, which allows it to avoid heap allocation when the actual number of
1145// elements is below that threshold. This allows normal "small" cases to be
1146// fast without losing generality for large inputs.
1147//
1148// Note that this does not attempt to be exception safe.
1149
1150template <typename T, unsigned int N>
1151class R__CLING_PTRCHECK(off) RVecN : public Detail::VecOps::RVecImpl<T>, Internal::VecOps::SmallVectorStorage<T, N> {
1152public:
1153 RVecN() : Detail::VecOps::RVecImpl<T>(N) {}
1154
1156 {
1157 if (this->Owns()) {
1158 // Destroy the constructed elements in the vector.
1159 this->destroy_range(this->begin(), this->end());
1160 }
1161 }
1162
1163 explicit RVecN(size_t Size, const T &Value) : Detail::VecOps::RVecImpl<T>(N) { this->assign(Size, Value); }
1164
1165 explicit RVecN(size_t Size) : Detail::VecOps::RVecImpl<T>(N)
1166 {
1167 if (Size > N)
1168 this->grow(Size);
1169 this->fSize = Size;
1170 ROOT::Internal::VecOps::UninitializedValueConstruct(this->begin(), this->end());
1171 }
1172
1173 template <typename ItTy,
1174 typename = typename std::enable_if<std::is_convertible<
1175 typename std::iterator_traits<ItTy>::iterator_category, std::input_iterator_tag>::value>::type>
1176 RVecN(ItTy S, ItTy E) : Detail::VecOps::RVecImpl<T>(N)
1177 {
1178 this->append(S, E);
1179 }
1180
1181 RVecN(std::initializer_list<T> IL) : Detail::VecOps::RVecImpl<T>(N) { this->assign(IL); }
1182
1183 RVecN(const RVecN &RHS) : Detail::VecOps::RVecImpl<T>(N)
1184 {
1185 if (!RHS.empty())
1187 }
1188
1189 RVecN &operator=(const RVecN &RHS)
1190 {
1192 return *this;
1193 }
1194
1195 RVecN(RVecN &&RHS) : Detail::VecOps::RVecImpl<T>(N)
1196 {
1197 if (!RHS.empty())
1199 }
1200
1201 RVecN(Detail::VecOps::RVecImpl<T> &&RHS) : Detail::VecOps::RVecImpl<T>(N)
1202 {
1203 if (!RHS.empty())
1205 }
1206
1207 RVecN(const std::vector<T> &RHS) : RVecN(RHS.begin(), RHS.end()) {}
1208
1210 {
1212 return *this;
1213 }
1214
1215 RVecN(T* p, size_t n) : Detail::VecOps::RVecImpl<T>(N)
1216 {
1217 this->fBeginX = p;
1218 this->fSize = n;
1219 this->fCapacity = -1;
1220 }
1221
1223 {
1225 return *this;
1226 }
1227
1228 RVecN &operator=(std::initializer_list<T> IL)
1229 {
1230 this->assign(IL);
1231 return *this;
1232 }
1233
1240
1242 {
1243 return begin()[idx];
1244 }
1245
1247 {
1248 return begin()[idx];
1249 }
1250
1251 template <typename V, unsigned M, typename = std::enable_if<std::is_convertible<V, bool>::value>>
1252 RVecN operator[](const RVecN<V, M> &conds) const
1253 {
1254 const size_type n = conds.size();
1255
1256 if (n != this->size()) {
1257 std::string msg = "Cannot index RVecN of size " + std::to_string(this->size()) +
1258 " with condition vector of different size (" + std::to_string(n) + ").";
1259 throw std::runtime_error(std::move(msg));
1260 }
1261
1262 size_type n_true = 0ull;
1263 for (auto c : conds)
1264 n_true += c; // relies on bool -> int conversion, faster than branching
1265
1266 RVecN ret;
1267 ret.reserve(n_true);
1268 size_type j = 0u;
1269 for (size_type i = 0u; i < n; ++i) {
1270 if (conds[i]) {
1271 ret.push_back(this->operator[](i));
1272 ++j;
1273 }
1274 }
1275 return ret;
1276 }
1277
1278 // conversion
1279 template <typename U, unsigned M, typename = std::enable_if<std::is_convertible<T, U>::value>>
1280 operator RVecN<U, M>() const
1281 {
1282 return RVecN<U, M>(this->begin(), this->end());
1283 }
1284
1286 {
1287 if (pos >= size_type(this->fSize)) {
1288 std::string msg = "RVecN::at: size is " + std::to_string(this->fSize) + " but out-of-bounds index " +
1289 std::to_string(pos) + " was requested.";
1290 throw std::out_of_range(std::move(msg));
1291 }
1292 return this->operator[](pos);
1293 }
1294
1296 {
1297 if (pos >= size_type(this->fSize)) {
1298 std::string msg = "RVecN::at: size is " + std::to_string(this->fSize) + " but out-of-bounds index " +
1299 std::to_string(pos) + " was requested.";
1300 throw std::out_of_range(std::move(msg));
1301 }
1302 return this->operator[](pos);
1303 }
1304
1305 /// No exception thrown. The user specifies the desired value in case the RVecN is shorter than `pos`.
1307 {
1308 if (pos >= size_type(this->fSize))
1309 return fallback;
1310 return this->operator[](pos);
1311 }
1312
1313 /// No exception thrown. The user specifies the desired value in case the RVecN is shorter than `pos`.
1314 value_type at(size_type pos, value_type fallback) const
1315 {
1316 if (pos >= size_type(this->fSize))
1317 return fallback;
1318 return this->operator[](pos);
1319 }
1320};
1321
1322// clang-format off
1323/**
1324\class ROOT::VecOps::RVec
1325\brief A "std::vector"-like collection of values implementing handy operation to analyse them
1326\tparam T The type of the contained objects
1327
1328A RVec is a container designed to make analysis of values' collections fast and easy.
1329Its storage is contiguous in memory and its interface is designed such to resemble to the one
1330of the stl vector. In addition the interface features methods and
1331[external functions](https://root.cern/doc/master/namespaceROOT_1_1VecOps.html) to ease the manipulation and analysis
1332of the data in the RVec.
1333
1334\note ROOT::VecOps::RVec can also be spelled simply ROOT::RVec. Shorthand aliases such as ROOT::RVecI or ROOT::RVecD
1335are also available as template instantiations of RVec of fundamental types. The full list of available aliases:
1336- RVecB (`bool`)
1337- RVecC (`char`)
1338- RVecD (`double`)
1339- RVecF (`float`)
1340- RVecI (`int`)
1341- RVecL (`long`)
1342- RVecLL (`long long`)
1343- RVecU (`unsigned`)
1344- RVecUL (`unsigned long`)
1345- RVecULL (`unsigned long long`)
1346
1347\note RVec does not attempt to be exception safe. Exceptions thrown by element constructors during insertions, swaps or
1348other operations will be propagated potentially leaving the RVec object in an invalid state.
1349
1350\note RVec methods (e.g. `at` or `size`) follow the STL naming convention instead of the ROOT naming convention in order
1351to make RVec a drop-in replacement for `std::vector`.
1352
1353\htmlonly
1354<a href="https://doi.org/10.5281/zenodo.1253756"><img src="https://zenodo.org/badge/DOI/10.5281/zenodo.1253756.svg" alt="DOI"></a>
1355\endhtmlonly
1356
1357## Table of Contents
1358- [Example](\ref example)
1359- [Owning and adopting memory](\ref owningandadoptingmemory)
1360- [Sorting and manipulation of indices](\ref sorting)
1361- [Usage in combination with RDataFrame](\ref usagetdataframe)
1362- [Reference for the RVec class](\ref RVecdoxyref)
1363- [Reference for RVec helper functions](https://root.cern/doc/master/namespaceROOT_1_1VecOps.html)
1364
1365\anchor example
1366## Example
1367Suppose to have an event featuring a collection of muons with a certain pseudorapidity,
1368momentum and charge, e.g.:
1369~~~{.cpp}
1370std::vector<short> mu_charge {1, 1, -1, -1, -1, 1, 1, -1};
1371std::vector<float> mu_pt {56, 45, 32, 24, 12, 8, 7, 6.2};
1372std::vector<float> mu_eta {3.1, -.2, -1.1, 1, 4.1, 1.6, 2.4, -.5};
1373~~~
1374Suppose you want to extract the transverse momenta of the muons satisfying certain
1375criteria, for example consider only negatively charged muons with a pseudorapidity
1376smaller or equal to 2 and with a transverse momentum greater than 10 GeV.
1377Such a selection would require, among the other things, the management of an explicit
1378loop, for example:
1379~~~{.cpp}
1380std::vector<float> goodMuons_pt;
1381const auto size = mu_charge.size();
1382for (size_t i=0; i < size; ++i) {
1383 if (mu_pt[i] > 10 && abs(mu_eta[i]) <= 2. && mu_charge[i] == -1) {
1384 goodMuons_pt.emplace_back(mu_pt[i]);
1385 }
1386}
1387~~~
1388These operations become straightforward with RVec - we just need to *write what
1389we mean*:
1390~~~{.cpp}
1391auto goodMuons_pt = mu_pt[ (mu_pt > 10.f && abs(mu_eta) <= 2.f && mu_charge == -1) ]
1392~~~
1393Now the clean collection of transverse momenta can be used within the rest of the data analysis, for
1394example to fill a histogram.
1395
1396\anchor owningandadoptingmemory
1397## Owning and adopting memory
1398RVec has contiguous memory associated to it. It can own it or simply adopt it. In the latter case,
1399it can be constructed with the address of the memory associated to it and its length. For example:
1400~~~{.cpp}
1401std::vector<int> myStlVec {1,2,3};
1402RVec<int> myRVec(myStlVec.data(), myStlVec.size());
1403~~~
1404In this case, the memory associated to myStlVec and myRVec is the same, myRVec simply "adopted it".
1405If any method which implies a re-allocation is called, e.g. *emplace_back* or *resize*, the adopted
1406memory is released and new one is allocated. The previous content is copied in the new memory and
1407preserved.
1408
1409\anchor sorting
1410## Sorting and manipulation of indices
1411
1412### Sorting
1413RVec complies to the STL interfaces when it comes to iterations. As a result, standard algorithms
1414can be used, for example sorting:
1415~~~{.cpp}
1416RVec<double> v{6., 4., 5.};
1417std::sort(v.begin(), v.end());
1418~~~
1419
1420For convenience, helpers are provided too:
1421~~~{.cpp}
1422auto sorted_v = Sort(v);
1423auto reversed_v = Reverse(v);
1424~~~
1425
1426### Manipulation of indices
1427
1428It is also possible to manipulated the RVecs acting on their indices. For example,
1429the following syntax
1430~~~{.cpp}
1431RVecD v0 {9., 7., 8.};
1432auto v1 = Take(v0, {1, 2, 0});
1433~~~
1434will yield a new RVec<double> the content of which is the first, second and zeroth element of
1435v0, i.e. `{7., 8., 9.}`.
1436
1437The `Argsort` and `StableArgsort` helper extracts the indices which order the content of a `RVec`.
1438For example, this snippet accomplishes in a more expressive way what we just achieved:
1439~~~{.cpp}
1440auto v1_indices = Argsort(v0); // The content of v1_indices is {1, 2, 0}.
1441v1 = Take(v0, v1_indices);
1442~~~
1443
1444The `Take` utility allows to extract portions of the `RVec`. The content to be *taken*
1445can be specified with an `RVec` of indices or an integer. If the integer is negative,
1446elements will be picked starting from the end of the container:
1447~~~{.cpp}
1448RVecF vf {1.f, 2.f, 3.f, 4.f};
1449auto vf_1 = Take(vf, {1, 3}); // The content is {2.f, 4.f}
1450auto vf_2 = Take(vf, 2); // The content is {1.f, 2.f}
1451auto vf_3 = Take(vf, -3); // The content is {2.f, 3.f, 4.f}
1452~~~
1453
1454\anchor usagetdataframe
1455## Usage in combination with RDataFrame
1456RDataFrame leverages internally RVecs. Suppose to have a dataset stored in a
1457TTree which holds these columns (here we choose C arrays to represent the
1458collections, they could be as well std::vector instances):
1459~~~{.bash}
1460 nPart "nPart/I" An integer representing the number of particles
1461 px "px[nPart]/D" The C array of the particles' x component of the momentum
1462 py "py[nPart]/D" The C array of the particles' y component of the momentum
1463 E "E[nPart]/D" The C array of the particles' Energy
1464~~~
1465Suppose you'd like to plot in a histogram the transverse momenta of all particles
1466for which the energy is greater than 200 MeV.
1467The code required would just be:
1468~~~{.cpp}
1469RDataFrame d("mytree", "myfile.root");
1470auto cutPt = [](RVecD &pxs, RVecD &pys, RVecD &Es) {
1471 auto all_pts = sqrt(pxs * pxs + pys * pys);
1472 auto good_pts = all_pts[Es > 200.];
1473 return good_pts;
1474 };
1475
1476auto hpt = d.Define("pt", cutPt, {"px", "py", "E"})
1477 .Histo1D("pt");
1478hpt->Draw();
1479~~~
1480And if you'd like to express your selection as a string:
1481~~~{.cpp}
1482RDataFrame d("mytree", "myfile.root");
1483auto hpt = d.Define("pt", "sqrt(pxs * pxs + pys * pys)[E>200]")
1484 .Histo1D("pt");
1485hpt->Draw();
1486~~~
1487\anchor RVecdoxyref
1488**/
1489// clang-format on
1490
1491template <typename T>
1492class R__CLING_PTRCHECK(off) RVec : public RVecN<T, Internal::VecOps::RVecInlineStorageSize<T>::value> {
1494
1495 friend void Internal::VecOps::ResetView<>(RVec<T> &v, T *addr, std::size_t sz);
1496
1497public:
1502 using SuperClass::begin;
1503 using SuperClass::size;
1504
1505 RVec() {}
1506
1507 explicit RVec(size_t Size, const T &Value) : SuperClass(Size, Value) {}
1508
1509 explicit RVec(size_t Size) : SuperClass(Size) {}
1510
1511 template <typename ItTy,
1512 typename = typename std::enable_if<std::is_convertible<
1513 typename std::iterator_traits<ItTy>::iterator_category, std::input_iterator_tag>::value>::type>
1514 RVec(ItTy S, ItTy E) : SuperClass(S, E)
1515 {
1516 }
1517
1518 RVec(std::initializer_list<T> IL) : SuperClass(IL) {}
1519
1520 RVec(const RVec &RHS) : SuperClass(RHS) {}
1521
1522 RVec &operator=(const RVec &RHS)
1523 {
1524 SuperClass::operator=(RHS);
1525 return *this;
1526 }
1527
1528 RVec(RVec &&RHS) : SuperClass(std::move(RHS)) {}
1529
1531 {
1532 SuperClass::operator=(std::move(RHS));
1533 return *this;
1534 }
1535
1536 RVec(Detail::VecOps::RVecImpl<T> &&RHS) : SuperClass(std::move(RHS)) {}
1537
1538 template <unsigned N>
1539 RVec(RVecN<T, N> &&RHS) : SuperClass(std::move(RHS)) {}
1540
1541 template <unsigned N>
1542 RVec(const RVecN<T, N> &RHS) : SuperClass(RHS) {}
1543
1544 RVec(const std::vector<T> &RHS) : SuperClass(RHS) {}
1545
1546 RVec(T* p, size_t n) : SuperClass(p, n) {}
1547
1548 // conversion
1549 template <typename U, typename = std::enable_if<std::is_convertible<T, U>::value>>
1550 operator RVec<U>() const
1551 {
1552 return RVec<U>(this->begin(), this->end());
1553 }
1554
1555 using SuperClass::operator[];
1556
1557 template <typename V, typename = std::enable_if<std::is_convertible<V, bool>::value>>
1558 RVec operator[](const RVec<V> &conds) const
1559 {
1560 return RVec(SuperClass::operator[](conds));
1561 }
1562
1563 using SuperClass::at;
1564
1565 friend bool ROOT::Detail::VecOps::IsSmall<T>(const RVec<T> &v);
1566
1567 friend bool ROOT::Detail::VecOps::IsAdopting<T>(const RVec<T> &v);
1568};
1569
1570template <typename T, unsigned N>
1571inline size_t CapacityInBytes(const RVecN<T, N> &X)
1572{
1573 return X.capacity_in_bytes();
1574}
1575
1576///@name RVec Unary Arithmetic Operators
1577///@{
1578
1579#define RVEC_UNARY_OPERATOR(OP) \
1580template <typename T> \
1581RVec<T> operator OP(const RVec<T> &v) \
1582{ \
1583 RVec<T> ret(v); \
1584 for (auto &x : ret) \
1585 x = OP x; \
1586return ret; \
1587} \
1588
1593#undef RVEC_UNARY_OPERATOR
1594
1595///@}
1596///@name RVec Binary Arithmetic Operators
1597///@{
1598
1599#define ERROR_MESSAGE(OP) \
1600 "Cannot call operator " #OP " on vectors of different sizes."
1601
1602#define RVEC_BINARY_OPERATOR(OP) \
1603template <typename T0, typename T1> \
1604auto operator OP(const RVec<T0> &v, const T1 &y) \
1605 -> RVec<decltype(v[0] OP y)> \
1606{ \
1607 RVec<decltype(v[0] OP y)> ret(v.size()); \
1608 auto op = [&y](const T0 &x) { return x OP y; }; \
1609 std::transform(v.begin(), v.end(), ret.begin(), op); \
1610 return ret; \
1611} \
1612 \
1613template <typename T0, typename T1> \
1614auto operator OP(const T0 &x, const RVec<T1> &v) \
1615 -> RVec<decltype(x OP v[0])> \
1616{ \
1617 RVec<decltype(x OP v[0])> ret(v.size()); \
1618 auto op = [&x](const T1 &y) { return x OP y; }; \
1619 std::transform(v.begin(), v.end(), ret.begin(), op); \
1620 return ret; \
1621} \
1622 \
1623template <typename T0, typename T1> \
1624auto operator OP(const RVec<T0> &v0, const RVec<T1> &v1) \
1625 -> RVec<decltype(v0[0] OP v1[0])> \
1626{ \
1627 if (v0.size() != v1.size()) \
1628 throw std::runtime_error(ERROR_MESSAGE(OP)); \
1629 \
1630 RVec<decltype(v0[0] OP v1[0])> ret(v0.size()); \
1631 auto op = [](const T0 &x, const T1 &y) { return x OP y; }; \
1632 std::transform(v0.begin(), v0.end(), v1.begin(), ret.begin(), op); \
1633 return ret; \
1634} \
1635
1644#undef RVEC_BINARY_OPERATOR
1645
1646///@}
1647///@name RVec Assignment Arithmetic Operators
1648///@{
1649
1650#define RVEC_ASSIGNMENT_OPERATOR(OP) \
1651template <typename T0, typename T1> \
1652RVec<T0>& operator OP(RVec<T0> &v, const T1 &y) \
1653{ \
1654 auto op = [&y](T0 &x) { return x OP y; }; \
1655 std::transform(v.begin(), v.end(), v.begin(), op); \
1656 return v; \
1657} \
1658 \
1659template <typename T0, typename T1> \
1660RVec<T0>& operator OP(RVec<T0> &v0, const RVec<T1> &v1) \
1661{ \
1662 if (v0.size() != v1.size()) \
1663 throw std::runtime_error(ERROR_MESSAGE(OP)); \
1664 \
1665 auto op = [](T0 &x, const T1 &y) { return x OP y; }; \
1666 std::transform(v0.begin(), v0.end(), v1.begin(), v0.begin(), op); \
1667 return v0; \
1668} \
1669
1680#undef RVEC_ASSIGNMENT_OPERATOR
1681
1682///@}
1683///@name RVec Comparison and Logical Operators
1684///@{
1685
1686#define RVEC_LOGICAL_OPERATOR(OP) \
1687template <typename T0, typename T1> \
1688auto operator OP(const RVec<T0> &v, const T1 &y) \
1689 -> RVec<int> /* avoid std::vector<bool> */ \
1690{ \
1691 RVec<int> ret(v.size()); \
1692 auto op = [y](const T0 &x) -> int { return x OP y; }; \
1693 std::transform(v.begin(), v.end(), ret.begin(), op); \
1694 return ret; \
1695} \
1696 \
1697template <typename T0, typename T1> \
1698auto operator OP(const T0 &x, const RVec<T1> &v) \
1699 -> RVec<int> /* avoid std::vector<bool> */ \
1700{ \
1701 RVec<int> ret(v.size()); \
1702 auto op = [x](const T1 &y) -> int { return x OP y; }; \
1703 std::transform(v.begin(), v.end(), ret.begin(), op); \
1704 return ret; \
1705} \
1706 \
1707template <typename T0, typename T1> \
1708auto operator OP(const RVec<T0> &v0, const RVec<T1> &v1) \
1709 -> RVec<int> /* avoid std::vector<bool> */ \
1710{ \
1711 if (v0.size() != v1.size()) \
1712 throw std::runtime_error(ERROR_MESSAGE(OP)); \
1713 \
1714 RVec<int> ret(v0.size()); \
1715 auto op = [](const T0 &x, const T1 &y) -> int { return x OP y; }; \
1716 std::transform(v0.begin(), v0.end(), v1.begin(), ret.begin(), op); \
1717 return ret; \
1718} \
1719
1728#undef RVEC_LOGICAL_OPERATOR
1729
1730///@}
1731///@name RVec Standard Mathematical Functions
1732///@{
1733
1734/// \cond
1735template <typename T> struct PromoteTypeImpl;
1736
1737template <> struct PromoteTypeImpl<float> { using Type = float; };
1738template <> struct PromoteTypeImpl<double> { using Type = double; };
1739template <> struct PromoteTypeImpl<long double> { using Type = long double; };
1740
1741template <typename T> struct PromoteTypeImpl { using Type = double; };
1742
1743template <typename T>
1744using PromoteType = typename PromoteTypeImpl<T>::Type;
1745
1746template <typename U, typename V>
1747using PromoteTypes = decltype(PromoteType<U>() + PromoteType<V>());
1748
1749/// \endcond
1750
1751#define RVEC_UNARY_FUNCTION(NAME, FUNC) \
1752 template <typename T> \
1753 RVec<PromoteType<T>> NAME(const RVec<T> &v) \
1754 { \
1755 RVec<PromoteType<T>> ret(v.size()); \
1756 auto f = [](const T &x) { return FUNC(x); }; \
1757 std::transform(v.begin(), v.end(), ret.begin(), f); \
1758 return ret; \
1759 }
1760
1761#define RVEC_BINARY_FUNCTION(NAME, FUNC) \
1762 template <typename T0, typename T1> \
1763 RVec<PromoteTypes<T0, T1>> NAME(const T0 &x, const RVec<T1> &v) \
1764 { \
1765 RVec<PromoteTypes<T0, T1>> ret(v.size()); \
1766 auto f = [&x](const T1 &y) { return FUNC(x, y); }; \
1767 std::transform(v.begin(), v.end(), ret.begin(), f); \
1768 return ret; \
1769 } \
1770 \
1771 template <typename T0, typename T1> \
1772 RVec<PromoteTypes<T0, T1>> NAME(const RVec<T0> &v, const T1 &y) \
1773 { \
1774 RVec<PromoteTypes<T0, T1>> ret(v.size()); \
1775 auto f = [&y](const T1 &x) { return FUNC(x, y); }; \
1776 std::transform(v.begin(), v.end(), ret.begin(), f); \
1777 return ret; \
1778 } \
1779 \
1780 template <typename T0, typename T1> \
1781 RVec<PromoteTypes<T0, T1>> NAME(const RVec<T0> &v0, const RVec<T1> &v1) \
1782 { \
1783 if (v0.size() != v1.size()) \
1784 throw std::runtime_error(ERROR_MESSAGE(NAME)); \
1785 \
1786 RVec<PromoteTypes<T0, T1>> ret(v0.size()); \
1787 auto f = [](const T0 &x, const T1 &y) { return FUNC(x, y); }; \
1788 std::transform(v0.begin(), v0.end(), v1.begin(), ret.begin(), f); \
1789 return ret; \
1790 } \
1791
1792#define RVEC_STD_UNARY_FUNCTION(F) RVEC_UNARY_FUNCTION(F, std::F)
1793#define RVEC_STD_BINARY_FUNCTION(F) RVEC_BINARY_FUNCTION(F, std::F)
1794
1799
1803
1808
1813
1821
1828
1835
1840#undef RVEC_STD_UNARY_FUNCTION
1841
1842///@}
1843///@name RVec Fast Mathematical Functions with Vdt
1844///@{
1845
1846#ifdef R__HAS_VDT
1847#define RVEC_VDT_UNARY_FUNCTION(F) RVEC_UNARY_FUNCTION(F, vdt::F)
1848
1849RVEC_VDT_UNARY_FUNCTION(fast_expf)
1850RVEC_VDT_UNARY_FUNCTION(fast_logf)
1851RVEC_VDT_UNARY_FUNCTION(fast_sinf)
1852RVEC_VDT_UNARY_FUNCTION(fast_cosf)
1853RVEC_VDT_UNARY_FUNCTION(fast_tanf)
1854RVEC_VDT_UNARY_FUNCTION(fast_asinf)
1855RVEC_VDT_UNARY_FUNCTION(fast_acosf)
1856RVEC_VDT_UNARY_FUNCTION(fast_atanf)
1857
1858RVEC_VDT_UNARY_FUNCTION(fast_exp)
1859RVEC_VDT_UNARY_FUNCTION(fast_log)
1860RVEC_VDT_UNARY_FUNCTION(fast_sin)
1861RVEC_VDT_UNARY_FUNCTION(fast_cos)
1862RVEC_VDT_UNARY_FUNCTION(fast_tan)
1863RVEC_VDT_UNARY_FUNCTION(fast_asin)
1864RVEC_VDT_UNARY_FUNCTION(fast_acos)
1865RVEC_VDT_UNARY_FUNCTION(fast_atan)
1866#undef RVEC_VDT_UNARY_FUNCTION
1867
1868#endif // R__HAS_VDT
1869
1870#undef RVEC_UNARY_FUNCTION
1871
1872///@}
1873
1874/// Inner product
1875///
1876/// Example code, at the ROOT prompt:
1877/// ~~~{.cpp}
1878/// using namespace ROOT::VecOps;
1879/// RVec<float> v1 {1., 2., 3.};
1880/// RVec<float> v2 {4., 5., 6.};
1881/// auto v1_dot_v2 = Dot(v1, v2);
1882/// v1_dot_v2
1883/// // (float) 32.0000f
1884/// ~~~
1885template <typename T, typename V>
1886auto Dot(const RVec<T> &v0, const RVec<V> &v1) -> decltype(v0[0] * v1[0])
1887{
1888 if (v0.size() != v1.size())
1889 throw std::runtime_error("Cannot compute inner product of vectors of different sizes");
1890 return std::inner_product(v0.begin(), v0.end(), v1.begin(), decltype(v0[0] * v1[0])(0));
1891}
1892
1893/// Sum elements of an RVec
1894///
1895/// Example code, at the ROOT prompt:
1896/// ~~~{.cpp}
1897/// using namespace ROOT::VecOps;
1898/// RVecF v {1.f, 2.f, 3.f};
1899/// auto v_sum = Sum(v);
1900/// v_sum
1901/// // (float) 6.f
1902/// auto v_sum_d = Sum(v, 0.);
1903/// v_sum_d
1904/// // (double) 6.0000000
1905/// ~~~
1906/// ~~~{.cpp}
1907/// using namespace ROOT::VecOps;
1908/// const ROOT::Math::PtEtaPhiMVector lv0 {15.5f, .3f, .1f, 105.65f},
1909/// lv1 {34.32f, 2.2f, 3.02f, 105.65f},
1910/// lv2 {12.95f, 1.32f, 2.2f, 105.65f};
1911/// RVec<ROOT::Math::PtEtaPhiMVector> v {lv0, lv1, lv2};
1912/// auto v_sum_lv = Sum(v, ROOT::Math::PtEtaPhiMVector());
1913/// v_sum_lv
1914/// // (ROOT::Math::LorentzVector<ROOT::Math::PtEtaPhiM4D<double> > &) (30.8489,2.46534,2.58947,361.084)
1915/// ~~~
1916template <typename T>
1917T Sum(const RVec<T> &v, const T zero = T(0))
1918{
1919 return std::accumulate(v.begin(), v.end(), zero);
1920}
1921
1922inline std::size_t Sum(const RVec<bool> &v, std::size_t zero = 0ul)
1923{
1924 return std::accumulate(v.begin(), v.end(), zero);
1925}
1926
1927/// Return the product of the elements of the RVec.
1928template <typename T>
1929T Product(const RVec<T> &v, const T init = T(1)) // initialize with identity
1930{
1931 return std::accumulate(v.begin(), v.end(), init, std::multiplies<T>());
1932}
1933
1934/// Get the mean of the elements of an RVec
1935///
1936/// The return type is a double precision floating point number.
1937///
1938/// Example code, at the ROOT prompt:
1939/// ~~~{.cpp}
1940/// using namespace ROOT::VecOps;
1941/// RVecF v {1.f, 2.f, 4.f};
1942/// auto v_mean = Mean(v);
1943/// v_mean
1944/// // (double) 2.3333333
1945/// ~~~
1946template <typename T>
1947double Mean(const RVec<T> &v)
1948{
1949 if (v.empty()) return 0.;
1950 return double(Sum(v)) / v.size();
1951}
1952
1953/// Get the mean of the elements of an RVec with custom initial value
1954///
1955/// The return type will be deduced from the `zero` parameter
1956///
1957/// Example code, at the ROOT prompt:
1958/// ~~~{.cpp}
1959/// using namespace ROOT::VecOps;
1960/// RVecF v {1.f, 2.f, 4.f};
1961/// auto v_mean_f = Mean(v, 0.f);
1962/// v_mean_f
1963/// // (float) 2.33333f
1964/// auto v_mean_d = Mean(v, 0.);
1965/// v_mean_d
1966/// // (double) 2.3333333
1967/// ~~~
1968/// ~~~{.cpp}
1969/// using namespace ROOT::VecOps;
1970/// const ROOT::Math::PtEtaPhiMVector lv0 {15.5f, .3f, .1f, 105.65f},
1971/// lv1 {34.32f, 2.2f, 3.02f, 105.65f},
1972/// lv2 {12.95f, 1.32f, 2.2f, 105.65f};
1973/// RVec<ROOT::Math::PtEtaPhiMVector> v {lv0, lv1, lv2};
1974/// auto v_mean_lv = Mean(v, ROOT::Math::PtEtaPhiMVector());
1975/// v_mean_lv
1976/// // (ROOT::Math::LorentzVector<ROOT::Math::PtEtaPhiM4D<double> > &) (10.283,2.46534,2.58947,120.361)
1977/// ~~~
1978template <typename T, typename R = T>
1979R Mean(const RVec<T> &v, const R zero)
1980{
1981 if (v.empty()) return zero;
1982 return Sum(v, zero) / v.size();
1983}
1984
1985/// Get the greatest element of an RVec
1986///
1987/// Example code, at the ROOT prompt:
1988/// ~~~{.cpp}
1989/// using namespace ROOT::VecOps;
1990/// RVecF v {1.f, 2.f, 4.f};
1991/// auto v_max = Max(v);
1992/// v_max
1993/// (float) 4.00000f
1994/// ~~~
1995template <typename T>
1996T Max(const RVec<T> &v)
1997{
1998 return *std::max_element(v.begin(), v.end());
1999}
2000
2001/// Get the smallest element of an RVec
2002///
2003/// Example code, at the ROOT prompt:
2004/// ~~~{.cpp}
2005/// using namespace ROOT::VecOps;
2006/// RVecF v {1.f, 2.f, 4.f};
2007/// auto v_min = Min(v);
2008/// v_min
2009/// (float) 1.00000f
2010/// ~~~
2011template <typename T>
2012T Min(const RVec<T> &v)
2013{
2014 return *std::min_element(v.begin(), v.end());
2015}
2016
2017/// Get the index of the greatest element of an RVec
2018/// In case of multiple occurrences of the maximum values,
2019/// the index corresponding to the first occurrence is returned.
2020///
2021/// Example code, at the ROOT prompt:
2022/// ~~~{.cpp}
2023/// using namespace ROOT::VecOps;
2024/// RVecF v {1.f, 2.f, 4.f};
2025/// auto v_argmax = ArgMax(v);
2026/// v_argmax
2027/// // (unsigned long) 2
2028/// ~~~
2029template <typename T>
2030std::size_t ArgMax(const RVec<T> &v)
2031{
2032 return std::distance(v.begin(), std::max_element(v.begin(), v.end()));
2033}
2034
2035/// Get the index of the smallest element of an RVec
2036/// In case of multiple occurrences of the minimum values,
2037/// the index corresponding to the first occurrence is returned.
2038///
2039/// Example code, at the ROOT prompt:
2040/// ~~~{.cpp}
2041/// using namespace ROOT::VecOps;
2042/// RVecF v {1.f, 2.f, 4.f};
2043/// auto v_argmin = ArgMin(v);
2044/// v_argmin
2045/// // (unsigned long) 0
2046/// ~~~
2047template <typename T>
2048std::size_t ArgMin(const RVec<T> &v)
2049{
2050 return std::distance(v.begin(), std::min_element(v.begin(), v.end()));
2051}
2052
2053/// Get the variance of the elements of an RVec
2054///
2055/// The return type is a double precision floating point number.
2056/// Example code, at the ROOT prompt:
2057/// ~~~{.cpp}
2058/// using namespace ROOT::VecOps;
2059/// RVecF v {1.f, 2.f, 4.f};
2060/// auto v_var = Var(v);
2061/// v_var
2062/// // (double) 2.3333333
2063/// ~~~
2064template <typename T>
2065double Var(const RVec<T> &v)
2066{
2067 const std::size_t size = v.size();
2068 if (size < std::size_t(2)) return 0.;
2069 T sum_squares(0), squared_sum(0);
2070 auto pred = [&sum_squares, &squared_sum](const T& x) {sum_squares+=x*x; squared_sum+=x;};
2071 std::for_each(v.begin(), v.end(), pred);
2072 squared_sum *= squared_sum;
2073 const auto dsize = (double) size;
2074 return 1. / (dsize - 1.) * (sum_squares - squared_sum / dsize );
2075}
2076
2077/// Get the standard deviation of the elements of an RVec
2078///
2079/// The return type is a double precision floating point number.
2080/// Example code, at the ROOT prompt:
2081/// ~~~{.cpp}
2082/// using namespace ROOT::VecOps;
2083/// RVecF v {1.f, 2.f, 4.f};
2084/// auto v_sd = StdDev(v);
2085/// v_sd
2086/// // (double) 1.5275252
2087/// ~~~
2088template <typename T>
2089double StdDev(const RVec<T> &v)
2090{
2091 return std::sqrt(Var(v));
2092}
2093
2094/// Create new collection applying a callable to the elements of the input collection
2095///
2096/// Example code, at the ROOT prompt:
2097/// ~~~{.cpp}
2098/// using namespace ROOT::VecOps;
2099/// RVecF v {1.f, 2.f, 4.f};
2100/// auto v_square = Map(v, [](float f){return f* 2.f;});
2101/// v_square
2102/// // (ROOT::VecOps::RVec<float> &) { 2.00000f, 4.00000f, 8.00000f }
2103///
2104/// RVecF x({1.f, 2.f, 3.f});
2105/// RVecF y({4.f, 5.f, 6.f});
2106/// RVecF z({7.f, 8.f, 9.f});
2107/// auto mod = [](float x, float y, float z) { return sqrt(x * x + y * y + z * z); };
2108/// auto v_mod = Map(x, y, z, mod);
2109/// v_mod
2110/// // (ROOT::VecOps::RVec<float> &) { 8.12404f, 9.64365f, 11.2250f }
2111/// ~~~
2112template <typename... Args>
2113auto Map(Args &&... args)
2114{
2115 /*
2116 Here the strategy in order to generalise the previous implementation of Map, i.e.
2117 `RVec Map(RVec, F)`, here we need to move the last parameter of the pack in first
2118 position in order to be able to invoke the Map function with automatic type deduction.
2119 This is achieved in two steps:
2120 1. Forward as tuple the pack to MapFromTuple
2121 2. Invoke the MapImpl helper which has the signature `template<...T, F> RVec MapImpl(F &&f, RVec<T>...)`
2122 */
2123
2124 // check the first N - 1 arguments are RVecs
2125 constexpr auto nArgs = sizeof...(Args);
2127 static_assert(ROOT::Internal::VecOps::All(isRVec, nArgs - 1),
2128 "Map: the first N-1 arguments must be RVecs or references to RVecs");
2129
2130 return ROOT::Internal::VecOps::MapFromTuple(std::forward_as_tuple(args...),
2131 std::make_index_sequence<sizeof...(args) - 1>());
2132}
2133
2134/// Create a new collection with the elements passing the filter expressed by the predicate
2135///
2136/// Example code, at the ROOT prompt:
2137/// ~~~{.cpp}
2138/// using namespace ROOT::VecOps;
2139/// RVecI v {1, 2, 4};
2140/// auto v_even = Filter(v, [](int i){return 0 == i%2;});
2141/// v_even
2142/// // (ROOT::VecOps::RVec<int> &) { 2, 4 }
2143/// ~~~
2144template <typename T, typename F>
2146{
2147 const auto thisSize = v.size();
2148 RVec<T> w;
2149 w.reserve(thisSize);
2150 for (auto &&val : v) {
2151 if (f(val))
2152 w.emplace_back(val);
2153 }
2154 return w;
2155}
2156
2157/// Return true if any of the elements equates to true, return false otherwise.
2158///
2159/// Example code, at the ROOT prompt:
2160/// ~~~{.cpp}
2161/// using namespace ROOT::VecOps;
2162/// RVecI v {0, 1, 0};
2163/// auto anyTrue = Any(v);
2164/// anyTrue
2165/// // (bool) true
2166/// ~~~
2167template <typename T>
2168auto Any(const RVec<T> &v) -> decltype(v[0] == true)
2169{
2170 for (auto &&e : v)
2171 if (static_cast<bool>(e) == true)
2172 return true;
2173 return false;
2174}
2175
2176/// Return true if all of the elements equate to true, return false otherwise.
2177///
2178/// Example code, at the ROOT prompt:
2179/// ~~~{.cpp}
2180/// using namespace ROOT::VecOps;
2181/// RVecI v {0, 1, 0};
2182/// auto allTrue = All(v);
2183/// allTrue
2184/// // (bool) false
2185/// ~~~
2186template <typename T>
2187auto All(const RVec<T> &v) -> decltype(v[0] == false)
2188{
2189 for (auto &&e : v)
2190 if (static_cast<bool>(e) == false)
2191 return false;
2192 return true;
2193}
2194
2195template <typename T>
2196void swap(RVec<T> &lhs, RVec<T> &rhs)
2197{
2198 lhs.swap(rhs);
2199}
2200
2201/// Return an RVec of indices that sort the input RVec
2202///
2203/// Example code, at the ROOT prompt:
2204/// ~~~{.cpp}
2205/// using namespace ROOT::VecOps;
2206/// RVecD v {2., 3., 1.};
2207/// auto sortIndices = Argsort(v)
2208/// // (ROOT::VecOps::RVec<unsigned long> &) { 2, 0, 1 }
2209/// auto values = Take(v, sortIndices)
2210/// // (ROOT::VecOps::RVec<double> &) { 1.0000000, 2.0000000, 3.0000000 }
2211/// ~~~
2212template <typename T>
2214{
2215 using size_type = typename RVec<T>::size_type;
2216 RVec<size_type> i(v.size());
2217 std::iota(i.begin(), i.end(), 0);
2218 std::sort(i.begin(), i.end(), [&v](size_type i1, size_type i2) { return v[i1] < v[i2]; });
2219 return i;
2220}
2221
2222/// Return an RVec of indices that sort the input RVec based on a comparison function.
2223///
2224/// Example code, at the ROOT prompt:
2225/// ~~~{.cpp}
2226/// using namespace ROOT::VecOps;
2227/// RVecD v {2., 3., 1.};
2228/// auto sortIndices = Argsort(v, [](double x, double y) {return x > y;})
2229/// // (ROOT::VecOps::RVec<unsigned long> &) { 1, 0, 2 }
2230/// auto values = Take(v, sortIndices)
2231/// // (ROOT::VecOps::RVec<double> &) { 3.0000000, 2.0000000, 1.0000000 }
2232/// ~~~
2233template <typename T, typename Compare>
2235{
2236 using size_type = typename RVec<T>::size_type;
2237 RVec<size_type> i(v.size());
2238 std::iota(i.begin(), i.end(), 0);
2239 std::sort(i.begin(), i.end(),
2240 [&v, &c](size_type i1, size_type i2) { return c(v[i1], v[i2]); });
2241 return i;
2242}
2243
2244/// Return an RVec of indices that sort the input RVec
2245/// while keeping the order of equal elements.
2246/// This is the stable variant of `Argsort`.
2247///
2248/// Example code, at the ROOT prompt:
2249/// ~~~{.cpp}
2250/// using namespace ROOT::VecOps;
2251/// RVecD v {2., 3., 2., 1.};
2252/// auto sortIndices = StableArgsort(v)
2253/// // (ROOT::VecOps::RVec<unsigned long> &) { 3, 0, 2, 1 }
2254/// auto values = Take(v, sortIndices)
2255/// // (ROOT::VecOps::RVec<double> &) { 1.0000000, 2.0000000, 2.0000000, 3.0000000 }
2256/// ~~~
2257template <typename T>
2259{
2260 using size_type = typename RVec<T>::size_type;
2261 RVec<size_type> i(v.size());
2262 std::iota(i.begin(), i.end(), 0);
2263 std::stable_sort(i.begin(), i.end(), [&v](size_type i1, size_type i2) { return v[i1] < v[i2]; });
2264 return i;
2265}
2266
2267/// Return an RVec of indices that sort the input RVec based on a comparison function
2268/// while keeping the order of equal elements.
2269/// This is the stable variant of `Argsort`.
2270///
2271/// Example code, at the ROOT prompt:
2272/// ~~~{.cpp}
2273/// using namespace ROOT::VecOps;
2274/// RVecD v {2., 3., 2., 1.};
2275/// auto sortIndices = StableArgsort(v, [](double x, double y) {return x > y;})
2276/// // (ROOT::VecOps::RVec<unsigned long> &) { 1, 0, 2, 3 }
2277/// auto values = Take(v, sortIndices)
2278/// // (ROOT::VecOps::RVec<double> &) { 3.0000000, 2.0000000, 2.0000000, 1.0000000 }
2279/// ~~~
2280template <typename T, typename Compare>
2282{
2283 using size_type = typename RVec<T>::size_type;
2284 RVec<size_type> i(v.size());
2285 std::iota(i.begin(), i.end(), 0);
2286 std::stable_sort(i.begin(), i.end(), [&v, &c](size_type i1, size_type i2) { return c(v[i1], v[i2]); });
2287 return i;
2288}
2289
2290/// Return elements of a vector at given indices
2291///
2292/// Example code, at the ROOT prompt:
2293/// ~~~{.cpp}
2294/// using namespace ROOT::VecOps;
2295/// RVecD v {2., 3., 1.};
2296/// auto vTaken = Take(v, {0,2});
2297/// vTaken
2298/// // (ROOT::VecOps::RVec<double>) { 2.0000000, 1.0000000 }
2299/// ~~~
2300
2301template <typename T>
2302RVec<T> Take(const RVec<T> &v, const RVec<typename RVec<T>::size_type> &i)
2303{
2304 using size_type = typename RVec<T>::size_type;
2305 const size_type isize = i.size();
2306 RVec<T> r(isize);
2307 for (size_type k = 0; k < isize; k++)
2308 r[k] = v[i[k]];
2309 return r;
2310}
2311
2312/// Take version that defaults to (user-specified) output value if some index is out of range
2313template <typename T>
2314RVec<T> Take(const RVec<T> &v, const RVec<typename RVec<T>::size_type> &i, const T default_val)
2315{
2316 using size_type = typename RVec<T>::size_type;
2317 const size_type isize = i.size();
2318 RVec<T> r(isize);
2319 for (size_type k = 0; k < isize; k++)
2320 {
2321 if (k < v.size()){
2322 r[k] = v[i[k]];
2323 }
2324 else {
2325 r[k] = default_val;
2326 }
2327 }
2328 return r;
2329}
2330
2331/// Return first `n` elements of an RVec if `n > 0` and last `n` elements if `n < 0`.
2332///
2333/// Example code, at the ROOT prompt:
2334/// ~~~{.cpp}
2335/// using namespace ROOT::VecOps;
2336/// RVecD v {2., 3., 1.};
2337/// auto firstTwo = Take(v, 2);
2338/// firstTwo
2339/// // (ROOT::VecOps::RVec<double>) { 2.0000000, 3.0000000 }
2340/// auto lastOne = Take(v, -1);
2341/// lastOne
2342/// // (ROOT::VecOps::RVec<double>) { 1.0000000 }
2343/// ~~~
2344template <typename T>
2345RVec<T> Take(const RVec<T> &v, const int n)
2346{
2347 using size_type = typename RVec<T>::size_type;
2348 const size_type size = v.size();
2349 const size_type absn = std::abs(n);
2350 if (absn > size) {
2351 const auto msg = std::to_string(absn) + " elements requested from Take but input contains only " +
2352 std::to_string(size) + " elements.";
2353 throw std::runtime_error(msg);
2354 }
2355 RVec<T> r(absn);
2356 if (n < 0) {
2357 for (size_type k = 0; k < absn; k++)
2358 r[k] = v[size - absn + k];
2359 } else {
2360 for (size_type k = 0; k < absn; k++)
2361 r[k] = v[k];
2362 }
2363 return r;
2364}
2365
2366/// Return first `n` elements of an RVec if `n > 0` and last `n` elements if `n < 0`.
2367///
2368/// This Take version defaults to a user-specified value
2369/// `default_val` if the absolute value of `n` is
2370/// greater than the size of the RVec `v`
2371///
2372/// Example code, at the ROOT prompt:
2373/// ~~~{.cpp}
2374/// using ROOT::VecOps::RVec;
2375/// RVec<int> x{1,2,3,4};
2376/// Take(x,-5,1)
2377/// // (ROOT::VecOps::RVec<int>) { 1, 1, 2, 3, 4 }
2378/// Take(x,5,20)
2379/// // (ROOT::VecOps::RVec<int>) { 1, 2, 3, 4, 20 }
2380/// Take(x,-1,1)
2381/// // (ROOT::VecOps::RVec<int>) { 4 }
2382/// Take(x,4,1)
2383/// // (ROOT::VecOps::RVec<int>) { 1, 2, 3, 4 }
2384/// ~~~
2385template <typename T>
2386RVec<T> Take(const RVec<T> &v, const int n, const T default_val)
2387{
2388 using size_type = typename RVec<T>::size_type;
2389 const size_type size = v.size();
2390 const size_type absn = std::abs(n);
2391 // Base case, can be handled by another overload of Take
2392 if (absn <= size) {
2393 return Take(v, n);
2394 }
2395 RVec<T> temp = v;
2396 // Case when n is positive and n > v.size()
2397 if (n > 0) {
2398 temp.resize(n, default_val);
2399 return temp;
2400 }
2401 // Case when n is negative and abs(n) > v.size()
2402 const auto num_to_fill = absn - size;
2403 ROOT::VecOps::RVec<T> fill_front(num_to_fill, default_val);
2404 return Concatenate(fill_front, temp);
2405}
2406
2407/// Return a copy of the container without the elements at the specified indices.
2408///
2409/// Duplicated and out-of-range indices in idxs are ignored.
2410template <typename T>
2412{
2413 // clean up input indices
2414 std::sort(idxs.begin(), idxs.end());
2415 idxs.erase(std::unique(idxs.begin(), idxs.end()), idxs.end());
2416
2417 RVec<T> r;
2418 if (v.size() > idxs.size())
2419 r.reserve(v.size() - idxs.size());
2420
2421 auto discardIt = idxs.begin();
2422 using sz_t = typename RVec<T>::size_type;
2423 for (sz_t i = 0u; i < v.size(); ++i) {
2424 if (discardIt != idxs.end() && i == *discardIt)
2425 ++discardIt;
2426 else
2427 r.emplace_back(v[i]);
2428 }
2429
2430 return r;
2431}
2432
2433/// Return copy of reversed vector
2434///
2435/// Example code, at the ROOT prompt:
2436/// ~~~{.cpp}
2437/// using namespace ROOT::VecOps;
2438/// RVecD v {2., 3., 1.};
2439/// auto v_reverse = Reverse(v);
2440/// v_reverse
2441/// // (ROOT::VecOps::RVec<double> &) { 1.0000000, 3.0000000, 2.0000000 }
2442/// ~~~
2443template <typename T>
2445{
2446 RVec<T> r(v);
2447 std::reverse(r.begin(), r.end());
2448 return r;
2449}
2450
2451/// Return copy of RVec with elements sorted in ascending order
2452///
2453/// This helper is different from Argsort since it does not return an RVec of indices,
2454/// but an RVec of values.
2455///
2456/// Example code, at the ROOT prompt:
2457/// ~~~{.cpp}
2458/// using namespace ROOT::VecOps;
2459/// RVecD v {2., 3., 1.};
2460/// auto v_sorted = Sort(v);
2461/// v_sorted
2462/// // (ROOT::VecOps::RVec<double> &) { 1.0000000, 2.0000000, 3.0000000 }
2463/// ~~~
2464template <typename T>
2465RVec<T> Sort(const RVec<T> &v)
2466{
2467 RVec<T> r(v);
2468 std::sort(r.begin(), r.end());
2469 return r;
2470}
2471
2472/// Return copy of RVec with elements sorted based on a comparison operator
2473///
2474/// The comparison operator has to fulfill the same requirements of the
2475/// predicate of by std::sort.
2476///
2477///
2478/// This helper is different from Argsort since it does not return an RVec of indices,
2479/// but an RVec of values.
2480///
2481/// Example code, at the ROOT prompt:
2482/// ~~~{.cpp}
2483/// using namespace ROOT::VecOps;
2484/// RVecD v {2., 3., 1.};
2485/// auto v_sorted = Sort(v, [](double x, double y) {return 1/x < 1/y;});
2486/// v_sorted
2487/// // (ROOT::VecOps::RVec<double> &) { 3.0000000, 2.0000000, 1.0000000 }
2488/// ~~~
2489template <typename T, typename Compare>
2490RVec<T> Sort(const RVec<T> &v, Compare &&c)
2491{
2492 RVec<T> r(v);
2493 std::sort(r.begin(), r.end(), std::forward<Compare>(c));
2494 return r;
2495}
2496
2497/// Return copy of RVec with elements sorted in ascending order
2498/// while keeping the order of equal elements.
2499///
2500/// This is the stable variant of `Sort`.
2501///
2502/// This helper is different from StableArgsort since it does not return an RVec of indices,
2503/// but an RVec of values.
2504///
2505/// Example code, at the ROOT prompt:
2506/// ~~~{.cpp}
2507/// using namespace ROOT::VecOps;
2508/// RVecD v {2., 3., 2, 1.};
2509/// auto v_sorted = StableSort(v);
2510/// v_sorted
2511/// // (ROOT::VecOps::RVec<double> &) { 1.0000000, 2.0000000, 2.0000000, 3.0000000 }
2512/// ~~~
2513template <typename T>
2515{
2516 RVec<T> r(v);
2517 std::stable_sort(r.begin(), r.end());
2518 return r;
2519}
2520
2521// clang-format off
2522/// Return copy of RVec with elements sorted based on a comparison operator
2523/// while keeping the order of equal elements.
2524///
2525/// The comparison operator has to fulfill the same requirements of the
2526/// predicate of std::stable_sort.
2527///
2528/// This helper is different from StableArgsort since it does not return an RVec of indices,
2529/// but an RVec of values.
2530///
2531/// This is the stable variant of `Sort`.
2532///
2533/// Example code, at the ROOT prompt:
2534/// ~~~{.cpp}
2535/// using namespace ROOT::VecOps;
2536/// RVecD v {2., 3., 2., 1.};
2537/// auto v_sorted = StableSort(v, [](double x, double y) {return 1/x < 1/y;});
2538/// v_sorted
2539/// // (ROOT::VecOps::RVec<double> &) { 3.0000000, 2.0000000, 2.0000000, 1.0000000 }
2540/// ~~~
2541/// ~~~{.cpp}
2542/// using namespace ROOT::VecOps;
2543/// RVec<RVecD> v {{2., 4.}, {3., 1.}, {2, 1.}, {1., 4.}};
2544/// auto v_sorted = StableSort(StableSort(v, [](const RVecD &x, const RVecD &y) {return x[1] < y[1];}), [](const RVecD &x, const RVecD &y) {return x[0] < y[0];});
2545/// v_sorted
2546/// // (ROOT::VecOps::RVec<ROOT::VecOps::RVec<double> > &) { { 1.0000000, 4.0000000 }, { 2.0000000, 1.0000000 }, { 2.0000000, 4.0000000 }, { 3.0000000, 1.0000000 } }
2547/// ~~~
2548// clang-format off
2549template <typename T, typename Compare>
2551{
2552 RVec<T> r(v);
2553 std::stable_sort(r.begin(), r.end(), std::forward<Compare>(c));
2554 return r;
2555}
2556
2557/// Return the indices that represent all combinations of the elements of two
2558/// RVecs.
2559///
2560/// The type of the return value is an RVec of two RVecs containing indices.
2561///
2562/// Example code, at the ROOT prompt:
2563/// ~~~{.cpp}
2564/// using namespace ROOT::VecOps;
2565/// auto comb_idx = Combinations(3, 2);
2566/// comb_idx
2567/// // (ROOT::VecOps::RVec<ROOT::VecOps::RVec<unsigned long> > &) { { 0, 0, 1, 1, 2, 2 }, { 0, 1, 0, 1, 0, 1 } }
2568/// ~~~
2569inline RVec<RVec<std::size_t>> Combinations(const std::size_t size1, const std::size_t size2)
2570{
2571 using size_type = std::size_t;
2573 r[0].resize(size1*size2);
2574 r[1].resize(size1*size2);
2575 size_type c = 0;
2576 for(size_type i=0; i<size1; i++) {
2577 for(size_type j=0; j<size2; j++) {
2578 r[0][c] = i;
2579 r[1][c] = j;
2580 c++;
2581 }
2582 }
2583 return r;
2584}
2585
2586/// Return the indices that represent all combinations of the elements of two
2587/// RVecs.
2588///
2589/// The type of the return value is an RVec of two RVecs containing indices.
2590///
2591/// Example code, at the ROOT prompt:
2592/// ~~~{.cpp}
2593/// using namespace ROOT::VecOps;
2594/// RVecD v1 {1., 2., 3.};
2595/// RVecD v2 {-4., -5.};
2596/// auto comb_idx = Combinations(v1, v2);
2597/// comb_idx
2598/// // (ROOT::VecOps::RVec<ROOT::VecOps::RVec<unsigned long> > &) { { 0, 0, 1, 1, 2, 2 }, { 0, 1, 0, 1, 0, 1 } }
2599/// ~~~
2600template <typename T1, typename T2>
2602{
2603 return Combinations(v1.size(), v2.size());
2604}
2605
2606/// Return the indices that represent all unique combinations of the
2607/// elements of a given RVec.
2608///
2609/// ~~~{.cpp}
2610/// using namespace ROOT::VecOps;
2611/// RVecD v {1., 2., 3., 4.};
2612/// auto v_1 = Combinations(v, 1);
2613/// v_1
2614/// (ROOT::VecOps::RVec<ROOT::VecOps::RVec<unsigned long> > &) { { 0, 1, 2, 3 } }
2615/// auto v_2 = Combinations(v, 2);
2616/// v_2
2617/// (ROOT::VecOps::RVec<ROOT::VecOps::RVec<unsigned long> > &) { { 0, 0, 0, 1, 1, 2 }, { 1, 2, 3, 2, 3, 3 } }
2618/// auto v_3 = Combinations(v, 3);
2619/// v_3
2620/// (ROOT::VecOps::RVec<ROOT::VecOps::RVec<unsigned long> > &) { { 0, 0, 0, 1 }, { 1, 1, 2, 2 }, { 2, 3, 3, 3 } }
2621/// auto v_4 = Combinations(v, 4);
2622/// v_4
2623/// (ROOT::VecOps::RVec<ROOT::VecOps::RVec<unsigned long> > &) { { 0 }, { 1 }, { 2 }, { 3 } }
2624/// ~~~
2625template <typename T>
2627{
2628 using size_type = typename RVec<T>::size_type;
2629 const size_type s = v.size();
2630 if (n > s) {
2631 throw std::runtime_error("Cannot make unique combinations of size " + std::to_string(n) +
2632 " from vector of size " + std::to_string(s) + ".");
2633 }
2634
2636 for(size_type k=0; k<s; k++)
2637 indices[k] = k;
2638
2639 const auto innersize = [=] {
2640 size_type inners = s - n + 1;
2641 for (size_type m = s - n + 2; m <= s; ++m)
2642 inners *= m;
2643
2644 size_type factn = 1;
2645 for (size_type i = 2; i <= n; ++i)
2646 factn *= i;
2647 inners /= factn;
2648
2649 return inners;
2650 }();
2651
2653 size_type inneridx = 0;
2654 for (size_type k = 0; k < n; k++)
2655 c[k][inneridx] = indices[k];
2656 ++inneridx;
2657
2658 while (true) {
2659 bool run_through = true;
2660 long i = n - 1;
2661 for (; i>=0; i--) {
2662 if (indices[i] != i + s - n){
2663 run_through = false;
2664 break;
2665 }
2666 }
2667 if (run_through) {
2668 return c;
2669 }
2670 indices[i]++;
2671 for (long j=i+1; j<(long)n; j++)
2672 indices[j] = indices[j-1] + 1;
2673 for (size_type k = 0; k < n; k++)
2674 c[k][inneridx] = indices[k];
2675 ++inneridx;
2676 }
2677}
2678
2679/// Return the indices of the elements which are not zero
2680///
2681/// Example code, at the ROOT prompt:
2682/// ~~~{.cpp}
2683/// using namespace ROOT::VecOps;
2684/// RVecD v {2., 0., 3., 0., 1.};
2685/// auto nonzero_idx = Nonzero(v);
2686/// nonzero_idx
2687/// // (ROOT::VecOps::RVec<unsigned long> &) { 0, 2, 4 }
2688/// ~~~
2689template <typename T>
2691{
2692 using size_type = typename RVec<T>::size_type;
2694 const auto size = v.size();
2695 r.reserve(size);
2696 for(size_type i=0; i<size; i++) {
2697 if(v[i] != 0) {
2698 r.emplace_back(i);
2699 }
2700 }
2701 return r;
2702}
2703
2704/// Return the intersection of elements of two RVecs.
2705///
2706/// Each element of v1 is looked up in v2 and added to the returned vector if
2707/// found. Following, the order of v1 is preserved. If v2 is already sorted, the
2708/// optional argument v2_is_sorted can be used to toggle of the internal sorting
2709/// step, therewith optimising runtime.
2710///
2711/// Example code, at the ROOT prompt:
2712/// ~~~{.cpp}
2713/// using namespace ROOT::VecOps;
2714/// RVecD v1 {1., 2., 3.};
2715/// RVecD v2 {-4., -5., 2., 1.};
2716/// auto v1_intersect_v2 = Intersect(v1, v2);
2717/// v1_intersect_v2
2718/// // (ROOT::VecOps::RVec<double> &) { 1.0000000, 2.0000000 }
2719/// ~~~
2720template <typename T>
2721RVec<T> Intersect(const RVec<T>& v1, const RVec<T>& v2, bool v2_is_sorted = false)
2722{
2723 RVec<T> v2_sorted;
2724 if (!v2_is_sorted) v2_sorted = Sort(v2);
2725 const auto v2_begin = v2_is_sorted ? v2.begin() : v2_sorted.begin();
2726 const auto v2_end = v2_is_sorted ? v2.end() : v2_sorted.end();
2727 RVec<T> r;
2728 const auto size = v1.size();
2729 r.reserve(size);
2730 using size_type = typename RVec<T>::size_type;
2731 for(size_type i=0; i<size; i++) {
2732 if (std::binary_search(v2_begin, v2_end, v1[i])) {
2733 r.emplace_back(v1[i]);
2734 }
2735 }
2736 return r;
2737}
2738
2739/// Return the elements of v1 if the condition c is true and v2 if the
2740/// condition c is false.
2741///
2742/// Example code, at the ROOT prompt:
2743/// ~~~{.cpp}
2744/// using namespace ROOT::VecOps;
2745/// RVecD v1 {1., 2., 3.};
2746/// RVecD v2 {-1., -2., -3.};
2747/// auto c = v1 > 1;
2748/// c
2749/// // (ROOT::VecOps::RVec<int> &) { 0, 1, 1 }
2750/// auto if_c_v1_else_v2 = Where(c, v1, v2);
2751/// if_c_v1_else_v2
2752/// // (ROOT::VecOps::RVec<double> &) { -1.0000000, 2.0000000, 3.0000000 }
2753/// ~~~
2754template <typename T>
2755RVec<T> Where(const RVec<int>& c, const RVec<T>& v1, const RVec<T>& v2)
2756{
2757 using size_type = typename RVec<T>::size_type;
2758 const size_type size = c.size();
2759 RVec<T> r;
2760 r.reserve(size);
2761 for (size_type i=0; i<size; i++) {
2762 r.emplace_back(c[i] != 0 ? v1[i] : v2[i]);
2763 }
2764 return r;
2765}
2766
2767/// Return the elements of v1 if the condition c is true and sets the value v2
2768/// if the condition c is false.
2769///
2770/// Example code, at the ROOT prompt:
2771/// ~~~{.cpp}
2772/// using namespace ROOT::VecOps;
2773/// RVecD v1 {1., 2., 3.};
2774/// double v2 = 4.;
2775/// auto c = v1 > 1;
2776/// c
2777/// // (ROOT::VecOps::RVec<int> &) { 0, 1, 1 }
2778/// auto if_c_v1_else_v2 = Where(c, v1, v2);
2779/// if_c_v1_else_v2
2780/// // (ROOT::VecOps::RVec<double>) { 4.0000000, 2.0000000, 3.0000000 }
2781/// ~~~
2782template <typename T>
2784{
2785 using size_type = typename RVec<T>::size_type;
2786 const size_type size = c.size();
2787 RVec<T> r;
2788 r.reserve(size);
2789 for (size_type i=0; i<size; i++) {
2790 r.emplace_back(c[i] != 0 ? v1[i] : v2);
2791 }
2792 return r;
2793}
2794
2795/// Return the elements of v2 if the condition c is false and sets the value v1
2796/// if the condition c is true.
2797///
2798/// Example code, at the ROOT prompt:
2799/// ~~~{.cpp}
2800/// using namespace ROOT::VecOps;
2801/// double v1 = 4.;
2802/// RVecD v2 {1., 2., 3.};
2803/// auto c = v2 > 1;
2804/// c
2805/// // (ROOT::VecOps::RVec<int> &) { 0, 1, 1 }
2806/// auto if_c_v1_else_v2 = Where(c, v1, v2);
2807/// if_c_v1_else_v2
2808/// // (ROOT::VecOps::RVec<double>) { 1.0000000, 4.0000000, 4.0000000 }
2809/// ~~~
2810template <typename T>
2812{
2813 using size_type = typename RVec<T>::size_type;
2814 const size_type size = c.size();
2815 RVec<T> r;
2816 r.reserve(size);
2817 for (size_type i=0; i<size; i++) {
2818 r.emplace_back(c[i] != 0 ? v1 : v2[i]);
2819 }
2820 return r;
2821}
2822
2823/// Return a vector with the value v2 if the condition c is false and sets the
2824/// value v1 if the condition c is true.
2825///
2826/// Example code, at the ROOT prompt:
2827/// ~~~{.cpp}
2828/// using namespace ROOT::VecOps;
2829/// double v1 = 4.;
2830/// double v2 = 2.;
2831/// RVecI c {0, 1, 1};
2832/// auto if_c_v1_else_v2 = Where(c, v1, v2);
2833/// if_c_v1_else_v2
2834/// // (ROOT::VecOps::RVec<double>) { 2.0000000, 4.0000000, 4.0000000 }
2835/// ~~~
2836template <typename T>
2838{
2839 using size_type = typename RVec<T>::size_type;
2840 const size_type size = c.size();
2841 RVec<T> r;
2842 r.reserve(size);
2843 for (size_type i=0; i<size; i++) {
2844 r.emplace_back(c[i] != 0 ? v1 : v2);
2845 }
2846 return r;
2847}
2848
2849/// Return the concatenation of two RVecs.
2850///
2851/// Example code, at the ROOT prompt:
2852/// ~~~{.cpp}
2853/// using namespace ROOT::VecOps;
2854/// RVecF rvf {0.f, 1.f, 2.f};
2855/// RVecI rvi {7, 8, 9};
2856/// Concatenate(rvf, rvi)
2857/// // (ROOT::VecOps::RVec<float>) { 0.00000f, 1.00000f, 2.00000f, 7.00000f, 8.00000f, 9.00000f }
2858/// ~~~
2859template <typename T0, typename T1, typename Common_t = typename std::common_type<T0, T1>::type>
2861{
2862 RVec<Common_t> res;
2863 res.reserve(v0.size() + v1.size());
2864 std::copy(v0.begin(), v0.end(), std::back_inserter(res));
2865 std::copy(v1.begin(), v1.end(), std::back_inserter(res));
2866 return res;
2867}
2868
2869/// Return the angle difference \f$\Delta \phi\f$ of two scalars.
2870///
2871/// The function computes the closest angle from v1 to v2 with sign and is
2872/// therefore in the range \f$[-\pi, \pi]\f$.
2873/// The computation is done per default in radians \f$c = \pi\f$ but can be switched
2874/// to degrees \f$c = 180\f$.
2875template <typename T>
2876T DeltaPhi(T v1, T v2, const T c = M_PI)
2877{
2878 static_assert(std::is_floating_point<T>::value,
2879 "DeltaPhi must be called with floating point values.");
2880 auto r = std::fmod(v2 - v1, 2.0 * c);
2881 if (r < -c) {
2882 r += 2.0 * c;
2883 }
2884 else if (r > c) {
2885 r -= 2.0 * c;
2886 }
2887 return r;
2888}
2889
2890/// Return the angle difference \f$\Delta \phi\f$ in radians of two vectors.
2891///
2892/// The function computes the closest angle from v1 to v2 with sign and is
2893/// therefore in the range \f$[-\pi, \pi]\f$.
2894/// The computation is done per default in radians \f$c = \pi\f$ but can be switched
2895/// to degrees \f$c = 180\f$.
2896template <typename T>
2897RVec<T> DeltaPhi(const RVec<T>& v1, const RVec<T>& v2, const T c = M_PI)
2898{
2899 using size_type = typename RVec<T>::size_type;
2900 const size_type size = v1.size();
2901 auto r = RVec<T>(size);
2902 for (size_type i = 0; i < size; i++) {
2903 r[i] = DeltaPhi(v1[i], v2[i], c);
2904 }
2905 return r;
2906}
2907
2908/// Return the angle difference \f$\Delta \phi\f$ in radians of a vector and a scalar.
2909///
2910/// The function computes the closest angle from v1 to v2 with sign and is
2911/// therefore in the range \f$[-\pi, \pi]\f$.
2912/// The computation is done per default in radians \f$c = \pi\f$ but can be switched
2913/// to degrees \f$c = 180\f$.
2914template <typename T>
2915RVec<T> DeltaPhi(const RVec<T>& v1, T v2, const T c = M_PI)
2916{
2917 using size_type = typename RVec<T>::size_type;
2918 const size_type size = v1.size();
2919 auto r = RVec<T>(size);
2920 for (size_type i = 0; i < size; i++) {
2921 r[i] = DeltaPhi(v1[i], v2, c);
2922 }
2923 return r;
2924}
2925
2926/// Return the angle difference \f$\Delta \phi\f$ in radians of a scalar and a vector.
2927///
2928/// The function computes the closest angle from v1 to v2 with sign and is
2929/// therefore in the range \f$[-\pi, \pi]\f$.
2930/// The computation is done per default in radians \f$c = \pi\f$ but can be switched
2931/// to degrees \f$c = 180\f$.
2932template <typename T>
2933RVec<T> DeltaPhi(T v1, const RVec<T>& v2, const T c = M_PI)
2934{
2935 using size_type = typename RVec<T>::size_type;
2936 const size_type size = v2.size();
2937 auto r = RVec<T>(size);
2938 for (size_type i = 0; i < size; i++) {
2939 r[i] = DeltaPhi(v1, v2[i], c);
2940 }
2941 return r;
2942}
2943
2944/// Return the square of the distance on the \f$\eta\f$-\f$\phi\f$ plane (\f$\Delta R\f$) from
2945/// the collections eta1, eta2, phi1 and phi2.
2946///
2947/// The function computes \f$\Delta R^2 = (\eta_1 - \eta_2)^2 + (\phi_1 - \phi_2)^2\f$
2948/// of the given collections eta1, eta2, phi1 and phi2. The angle \f$\phi\f$ can
2949/// be set to radian or degrees using the optional argument c, see the documentation
2950/// of the DeltaPhi helper.
2951template <typename T>
2952RVec<T> DeltaR2(const RVec<T>& eta1, const RVec<T>& eta2, const RVec<T>& phi1, const RVec<T>& phi2, const T c = M_PI)
2953{
2954 const auto dphi = DeltaPhi(phi1, phi2, c);
2955 return (eta1 - eta2) * (eta1 - eta2) + dphi * dphi;
2956}
2957
2958/// Return the distance on the \f$\eta\f$-\f$\phi\f$ plane (\f$\Delta R\f$) from
2959/// the collections eta1, eta2, phi1 and phi2.
2960///
2961/// The function computes \f$\Delta R = \sqrt{(\eta_1 - \eta_2)^2 + (\phi_1 - \phi_2)^2}\f$
2962/// of the given collections eta1, eta2, phi1 and phi2. The angle \f$\phi\f$ can
2963/// be set to radian or degrees using the optional argument c, see the documentation
2964/// of the DeltaPhi helper.
2965template <typename T>
2966RVec<T> DeltaR(const RVec<T>& eta1, const RVec<T>& eta2, const RVec<T>& phi1, const RVec<T>& phi2, const T c = M_PI)
2967{
2968 return sqrt(DeltaR2(eta1, eta2, phi1, phi2, c));
2969}
2970
2971/// Return the distance on the \f$\eta\f$-\f$\phi\f$ plane (\f$\Delta R\f$) from
2972/// the scalars eta1, eta2, phi1 and phi2.
2973///
2974/// The function computes \f$\Delta R = \sqrt{(\eta_1 - \eta_2)^2 + (\phi_1 - \phi_2)^2}\f$
2975/// of the given scalars eta1, eta2, phi1 and phi2. The angle \f$\phi\f$ can
2976/// be set to radian or degrees using the optional argument c, see the documentation
2977/// of the DeltaPhi helper.
2978template <typename T>
2979T DeltaR(T eta1, T eta2, T phi1, T phi2, const T c = M_PI)
2980{
2981 const auto dphi = DeltaPhi(phi1, phi2, c);
2982 return std::sqrt((eta1 - eta2) * (eta1 - eta2) + dphi * dphi);
2983}
2984
2985/// Return the invariant mass of two particles given the collections of the quantities
2986/// transverse momentum (pt), rapidity (eta), azimuth (phi) and mass.
2987///
2988/// The function computes the invariant mass of two particles with the four-vectors
2989/// (pt1, eta2, phi1, mass1) and (pt2, eta2, phi2, mass2).
2990template <typename T>
2992 const RVec<T>& pt1, const RVec<T>& eta1, const RVec<T>& phi1, const RVec<T>& mass1,
2993 const RVec<T>& pt2, const RVec<T>& eta2, const RVec<T>& phi2, const RVec<T>& mass2)
2994{
2995 std::size_t size = pt1.size();
2996
2997 R__ASSERT(eta1.size() == size && phi1.size() == size && mass1.size() == size);
2998 R__ASSERT(pt2.size() == size && phi2.size() == size && mass2.size() == size);
2999
3000 RVec<T> inv_masses(size);
3001
3002 for (std::size_t i = 0u; i < size; ++i) {
3003 // Conversion from (pt, eta, phi, mass) to (x, y, z, e) coordinate system
3004 const auto x1 = pt1[i] * std::cos(phi1[i]);
3005 const auto y1 = pt1[i] * std::sin(phi1[i]);
3006 const auto z1 = pt1[i] * std::sinh(eta1[i]);
3007 const auto e1 = std::sqrt(x1 * x1 + y1 * y1 + z1 * z1 + mass1[i] * mass1[i]);
3008
3009 const auto x2 = pt2[i] * std::cos(phi2[i]);
3010 const auto y2 = pt2[i] * std::sin(phi2[i]);
3011 const auto z2 = pt2[i] * std::sinh(eta2[i]);
3012 const auto e2 = std::sqrt(x2 * x2 + y2 * y2 + z2 * z2 + mass2[i] * mass2[i]);
3013
3014 // Addition of particle four-vector elements
3015 const auto e = e1 + e2;
3016 const auto x = x1 + x2;
3017 const auto y = y1 + y2;
3018 const auto z = z1 + z2;
3019
3020 inv_masses[i] = std::sqrt(e * e - x * x - y * y - z * z);
3021 }
3022
3023 // Return invariant mass with (+, -, -, -) metric
3024 return inv_masses;
3025}
3026
3027/// Return the invariant mass of multiple particles given the collections of the
3028/// quantities transverse momentum (pt), rapidity (eta), azimuth (phi) and mass.
3029///
3030/// The function computes the invariant mass of multiple particles with the
3031/// four-vectors (pt, eta, phi, mass).
3032template <typename T>
3033T InvariantMass(const RVec<T>& pt, const RVec<T>& eta, const RVec<T>& phi, const RVec<T>& mass)
3034{
3035 const std::size_t size = pt.size();
3036
3037 R__ASSERT(eta.size() == size && phi.size() == size && mass.size() == size);
3038
3039 T x_sum = 0.;
3040 T y_sum = 0.;
3041 T z_sum = 0.;
3042 T e_sum = 0.;
3043
3044 for (std::size_t i = 0u; i < size; ++ i) {
3045 // Convert to (e, x, y, z) coordinate system and update sums
3046 const auto x = pt[i] * std::cos(phi[i]);
3047 x_sum += x;
3048 const auto y = pt[i] * std::sin(phi[i]);
3049 y_sum += y;
3050 const auto z = pt[i] * std::sinh(eta[i]);
3051 z_sum += z;
3052 const auto e = std::sqrt(x * x + y * y + z * z + mass[i] * mass[i]);
3053 e_sum += e;
3054 }
3055
3056 // Return invariant mass with (+, -, -, -) metric
3057 return std::sqrt(e_sum * e_sum - x_sum * x_sum - y_sum * y_sum - z_sum * z_sum);
3058}
3059
3060////////////////////////////////////////////////////////////////////////////
3061/// \brief Build an RVec of objects starting from RVecs of input to their constructors.
3062/// \tparam T Type of the objects contained in the created RVec.
3063/// \tparam Args_t Pack of types templating the input RVecs.
3064/// \param[in] args The RVecs containing the values used to initialise the output objects.
3065/// \return The RVec of objects initialised with the input parameters.
3066///
3067/// Example code, at the ROOT prompt:
3068/// ~~~{.cpp}
3069/// using namespace ROOT::VecOps;
3070/// RVecF pts = {15.5, 34.32, 12.95};
3071/// RVecF etas = {0.3, 2.2, 1.32};
3072/// RVecF phis = {0.1, 3.02, 2.2};
3073/// RVecF masses = {105.65, 105.65, 105.65};
3074/// auto fourVecs = Construct<ROOT::Math::PtEtaPhiMVector>(pts, etas, phis, masses);
3075/// cout << fourVecs << endl;
3076/// // { (15.5,0.3,0.1,105.65), (34.32,2.2,3.02,105.65), (12.95,1.32,2.2,105.65) }
3077/// ~~~
3078template <typename T, typename... Args_t>
3080{
3081 const auto size = ::ROOT::Internal::VecOps::GetVectorsSize("Construct", args...);
3082 RVec<T> ret;
3083 ret.reserve(size);
3084 for (auto i = 0UL; i < size; ++i) {
3085 ret.emplace_back(args[i]...);
3086 }
3087 return ret;
3088}
3089
3090/// For any Rvec v produce another RVec with entries starting from 0, and incrementing by 1 until a N = v.size() is reached.
3091/// Example code, at the ROOT prompt:
3092/// ~~~{.cpp}
3093/// using namespace ROOT::VecOps;
3094/// RVecF v = {1., 2., 3.};
3095/// cout << Enumerate(v1) << "\n";
3096/// // { 0, 1, 2 }
3097/// ~~~
3098template <typename T>
3100{
3101 const auto size = v.size();
3102 RVec<T> ret;
3103 ret.reserve(size);
3104 for (auto i = 0UL; i < size; ++i) {
3105 ret.emplace_back(i);
3106 }
3107 return ret;
3108}
3109
3110/// Produce RVec with entries starting from 0, and incrementing by 1 until a user-specified N is reached.
3111/// Example code, at the ROOT prompt:
3112/// ~~~{.cpp}
3113/// using namespace ROOT::VecOps;
3114/// cout << Range(3) << "\n";
3115/// // { 0, 1, 2 }
3116/// ~~~
3117inline RVec<std::size_t> Range(std::size_t length)
3118{
3120 ret.reserve(length);
3121 for (auto i = 0UL; i < length; ++i) {
3122 ret.emplace_back(i);
3123 }
3124 return ret;
3125}
3126
3127/// Produce RVec with entries equal to begin, begin+1, ..., end-1.
3128/// An empty RVec is returned if begin >= end.
3129inline RVec<std::size_t> Range(std::size_t begin, std::size_t end)
3130{
3132 ret.reserve(begin < end ? end - begin : 0u);
3133 for (auto i = begin; i < end; ++i)
3134 ret.push_back(i);
3135 return ret;
3136}
3137
3138/// Allows for negative begin, end, and/or stride. Produce RVec<int> with entries equal to begin, begin+stride, ... , N,
3139/// where N is the first integer such that N+stride exceeds or equals N in the positive or negative direction (same as in Python).
3140/// An empty RVec is returned if begin >= end and stride > 0 or if
3141/// begin < end and stride < 0. Throws a runtime_error if stride==0
3142/// Example code, at the ROOT prompt:
3143/// ~~~{.cpp}
3144/// using namespace ROOT::VecOps;
3145/// cout << Range(1, 5, 2) << "\n";
3146/// // { 1, 3 }
3147/// cout << Range(-1, -11, -4) << "\n";
3148/// // { -1, -5, -9 }
3149/// ~~~
3150inline RVec<long long int> Range(long long int begin, long long int end, long long int stride)
3151{
3152 if (stride==0ll)
3153 {
3154 throw std::runtime_error("Range: the stride must not be zero");
3155 }
3157 float ret_cap = std::ceil(static_cast<float>(end-begin) / stride); //the capacity to reserve
3158 //ret_cap < 0 if either begin > end & stride > 0, or begin < end & stride < 0. In both cases, an empty RVec should be returned
3159 if (ret_cap < 0)
3160 {
3161 return ret;
3162 }
3163 ret.reserve(static_cast<size_t>(ret_cap));
3164 if (stride > 0)
3165 {
3166 for (auto i = begin; i < end; i+=stride)
3167 ret.push_back(i);
3168 }
3169 else
3170 {
3171 for (auto i = begin; i > end; i+=stride)
3172 ret.push_back(i);
3173 }
3174 return ret;
3175}
3176
3177////////////////////////////////////////////////////////////////////////////////
3178/// Print a RVec at the prompt:
3179template <class T>
3180std::ostream &operator<<(std::ostream &os, const RVec<T> &v)
3181{
3182 // In order to print properly, convert to 64 bit int if this is a char
3183 constexpr bool mustConvert = std::is_same<char, T>::value || std::is_same<signed char, T>::value ||
3184 std::is_same<unsigned char, T>::value || std::is_same<wchar_t, T>::value ||
3185 std::is_same<char16_t, T>::value || std::is_same<char32_t, T>::value;
3186 using Print_t = typename std::conditional<mustConvert, long long int, T>::type;
3187 os << "{ ";
3188 auto size = v.size();
3189 if (size) {
3190 for (std::size_t i = 0; i < size - 1; ++i) {
3191 os << (Print_t)v[i] << ", ";
3192 }
3193 os << (Print_t)v[size - 1];
3194 }
3195 os << " }";
3196 return os;
3197}
3198
3199#if (_VECOPS_USE_EXTERN_TEMPLATES)
3200
3201#define RVEC_EXTERN_UNARY_OPERATOR(T, OP) \
3202 extern template RVec<T> operator OP<T>(const RVec<T> &);
3203
3204#define RVEC_EXTERN_BINARY_OPERATOR(T, OP) \
3205 extern template auto operator OP<T, T>(const T &x, const RVec<T> &v) \
3206 -> RVec<decltype(x OP v[0])>; \
3207 extern template auto operator OP<T, T>(const RVec<T> &v, const T &y) \
3208 -> RVec<decltype(v[0] OP y)>; \
3209 extern template auto operator OP<T, T>(const RVec<T> &v0, const RVec<T> &v1)\
3210 -> RVec<decltype(v0[0] OP v1[0])>;
3211
3212#define RVEC_EXTERN_ASSIGN_OPERATOR(T, OP) \
3213 extern template RVec<T> &operator OP<T, T>(RVec<T> &, const T &); \
3214 extern template RVec<T> &operator OP<T, T>(RVec<T> &, const RVec<T> &);
3215
3216#define RVEC_EXTERN_LOGICAL_OPERATOR(T, OP) \
3217 extern template RVec<int> operator OP<T, T>(const RVec<T> &, const T &); \
3218 extern template RVec<int> operator OP<T, T>(const T &, const RVec<T> &); \
3219 extern template RVec<int> operator OP<T, T>(const RVec<T> &, const RVec<T> &);
3220
3221#define RVEC_EXTERN_FLOAT_TEMPLATE(T) \
3222 extern template class RVec<T>; \
3223 RVEC_EXTERN_UNARY_OPERATOR(T, +) \
3224 RVEC_EXTERN_UNARY_OPERATOR(T, -) \
3225 RVEC_EXTERN_UNARY_OPERATOR(T, !) \
3226 RVEC_EXTERN_BINARY_OPERATOR(T, +) \
3227 RVEC_EXTERN_BINARY_OPERATOR(T, -) \
3228 RVEC_EXTERN_BINARY_OPERATOR(T, *) \
3229 RVEC_EXTERN_BINARY_OPERATOR(T, /) \
3230 RVEC_EXTERN_ASSIGN_OPERATOR(T, +=) \
3231 RVEC_EXTERN_ASSIGN_OPERATOR(T, -=) \
3232 RVEC_EXTERN_ASSIGN_OPERATOR(T, *=) \
3233 RVEC_EXTERN_ASSIGN_OPERATOR(T, /=) \
3234 RVEC_EXTERN_LOGICAL_OPERATOR(T, <) \
3235 RVEC_EXTERN_LOGICAL_OPERATOR(T, >) \
3236 RVEC_EXTERN_LOGICAL_OPERATOR(T, ==) \
3237 RVEC_EXTERN_LOGICAL_OPERATOR(T, !=) \
3238 RVEC_EXTERN_LOGICAL_OPERATOR(T, <=) \
3239 RVEC_EXTERN_LOGICAL_OPERATOR(T, >=) \
3240 RVEC_EXTERN_LOGICAL_OPERATOR(T, &&) \
3241 RVEC_EXTERN_LOGICAL_OPERATOR(T, ||)
3242
3243#define RVEC_EXTERN_INTEGER_TEMPLATE(T) \
3244 extern template class RVec<T>; \
3245 RVEC_EXTERN_UNARY_OPERATOR(T, +) \
3246 RVEC_EXTERN_UNARY_OPERATOR(T, -) \
3247 RVEC_EXTERN_UNARY_OPERATOR(T, ~) \
3248 RVEC_EXTERN_UNARY_OPERATOR(T, !) \
3249 RVEC_EXTERN_BINARY_OPERATOR(T, +) \
3250 RVEC_EXTERN_BINARY_OPERATOR(T, -) \
3251 RVEC_EXTERN_BINARY_OPERATOR(T, *) \
3252 RVEC_EXTERN_BINARY_OPERATOR(T, /) \
3253 RVEC_EXTERN_BINARY_OPERATOR(T, %) \
3254 RVEC_EXTERN_BINARY_OPERATOR(T, &) \
3255 RVEC_EXTERN_BINARY_OPERATOR(T, |) \
3256 RVEC_EXTERN_BINARY_OPERATOR(T, ^) \
3257 RVEC_EXTERN_ASSIGN_OPERATOR(T, +=) \
3258 RVEC_EXTERN_ASSIGN_OPERATOR(T, -=) \
3259 RVEC_EXTERN_ASSIGN_OPERATOR(T, *=) \
3260 RVEC_EXTERN_ASSIGN_OPERATOR(T, /=) \
3261 RVEC_EXTERN_ASSIGN_OPERATOR(T, %=) \
3262 RVEC_EXTERN_ASSIGN_OPERATOR(T, &=) \
3263 RVEC_EXTERN_ASSIGN_OPERATOR(T, |=) \
3264 RVEC_EXTERN_ASSIGN_OPERATOR(T, ^=) \
3265 RVEC_EXTERN_ASSIGN_OPERATOR(T, >>=) \
3266 RVEC_EXTERN_ASSIGN_OPERATOR(T, <<=) \
3267 RVEC_EXTERN_LOGICAL_OPERATOR(T, <) \
3268 RVEC_EXTERN_LOGICAL_OPERATOR(T, >) \
3269 RVEC_EXTERN_LOGICAL_OPERATOR(T, ==) \
3270 RVEC_EXTERN_LOGICAL_OPERATOR(T, !=) \
3271 RVEC_EXTERN_LOGICAL_OPERATOR(T, <=) \
3272 RVEC_EXTERN_LOGICAL_OPERATOR(T, >=) \
3273 RVEC_EXTERN_LOGICAL_OPERATOR(T, &&) \
3274 RVEC_EXTERN_LOGICAL_OPERATOR(T, ||)
3275
3276RVEC_EXTERN_INTEGER_TEMPLATE(char)
3277RVEC_EXTERN_INTEGER_TEMPLATE(short)
3278RVEC_EXTERN_INTEGER_TEMPLATE(int)
3279RVEC_EXTERN_INTEGER_TEMPLATE(long)
3280//RVEC_EXTERN_INTEGER_TEMPLATE(long long)
3281
3282RVEC_EXTERN_INTEGER_TEMPLATE(unsigned char)
3283RVEC_EXTERN_INTEGER_TEMPLATE(unsigned short)
3284RVEC_EXTERN_INTEGER_TEMPLATE(unsigned int)
3285RVEC_EXTERN_INTEGER_TEMPLATE(unsigned long)
3286//RVEC_EXTERN_INTEGER_TEMPLATE(unsigned long long)
3287
3288RVEC_EXTERN_FLOAT_TEMPLATE(float)
3289RVEC_EXTERN_FLOAT_TEMPLATE(double)
3290
3291#undef RVEC_EXTERN_UNARY_OPERATOR
3292#undef RVEC_EXTERN_BINARY_OPERATOR
3293#undef RVEC_EXTERN_ASSIGN_OPERATOR
3294#undef RVEC_EXTERN_LOGICAL_OPERATOR
3295#undef RVEC_EXTERN_INTEGER_TEMPLATE
3296#undef RVEC_EXTERN_FLOAT_TEMPLATE
3297
3298#define RVEC_EXTERN_UNARY_FUNCTION(T, NAME, FUNC) \
3299 extern template RVec<PromoteType<T>> NAME(const RVec<T> &);
3300
3301#define RVEC_EXTERN_STD_UNARY_FUNCTION(T, F) RVEC_EXTERN_UNARY_FUNCTION(T, F, std::F)
3302
3303#define RVEC_EXTERN_BINARY_FUNCTION(T0, T1, NAME, FUNC) \
3304 extern template RVec<PromoteTypes<T0, T1>> NAME(const RVec<T0> &, const T1 &); \
3305 extern template RVec<PromoteTypes<T0, T1>> NAME(const T0 &, const RVec<T1> &); \
3306 extern template RVec<PromoteTypes<T0, T1>> NAME(const RVec<T0> &, const RVec<T1> &);
3307
3308#define RVEC_EXTERN_STD_BINARY_FUNCTION(T, F) RVEC_EXTERN_BINARY_FUNCTION(T, T, F, std::F)
3309
3310#define RVEC_EXTERN_STD_FUNCTIONS(T) \
3311 RVEC_EXTERN_STD_UNARY_FUNCTION(T, abs) \
3312 RVEC_EXTERN_STD_BINARY_FUNCTION(T, fdim) \
3313 RVEC_EXTERN_STD_BINARY_FUNCTION(T, fmod) \
3314 RVEC_EXTERN_STD_BINARY_FUNCTION(T, remainder) \
3315 RVEC_EXTERN_STD_UNARY_FUNCTION(T, exp) \
3316 RVEC_EXTERN_STD_UNARY_FUNCTION(T, exp2) \
3317 RVEC_EXTERN_STD_UNARY_FUNCTION(T, expm1) \
3318 RVEC_EXTERN_STD_UNARY_FUNCTION(T, log) \
3319 RVEC_EXTERN_STD_UNARY_FUNCTION(T, log10) \
3320 RVEC_EXTERN_STD_UNARY_FUNCTION(T, log2) \
3321 RVEC_EXTERN_STD_UNARY_FUNCTION(T, log1p) \
3322 RVEC_EXTERN_STD_BINARY_FUNCTION(T, pow) \
3323 RVEC_EXTERN_STD_UNARY_FUNCTION(T, sqrt) \
3324 RVEC_EXTERN_STD_UNARY_FUNCTION(T, cbrt) \
3325 RVEC_EXTERN_STD_BINARY_FUNCTION(T, hypot) \
3326 RVEC_EXTERN_STD_UNARY_FUNCTION(T, sin) \
3327 RVEC_EXTERN_STD_UNARY_FUNCTION(T, cos) \
3328 RVEC_EXTERN_STD_UNARY_FUNCTION(T, tan) \
3329 RVEC_EXTERN_STD_UNARY_FUNCTION(T, asin) \
3330 RVEC_EXTERN_STD_UNARY_FUNCTION(T, acos) \
3331 RVEC_EXTERN_STD_UNARY_FUNCTION(T, atan) \
3332 RVEC_EXTERN_STD_BINARY_FUNCTION(T, atan2) \
3333 RVEC_EXTERN_STD_UNARY_FUNCTION(T, sinh) \
3334 RVEC_EXTERN_STD_UNARY_FUNCTION(T, cosh) \
3335 RVEC_EXTERN_STD_UNARY_FUNCTION(T, tanh) \
3336 RVEC_EXTERN_STD_UNARY_FUNCTION(T, asinh) \
3337 RVEC_EXTERN_STD_UNARY_FUNCTION(T, acosh) \
3338 RVEC_EXTERN_STD_UNARY_FUNCTION(T, atanh) \
3339 RVEC_EXTERN_STD_UNARY_FUNCTION(T, floor) \
3340 RVEC_EXTERN_STD_UNARY_FUNCTION(T, ceil) \
3341 RVEC_EXTERN_STD_UNARY_FUNCTION(T, trunc) \
3342 RVEC_EXTERN_STD_UNARY_FUNCTION(T, round) \
3343 RVEC_EXTERN_STD_UNARY_FUNCTION(T, erf) \
3344 RVEC_EXTERN_STD_UNARY_FUNCTION(T, erfc) \
3345 RVEC_EXTERN_STD_UNARY_FUNCTION(T, lgamma) \
3346 RVEC_EXTERN_STD_UNARY_FUNCTION(T, tgamma) \
3347
3348RVEC_EXTERN_STD_FUNCTIONS(float)
3349RVEC_EXTERN_STD_FUNCTIONS(double)
3350#undef RVEC_EXTERN_STD_UNARY_FUNCTION
3351#undef RVEC_EXTERN_STD_BINARY_FUNCTION
3352#undef RVEC_EXTERN_STD_UNARY_FUNCTIONS
3353
3354#ifdef R__HAS_VDT
3355
3356#define RVEC_EXTERN_VDT_UNARY_FUNCTION(T, F) RVEC_EXTERN_UNARY_FUNCTION(T, F, vdt::F)
3357
3358RVEC_EXTERN_VDT_UNARY_FUNCTION(float, fast_expf)
3359RVEC_EXTERN_VDT_UNARY_FUNCTION(float, fast_logf)
3360RVEC_EXTERN_VDT_UNARY_FUNCTION(float, fast_sinf)
3361RVEC_EXTERN_VDT_UNARY_FUNCTION(float, fast_cosf)
3362RVEC_EXTERN_VDT_UNARY_FUNCTION(float, fast_tanf)
3363RVEC_EXTERN_VDT_UNARY_FUNCTION(float, fast_asinf)
3364RVEC_EXTERN_VDT_UNARY_FUNCTION(float, fast_acosf)
3365RVEC_EXTERN_VDT_UNARY_FUNCTION(float, fast_atanf)
3366
3367RVEC_EXTERN_VDT_UNARY_FUNCTION(double, fast_exp)
3368RVEC_EXTERN_VDT_UNARY_FUNCTION(double, fast_log)
3369RVEC_EXTERN_VDT_UNARY_FUNCTION(double, fast_sin)
3370RVEC_EXTERN_VDT_UNARY_FUNCTION(double, fast_cos)
3371RVEC_EXTERN_VDT_UNARY_FUNCTION(double, fast_tan)
3372RVEC_EXTERN_VDT_UNARY_FUNCTION(double, fast_asin)
3373RVEC_EXTERN_VDT_UNARY_FUNCTION(double, fast_acos)
3374RVEC_EXTERN_VDT_UNARY_FUNCTION(double, fast_atan)
3375
3376#endif // R__HAS_VDT
3377
3378#endif // _VECOPS_USE_EXTERN_TEMPLATES
3379
3380/** @} */ // end of Doxygen group vecops
3381
3382} // End of VecOps NS
3383
3384// Allow to use RVec as ROOT::RVec
3385using ROOT::VecOps::RVec;
3386
3397
3398} // End of ROOT NS
3399
3400#endif // ROOT_RVEC
size_t fSize
#define R__unlikely(expr)
Definition RConfig.hxx:586
#define f(i)
Definition RSha256.hxx:104
#define c(i)
Definition RSha256.hxx:101
#define e(i)
Definition RSha256.hxx:103
#define R__RVEC_NODISCARD
Definition RVec.hxx:19
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define M_PI
Definition Rotated.cxx:105
TBuffer & operator<<(TBuffer &buf, const Tmpl *obj)
Definition TBuffer.h:399
#define R__CLING_PTRCHECK(ONOFF)
Definition Rtypes.h:500
static Double_t Product(const Double_t *x, const Float_t *y)
Product.
Definition TCTUB.cxx:101
#define X(type, name)
#define R__ASSERT(e)
Definition TError.h:118
#define N
Int_t Compare(const void *item1, const void *item2)
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h length
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
Option_t Option_t TPoint TPoint const char y1
#define free
Definition civetweb.c:1539
#define malloc
Definition civetweb.c:1536
This class consists of common code factored out of the SmallVector class to reduce code duplication b...
Definition RVec.hxx:561
void assign(size_type NumElts, const T &Elt)
Definition RVec.hxx:680
typename SuperClass::iterator iterator
Definition RVec.hxx:565
typename SuperClass::size_type size_type
Definition RVec.hxx:568
void append(in_iter in_start, in_iter in_end)
Add the specified range to the end of the SmallVector.
Definition RVec.hxx:654
iterator insert(iterator I, T &&Elt)
Definition RVec.hxx:741
void resize(size_type N)
Definition RVec.hxx:596
void assign(std::initializer_list< T > IL)
Definition RVec.hxx:698
typename SuperClass::const_iterator const_iterator
Definition RVec.hxx:566
void resize(size_type N, const T &NV)
Definition RVec.hxx:611
void reserve(size_type N)
Definition RVec.hxx:625
iterator insert(iterator I, ItTy From, ItTy To)
Definition RVec.hxx:859
reference emplace_back(ArgTypes &&...Args)
Definition RVec.hxx:920
void assign(in_iter in_start, in_iter in_end)
Definition RVec.hxx:692
iterator insert(iterator I, const T &Elt)
Definition RVec.hxx:773
void swap(RVecImpl &RHS)
Definition RVec.hxx:935
iterator insert(iterator I, size_type NumToInsert, const T &Elt)
Definition RVec.hxx:804
RVecImpl & operator=(const RVecImpl &RHS)
Definition RVec.hxx:995
iterator erase(const_iterator CS, const_iterator CE)
Definition RVec.hxx:721
typename SuperClass::reference reference
Definition RVec.hxx:567
void append(size_type NumInputs, const T &Elt)
Append NumInputs copies of Elt to the end.
Definition RVec.hxx:665
iterator erase(const_iterator CI)
Definition RVec.hxx:704
RVecImpl & operator=(RVecImpl &&RHS)
Definition RVec.hxx:1048
void pop_back_n(size_type NumItems)
Definition RVec.hxx:631
RVecImpl(const RVecImpl &)=delete
void append(std::initializer_list< T > IL)
Definition RVec.hxx:674
void insert(iterator I, std::initializer_list< T > IL)
Definition RVec.hxx:917
This is all the stuff common to all SmallVectors.
Definition RVec.hxx:138
SmallVectorBase(void *FirstEl, size_t TotalCapacity)
Definition RVec.hxx:156
static constexpr size_t SizeTypeMax()
The maximum value of the Size_T used.
Definition RVec.hxx:153
Size_T fCapacity
Always >= -1. fCapacity == -1 indicates the RVec is in "memory adoption" mode.
Definition RVec.hxx:150
bool Owns() const
If false, the RVec is in "memory adoption" mode, i.e. it is acting as a view on a memory buffer it do...
Definition RVec.hxx:171
size_t capacity() const noexcept
Definition RVec.hxx:175
void set_size(size_t N)
Set the array size to N, which the current array must have enough capacity for.
Definition RVec.hxx:188
void grow(size_t MinSize=0)
Double the size of the allocated memory, guaranteeing space for at least one more element or MinSize ...
Definition RVec.hxx:472
static void uninitialized_move(It1 I, It1 E, It2 Dest)
Move the range [I, E) onto the uninitialized memory starting with "Dest", constructing elements into ...
Definition RVec.hxx:440
typename SuperClass::const_iterator const_iterator
Definition RVec.hxx:479
static void uninitialized_copy(T1 *I, T1 *E, T2 *Dest, typename std::enable_if< std::is_same< typename std::remove_const< T1 >::type, T2 >::value >::type *=nullptr)
Copy the range [I, E) onto the uninitialized memory starting with "Dest", constructing elements into ...
Definition RVec.hxx:458
static void uninitialized_copy(It1 I, It1 E, It2 Dest)
Copy the range [I, E) onto the uninitialized memory starting with "Dest", constructing elements into ...
Definition RVec.hxx:449
SmallVectorTemplateBase<TriviallyCopyable = false> - This is where we put method implementations that...
Definition RVec.hxx:328
void grow(size_t MinSize=0)
Grow the allocated memory (without initializing new elements), doubling the size of the allocated mem...
static void uninitialized_move(It1 I, It1 E, It2 Dest)
Move the range [I, E) into the uninitialized memory starting with "Dest", constructing elements as ne...
Definition RVec.hxx:343
static void uninitialized_copy(It1 I, It1 E, It2 Dest)
Copy the range [I, E) onto the uninitialized memory starting with "Dest", constructing elements as ne...
Definition RVec.hxx:351
This is the part of SmallVectorTemplateBase which does not depend on whether the type T is a POD.
Definition RVec.hxx:206
const_iterator cbegin() const noexcept
Definition RVec.hxx:261
void grow_pod(size_t MinSize, size_t TSize)
Definition RVec.hxx:222
const_iterator cend() const noexcept
Definition RVec.hxx:264
void resetToSmall()
Put this vector in a state of being small.
Definition RVec.hxx:229
std::reverse_iterator< iterator > reverse_iterator
Definition RVec.hxx:247
bool isSmall() const
Return true if this is a smallvector which has not had dynamic memory allocated for it.
Definition RVec.hxx:226
const_reverse_iterator crend() const noexcept
Definition RVec.hxx:272
const_iterator end() const noexcept
Definition RVec.hxx:263
const_reverse_iterator crbegin() const noexcept
Definition RVec.hxx:269
pointer data() noexcept
Return a pointer to the vector's buffer, even if empty().
Definition RVec.hxx:280
const_reverse_iterator rbegin() const noexcept
Definition RVec.hxx:268
std::reverse_iterator< const_iterator > const_reverse_iterator
Definition RVec.hxx:246
const_iterator begin() const noexcept
Definition RVec.hxx:260
const_pointer data() const noexcept
Return a pointer to the vector's buffer, even if empty().
Definition RVec.hxx:282
void * getFirstEl() const
Find the address of the first element.
Definition RVec.hxx:212
const_reverse_iterator rend() const noexcept
Definition RVec.hxx:271
RVecN(size_t Size)
Definition RVec.hxx:1165
RVecN(Detail::VecOps::RVecImpl< T > &&RHS)
Definition RVec.hxx:1201
reference operator[](size_type idx)
Definition RVec.hxx:1241
typename Internal::VecOps::SmallVectorTemplateCommon< T >::const_reference const_reference
Definition RVec.hxx:1235
RVecN operator[](const RVecN< V, M > &conds) const
Definition RVec.hxx:1252
RVecN(std::initializer_list< T > IL)
Definition RVec.hxx:1181
const_reference at(size_type pos) const
Definition RVec.hxx:1295
RVecN(const RVecN &RHS)
Definition RVec.hxx:1183
RVecN & operator=(Detail::VecOps::RVecImpl< T > &&RHS)
Definition RVec.hxx:1222
RVecN & operator=(RVecN &&RHS)
Definition RVec.hxx:1209
typename Internal::VecOps::SmallVectorTemplateCommon< T >::size_type size_type
Definition RVec.hxx:1236
value_type at(size_type pos, value_type fallback) const
No exception thrown. The user specifies the desired value in case the RVecN is shorter than pos.
Definition RVec.hxx:1314
RVecN & operator=(std::initializer_list< T > IL)
Definition RVec.hxx:1228
RVecN & operator=(const RVecN &RHS)
Definition RVec.hxx:1189
RVecN(const std::vector< T > &RHS)
Definition RVec.hxx:1207
RVecN(size_t Size, const T &Value)
Definition RVec.hxx:1163
RVecN(RVecN &&RHS)
Definition RVec.hxx:1195
RVecN(ItTy S, ItTy E)
Definition RVec.hxx:1176
reference at(size_type pos)
Definition RVec.hxx:1285
value_type at(size_type pos, value_type fallback)
No exception thrown. The user specifies the desired value in case the RVecN is shorter than pos.
Definition RVec.hxx:1306
RVecN(T *p, size_t n)
Definition RVec.hxx:1215
typename Internal::VecOps::SmallVectorTemplateCommon< T >::reference reference
Definition RVec.hxx:1234
typename Internal::VecOps::SmallVectorTemplateCommon< T >::value_type value_type
Definition RVec.hxx:1237
const_reference operator[](size_type idx) const
Definition RVec.hxx:1246
A "std::vector"-like collection of values implementing handy operation to analyse them.
Definition RVec.hxx:1492
RVec(RVecN< T, N > &&RHS)
Definition RVec.hxx:1539
typename SuperClass::reference reference
Definition RVec.hxx:1498
RVec(const RVecN< T, N > &RHS)
Definition RVec.hxx:1542
RVec(size_t Size, const T &Value)
Definition RVec.hxx:1507
RVec(const RVec &RHS)
Definition RVec.hxx:1520
RVec & operator=(RVec &&RHS)
Definition RVec.hxx:1530
RVec(T *p, size_t n)
Definition RVec.hxx:1546
RVec operator[](const RVec< V > &conds) const
Definition RVec.hxx:1558
RVec(std::initializer_list< T > IL)
Definition RVec.hxx:1518
typename SuperClass::const_reference const_reference
Definition RVec.hxx:1499
RVec(size_t Size)
Definition RVec.hxx:1509
RVec(ItTy S, ItTy E)
Definition RVec.hxx:1514
RVec(const std::vector< T > &RHS)
Definition RVec.hxx:1544
typename SuperClass::size_type size_type
Definition RVec.hxx:1500
RVec(Detail::VecOps::RVecImpl< T > &&RHS)
Definition RVec.hxx:1536
RVec(RVec &&RHS)
Definition RVec.hxx:1528
typename SuperClass::value_type value_type
Definition RVec.hxx:1501
RVec & operator=(const RVec &RHS)
Definition RVec.hxx:1522
TPaveText * pt
Type
enumeration specifying the integration types.
RVec< T > Reverse(const RVec< T > &v)
Return copy of reversed vector.
Definition RVec.hxx:2444
RVec< T > Intersect(const RVec< T > &v1, const RVec< T > &v2, bool v2_is_sorted=false)
Return the intersection of elements of two RVecs.
Definition RVec.hxx:2721
RVec< typename RVec< T >::size_type > Nonzero(const RVec< T > &v)
Return the indices of the elements which are not zero.
Definition RVec.hxx:2690
#define RVEC_UNARY_OPERATOR(OP)
Definition RVec.hxx:1579
#define RVEC_ASSIGNMENT_OPERATOR(OP)
Definition RVec.hxx:1650
RVec< typename RVec< T >::size_type > StableArgsort(const RVec< T > &v)
Return an RVec of indices that sort the input RVec while keeping the order of equal elements.
Definition RVec.hxx:2258
RVec< Common_t > Concatenate(const RVec< T0 > &v0, const RVec< T1 > &v1)
Return the concatenation of two RVecs.
Definition RVec.hxx:2860
RVec< T > InvariantMasses(const RVec< T > &pt1, const RVec< T > &eta1, const RVec< T > &phi1, const RVec< T > &mass1, const RVec< T > &pt2, const RVec< T > &eta2, const RVec< T > &phi2, const RVec< T > &mass2)
Return the invariant mass of two particles given the collections of the quantities transverse momentu...
Definition RVec.hxx:2991
T Sum(const RVec< T > &v, const T zero=T(0))
Sum elements of an RVec.
Definition RVec.hxx:1917
RVec< T > Take(const RVec< T > &v, const RVec< typename RVec< T >::size_type > &i)
Return elements of a vector at given indices.
Definition RVec.hxx:2302
RVec< T > Construct(const RVec< Args_t > &... args)
Build an RVec of objects starting from RVecs of input to their constructors.
Definition RVec.hxx:3079
#define RVEC_STD_BINARY_FUNCTION(F)
Definition RVec.hxx:1793
#define RVEC_BINARY_OPERATOR(OP)
Definition RVec.hxx:1602
RVec< T > Drop(const RVec< T > &v, RVec< typename RVec< T >::size_type > idxs)
Return a copy of the container without the elements at the specified indices.
Definition RVec.hxx:2411
size_t CapacityInBytes(const RVecN< T, N > &X)
Definition RVec.hxx:1571
#define RVEC_LOGICAL_OPERATOR(OP)
Definition RVec.hxx:1686
RVec< RVec< std::size_t > > Combinations(const std::size_t size1, const std::size_t size2)
Return the indices that represent all combinations of the elements of two RVecs.
Definition RVec.hxx:2569
#define RVEC_STD_UNARY_FUNCTION(F)
Definition RVec.hxx:1792
RVec< typename RVec< T >::size_type > Enumerate(const RVec< T > &v)
For any Rvec v produce another RVec with entries starting from 0, and incrementing by 1 until a N = v...
Definition RVec.hxx:3099
auto Map(Args &&... args)
Create new collection applying a callable to the elements of the input collection.
Definition RVec.hxx:2113
RVec< T > Where(const RVec< int > &c, const RVec< T > &v1, const RVec< T > &v2)
Return the elements of v1 if the condition c is true and v2 if the condition c is false.
Definition RVec.hxx:2755
auto Any(const RVec< T > &v) -> decltype(v[0]==true)
Return true if any of the elements equates to true, return false otherwise.
Definition RVec.hxx:2168
RVec< typename RVec< T >::size_type > Argsort(const RVec< T > &v)
Return an RVec of indices that sort the input RVec.
Definition RVec.hxx:2213
std::size_t ArgMin(const RVec< T > &v)
Get the index of the smallest element of an RVec In case of multiple occurrences of the minimum value...
Definition RVec.hxx:2048
RVec< T > StableSort(const RVec< T > &v)
Return copy of RVec with elements sorted in ascending order while keeping the order of equal elements...
Definition RVec.hxx:2514
double Var(const RVec< T > &v)
Get the variance of the elements of an RVec.
Definition RVec.hxx:2065
RVec< T > Filter(const RVec< T > &v, F &&f)
Create a new collection with the elements passing the filter expressed by the predicate.
Definition RVec.hxx:2145
std::size_t ArgMax(const RVec< T > &v)
Get the index of the greatest element of an RVec In case of multiple occurrences of the maximum value...
Definition RVec.hxx:2030
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
#define T2
Definition md5.inl:147
#define F(x, y, z)
#define I(x, y, z)
#define T1
Definition md5.inl:146
bool IsSmall(const ROOT::VecOps::RVec< T > &v)
Definition RVec.hxx:1118
bool IsAdopting(const ROOT::VecOps::RVec< T > &v)
Definition RVec.hxx:1124
auto MapImpl(F &&f, RVecs &&... vs) -> RVec< decltype(f(vs[0]...))>
Definition RVec.hxx:105
void ResetView(RVec< T > &v, T *addr, std::size_t sz)
An unsafe function to reset the buffer for which this RVec is acting as a view.
Definition RVec.hxx:546
uint64_t NextPowerOf2(uint64_t A)
Return the next power of two (in 64-bits) that is strictly greater than A.
Definition RVec.hxx:126
constexpr bool All(const bool *vals, std::size_t size)
Definition RVec.hxx:79
std::size_t GetVectorsSize(const std::string &id, const RVec< T > &... vs)
Definition RVec.hxx:88
void UninitializedValueConstruct(ForwardIt first, ForwardIt last)
Definition RVec.hxx:530
auto MapFromTuple(Tuple_t &&t, std::index_sequence< Is... >) -> decltype(MapImpl(std::get< std::tuple_size< Tuple_t >::value - 1 >(t), std::get< Is >(t)...))
Definition RVec.hxx:117
This file contains a specialised ROOT message handler to test for diagnostic in unit tests.
Definition first.py:1
#define Dot(u, v)
Definition normal.c:49
The size of the inline storage of an RVec.
Definition RVec.hxx:512
Used to figure out the offset of the first element of an RVec.
Definition RVec.hxx:199
Storage for the SmallVector elements.
Definition RVec.hxx:497
Ta Range(0, 0, 1, 1)
TMarker m
Definition textangle.C:8