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