GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TEpicsFrag.cxx
Go to the documentation of this file.
1#include "TEpicsFrag.h"
2
3#include <iostream>
4#include <iomanip>
5#include <stdexcept>
6#include <limits>
7
8#include <ctime>
9
10#include "TChannel.h"
11#include "TRunInfo.h"
12
13///////////////////////////////////////////////////////////////
14///
15/// \class TEpicsFrag
16///
17/// This Class should contain all the information found in
18/// NOT typeid 1 midas events. aka Epics (scaler) Events.
19///
20///
21///////////////////////////////////////////////////////////////
22
23std::vector<std::string> TEpicsFrag::fNameList;
24std::map<Long64_t, TEpicsFrag> TEpicsFrag::fScalerMap;
25Long64_t TEpicsFrag::fSmallestTime = std::numeric_limits<Long64_t>::max();
26
27void TEpicsFrag::Clear(Option_t*)
28{
29 // Clears the TEpicsFrag.
30 fDaqTimeStamp = 0;
31 fDaqId = -1;
32
33 fName.clear();
34 fData.clear();
35}
36
37void TEpicsFrag::Print(Option_t*) const
38{
39 // Prints the TEpicsFrag. This includes Midas information as well the data
40 // kep inside of the scaler.
41 size_t largest = fData.size();
42 std::cout << "------ EPICS " << largest << " Varibles Found ------" << std::endl;
43
44 // TODO maybe we can change this to std::array?
45 std::array<char, 20> buff;
46 ctime(&fDaqTimeStamp);
47 struct tm* timeInfo = localtime(&fDaqTimeStamp);
48 strftime(buff.data(), buff.size(), "%b %d %H:%M:%S", timeInfo);
49
50 std::cout << " DaqTimeStamp: " << buff.data() << std::endl;
51 std::cout << " DaqId: " << fDaqId << std::endl;
52 for(size_t i = 0; i < largest; i++) {
53 std::cout << std::setw(3) << i << ": ";
54 std::cout << std::setw(30) << fName.at(i) << " --- ";
55 std::cout << fData.at(i);
56 std::cout << std::endl;
57 }
58}
59
60void TEpicsFrag::AddEpicsVariable(const char* name)
61{
62 fNameList.emplace_back(name);
63}
64
65std::string TEpicsFrag::GetEpicsVariableName(const int& index)
66{
67 try {
68 return fNameList.at(index);
69 } catch(const std::out_of_range& oor) {
70 std::cout << DRED << "Could not find variable at position " << index << ", returning nothing" << std::endl;
71 return "";
72 }
73}
74
76{
77 int idx = 0;
78 for(const auto& name : fNameList) {
79 std::cout << idx++ << ": " << name << std::endl;
80 }
81}
82
83void TEpicsFrag::SetEpicsNameList(const std::vector<std::string>& names)
84{
85 fNameList.clear();
86 for(const auto& name : names) {
87 fNameList.push_back(name);
88 }
89}
90
92{
93 if(tree == nullptr) {
94 std::cout << DRED << "Could not build map from tree" << RESET_COLOR << std::endl;
95 }
96 // Loop through the tree and insert the scalers into the map
97 fScalerMap.clear();
98 TEpicsFrag* my_frag = nullptr;
99 if(tree->SetBranchAddress("TEpicsFrag", &my_frag) == 0) {
100 for(int i = 0; i < tree->GetEntries(); ++i) {
101 tree->GetEntry(i);
102 if((static_cast<Long64_t>(my_frag->fDaqTimeStamp) - static_cast<Long64_t>(TRunInfo::RunStart())) < fSmallestTime) {
103 fSmallestTime = static_cast<Long64_t>(my_frag->fDaqTimeStamp) - static_cast<Long64_t>(TRunInfo::RunStart());
104 }
105 fScalerMap[static_cast<Long64_t>(my_frag->fDaqTimeStamp) -
106 static_cast<Long64_t>(TRunInfo::RunStart())] = *my_frag;
107 }
108 } else {
109 std::cout << DRED << "Could not build map from tree" << RESET_COLOR << std::endl;
110 }
111}
112
114{
115 TTree* scaler_tree = static_cast<TTree*>(gDirectory->Get("EpicsTree"));
116 if(scaler_tree == nullptr) {
117 return;
118 }
119
120 BuildScalerMap(scaler_tree);
121}
122
124{
125 if(fScalerMap.empty()) {
127 if(fScalerMap.empty()) {
128 std::cout << DRED << "Could not build the epics map" << RESET_COLOR << std::endl;
129 return nullptr;
130 }
131 }
132 if(time < fSmallestTime) {
133 time = fSmallestTime;
134 }
135 return &((--(fScalerMap.upper_bound(time)))->second);
136}
137
139{
140 if(fScalerMap.empty()) {
142 if(fScalerMap.empty()) {
143 std::cout << DRED << "Could not build the epics map" << RESET_COLOR << std::endl;
144 return;
145 }
146 }
147 for(const auto& item : fScalerMap) {
148 std::cout << item.first << " " << item.second.fDaqTimeStamp << std::endl;
149 }
150}
#define DRED
Definition Globals.h:18
#define RESET_COLOR
Definition Globals.h:5
static void SetEpicsNameList(const std::vector< std::string > &names)
static std::map< Long64_t, TEpicsFrag > fScalerMap
Definition TEpicsFrag.h:80
void Print(Option_t *opt="") const override
!
static void BuildScalerMap()
static void PrintVariableNames()
void Clear(Option_t *opt="") override
!
static std::vector< std::string > fNameList
Definition TEpicsFrag.h:79
Int_t fDaqId
Definition TEpicsFrag.h:84
std::vector< std::string > fName
The name of the scaler.
Definition TEpicsFrag.h:87
time_t fDaqTimeStamp
Definition TEpicsFrag.h:83
static Long64_t fSmallestTime
Definition TEpicsFrag.h:81
static void PrintScalerMap()
std::vector< float > fData
The data in the scaler.
Definition TEpicsFrag.h:86
static std::string GetEpicsVariableName(const int &index)
static void AddEpicsVariable(const char *name)
static TEpicsFrag * GetScalerAtTime(Long64_t time)
static double RunStart()
Definition TRunInfo.h:227