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 "TChannel.h"
6
8
9TParsingDiagnosticsData::TParsingDiagnosticsData(const std::shared_ptr<const TFragment>& frag)
10 : fMinChannelId(frag->GetChannelId()), fMaxChannelId(frag->GetChannelId()), fNumberOfHits(1),
11 fDeadTime(frag->GetDeadTime()), fMinTimeStamp(frag->GetTimeStampNs()), fMaxTimeStamp(frag->GetTimeStampNs())
12{
13}
14
15void TParsingDiagnosticsData::Update(const std::shared_ptr<const TFragment>& frag)
16{
17 UInt_t channelId = frag->GetChannelId();
18 auto timeStamp = frag->GetTimeStampNs();
19
20 // update minimum and maximum channel id if necessary
21 if(channelId < fMinChannelId) {
22 fMinChannelId = channelId;
23 }
24 if(channelId > fMaxChannelId) {
25 fMaxChannelId = channelId;
26 }
27
29
30 // increment the dead time and set per channel min/max timestamps
31 fDeadTime += frag->GetDeadTime();
32 if(timeStamp < fMinTimeStamp) {
33 fMinTimeStamp = timeStamp;
34 }
35 if(timeStamp > fMaxTimeStamp) {
36 fMaxTimeStamp = timeStamp;
37 }
38}
39
40void TParsingDiagnosticsData::Print(UInt_t address) const
41{
42 std::cout << "channel " << hex(address, 4) << ": " << fDeadTime / 100000 << " ms deadtime out of ";
44 std::cout << "empty channel" << std::endl;
45 } else {
46 std::cout << std::setw(12) << (fMaxTimeStamp - fMinTimeStamp) / 100000
47 << "ms = " << 100. * static_cast<double>(fDeadTime) / static_cast<double>(fMaxTimeStamp - fMinTimeStamp) << " %"
48 << std::endl;
49 }
50}
51
56
61
66
77
79{
80 delete fIdHist;
81 fIdHist = nullptr;
87 fMinNetworkPacketNumber = 0x7fffffff; // just a large number
90}
91
92void TParsingDiagnostics::Print(Option_t*) const
93{
94 std::cout << "Total run time of this (sub-)run is " << fMaxDaqTimeStamp - fMinDaqTimeStamp << " s" << std::endl
95 << "PPG cycle is " << fPPGCycleLength / 100000 << " ms long." << std::endl
96 << "Found " << fNumberOfNetworkPackets << " network packets in range " << fMinNetworkPacketNumber << " - "
97 << fMaxNetworkPacketNumber << " => "
98 << 100. * static_cast<double>(fNumberOfNetworkPackets) / static_cast<double>(fMaxNetworkPacketNumber - fMinNetworkPacketNumber + 1)
99 << " % packet survival." << std::endl;
100 // loop over number of good fragments per detector type
101 for(const auto& fNumberOfGoodFragment : fNumberOfGoodFragments) {
102 std::cout << "detector type " << std::setw(2) << fNumberOfGoodFragment.first << ": " << std::setw(12)
103 << fNumberOfGoodFragment.second << " good, ";
104 // check if we have corresponding bad fragment for this detector type
105 if(fNumberOfBadFragments.find(fNumberOfGoodFragment.first) == fNumberOfBadFragments.end()) {
106 std::cout << " no";
107 } else {
108 std::cout << std::setw(12) << fNumberOfBadFragments.at(fNumberOfGoodFragment.first) << " ("
109 << 100. * static_cast<double>(fNumberOfBadFragments.at(fNumberOfGoodFragment.first)) / static_cast<double>(fNumberOfGoodFragment.second)
110 << " %)";
111 }
112 std::cout << " bad fragments." << std::endl;
113 }
114 for(const auto& iter : fChannelAddressData) {
115 iter.second.Print(iter.first);
116 }
117}
118
119void TParsingDiagnostics::GoodFragment(const std::shared_ptr<const TFragment>& frag)
120{
121 /// increment the counter of good fragments for this detector type and check if any trigger ids have been lost
122 fNumberOfGoodFragments[frag->GetDetectorType()]++;
123
124 UInt_t channelAddress = frag->GetAddress();
125 if(fChannelAddressData.find(channelAddress) == fChannelAddressData.end()) {
126 fChannelAddressData[channelAddress] = TParsingDiagnosticsData(frag);
127 } else {
128 fChannelAddressData[channelAddress].Update(frag);
129 }
130
131 // check if this is a new minimum/maximum network packet id
132 if(frag->GetNetworkPacketNumber() > 0) {
134 if(frag->GetNetworkPacketNumber() < fMinNetworkPacketNumber) {
135 fMinNetworkPacketNumber = frag->GetNetworkPacketNumber();
136 }
137 if(frag->GetNetworkPacketNumber() > fMaxNetworkPacketNumber) {
138 fMaxNetworkPacketNumber = frag->GetNetworkPacketNumber();
139 }
140 }
141
142 if(fMinDaqTimeStamp == 0 || frag->GetDaqTimeStamp() < fMinDaqTimeStamp) {
143 fMinDaqTimeStamp = frag->GetDaqTimeStamp();
144 }
145 if(fMaxDaqTimeStamp == 0 || frag->GetDaqTimeStamp() > fMaxDaqTimeStamp) {
146 fMaxDaqTimeStamp = frag->GetDaqTimeStamp();
147 }
148}
149
151{
152 /// store different TPPG diagnostics like cycle length, length of each state, offset, how often each state was found
153 if(ppg == nullptr) {
154 return;
155 }
157}
158
159void TParsingDiagnostics::Draw(Option_t* opt)
160{
161 UInt_t minChannel = fChannelAddressData.begin()->first;
162 UInt_t maxChannel = fChannelAddressData.begin()->first;
163
164 for(const auto& iter : fChannelAddressData) {
165 if(iter.first < minChannel) { minChannel = iter.first; }
166 if(iter.first > maxChannel) { maxChannel = iter.first; }
167 }
168
169 // check that the histogram (if it already exists) has the right number of bins
170 if(fIdHist != nullptr && fIdHist->GetNbinsX() != static_cast<Int_t>(maxChannel - minChannel + 1)) {
171 delete fIdHist;
172 fIdHist = nullptr;
173 }
174 if(fIdHist == nullptr) {
175 fIdHist = new TH1F("IdHist", "Event survival;channel number;survival rate [%]", maxChannel - minChannel + 1,
176 minChannel, maxChannel + 1);
177 } else {
178 // the histogram already had the right number of bins, but to be save we set the range
179 fIdHist->SetAxisRange(minChannel, maxChannel + 1);
180 }
181
182 for(const auto& iter : fChannelAddressData) {
183 if(iter.second.MinChannelId() != 0 || iter.second.MinChannelId() != iter.second.MaxChannelId()) {
184 fIdHist->SetBinContent(fIdHist->GetXaxis()->FindBin(iter.first),
185 100. * static_cast<double>(iter.second.NumberOfHits()) / static_cast<double>(iter.second.MaxChannelId() - iter.second.MinChannelId() + 1));
186 }
187 }
188
189 fIdHist->Draw(opt);
190}
191
192void TParsingDiagnostics::WriteToFile(const char* fileName) const
193{
194 std::ofstream statsOut(fileName);
195 statsOut << std::endl
196 << "Run time to the nearest second = " << fMaxDaqTimeStamp - fMinDaqTimeStamp << std::endl
197 << std::endl;
198
199 statsOut << "Good fragments:";
200 for(const auto& iter : fNumberOfGoodFragments) {
201 statsOut << " " << iter.second << " of type " << iter.first;
202 }
203 statsOut << std::endl;
204
205 statsOut << "Bad fragments:";
206 for(const auto& iter : fNumberOfBadFragments) {
207 statsOut << " " << iter.second << " of type " << iter.first;
208 }
209 statsOut << std::endl;
210
211 for(const auto& iter : fChannelAddressData) {
212 TChannel* chan = TChannel::GetChannel(iter.first, false);
213 if(chan == nullptr) {
214 continue;
215 }
216 statsOut << hex(iter.first, 4) << ":\t" << chan->GetName()
217 << "\tdead time: " << static_cast<float>(iter.second.DeadTime()) / 1e9 << " seconds." << std::endl;
218 }
219 statsOut << std::endl;
220
221 statsOut.close();
222}
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:459
Definition TPPG.h:132
ULong64_t GetCycleLength()
Definition TPPG.cxx:629
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