![]() |
GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
|
![]() |
The TScaler is designed to hold all of the information about the scaler status.
Public Member Functions | |
TScaler (bool loadIntoMap=false) | |
TScaler (const TScaler &) | |
TScaler (TScaler &&) noexcept=default | |
TScaler (TTree *, bool loadIntoMap=false) | |
~TScaler () | |
void | Clear (Option_t *opt="") override |
void | Copy (TObject &obj) const override |
TH1D * | Draw (UInt_t address, size_t index=0, Option_t *option="") |
TH1D * | Draw (UInt_t lowAddress, UInt_t highAddress, size_t index=0, Option_t *option="") |
TH1D * | DrawRawTimes (UInt_t address, Double_t lowtime, Double_t hightime, size_t index=0, Option_t *option="") |
std::map< UInt_t, std::map< ULong64_t, int > > | GetNumberOfTimePeriods () |
std::vector< UInt_t > | GetScaler (UInt_t address, ULong64_t time) const |
UInt_t | GetScaler (UInt_t address, ULong64_t time, size_t index) const |
UInt_t | GetScalerDifference (UInt_t address, ULong64_t time, size_t index) const |
ULong64_t | GetTimePeriod (UInt_t address) |
std::map< UInt_t, ULong64_t > | GetTimePeriodMap () |
void | ListHistograms () |
ULong64_t | NumberOfScalerReadouts () const |
TScaler & | operator= (const TScaler &)=default |
TScaler & | operator= (TScaler &&) noexcept=default |
Static Public Member Functions | |
static double | GetLastScaler (TTree *tree=nullptr, UInt_t address=0xffff, size_t nominator=2, size_t denominator=1) |
Private Member Functions | |
void | ReadTree (bool loadIntoMap) |
Private Attributes | |
Long64_t | fEntries |
std::map< UInt_t, TH1D * > | fHist |
std::map< std::pair< UInt_t, UInt_t >, TH1D * > | fHistRange |
std::map< UInt_t, std::map< ULong64_t, int > > | fNumberOfTimePeriods |
TPPG * | fPPG |
TScalerData * | fScalerData |
std::map< UInt_t, std::map< ULong64_t, std::vector< UInt_t > > > | fScalerMap |
std::map< UInt_t, ULong64_t > | fTimePeriod |
std::map< ULong64_t, int > | fTotalNumberOfTimePeriods |
ULong64_t | fTotalTimePeriod |
TTree * | fTree |
|
explicit |
This constructor tries to find the "ScalerTree" and uses it (if requested) to load the scaler data into the map.
[in] | loadIntoMap | Flag telling TScaler to load all scaler data into fScalerMap. |
Definition at line 46 of file TScaler.cxx.
References ReadTree().
|
explicit |
Definition at line 54 of file TScaler.cxx.
References ReadTree().
TScaler::TScaler | ( | const TScaler & | rhs | ) |
Definition at line 74 of file TScaler.cxx.
References Copy().
|
defaultnoexcept |
TScaler::~TScaler | ( | ) |
Definition at line 79 of file TScaler.cxx.
References fEntries, fHist, fHistRange, fNumberOfTimePeriods, fPPG, fScalerData, fScalerMap, fTimePeriod, fTotalNumberOfTimePeriods, fTotalTimePeriod, and fTree.
|
override |
Definition at line 228 of file TScaler.cxx.
References fEntries, fHist, fHistRange, fNumberOfTimePeriods, fPPG, fScalerData, fScalerMap, fTimePeriod, fTotalNumberOfTimePeriods, fTotalTimePeriod, and fTree.
Referenced by Copy(), and MakeSpectra().
|
override |
Definition at line 95 of file TScaler.cxx.
References Clear(), fEntries, fNumberOfTimePeriods, fScalerData, fScalerMap, fTimePeriod, fTotalNumberOfTimePeriods, fTotalTimePeriod, and fTree.
Referenced by TScaler().
TH1D * TScaler::Draw | ( | UInt_t | address, |
size_t | index = 0, | ||
Option_t * | option = "" ) |
Draw scaler differences (i.e. current scaler minus last scaler) vs. time in cycle. Passing "redraw" as option forces re-drawing of the histogram (e.g. for a different index).
Definition at line 255 of file TScaler.cxx.
References fEntries, fHist, fPPG, fScalerData, fTree, TSingleton< TPPG >::Get(), TScalerData::GetAddress(), TPPG::GetCycleLength(), GetScaler(), TScalerData::GetScaler(), TPPG::GetTimeInCycle(), GetTimePeriod(), and TScalerData::GetTimeStamp().
TH1D * TScaler::Draw | ( | UInt_t | lowAddress, |
UInt_t | highAddress, | ||
size_t | index = 0, | ||
Option_t * | option = "" ) |
Draw scaler differences (i.e. current scaler minus last scaler) vs. time in cycle. Passing "redraw" as option forces re-drawing of the histogram (e.g. for a different index). The "single" options creates spectra for all single addresses in the given range (instead of one accumulative spectrum).
Definition at line 322 of file TScaler.cxx.
References fEntries, fHist, fHistRange, fPPG, fScalerData, fTree, TSingleton< TPPG >::Get(), TScalerData::GetAddress(), TPPG::GetCycleLength(), TScalerData::GetScaler(), TPPG::GetTimeInCycle(), GetTimePeriod(), and TScalerData::GetTimeStamp().
TH1D * TScaler::DrawRawTimes | ( | UInt_t | address, |
Double_t | lowtime, | ||
Double_t | hightime, | ||
size_t | index = 0, | ||
Option_t * | option = "" ) |
Draw scaler differences (i.e. current scaler minus last scaler) vs. time.
Definition at line 450 of file TScaler.cxx.
References fEntries, fScalerData, fTree, TScalerData::GetAddress(), GetScaler(), TScalerData::GetScaler(), GetTimePeriod(), and TScalerData::GetTimeStamp().
|
inlinestatic |
This function returns the ratio of the two scalers nominator and denominator for the last entry with a matching address for a given tree. If no tree is provided it tries to get the "RateScaler" tree from the current directory.
Definition at line 134 of file TScaler.h.
References TScalerData::GetAddress(), TScalerData::GetScaler(), and hex().
|
inline |
Definition at line 124 of file TScaler.h.
References fNumberOfTimePeriods.
std::vector< UInt_t > TScaler::GetScaler | ( | UInt_t | address, |
ULong64_t | time ) const |
Returns the vector of scaler values for address "address" at the time "time". If the time is after the last entry, we return the last entry, otherwise we return the following entry (or the current entry if the time is an exact match).
Definition at line 108 of file TScaler.cxx.
References fEntries, fScalerData, fScalerMap, fTree, TScalerData::GetAddress(), TScalerData::GetScaler(), and TScalerData::GetTimeStamp().
Referenced by Draw(), DrawRawTimes(), and GetScaler().
UInt_t TScaler::GetScaler | ( | UInt_t | address, |
ULong64_t | time, | ||
size_t | index ) const |
Definition at line 154 of file TScaler.cxx.
References GetScaler().
UInt_t TScaler::GetScalerDifference | ( | UInt_t | address, |
ULong64_t | time, | ||
size_t | index ) const |
Returns the difference between "index"th scaler value for address "address" at the time "time" and the previous scaler value. If the time is after our last entry, we return the last entry divided by the number of entries, if this is before the first scaler, we just return the first scaler, otherwise we return the scaler minus the previous scaler.
Definition at line 163 of file TScaler.cxx.
References fEntries, fScalerData, fScalerMap, fTree, TScalerData::GetAddress(), TScalerData::GetScaler(), and TScalerData::GetTimeStamp().
ULong64_t TScaler::GetTimePeriod | ( | UInt_t | address | ) |
Get time period of scaler readouts for address "address" by calculating all time differences and choosing the one that occurs most often. Returns 0 if the address doesn't exist in the map.
Definition at line 491 of file TScaler.cxx.
References fEntries, fNumberOfTimePeriods, fScalerData, fTimePeriod, fTree, TScalerData::GetAddress(), and TScalerData::GetTimeStamp().
Referenced by Draw(), Draw(), and DrawRawTimes().
|
inline |
Definition at line 123 of file TScaler.h.
References fTimePeriod.
void TScaler::ListHistograms | ( | ) |
Definition at line 523 of file TScaler.cxx.
References fHist, fHistRange, and hex().
|
inline |
|
private |
Definition at line 60 of file TScaler.cxx.
References fEntries, fScalerData, fScalerMap, fTree, TScalerData::GetAddress(), TScalerData::GetScaler(), and TScalerData::GetTimeStamp().
|
private |
Definition at line 167 of file TScaler.h.
Referenced by ~TScaler(), Clear(), Copy(), Draw(), Draw(), DrawRawTimes(), GetScaler(), GetScalerDifference(), GetTimePeriod(), NumberOfScalerReadouts(), and ReadTree().
|
private |
! map to store histograms that have already been drawn
Definition at line 174 of file TScaler.h.
Referenced by ~TScaler(), Clear(), Draw(), Draw(), and ListHistograms().
|
private |
! map to store histograms for address-ranges
Definition at line 175 of file TScaler.h.
Referenced by ~TScaler(), Clear(), Draw(), and ListHistograms().
|
private |
!
Definition at line 170 of file TScaler.h.
Referenced by ~TScaler(), Clear(), Copy(), GetNumberOfTimePeriods(), and GetTimePeriod().
|
private |
|
private |
Definition at line 166 of file TScaler.h.
Referenced by ~TScaler(), Clear(), Copy(), Draw(), Draw(), DrawRawTimes(), GetScaler(), GetScalerDifference(), GetTimePeriod(), and ReadTree().
|
private |
! an address-map of timestamp mapped scaler values
Definition at line 168 of file TScaler.h.
Referenced by ~TScaler(), Clear(), Copy(), GetScaler(), GetScalerDifference(), and ReadTree().
|
private |
! a map between addresses and time differences (used to calculate the time period)
Definition at line 169 of file TScaler.h.
Referenced by ~TScaler(), Clear(), Copy(), GetTimePeriod(), and GetTimePeriodMap().
|
private |
|
private |
|
private |
Definition at line 165 of file TScaler.h.
Referenced by ~TScaler(), Clear(), Copy(), Draw(), Draw(), DrawRawTimes(), GetScaler(), GetScalerDifference(), GetTimePeriod(), and ReadTree().