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
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
|
/*
* See header file for a description of this class.
*
* \author G. Mila - INFN Torino
*/
#include "DQMOffline/Muon/interface/DTSegmentsTask.h"
// Framework
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ServiceRegistry/interface/Service.h"
#include "DQMServices/Core/interface/DQMStore.h"
//Geometry
#include "DataFormats/GeometryVector/interface/Pi.h"
#include "DataFormats/MuonDetId/interface/DTChamberId.h"
#include <iterator>
using namespace edm;
using namespace std;
DTSegmentsTask::DTSegmentsTask(const edm::ParameterSet& pset)
: statusMapToken_(esConsumes<DTStatusFlag, DTStatusFlagRcd>()) {
debug = pset.getUntrackedParameter<bool>("debug", false);
parameters = pset;
// should be init from pset, but it was commented out...
// better false than undefined
checkNoisyChannels = false;
// the name of the 4D rec hits collection
theRecHits4DLabel_ = consumes<DTRecSegment4DCollection>(parameters.getParameter<std::string>("recHits4DLabel"));
}
DTSegmentsTask::~DTSegmentsTask() {}
void DTSegmentsTask::bookHistograms(DQMStore::IBooker& ibooker,
edm::Run const& /*iRun*/,
edm::EventSetup const& /* iSetup */) {
ibooker.cd();
ibooker.setCurrentFolder("Muons/DTSegmentsMonitor");
// histos for phi segments
phiHistos.push_back(
ibooker.book2D("phiSegments_numHitsVsWheel", "phiSegments_numHitsVsWheel", 5, -2.5, 2.5, 20, 0, 20));
phiHistos[0]->setBinLabel(1, "W-2", 1);
phiHistos[0]->setBinLabel(2, "W-1", 1);
phiHistos[0]->setBinLabel(3, "W0", 1);
phiHistos[0]->setBinLabel(4, "W1", 1);
phiHistos[0]->setBinLabel(5, "W2", 1);
phiHistos.push_back(
ibooker.book2D("phiSegments_numHitsVsSector", "phiSegments_numHitsVsSector", 14, 0.5, 14.5, 20, 0, 20));
phiHistos[1]->setBinLabel(1, "Sec1", 1);
phiHistos[1]->setBinLabel(2, "Sec2", 1);
phiHistos[1]->setBinLabel(3, "Sec3", 1);
phiHistos[1]->setBinLabel(4, "Sec4", 1);
phiHistos[1]->setBinLabel(5, "Sec5", 1);
phiHistos[1]->setBinLabel(6, "Sec6", 1);
phiHistos[1]->setBinLabel(7, "Sec7", 1);
phiHistos[1]->setBinLabel(8, "Sec8", 1);
phiHistos[1]->setBinLabel(9, "Sec9", 1);
phiHistos[1]->setBinLabel(10, "Sec10", 1);
phiHistos[1]->setBinLabel(11, "Sec11", 1);
phiHistos[1]->setBinLabel(12, "Sec12", 1);
phiHistos[1]->setBinLabel(13, "Sec13", 1);
phiHistos[1]->setBinLabel(14, "Sec14", 1);
phiHistos.push_back(
ibooker.book2D("phiSegments_numHitsVsStation", "phiSegments_numHitsVsStation", 4, 0.5, 4.5, 20, 0, 20));
phiHistos[2]->setBinLabel(1, "St1", 1);
phiHistos[2]->setBinLabel(2, "St2", 1);
phiHistos[2]->setBinLabel(3, "St3", 1);
phiHistos[2]->setBinLabel(4, "St4", 1);
// histos for theta segments
thetaHistos.push_back(
ibooker.book2D("thetaSegments_numHitsVsWheel", "thetaSegments_numHitsVsWheel", 5, -2.5, 2.5, 20, 0, 20));
thetaHistos[0]->setBinLabel(1, "W-2", 1);
thetaHistos[0]->setBinLabel(2, "W-1", 1);
thetaHistos[0]->setBinLabel(3, "W0", 1);
thetaHistos[0]->setBinLabel(4, "W1", 1);
thetaHistos[0]->setBinLabel(5, "W2", 1);
thetaHistos.push_back(
ibooker.book2D("thetaSegments_numHitsVsSector", "thetaSegments_numHitsVsSector", 14, 0.5, 14.5, 20, 0, 20));
thetaHistos[1]->setBinLabel(1, "Sec1", 1);
thetaHistos[1]->setBinLabel(2, "Sec2", 1);
thetaHistos[1]->setBinLabel(3, "Sec3", 1);
thetaHistos[1]->setBinLabel(4, "Sec4", 1);
thetaHistos[1]->setBinLabel(5, "Sec5", 1);
thetaHistos[1]->setBinLabel(6, "Sec6", 1);
thetaHistos[1]->setBinLabel(7, "Sec7", 1);
thetaHistos[1]->setBinLabel(8, "Sec8", 1);
thetaHistos[1]->setBinLabel(9, "Sec9", 1);
thetaHistos[1]->setBinLabel(10, "Sec10", 1);
thetaHistos[1]->setBinLabel(11, "Sec11", 1);
thetaHistos[1]->setBinLabel(12, "Sec12", 1);
thetaHistos[1]->setBinLabel(13, "Sec13", 1);
thetaHistos[1]->setBinLabel(14, "Sec14", 1);
thetaHistos.push_back(
ibooker.book2D("thetaSegments_numHitsVsStation", "thetaSegments_numHitsVsStation", 4, 0.5, 4.5, 20, 0, 20));
thetaHistos[2]->setBinLabel(1, "St1", 1);
thetaHistos[2]->setBinLabel(2, "St2", 1);
thetaHistos[2]->setBinLabel(3, "St3", 1);
thetaHistos[2]->setBinLabel(4, "St4", 1);
}
void DTSegmentsTask::analyze(const edm::Event& event, const edm::EventSetup& setup) {
// Get the map of noisy channels
// bool checkNoisyChannels = parameters.getUntrackedParameter<bool>("checkNoisyChannels",false);
ESHandle<DTStatusFlag> statusMap;
if (checkNoisyChannels) {
statusMap = setup.getHandle(statusMapToken_);
}
// Get the 4D segment collection from the event
edm::Handle<DTRecSegment4DCollection> all4DSegments;
event.getByToken(theRecHits4DLabel_, all4DSegments);
// Loop over all chambers containing a segment
DTRecSegment4DCollection::id_iterator chamberId;
for (chamberId = all4DSegments->id_begin(); chamberId != all4DSegments->id_end(); ++chamberId) {
// Get the range for the corresponding ChamerId
DTRecSegment4DCollection::range range = all4DSegments->get(*chamberId);
int nsegm = distance(range.first, range.second);
if (debug)
cout << " Chamber: " << *chamberId << " has " << nsegm << " 4D segments" << endl;
// Loop over the rechits of this ChamerId
for (DTRecSegment4DCollection::const_iterator segment4D = range.first; segment4D != range.second; ++segment4D) {
//FOR NOISY CHANNELS////////////////////////////////
bool segmNoisy = false;
if ((*segment4D).hasPhi()) {
const DTChamberRecSegment2D* phiSeg = (*segment4D).phiSegment();
vector<DTRecHit1D> phiHits = phiSeg->specificRecHits();
map<DTSuperLayerId, vector<DTRecHit1D> > hitsBySLMap;
for (vector<DTRecHit1D>::const_iterator hit = phiHits.begin(); hit != phiHits.end(); ++hit) {
DTWireId wireId = (*hit).wireId();
// Check for noisy channels to skip them
if (checkNoisyChannels) {
bool isNoisy = false;
bool isFEMasked = false;
bool isTDCMasked = false;
bool isTrigMask = false;
bool isDead = false;
bool isNohv = false;
statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
if (isNoisy) {
if (debug)
cout << "Wire: " << wireId << " is noisy, skipping!" << endl;
segmNoisy = true;
}
}
}
}
if ((*segment4D).hasZed()) {
const DTSLRecSegment2D* zSeg = (*segment4D).zSegment(); // zSeg lives in the SL RF
// Check for noisy channels to skip them
vector<DTRecHit1D> zHits = zSeg->specificRecHits();
for (vector<DTRecHit1D>::const_iterator hit = zHits.begin(); hit != zHits.end(); ++hit) {
DTWireId wireId = (*hit).wireId();
if (checkNoisyChannels) {
bool isNoisy = false;
bool isFEMasked = false;
bool isTDCMasked = false;
bool isTrigMask = false;
bool isDead = false;
bool isNohv = false;
//cout<<"wire id "<<wireId<<endl;
statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
if (isNoisy) {
if (debug)
cout << "Wire: " << wireId << " is noisy, skipping!" << endl;
segmNoisy = true;
}
}
}
}
if (segmNoisy) {
if (debug)
cout << "skipping the segment: it contains noisy cells" << endl;
continue;
}
//END FOR NOISY CHANNELS////////////////////////////////
// Fill the histos
int nHits = 0;
if ((*segment4D).hasPhi()) {
nHits = (((*segment4D).phiSegment())->specificRecHits()).size();
if (debug)
cout << "Phi segment with number of hits: " << nHits << endl;
phiHistos[0]->Fill((*chamberId).wheel(), nHits);
phiHistos[1]->Fill((*chamberId).sector(), nHits);
phiHistos[2]->Fill((*chamberId).station(), nHits);
}
if ((*segment4D).hasZed()) {
nHits = (((*segment4D).zSegment())->specificRecHits()).size();
if (debug)
cout << "Zed segment with number of hits: " << nHits << endl;
thetaHistos[0]->Fill((*chamberId).wheel(), nHits);
thetaHistos[1]->Fill((*chamberId).sector(), nHits);
thetaHistos[2]->Fill((*chamberId).station(), nHits);
}
} //loop over segments
} // loop over chambers
}
|