GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TScaler.h
Go to the documentation of this file.
1#ifndef TSCALER_H
2#define TSCALER_H
3
4/** \addtogroup Sorting
5 * @{
6 */
7
8/*
9 * Author: V. Bildstein, <vbildste@uoguelph.ca>
10 *
11 * Based on the TPPG class
12 *
13 * Please indicate changes with your initials.
14 *
15 */
16
17//////////////////////////////////////////////////////////////////////////
18///
19/// \class TScaler
20///
21/// The TScaler is designed to hold all of the information about the
22/// scaler status.
23///
24//////////////////////////////////////////////////////////////////////////
25
26#include <iostream>
27#include <map>
28
29#include "TObject.h"
30#include "Globals.h"
31#include "TTree.h"
32#include "TH1.h"
33
34#include "TPPG.h"
35
36class TScalerData : public TObject {
37public:
40 TScalerData(TScalerData&&) noexcept = default;
41 TScalerData& operator=(const TScalerData&) = default;
42 TScalerData& operator=(TScalerData&&) noexcept = default;
43 ~TScalerData() = default;
44
45 void Copy(TObject& rhs) const override;
46
47 void SetAddress(UInt_t address) { fAddress = address; }
48 void SetNetworkPacketId(UInt_t networkId) { fNetworkPacketId = networkId; }
49 void SetLowTimeStamp(UInt_t lowTime) { fLowTimeStamp = lowTime; }
50 void SetHighTimeStamp(UInt_t highTime) { fHighTimeStamp = highTime; }
51 void SetScaler(size_t index, UInt_t scaler)
52 {
53 if(index < fScaler.size()) {
54 fScaler[index] = scaler;
55 } else {
56 std::cout << "Failed to set scaler " << scaler << ", index " << index << " is out of range 0 - "
57 << fScaler.size() << std::endl;
58 }
59 }
60 void SetScaler(UInt_t* data, int size)
61 {
62 fScaler = std::vector<uint32_t>(data, data + size);
63 }
64
65 UInt_t GetAddress() const { return fAddress; }
66 UInt_t GetNetworkPacketId() const { return fNetworkPacketId; }
67 UInt_t GetLowTimeStamp() const { return fLowTimeStamp; }
68 UInt_t GetHighTimeStamp() const { return fHighTimeStamp; }
69 std::vector<UInt_t> GetScaler() const { return fScaler; }
70 UInt_t GetScaler(size_t index) const
71 {
72 if(index < fScaler.size()) {
73 return fScaler[index];
74 }
75 return 0;
76 }
77
78 ULong64_t GetTimeStamp() const
79 {
80 ULong64_t time = GetHighTimeStamp();
81 time = time << 28;
82 time |= GetLowTimeStamp() & 0x0fffffff;
83 return 10 * time; // convert from raw 10 ns units to ns
84 }
85
86 void ResizeScaler(size_t newSize = 1) { fScaler.resize(newSize); }
87
88 void Print(Option_t* opt = "") const override;
89 void Clear(Option_t* opt = "") override;
90
91private:
93 UInt_t fAddress{0};
94 std::vector<UInt_t> fScaler;
95 UInt_t fLowTimeStamp{0};
96 UInt_t fHighTimeStamp{0};
97
98 /// \cond CLASSIMP
99 ClassDefOverride(TScalerData, 2) // NOLINT(readability-else-after-return)
100 /// \endcond
101};
102
103class TScaler : public TObject {
104public:
105 explicit TScaler(bool loadIntoMap = false);
106 explicit TScaler(TTree*, bool loadIntoMap = false);
107 TScaler(const TScaler&);
108 TScaler(TScaler&&) noexcept = default;
109 TScaler& operator=(const TScaler&) = default;
110 TScaler& operator=(TScaler&&) noexcept = default;
111 ~TScaler();
112
113 void Copy(TObject& obj) const override;
114
115 std::vector<UInt_t> GetScaler(UInt_t address, ULong64_t time) const;
116 UInt_t GetScaler(UInt_t address, ULong64_t time, size_t index) const;
117 UInt_t GetScalerDifference(UInt_t address, ULong64_t time, size_t index) const;
118 ULong64_t NumberOfScalerReadouts() const { return fEntries; };
119
120 ULong64_t GetTimePeriod(UInt_t address);
121
122 std::map<UInt_t, ULong64_t> GetTimePeriodMap() { return fTimePeriod; }
123 std::map<UInt_t, std::map<ULong64_t, int>> GetNumberOfTimePeriods() { return fNumberOfTimePeriods; }
124
125 void Clear(Option_t* opt = "") override;
126 using TObject::Draw; // This is to remove hidden overload
127 TH1D* Draw(UInt_t address, size_t index = 0, Option_t* option = "");
128 TH1D* Draw(UInt_t lowAddress, UInt_t highAddress, size_t index = 0, Option_t* option = "");
129 TH1D* DrawRawTimes(UInt_t address, Double_t lowtime, Double_t hightime, size_t index = 0, Option_t* option = "");
130
131 void ListHistograms();
132
133 static double GetLastScaler(TTree* tree = nullptr, UInt_t address = 0xffff, size_t nominator = 2, size_t denominator = 1)
134 {
135 /// This function returns the ratio of the two scalers nominator and denominator for the last entry with a matching address for a given tree.
136 /// If no tree is provided it tries to get the "RateScaler" tree from the current directory.
137 if(tree == nullptr) {
138 tree = static_cast<TTree*>(gDirectory->Get("RateScaler"));
139 if(tree == nullptr) {
140 std::cerr << __PRETTY_FUNCTION__ << ": no tree provided and failed to find \"RateScaler\" in " << gDirectory->GetName() << std::endl; // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
141 return -1.;
142 }
143 }
144 TScalerData* scalerData = nullptr;
145 tree->SetBranchAddress("ScalerData", &scalerData);
146
147 for(Long64_t entry = tree->GetEntries() - 1; entry >= 0; --entry) {
148 tree->GetEntry(entry);
149 if(scalerData->GetAddress() != address) { continue; }
150 auto scalers = scalerData->GetScaler();
151 if(nominator >= scalers.size() || denominator >= scalers.size()) {
152 std::cerr << __PRETTY_FUNCTION__ << ": trying to get nominator " << nominator << " and denominator " << denominator << " from vector of size " << scalers.size() << std::endl; // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
153 return -1.;
154 }
155 return static_cast<double>(scalers[nominator]) / scalers[denominator];
156 }
157 std::cerr << __PRETTY_FUNCTION__ << ": failed to find any entry for address " << hex(address, 4) << std::endl; // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
158 return -1.;
159 }
160
161private:
162 void ReadTree(bool loadIntoMap);
163
164 TTree* fTree;
166 Long64_t fEntries;
167 std::map<UInt_t, std::map<ULong64_t, std::vector<UInt_t>>> fScalerMap; //!<! an address-map of timestamp mapped scaler values
168 std::map<UInt_t, ULong64_t> fTimePeriod; //!<! a map between addresses and time differences (used to calculate the time period)
169 std::map<UInt_t, std::map<ULong64_t, int>> fNumberOfTimePeriods; //!<!
170 ULong64_t fTotalTimePeriod; //!<!
171 std::map<ULong64_t, int> fTotalNumberOfTimePeriods; //!<!
172 TPPG* fPPG; //!<!
173 std::map<UInt_t, TH1D*> fHist; //!<! map to store histograms that have already been drawn
174 std::map<std::pair<UInt_t, UInt_t>, TH1D*> fHistRange; //!<! map to store histograms for address-ranges
175
176 /// \cond CLASSIMP
177 ClassDefOverride(TScaler, 2) // NOLINT(readability-else-after-return)
178 /// \endcond
179};
180/*! @} */
181#endif
std::string hex(T val, int width=-1)
Definition Globals.h:129
Definition TPPG.h:130
void Clear(Option_t *opt="") override
Definition TScaler.cxx:27
std::vector< UInt_t > fScaler
Definition TScaler.h:94
UInt_t GetScaler(size_t index) const
Definition TScaler.h:70
ULong64_t GetTimeStamp() const
Definition TScaler.h:78
UInt_t fAddress
Definition TScaler.h:93
UInt_t fNetworkPacketId
Definition TScaler.h:92
void Copy(TObject &rhs) const override
Definition TScaler.cxx:18
UInt_t fLowTimeStamp
Definition TScaler.h:95
UInt_t GetNetworkPacketId() const
Definition TScaler.h:66
UInt_t fHighTimeStamp
Definition TScaler.h:96
std::vector< UInt_t > GetScaler() const
Definition TScaler.h:69
void Print(Option_t *opt="") const override
Definition TScaler.cxx:37
void SetScaler(size_t index, UInt_t scaler)
Definition TScaler.h:51
void SetAddress(UInt_t address)
Definition TScaler.h:47
UInt_t GetAddress() const
Definition TScaler.h:65
void SetScaler(UInt_t *data, int size)
Definition TScaler.h:60
void SetNetworkPacketId(UInt_t networkId)
Definition TScaler.h:48
void ResizeScaler(size_t newSize=1)
Definition TScaler.h:86
UInt_t GetLowTimeStamp() const
Definition TScaler.h:67
void SetLowTimeStamp(UInt_t lowTime)
Definition TScaler.h:49
UInt_t GetHighTimeStamp() const
Definition TScaler.h:68
TScalerData(TScalerData &&) noexcept=default
void SetHighTimeStamp(UInt_t highTime)
Definition TScaler.h:50
std::map< UInt_t, std::map< ULong64_t, int > > fNumberOfTimePeriods
!
Definition TScaler.h:169
std::map< UInt_t, std::map< ULong64_t, int > > GetNumberOfTimePeriods()
Definition TScaler.h:123
ULong64_t fTotalTimePeriod
!
Definition TScaler.h:170
std::map< UInt_t, ULong64_t > GetTimePeriodMap()
Definition TScaler.h:122
TTree * fTree
Definition TScaler.h:164
void Clear(Option_t *opt="") override
Definition TScaler.cxx:228
TH1D * DrawRawTimes(UInt_t address, Double_t lowtime, Double_t hightime, size_t index=0, Option_t *option="")
Definition TScaler.cxx:450
static double GetLastScaler(TTree *tree=nullptr, UInt_t address=0xffff, size_t nominator=2, size_t denominator=1)
Definition TScaler.h:133
void ListHistograms()
Definition TScaler.cxx:523
std::map< ULong64_t, int > fTotalNumberOfTimePeriods
!
Definition TScaler.h:171
std::map< UInt_t, ULong64_t > fTimePeriod
! a map between addresses and time differences (used to calculate the time period)
Definition TScaler.h:168
std::vector< UInt_t > GetScaler(UInt_t address, ULong64_t time) const
Definition TScaler.cxx:108
std::map< std::pair< UInt_t, UInt_t >, TH1D * > fHistRange
! map to store histograms for address-ranges
Definition TScaler.h:174
TScalerData * fScalerData
Definition TScaler.h:165
TScaler(bool loadIntoMap=false)
Definition TScaler.cxx:46
std::map< UInt_t, TH1D * > fHist
! map to store histograms that have already been drawn
Definition TScaler.h:173
UInt_t GetScalerDifference(UInt_t address, ULong64_t time, size_t index) const
Definition TScaler.cxx:163
Long64_t fEntries
Definition TScaler.h:166
std::map< UInt_t, std::map< ULong64_t, std::vector< UInt_t > > > fScalerMap
! an address-map of timestamp mapped scaler values
Definition TScaler.h:167
void Copy(TObject &obj) const override
Definition TScaler.cxx:95
ULong64_t GetTimePeriod(UInt_t address)
Definition TScaler.cxx:491
TScaler(TScaler &&) noexcept=default
TPPG * fPPG
!
Definition TScaler.h:172
void ReadTree(bool loadIntoMap)
Definition TScaler.cxx:60
TH1D * Draw(UInt_t address, size_t index=0, Option_t *option="")
Definition TScaler.cxx:255
ULong64_t NumberOfScalerReadouts() const
Definition TScaler.h:118