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