GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TSortingDiagnostics.cxx
Go to the documentation of this file.
2
3#include <fstream>
4#include <string>
5
6#include "TChannel.h"
7#include "TGRSIOptions.h"
8
13
20
22{
26}
27
28void TSortingDiagnostics::OutOfTimeOrder(double newFragTime, double oldFragTime, int64_t newEntry)
29{
30 fFragmentsOutOfTimeOrder[oldFragTime] = std::make_pair(oldFragTime - newFragTime, newEntry);
31 // try and find a time before newFragTime
32 size_t entry = 0;
33 if(!fPreviousTimes.empty()) {
34 for(entry = fPreviousTimes.size() - 1; entry > 0; --entry) {
35 if(fPreviousTimes[entry] < newFragTime) {
36 break;
37 }
38 }
39 }
40 int64_t entryDiff = newEntry - (entry * TGRSIOptions::Get()->SortDepth());
41 if(entryDiff > fMaxEntryDiff) {
42 fMaxEntryDiff = entryDiff;
43 }
44}
45
46void TSortingDiagnostics::OutOfOrder(int64_t newFragTS, int64_t oldFragTS, int64_t newEntry)
47{
48 fFragmentsOutOfOrder[oldFragTS] = std::make_pair(oldFragTS - newFragTS, newEntry);
49 // try and find a timestamp before newFragTS
50 size_t entry = 0;
51 if(!fPreviousTimeStamps.empty()) {
52 for(entry = fPreviousTimeStamps.size() - 1; entry > 0; --entry) {
53 if(fPreviousTimeStamps[entry] < newFragTS) {
54 break;
55 }
56 }
57 }
58 int64_t entryDiff = newEntry - (entry * TGRSIOptions::Get()->SortDepth());
59 if(entryDiff > fMaxEntryDiff) {
60 fMaxEntryDiff = entryDiff;
61 }
62}
63
64void TSortingDiagnostics::MissingChannel(const UInt_t& address)
65{
66 if(fMissingChannels.find(address) != fMissingChannels.end()) {
67 ++(fMissingChannels[address]);
68 } else {
69 fMissingChannels[address] = 0;
70 }
71}
72
74{
77 } else {
79 std::cout << "Failed to find detector class " << channel->GetClassType() << " for channel:" << std::endl;
80 channel->Print();
81 }
82}
83
84void TSortingDiagnostics::Print(Option_t* opt) const
85{
86 TString option = opt;
87 option.ToUpper();
88 if(!fMissingChannels.empty()) {
89 std::cout << "Missing channels:" << std::endl;
90 for(auto iter : fMissingChannels) {
91 std::cout << hex(iter.first, 4) << ": " << iter.second << std::endl;
92 }
93 }
94 if(!fMissingDetectorClasses.empty()) {
95 std::cout << "Missing detector classes:" << std::endl;
96 for(auto iter : fMissingDetectorClasses) {
97 std::cout << (iter.first == nullptr ? "nullptr" : iter.first->GetName()) << ": " << iter.second << std::endl;
98 }
99 }
100 std::string color;
101 if(!fHitsRemoved.empty()) {
102 if(option.EqualTo("ERROR")) {
103 color = DRED;
104 }
105 std::cout << color << "Removed hits per detector class:" << RESET_COLOR << std::endl;
106 for(auto iter : fHitsRemoved) {
107 std::cout << iter.first->GetName() << ": " << iter.second.first << "/" << iter.second.second << " = " << 100. * static_cast<double>(iter.second.first) / static_cast<double>(iter.second.second) << "%" << std::endl;
108 }
109 } else {
110 if(option.EqualTo("ERROR")) {
111 color = DGREEN;
112 }
113 std::cout << color << "No hits were removed!" << RESET_COLOR << std::endl;
114 }
115 color = ""; // reset color string
116 if(fFragmentsOutOfOrder.empty() && fFragmentsOutOfTimeOrder.empty()) {
117 if(option.EqualTo("ERROR")) {
118 color = DGREEN;
119 }
120 std::cout << color << "No fragments out of order!" << RESET_COLOR << std::endl;
121 return;
122 }
123 if(option.EqualTo("ERROR")) {
124 color = DRED;
125 }
126 if(!fFragmentsOutOfOrder.empty()) {
127 std::cerr << color << NumberOfFragmentsOutOfOrder() << " fragments were out of order, maximum entry difference was "
128 << fMaxEntryDiff << "!" << std::endl
129 << "Please consider increasing the sort depth with --sort-depth=" << fMaxEntryDiff << RESET_COLOR
130 << std::endl;
131 }
132 if(!fFragmentsOutOfTimeOrder.empty()) {
133 std::cerr << color << NumberOfFragmentsOutOfTimeOrder() << " fragments were out of order, maximum entry difference was "
134 << fMaxEntryDiff << "!" << std::endl
135 << "Please consider increasing the sort depth with --sort-depth=" << fMaxEntryDiff << RESET_COLOR
136 << std::endl;
137 }
138}
139
141{
142}
143
144void TSortingDiagnostics::WriteToFile(const char* fileName) const
145{
146 std::ofstream statsOut(fileName);
147 statsOut << std::endl
148 << "Number of fragments out of order = " << NumberOfFragmentsOutOfOrder() << std::endl
149 << "Number of fragments out of time order = " << NumberOfFragmentsOutOfTimeOrder() << std::endl
150 << "Maximum entry difference = " << fMaxEntryDiff << std::endl
151 << std::endl;
152}
153
154void TSortingDiagnostics::RemovedHits(TClass* detClass, int64_t removed, int64_t total)
155{
156 if(fHitsRemoved.find(detClass) == fHitsRemoved.end()) {
157 fHitsRemoved[detClass] = std::make_pair(removed, total);
158 } else {
159 fHitsRemoved[detClass].first += removed;
160 fHitsRemoved[detClass].second += total;
161 }
162}
#define DRED
Definition Globals.h:18
#define DGREEN
Definition Globals.h:17
std::string hex(T val, int width=-1)
Definition Globals.h:129
#define RESET_COLOR
Definition Globals.h:5
TClass * GetClassType() const
void Print(Option_t *opt="") const override
Definition TChannel.cxx:833
static TGRSIOptions * Get(int argc=0, char **argv=nullptr)
Do not use!
int SortDepth() const
void OutOfTimeOrder(double newFragTime, double oldFragTime, int64_t newEntry)
void Copy(TObject &) const override
std::unordered_map< UInt_t, int64_t > fMissingChannels
counts of missing channels
void RemovedHits(TClass *detClass, int64_t removed, int64_t total)
void Draw(Option_t *opt="") override
void Clear(Option_t *opt="all") override
std::vector< double > fPreviousTimes
times of previous fragments, saved every 'BuildWindow' entries
std::unordered_map< TClass *, int64_t > fMissingDetectorClasses
counts of missing detector classes
std::unordered_map< int64_t, std::pair< int64_t, int64_t > > fFragmentsOutOfOrder
void MissingChannel(const UInt_t &address)
std::unordered_map< double, std::pair< double, double > > fFragmentsOutOfTimeOrder
void OutOfOrder(int64_t newFragTS, int64_t oldFragTS, int64_t newEntry)
void Print(Option_t *opt="") const override
void WriteToFile(const char *) const
size_t NumberOfFragmentsOutOfOrder() const
std::unordered_map< TClass *, std::pair< int64_t, int64_t > > fHitsRemoved
removed hits and total hits per detector class
size_t NumberOfFragmentsOutOfTimeOrder() const
void AddDetectorClass(TChannel *)
std::vector< Long_t > fPreviousTimeStamps
timestamps of previous fragments, saved every 'BuildWindow' entries