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