1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
|
#include <string>
#include <iostream>
#include <vector>
#include <array>
#include <cstdint>
#include "FWCore/Framework/interface/stream/EDAnalyzer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
#include "DataFormats/HcalRecHit/interface/CaloRecHitAuxSetter.h"
namespace {
std::ostream& operator<<(std::ostream& s, const HFQIE10Info& i) {
s << i.id() << ": " << i.energy() << " GeV"
<< ", t= " << i.timeRising() << " to " << i.timeFalling() << " ns";
return s;
}
std::ostream& operator<<(std::ostream& s, const HFPreRecHit& hit) {
s << "{ ";
const HFQIE10Info* i = hit.getHFQIE10Info(0);
if (i) {
s << *i;
}
s << " }, ";
s << "{ ";
i = hit.getHFQIE10Info(1);
if (i) {
s << *i;
}
s << " }";
return s;
}
template <std::size_t N>
void printBits(std::ostream& s, const std::array<uint32_t, N>& allbits, const std::vector<int>& bits) {
const int maxbit = N * 32;
const unsigned len = bits.size();
for (unsigned i = 0; i < len; ++i) {
const int bitnum = bits[i];
if (bitnum >= 0 && bitnum < maxbit) {
const unsigned ibit = bitnum % 32;
const bool bit = (allbits[bitnum / 32] & (1U << ibit)) >> ibit;
s << bit;
} else
s << '-';
}
}
void printRecHitAuxInfo(std::ostream& s, const HFPreRecHit& i, const std::vector<int>& bits, bool) {}
void printRecHitAuxInfo(std::ostream& s, const HBHERecHit& i, const std::vector<int>& bits, const bool plan1) {
if (plan1 && i.isMerged()) {
// This is a "Plan 1" combined rechit
std::vector<HcalDetId> ids;
i.getMergedIds(&ids);
const unsigned n = ids.size();
s << "; merged " << n << ": ";
for (unsigned j = 0; j < n; ++j) {
if (j)
s << ", ";
s << ids[j];
}
}
if (!bits.empty()) {
std::array<uint32_t, 4> allbits;
allbits[0] = i.flags();
allbits[1] = i.aux();
allbits[2] = i.auxHBHE();
allbits[3] = i.auxPhase1();
s << "; bits: ";
printBits(s, allbits, bits);
}
// Dump TDC data
s << "; tdc:";
const uint32_t auxTDC = i.auxTDC();
if (auxTDC) {
const unsigned six_bits_mask = 0x3f;
for (unsigned ts = 0; ts < 5; ++ts)
s << ' ' << CaloRecHitAuxSetter::getField(auxTDC, six_bits_mask, ts * 6);
} else {
s << " none";
}
}
void printRecHitAuxInfo(std::ostream& s, const HFRecHit& i, const std::vector<int>& bits, bool) {
if (!bits.empty()) {
std::array<uint32_t, 3> allbits;
allbits[0] = i.flags();
allbits[1] = i.aux();
allbits[2] = i.getAuxHF();
s << "; bits: ";
printBits(s, allbits, bits);
}
}
} // namespace
using namespace std;
class HcalRecHitDump : public edm::stream::EDAnalyzer<> {
public:
explicit HcalRecHitDump(edm::ParameterSet const& conf);
virtual void analyze(edm::Event const& e, edm::EventSetup const& c) override;
private:
string hbhePrefix_;
string hfPrefix_;
string hfprePrefix_;
std::vector<int> bits_;
bool printPlan1Info_;
edm::EDGetTokenT<HBHERecHitCollection> tok_hbhe_;
edm::EDGetTokenT<HFRecHitCollection> tok_hf_;
edm::EDGetTokenT<HFPreRecHitCollection> tok_prehf_;
unsigned long long counter_;
template <class Collection, class Token>
void analyzeT(edm::Event const& e,
const Token& tok,
const char* name,
const string& prefix,
const bool printPlan1Info = false) const {
cout << prefix << " rechit dump " << counter_ << endl;
edm::Handle<Collection> coll;
bool found = false;
try {
e.getByToken(tok, coll);
found = true;
} catch (...) {
cout << prefix << " Error: no " << name << " rechit data" << endl;
}
if (found) {
for (typename Collection::const_iterator j = coll->begin(); j != coll->end(); ++j) {
cout << prefix << *j;
printRecHitAuxInfo(cout, *j, bits_, printPlan1Info);
cout << endl;
}
}
}
};
HcalRecHitDump::HcalRecHitDump(edm::ParameterSet const& conf)
: hbhePrefix_(conf.getUntrackedParameter<string>("hbhePrefix", "")),
hfPrefix_(conf.getUntrackedParameter<string>("hfPrefix", "")),
hfprePrefix_(conf.getUntrackedParameter<string>("hfprePrefix", "")),
bits_(conf.getUntrackedParameter<std::vector<int> >("bits")),
printPlan1Info_(conf.getUntrackedParameter<bool>("printPlan1Info", false)),
counter_(0) {
if (!hbhePrefix_.empty())
tok_hbhe_ = consumes<HBHERecHitCollection>(conf.getParameter<edm::InputTag>("tagHBHE"));
if (!hfPrefix_.empty())
tok_hf_ = consumes<HFRecHitCollection>(conf.getParameter<edm::InputTag>("tagHF"));
if (!hfprePrefix_.empty())
tok_prehf_ = consumes<HFPreRecHitCollection>(conf.getParameter<edm::InputTag>("tagPreHF"));
}
void HcalRecHitDump::analyze(edm::Event const& e, edm::EventSetup const& c) {
if (!hbhePrefix_.empty())
analyzeT<HBHERecHitCollection>(e, tok_hbhe_, "HBHE", hbhePrefix_, printPlan1Info_);
if (!hfPrefix_.empty())
analyzeT<HFRecHitCollection>(e, tok_hf_, "HF", hfPrefix_);
if (!hfprePrefix_.empty())
analyzeT<HFPreRecHitCollection>(e, tok_prehf_, "PreHF", hfprePrefix_);
++counter_;
}
DEFINE_FWK_MODULE(HcalRecHitDump);
|