CMS 3D CMS Logo

DTResolutionAnalysisTask.cc
Go to the documentation of this file.
1 
2 /*
3  * See header file for a description of this class.
4  *
5  * \author G. Cerminara - INFN Torino
6  */
7 
9 
10 // Framework
14 //#include "FWCore/Framework/interface/LuminosityBlock.h"
16 
18 
20 
21 //Geometry
23 
24 //RecHit
26 
27 #include <iterator>
28 
29 using namespace edm;
30 using namespace std;
31 
33  : muonGeomToken_(esConsumes<edm::Transition::BeginRun>()) {
34  edm::LogVerbatim("DTDQM|DTMonitorModule|DTResolutionAnalysisTask")
35  << "[DTResolutionAnalysisTask] Constructor called!" << endl;
36 
37  // the name of the 4D rec hits collection
39  consumes<DTRecSegment4DCollection>(edm::InputTag(pset.getUntrackedParameter<string>("recHits4DLabel")));
40 
41  prescaleFactor = pset.getUntrackedParameter<int>("diagnosticPrescale", 1);
42  resetCycle = pset.getUntrackedParameter<int>("ResetCycle", -1);
43  // top folder for the histograms in DQMStore
44  topHistoFolder = pset.getUntrackedParameter<string>("topHistoFolder", "DT/02-Segments");
45 
46  thePhiHitsCut = pset.getUntrackedParameter<u_int32_t>("phiHitsCut", 8);
47  theZHitsCut = pset.getUntrackedParameter<u_int32_t>("zHitsCut", 4);
48 }
49 
51  edm::LogVerbatim("DTDQM|DTMonitorModule|DTResolutionAnalysisTask")
52  << "[DTResolutionAnalysisTask] Destructor called!" << endl;
53 }
54 
56  // Get the DT Geometry
57  dtGeom = &setup.getData(muonGeomToken_);
58 }
59 
61  edm::Run const& iRun,
62  edm::EventSetup const& /* iSetup */) {
63  // Book the histograms
64  vector<const DTChamber*> chambers = dtGeom->chambers();
65  for (vector<const DTChamber*>::const_iterator chamber = chambers.begin(); chamber != chambers.end();
66  ++chamber) { // Loop over all chambers
67  DTChamberId dtChId = (*chamber)->id();
68  for (int sl = 1; sl <= 3; ++sl) { // Loop over SLs
69  if (dtChId.station() == 4 && sl == 2)
70  continue;
71  const DTSuperLayerId dtSLId(dtChId, sl);
72  bookHistos(ibooker, dtSLId);
73  }
74  }
75 }
76 /*
77 void DTResolutionAnalysisTask::beginLuminosityBlock(const LuminosityBlock& lumiSeg, const EventSetup& context) {
78  edm::LogVerbatim("DTDQM|DTMonitorModule|DTResolutionAnalysisTask")
79  << "[DTResolutionTask]: Begin of LS transition" << endl;
80 
81  if (resetCycle != -1 && lumiSeg.id().luminosityBlock() % resetCycle == 0) {
82  for (map<DTSuperLayerId, vector<MonitorElement*> >::const_iterator histo = histosPerSL.begin();
83  histo != histosPerSL.end();
84  histo++) {
85  int size = (*histo).second.size();
86  for (int i = 0; i < size; i++) {
87  (*histo).second[i]->Reset();
88  }
89  }
90  }
91 }
92 */
94  edm::LogVerbatim("DTDQM|DTMonitorModule|DTResolutionAnalysisTask")
95  << "[DTResolutionAnalysisTask] Analyze #Run: " << event.id().run() << " #Event: " << event.id().event() << endl;
96 
97  // Get the 4D segment collection from the event
99  event.getByToken(recHits4DToken_, all4DSegments);
100 
101  // check the validity of the collection
102  if (!all4DSegments.isValid())
103  return;
104 
105  // Loop over all chambers containing a segment
107  for (chamberId = all4DSegments->id_begin(); chamberId != all4DSegments->id_end(); ++chamberId) {
108  // Get the range for the corresponding ChamerId
109  DTRecSegment4DCollection::range range = all4DSegments->get(*chamberId);
110 
111  // Get the chamber
113 
114  // Loop over the rechits of this ChamerId
115  for (DTRecSegment4DCollection::const_iterator segment4D = range.first; segment4D != range.second; ++segment4D) {
116  // If Statio != 4 skip RecHits with dimension != 4
117  // For the Station 4 consider 2D RecHits
118  if ((*chamberId).station() != 4 && (*segment4D).dimension() != 4) {
119  edm::LogVerbatim("DTDQM|DTMonitorModule|DTResolutionAnalysisTask")
120  << "[DTResolutionAnalysisTask]***Warning: RecSegment dimension is not 4 but " << (*segment4D).dimension()
121  << "!" << endl;
122  continue;
123  } else if ((*chamberId).station() == 4 && (*segment4D).dimension() != 2) {
124  edm::LogVerbatim("DTDQM|DTMonitorModule|DTResolutionAnalysisTask")
125  << "[DTResolutionAnalysisTask]***Warning: RecSegment dimension is not 2 but " << (*segment4D).dimension()
126  << "!" << endl;
127  continue;
128  }
129 
130  // Get all 1D RecHits at step 3 within the 4D segment
131  vector<DTRecHit1D> recHits1D_S3;
132 
133  // Get 1D RecHits at Step 3 and select only events with
134  // 8 hits in phi and 4 hits in theta (if any)
135 
136  if ((*segment4D).hasPhi()) { // has phi component
137  const DTChamberRecSegment2D* phiSeg = (*segment4D).phiSegment();
138  vector<DTRecHit1D> phiRecHits = phiSeg->specificRecHits();
139 
140  if (phiRecHits.size() < thePhiHitsCut) {
141  continue;
142  }
143  copy(phiRecHits.begin(), phiRecHits.end(), back_inserter(recHits1D_S3));
144  } else {
145  }
146 
147  if ((*segment4D).hasZed()) {
148  const DTSLRecSegment2D* zSeg = (*segment4D).zSegment();
149  vector<DTRecHit1D> zRecHits = zSeg->specificRecHits();
150  if (zRecHits.size() < theZHitsCut) {
151  continue;
152  }
153  copy(zRecHits.begin(), zRecHits.end(), back_inserter(recHits1D_S3));
154  }
155 
156  // Loop over 1D RecHit inside 4D segment
157  for (vector<DTRecHit1D>::const_iterator recHit1D = recHits1D_S3.begin(); recHit1D != recHits1D_S3.end();
158  recHit1D++) {
159  const DTWireId wireId = (*recHit1D).wireId();
160 
161  // Get the layer and the wire position
162  const DTLayer* layer = chamber->superLayer(wireId.superlayerId())->layer(wireId.layerId());
163  float wireX = layer->specificTopology().wirePosition(wireId.wire());
164 
165  // Distance of the 1D rechit from the wire
166  float distRecHitToWire = fabs(wireX - (*recHit1D).localPosition().x());
167 
168  // Extrapolate the segment to the z of the wire
169 
170  // Get wire position in chamber RF
171  LocalPoint wirePosInLay(wireX, (*recHit1D).localPosition().y(), (*recHit1D).localPosition().z());
172  GlobalPoint wirePosGlob = layer->toGlobal(wirePosInLay);
173  LocalPoint wirePosInChamber = chamber->toLocal(wirePosGlob);
174 
175  // Segment position at Wire z in chamber local frame
176  LocalPoint segPosAtZWire = (*segment4D).localPosition() + (*segment4D).localDirection() * wirePosInChamber.z() /
177  cos((*segment4D).localDirection().theta());
178 
179  // Compute the distance of the segment from the wire
180  int sl = wireId.superlayer();
181 
182  double distSegmToWire = -1;
183  if (sl == 1 || sl == 3) {
184  // RPhi SL
185  distSegmToWire = fabs(wirePosInChamber.x() - segPosAtZWire.x());
186  } else if (sl == 2) {
187  // RZ SL
188  distSegmToWire = fabs(wirePosInChamber.y() - segPosAtZWire.y());
189  }
190 
191  if (distSegmToWire > 2.1)
192  edm::LogVerbatim("DTDQM|DTMonitorModule|DTResolutionAnalysisTask")
193  << " Warning: dist segment-wire: " << distSegmToWire << endl;
194 
195  double residual = distRecHitToWire - distSegmToWire;
196  // FIXME: Fill the histos
197  fillHistos(wireId.superlayerId(), distSegmToWire, residual);
198 
199  } // End of loop over 1D RecHit inside 4D segment
200  } // End of loop over the rechits of this ChamerId
201  }
202  // -----------------------------------------------------------------------------
203 }
204 
205 // Book a set of histograms for a given SL
207  edm::LogVerbatim("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << " Booking histos for SL: " << slId << endl;
208 
209  // Compose the chamber name
210  stringstream wheel;
211  wheel << slId.wheel();
212  stringstream station;
213  station << slId.station();
214  stringstream sector;
215  sector << slId.sector();
216  stringstream superLayer;
217  superLayer << slId.superlayer();
218 
219  string slHistoName = "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str() + "_SL" + superLayer.str();
220 
221  ibooker.setCurrentFolder(topHistoFolder + "/Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" +
222  station.str());
223  // Create the monitor elements
224  vector<MonitorElement*> histos;
225  // Note the order matters
226  histos.push_back(ibooker.book1D(
227  "hResDist" + slHistoName, "Residuals on the distance from wire (rec_hit - segm_extr) (cm)", 200, -0.4, 0.4));
228  histosPerSL[slId] = histos;
229 }
230 
231 // Fill a set of histograms for a given SL
232 void DTResolutionAnalysisTask::fillHistos(DTSuperLayerId slId, float distExtr, float residual) {
233  vector<MonitorElement*> histos = histosPerSL[slId];
234  histos[0]->Fill(residual);
235 }
236 
237 // Local Variables:
238 // show-trailing-whitespace: t
239 // truncate-lines: t
240 // End:
Log< level::Info, true > LogVerbatim
int station() const
Return the station number.
Definition: DTChamberId.h:45
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
int wire() const
Return the wire number.
Definition: DTWireId.h:45
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
T z() const
Definition: PV3DBase.h:61
identifier iterator
Definition: RangeMap.h:130
~DTResolutionAnalysisTask() override
Destructor.
edm::ESGetToken< DTGeometry, MuonGeometryRecord > muonGeomToken_
void analyze(const edm::Event &event, const edm::EventSetup &setup) override
To reset the MEs.
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Transition
Definition: Transition.h:12
edm::EDGetTokenT< DTRecSegment4DCollection > recHits4DToken_
std::map< DTSuperLayerId, std::vector< MonitorElement * > > histosPerSL
int superlayer() const
Return the superlayer number (deprecated method name)
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
BookHistograms.
bool isValid() const
Definition: HandleBase.h:70
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:42
HLT enums.
int sector() const
Definition: DTChamberId.h:52
void bookHistos(DQMStore::IBooker &ibooker, DTSuperLayerId slId)
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
BeginRun.
void fillHistos(DTSuperLayerId slId, float distExtr, float residual)
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:48
const std::vector< const DTChamber * > & chambers() const
Return a vector of all Chamber.
Definition: DTGeometry.cc:84
DTLayerId layerId() const
Return the corresponding LayerId.
Definition: DTWireId.h:48
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
static char chambers[264][20]
Definition: ReadPGInfo.cc:243
DTResolutionAnalysisTask(const edm::ParameterSet &pset)
Constructor.
const DTChamber * chamber(const DTChamberId &id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:90
Definition: event.py:1
Definition: Run.h:45