GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TParsingDiagnostics.cxx
Go to the documentation of this file.
2
3#include <fstream>
4
5#include "TH1.h"
6
7#include "TChannel.h"
8
10
11TParsingDiagnosticsData::TParsingDiagnosticsData(const std::shared_ptr<const TFragment>& frag)
12 : fMinChannelId(frag->GetChannelId()), fMaxChannelId(frag->GetChannelId()), fNumberOfHits(1),
13 fDeadTime(frag->GetDeadTime()), fMinTimeStamp(frag->GetTimeStampNs()), fMaxTimeStamp(frag->GetTimeStampNs())
14{
15}
16
17void TParsingDiagnosticsData::Update(const std::shared_ptr<const TFragment>& frag)
18{
19 UInt_t channelId = frag->GetChannelId();
20 auto timeStamp = frag->GetTimeStampNs();
21
22 // update minimum and maximum channel id if necessary
23 fMinChannelId = std::min(fMinChannelId, channelId);
24 fMaxChannelId = std::max(fMaxChannelId, channelId);
25
27
28 // increment the dead time and set per channel min/max timestamps
29 fDeadTime += frag->GetDeadTime();
30 fMinTimeStamp = std::min(fMinTimeStamp, static_cast<int64_t>(timeStamp));
31 fMaxTimeStamp = std::max(fMaxTimeStamp, static_cast<int64_t>(timeStamp));
32}
33
34void TParsingDiagnosticsData::Print(UInt_t address) const
35{
36 std::cout << "channel " << hex(address, 4) << ": " << fDeadTime / 100000 << " ms deadtime out of ";
38 std::cout << "empty channel" << std::endl;
39 } else {
40 std::cout << std::setw(12) << (fMaxTimeStamp - fMinTimeStamp) / 100000
41 << "ms = " << 100. * static_cast<double>(fDeadTime) / static_cast<double>(fMaxTimeStamp - fMinTimeStamp) << " %"
42 << std::endl;
43 }
44}
45
50
55
60
71
73{
74 delete fIdHist;
75 fIdHist = nullptr;
81 fMinNetworkPacketNumber = 0x7fffffff; // just a large number
84}
85
86void TParsingDiagnostics::Print(Option_t*) const
87{
88 std::cout << "Total run time of this (sub-)run is " << fMaxDaqTimeStamp - fMinDaqTimeStamp << " s" << std::endl
89 << "PPG cycle is " << fPPGCycleLength / 100000 << " ms long." << std::endl
90 << "Found " << fNumberOfNetworkPackets << " network packets in range " << fMinNetworkPacketNumber << " - "
91 << fMaxNetworkPacketNumber << " => "
92 << 100. * static_cast<double>(fNumberOfNetworkPackets) / static_cast<double>(fMaxNetworkPacketNumber - fMinNetworkPacketNumber + 1)
93 << " % packet survival." << std::endl;
94 // loop over number of good fragments per detector type
95 for(const auto& fNumberOfGoodFragment : fNumberOfGoodFragments) {
96 std::cout << "detector type " << std::setw(2) << fNumberOfGoodFragment.first << ": " << std::setw(12)
97 << fNumberOfGoodFragment.second << " good, ";
98 // check if we have corresponding bad fragment for this detector type
99 if(fNumberOfBadFragments.find(fNumberOfGoodFragment.first) == fNumberOfBadFragments.end()) {
100 std::cout << " no";
101 } else {
102 std::cout << std::setw(12) << fNumberOfBadFragments.at(fNumberOfGoodFragment.first) << " ("
103 << 100. * static_cast<double>(fNumberOfBadFragments.at(fNumberOfGoodFragment.first)) / static_cast<double>(fNumberOfGoodFragment.second)
104 << " %)";
105 }
106 std::cout << " bad fragments." << std::endl;
107 }
108 for(const auto& iter : fChannelAddressData) {
109 iter.second.Print(iter.first);
110 }
111}
112
113void TParsingDiagnostics::GoodFragment(const std::shared_ptr<const TFragment>& frag)
114{
115 /// increment the counter of good fragments for this detector type and check if any trigger ids have been lost
116 fNumberOfGoodFragments[frag->GetDetectorType()]++;
117
118 UInt_t channelAddress = frag->GetAddress();
119 if(fChannelAddressData.find(channelAddress) == fChannelAddressData.end()) {
120 fChannelAddressData[channelAddress] = TParsingDiagnosticsData(frag);
121 } else {
122 fChannelAddressData[channelAddress].Update(frag);
123 }
124
125 // check if this is a new minimum/maximum network packet id
126 if(frag->GetNetworkPacketNumber() > 0) {
128 fMinNetworkPacketNumber = std::min(fMinNetworkPacketNumber, frag->GetNetworkPacketNumber());
129 fMaxNetworkPacketNumber = std::max(fMaxNetworkPacketNumber, frag->GetNetworkPacketNumber());
130 }
131
132 if(fMinDaqTimeStamp == 0 || frag->GetDaqTimeStamp() < fMinDaqTimeStamp) {
133 fMinDaqTimeStamp = frag->GetDaqTimeStamp();
134 }
135 if(fMaxDaqTimeStamp == 0 || frag->GetDaqTimeStamp() > fMaxDaqTimeStamp) {
136 fMaxDaqTimeStamp = frag->GetDaqTimeStamp();
137 }
138}
139
141{
142 /// store different TPPG diagnostics like cycle length, length of each state, offset, how often each state was found
143 if(ppg == nullptr) {
144 return;
145 }
147}
148
149void TParsingDiagnostics::Draw(Option_t* opt)
150{
151 UInt_t minChannel = fChannelAddressData.begin()->first;
152 UInt_t maxChannel = fChannelAddressData.begin()->first;
153
154 for(const auto& iter : fChannelAddressData) {
155 minChannel = std::min(minChannel, iter.first);
156 maxChannel = std::max(maxChannel, iter.first);
157 }
158
159 // check that the histogram (if it already exists) has the right number of bins
160 if(fIdHist != nullptr && fIdHist->GetNbinsX() != static_cast<Int_t>(maxChannel - minChannel + 1)) {
161 delete fIdHist;
162 fIdHist = nullptr;
163 }
164 if(fIdHist == nullptr) {
165 fIdHist = new TH1F("IdHist", "Event survival;channel number;survival rate [%]", maxChannel - minChannel + 1,
166 minChannel, maxChannel + 1);
167 } else {
168 // the histogram already had the right number of bins, but to be save we set the range
169 fIdHist->SetAxisRange(minChannel, maxChannel + 1);
170 }
171
172 for(const auto& iter : fChannelAddressData) {
173 if(iter.second.MinChannelId() != 0 || iter.second.MinChannelId() != iter.second.MaxChannelId()) {
174 fIdHist->SetBinContent(fIdHist->GetXaxis()->FindBin(iter.first),
175 100. * static_cast<double>(iter.second.NumberOfHits()) / static_cast<double>(iter.second.MaxChannelId() - iter.second.MinChannelId() + 1));
176 }
177 }
178
179 fIdHist->Draw(opt);
180}
181
182void TParsingDiagnostics::WriteToFile(const char* fileName) const
183{
184 std::ofstream statsOut(fileName);
185 statsOut << std::endl
186 << "Run time to the nearest second = " << fMaxDaqTimeStamp - fMinDaqTimeStamp << std::endl
187 << std::endl;
188
189 statsOut << "Good fragments:";
190 for(const auto& iter : fNumberOfGoodFragments) {
191 statsOut << " " << iter.second << " of type " << iter.first;
192 }
193 statsOut << std::endl;
194
195 statsOut << "Bad fragments:";
196 for(const auto& iter : fNumberOfBadFragments) {
197 statsOut << " " << iter.second << " of type " << iter.first;
198 }
199 statsOut << std::endl;
200
201 for(const auto& iter : fChannelAddressData) {
202 TChannel* chan = TChannel::GetChannel(iter.first, false);
203 if(chan == nullptr) {
204 continue;
205 }
206 statsOut << hex(iter.first, 4) << ":\t" << chan->GetName()
207 << "\tdead time: " << static_cast<float>(iter.second.DeadTime()) / 1e9 << " seconds." << std::endl;
208 }
209 statsOut << std::endl;
210
211 statsOut.close();
212}
std::string hex(T val, int width=-1)
Definition Globals.h:129
static TChannel * GetChannel(unsigned int temp_address, bool warn=false)
Definition TChannel.cxx:476
Definition TPPG.h:130
ULong64_t GetCycleLength()
Definition TPPG.cxx:631
int64_t fMinTimeStamp
minimum timestamp per channel address
UInt_t fMaxChannelId
maximum channel id per channel address
void Update(const std::shared_ptr< const TFragment > &frag)
int64_t fDeadTime
deadtime per channel address
UInt_t fMinChannelId
minimum channel id per channel address
void Print(UInt_t address) const
Long_t fNumberOfHits
number of hits per channel address
int64_t fMaxTimeStamp
maximum timestamp per channel address
void Copy(TObject &) const override
time_t fMaxDaqTimeStamp
maximum daq timestamp
Int_t fMinNetworkPacketNumber
minimum network packet id
time_t fMinDaqTimeStamp
minimum daq timestamp
void Print(Option_t *opt="") const override
void Clear(Option_t *opt="all") override
void Draw(Option_t *opt="") override
TH1F * fIdHist
histogram of event survival
void WriteToFile(const char *) const
std::unordered_map< UInt_t, TParsingDiagnosticsData > fChannelAddressData
unordered_map of data per channel address
Int_t fMaxNetworkPacketNumber
maximum network packet id
std::unordered_map< Short_t, Long_t > fNumberOfBadFragments
unordered_map of number of bad fragments per detector type
void GoodFragment(const std::shared_ptr< const TFragment > &)
std::unordered_map< Short_t, Long_t > fNumberOfGoodFragments
unordered_map of number of good fragments per detector type