CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTCalibValidation.cc
Go to the documentation of this file.
1 
2 /*
3  * See header file for a description of this class.
4  *
5  * \author G. Mila - INFN Torino
6  */
7 
9 
10 // Framework
18 
19 //Geometry
22 
23 //RecHit
26 
27 
28 #include <iterator>
29 
30 using namespace edm;
31 using namespace std;
32 
33 
35 
36  //debug = pset.getUntrackedParameter<bool>("debug",false);
37 
38  // Get the DQM needed services
39 
40  // To remove into CMSSW versions before 20X
41  theDbe = edm::Service<DQMStore>().operator->();
42  // To add into CMSSW versions before 20X
43  /*theDbe = edm::Service<DaqMonitorBEInterface>().operator->();
44  theDbe->setVerbose(1);
45  edm::Service<MonitorDaemon>().operator->();*/
46 
47  theDbe->setCurrentFolder("DT/DTCalibValidation");
48 
49  parameters = pset;
50 }
51 
52 
54 }
55 
56 
58 
59  // the name of the rechits collection at step 1
60  recHits1DToken_ = consumes<DTRecHitCollection>(
61  edm::InputTag(parameters.getUntrackedParameter<string>("recHits1DLabel")));
62  // the name of the 2D segments
63  segment2DToken_ = consumes<DTRecSegment2DCollection>(
64  edm::InputTag(parameters.getUntrackedParameter<string>("segment2DLabel")));
65  // the name of the 4D segments
66  segment4DToken_ = consumes<DTRecSegment4DCollection>(
67  edm::InputTag(parameters.getUntrackedParameter<string>("segment4DLabel")));
68  // the counter of segments not used to compute residuals
69  wrongSegment = 0;
70  // the counter of segments used to compute residuals
71  rightSegment = 0;
72  // the analysis type
73  detailedAnalysis = parameters.getUntrackedParameter<bool>("detailedAnalysis",false);
74 
75  nevent=0;
76 
77 }
78 
80 
81  // get the geometry
82  setup.get<MuonGeometryRecord>().get(dtGeom);
83 
84  // Loop over all the chambers
85  vector<const DTChamber*>::const_iterator ch_it = dtGeom->chambers().begin();
86  vector<const DTChamber*>::const_iterator ch_end = dtGeom->chambers().end();
87  for (; ch_it != ch_end; ++ch_it) {
88  vector<const DTSuperLayer*>::const_iterator sl_it = (*ch_it)->superLayers().begin();
89  vector<const DTSuperLayer*>::const_iterator sl_end = (*ch_it)->superLayers().end();
90  // Loop over the SLs
91  for(; sl_it != sl_end; ++sl_it) {
92  DTSuperLayerId slId = (*sl_it)->id();
93  if(detailedAnalysis){
94  // Loop over the 3 steps
95  for(int step = 1; step <= 3; ++step) bookHistos(slId,step);
96  } else {
97  // Only step 3
98  bookHistos(slId,3);
99  }
100  }
101  }
102 
103 }
104 
105 
107 
108  LogVerbatim("DTCalibValidation") << "Segments used to compute residuals: " << rightSegment;
109  LogVerbatim("DTCalibValidation") << "Segments not used to compute residuals: " << wrongSegment;
110 
111  //theDbe->showDirStructure();
112  bool outputMEsInRootFile = parameters.getParameter<bool>("OutputMEsInRootFile");
113  std::string outputFileName = parameters.getParameter<std::string>("OutputFileName");
114  if(outputMEsInRootFile){
115  theDbe->showDirStructure();
116  theDbe->save(outputFileName);
117  }
118 
119  theDbe->rmdir("DT/DTCalibValidation");
120 }
121 
122 
123 
125 
126  ++nevent;
127  LogTrace("DTCalibValidation") << "[DTCalibValidation] Analyze #Run: " << event.id().run()
128  << " #Event: " << nevent;
129 
130  // RecHit mapping at Step 1 -------------------------------
131  map<DTWireId,vector<DTRecHit1DPair> > recHitsPerWire_1S;
132 
133  // RecHit mapping at Step 2 ------------------------------
134  map<DTWireId,vector<DTRecHit1D> > recHitsPerWire_2S;
135 
136  if(detailedAnalysis){
137  LogTrace("DTCalibValidation") << " -- DTRecHit S1: begin analysis:";
138  // Get the rechit collection from the event
139  Handle<DTRecHitCollection> dtRecHits;
140  event.getByToken(recHits1DToken_, dtRecHits);
141  recHitsPerWire_1S = map1DRecHitsPerWire(dtRecHits.product());
142 
143  LogTrace("DTCalibValidation") << " -- DTRecHit S2: begin analysis:";
144  // Get the 2D rechits from the event
146  event.getByToken(segment2DToken_, segment2Ds);
147  recHitsPerWire_2S = map1DRecHitsPerWire(segment2Ds.product());
148  }
149 
150  // RecHit mapping at Step 3 ---------------------------------
151  LogTrace("DTCalibValidation") << " -- DTRecHit S3: begin analysis:";
152  // Get the 4D rechits from the event
154  event.getByToken(segment4DToken_, segment4Ds);
155  map<DTWireId,vector<DTRecHit1D> > recHitsPerWire_3S = map1DRecHitsPerWire(segment4Ds.product());
156 
157 
158  // Loop over all 4D segments
159  for(DTRecSegment4DCollection::const_iterator segment = segment4Ds->begin();
160  segment != segment4Ds->end();
161  ++segment) {
162 
163  if(detailedAnalysis){
164  LogTrace("DTCalibValidation") << "Anlysis on recHit at step 1";
165  compute(dtGeom.product(), (*segment), recHitsPerWire_1S, 1);
166 
167  LogTrace("DTCalibValidation") << "Anlysis on recHit at step 2";
168  compute(dtGeom.product(), (*segment), recHitsPerWire_2S, 2);
169  }
170 
171  LogTrace("DTCalibValidation") << "Anlysis on recHit at step 3";
172  compute(dtGeom.product(), (*segment), recHitsPerWire_3S, 3);
173  }
174 
175 }
176 
177 
178 // Return a map between DTRecHit1DPair and wireId
179 map<DTWireId, vector<DTRecHit1DPair> >
181  map<DTWireId, vector<DTRecHit1DPair> > ret;
182 
183  for(DTRecHitCollection::const_iterator rechit = dt1DRecHitPairs->begin();
184  rechit != dt1DRecHitPairs->end(); ++rechit) {
185  ret[(*rechit).wireId()].push_back(*rechit);
186  }
187 
188  return ret;
189 }
190 
191 
192 // Return a map between DTRecHit1D at S2 and wireId
193 map<DTWireId, vector<DTRecHit1D> >
195  map<DTWireId, vector<DTRecHit1D> > ret;
196 
197  // Loop over all 2D segments
198  for(DTRecSegment2DCollection::const_iterator segment = segment2Ds->begin();
199  segment != segment2Ds->end();
200  ++segment) {
201  vector<DTRecHit1D> component1DHits= (*segment).specificRecHits();
202  // Loop over all component 1D hits
203  for(vector<DTRecHit1D>::const_iterator hit = component1DHits.begin();
204  hit != component1DHits.end(); ++hit) {
205  ret[(*hit).wireId()].push_back(*hit);
206  }
207  }
208  return ret;
209 }
210 
211 // Return a map between DTRecHit1D at S3 and wireId
212 map<DTWireId, std::vector<DTRecHit1D> >
214  map<DTWireId, vector<DTRecHit1D> > ret;
215  // Loop over all 4D segments
216  for(DTRecSegment4DCollection::const_iterator segment = segment4Ds->begin();
217  segment != segment4Ds->end();
218  ++segment) {
219  // Get component 2D segments
220  vector<const TrackingRecHit*> segment2Ds = (*segment).recHits();
221  // Loop over 2D segments:
222  for(vector<const TrackingRecHit*>::const_iterator segment2D = segment2Ds.begin();
223  segment2D != segment2Ds.end();
224  ++segment2D) {
225  // Get 1D component rechits
226  vector<const TrackingRecHit*> hits = (*segment2D)->recHits();
227  // Loop over them
228  for(vector<const TrackingRecHit*>::const_iterator hit = hits.begin();
229  hit != hits.end(); ++hit) {
230  const DTRecHit1D* hit1D = dynamic_cast<const DTRecHit1D*>(*hit);
231  ret[hit1D->wireId()].push_back(*hit1D);
232  }
233  }
234  }
235 
236  return ret;
237 }
238 
239 
240 
241 // Find the RecHit closest to the segment4D
242 template <typename type>
243 const type*
245  DTWireId wireId,
246  const vector<type>& recHits,
247  const float segmDist) {
248  float res = 99999;
249  const type* theBestRecHit = 0;
250  // Loop over RecHits within the cell
251  for(typename vector<type>::const_iterator recHit = recHits.begin();
252  recHit != recHits.end();
253  ++recHit) {
254  float distTmp = recHitDistFromWire(*recHit, layer);
255  if(fabs(distTmp-segmDist) < res) {
256  res = fabs(distTmp-segmDist);
257  theBestRecHit = &(*recHit);
258  }
259  } // End of loop over RecHits within the cell
260 
261  return theBestRecHit;
262 }
263 
264 
265 // Compute the distance from wire (cm) of a hits in a DTRecHit1DPair
266 float
268  return fabs(hitPair.localPosition(DTEnums::Left).x() -
269  hitPair.localPosition(DTEnums::Right).x())/2.;
270 }
271 
272 
273 // Compute the distance from wire (cm) of a hits in a DTRecHit1D
274 float
276  return fabs(recHit.localPosition().x() - layer->specificTopology().wirePosition(recHit.wireId().wire()));
277 
278  //to check the compatibility position / distance
279  /* // Get the recHit and the wire position
280  const DTChamber* chamber = (*layer).chamber();
281  GlobalPoint recHitPosGlob = layer->toGlobal(recHit.localPosition());
282  LocalPoint recHitPosInChamber = chamber->toLocal(recHitPosGlob);
283  float wireX = layer->specificTopology().wirePosition(recHit.wireId().wire());
284  LocalPoint wirePosInLay(wireX,recHit.localPosition().y(),0);
285  GlobalPoint wirePosGlob = layer->toGlobal(wirePosInLay);
286  LocalPoint wirePosInChamber = chamber->toLocal(wirePosGlob);
287 
288  float recHitDist = -1;
289  if(recHit.wireId().layerId().superlayerId().superlayer() != 2)
290  recHitDist = fabs(recHitPosInChamber.x()-wirePosInChamber.x());
291  else
292  recHitDist = fabs(recHitPosInChamber.y()-wirePosInChamber.y());
293 
294  return recHitDist; */
295 
296 }
297 
298 
299 
300 
301 // Compute the position (cm) of a hits in a DTRecHit1DPair
302 float
303 DTCalibValidation::recHitPosition(const DTRecHit1DPair& hitPair, const DTLayer* layer, const DTChamber* chamber, float segmentPos, int sl) {
304 
305  // Get the layer and the wire position
306  GlobalPoint hitPosGlob_right = layer->toGlobal(hitPair.localPosition(DTEnums::Right));
307  LocalPoint hitPosInChamber_right = chamber->toLocal(hitPosGlob_right);
308  GlobalPoint hitPosGlob_left = layer->toGlobal(hitPair.localPosition(DTEnums::Left));
309  LocalPoint hitPosInChamber_left = chamber->toLocal(hitPosGlob_left);
310 
311  float recHitPos=-1;
312  if(sl != 2){
313  if(fabs(hitPosInChamber_left.x()-segmentPos)<fabs(hitPosInChamber_right.x()-segmentPos))
314  recHitPos = hitPosInChamber_left.x();
315  else
316  recHitPos = hitPosInChamber_right.x();
317  }
318  else{
319  if(fabs(hitPosInChamber_left.y()-segmentPos)<fabs(hitPosInChamber_right.y()-segmentPos))
320  recHitPos = hitPosInChamber_left.y();
321  else
322  recHitPos = hitPosInChamber_right.y();
323  }
324 
325  return recHitPos;
326 }
327 
328 
329 // Compute the position (cm) of a hits in a DTRecHit1D
330 float
331 DTCalibValidation::recHitPosition(const DTRecHit1D& recHit, const DTLayer* layer, const DTChamber* chamber, float segmentPos, int sl) {
332 
333  // Get the layer and the wire position
334  GlobalPoint recHitPosGlob = layer->toGlobal(recHit.localPosition());
335  LocalPoint recHitPosInChamber = chamber->toLocal(recHitPosGlob);
336 
337  float recHitPos = -1;
338  if(sl != 2)
339  recHitPos = recHitPosInChamber.x();
340  else
341  recHitPos = recHitPosInChamber.y();
342 
343  return recHitPos;
344 
345 }
346 
347 
348 
349 // Compute the residuals
350 template <typename type>
352  const DTRecSegment4D& segment,
353  const std::map<DTWireId, std::vector<type> >& recHitsPerWire,
354  int step) {
355  bool computeResidual = true;
356 
357  // Get all 1D RecHits at step 3 within the 4D segment
358  vector<DTRecHit1D> recHits1D_S3;
359 
360  // Get 1D RecHits at Step 3 and select only events with
361  // 8 hits in phi and 4 hits in theta (if any)
362  const DTChamberRecSegment2D* phiSeg = segment.phiSegment();
363  if(phiSeg){
364  vector<DTRecHit1D> phiRecHits = phiSeg->specificRecHits();
365  if(phiRecHits.size() != 8) {
366  LogTrace("DTCalibValidation") << "[DTCalibValidation] Phi segments has: " << phiRecHits.size()
367  << " hits, skipping"; // FIXME: info output
368  computeResidual = false;
369  }
370  copy(phiRecHits.begin(), phiRecHits.end(), back_inserter(recHits1D_S3));
371  }
372  if(!phiSeg){
373  LogTrace("DTCalibValidation") << " [DTCalibValidation] 4D segment has not the phi segment! ";
374  computeResidual = false;
375  }
376 
377  if(segment.dimension() == 4) {
378  const DTSLRecSegment2D* zSeg = segment.zSegment();
379  if(zSeg){
380  vector<DTRecHit1D> zRecHits = zSeg->specificRecHits();
381  if(zRecHits.size() != 4) {
382  LogTrace("DTCalibValidation") << "[DTCalibValidation] Theta segments has: " << zRecHits.size()
383  << " hits, skipping"; // FIXME: info output
384  computeResidual = false;
385  }
386  copy(zRecHits.begin(), zRecHits.end(), back_inserter(recHits1D_S3));
387  }
388  if(!zSeg){
389  LogTrace("DTCalibValidation") << " [DTCalibValidation] 4D segment has not the z segment! ";
390  computeResidual = false;
391  }
392  }
393 
394  if(!computeResidual)
395  ++wrongSegment;
396  if(computeResidual){
397  ++rightSegment;
398  // Loop over 1D RecHit inside 4D segment
399  for(vector<DTRecHit1D>::const_iterator recHit1D = recHits1D_S3.begin();
400  recHit1D != recHits1D_S3.end();
401  ++recHit1D) {
402  const DTWireId wireId = (*recHit1D).wireId();
403 
404  // Get the layer and the wire position
405  const DTLayer* layer = dtGeom->layer(wireId);
406  float wireX = layer->specificTopology().wirePosition(wireId.wire());
407 
408  // Extrapolate the segment to the z of the wire
409  // Get wire position in chamber RF
410  // (y and z must be those of the hit to be coherent in the transf. of RF in case of rotations of the layer alignment)
411  LocalPoint wirePosInLay(wireX,(*recHit1D).localPosition().y(),(*recHit1D).localPosition().z());
412  GlobalPoint wirePosGlob = layer->toGlobal(wirePosInLay);
413  const DTChamber* chamber = dtGeom->chamber((*recHit1D).wireId().layerId().chamberId());
414  LocalPoint wirePosInChamber = chamber->toLocal(wirePosGlob);
415 
416  // Segment position at Wire z in chamber local frame
417  LocalPoint segPosAtZWire = segment.localPosition()
418  + segment.localDirection()*wirePosInChamber.z()/cos(segment.localDirection().theta());
419 
420  // Compute the distance of the segment from the wire
421  int sl = wireId.superlayer();
422  float SegmDistance = -1;
423  if(sl == 1 || sl == 3) {
424  // RPhi SL
425  SegmDistance = fabs(wirePosInChamber.x() - segPosAtZWire.x());
426  LogTrace("DTCalibValidation") << "SegmDistance: " << SegmDistance;
427  } else if(sl == 2) {
428  // RZ SL
429  SegmDistance = fabs(segPosAtZWire.y() - wirePosInChamber.y());
430  LogTrace("DTCalibValidation") << "SegmDistance: " << SegmDistance;
431  }
432  if(SegmDistance > 2.1)
433  LogTrace("DTCalibValidation") << " Warning: dist segment-wire: " << SegmDistance;
434 
435  // Look for RecHits in the same cell
436  if(recHitsPerWire.find(wireId) == recHitsPerWire.end()) {
437  LogTrace("DTCalibValidation") << " No RecHit found at Step: " << step << " in cell: " << wireId;
438  } else {
439  vector<type> recHits = recHitsPerWire.at(wireId);
440  LogTrace("DTCalibValidation") << " " << recHits.size() << " RecHits, Step " << step << " in channel: " << wireId;
441 
442  // Get the layer
443  const DTLayer* layer = dtGeom->layer(wireId);
444  // Find the best RecHits
445  const type* theBestRecHit = findBestRecHit(layer, wireId, recHits, SegmDistance);
446  // Compute the distance of the recHit from the wire
447  float recHitWireDist = recHitDistFromWire(*theBestRecHit, layer);
448  LogTrace("DTCalibValidation") << "recHitWireDist: " << recHitWireDist;
449 
450  // Compute the residuals
451  float residualOnDistance = recHitWireDist - SegmDistance;
452  LogTrace("DTCalibValidation") << "WireId: " << wireId << " ResidualOnDistance: " << residualOnDistance;
453  float residualOnPosition = -1;
454  float recHitPos = -1;
455  if(sl == 1 || sl == 3) {
456  recHitPos = recHitPosition(*theBestRecHit, layer, chamber, segPosAtZWire.x(), sl);
457  residualOnPosition = recHitPos - segPosAtZWire.x();
458  }
459  else{
460  recHitPos = recHitPosition(*theBestRecHit, layer, chamber, segPosAtZWire.y(), sl);
461  residualOnPosition = recHitPos - segPosAtZWire.y();
462  }
463  LogTrace("DTCalibValidation") << "WireId: " << wireId << " ResidualOnPosition: " << residualOnPosition;
464 
465  // Fill the histos
466  if(sl == 1 || sl == 3)
467  fillHistos(wireId.superlayerId(), SegmDistance, residualOnDistance, (wirePosInChamber.x() - segPosAtZWire.x()), residualOnPosition, step);
468  else
469  fillHistos(wireId.superlayerId(), SegmDistance, residualOnDistance, (wirePosInChamber.y() - segPosAtZWire.y()), residualOnPosition, step);
470 
471 
472  }
473  }
474  }
475 
476 }
477 
478 
479 
480 // Book a set of histograms for a given SL
482  LogTrace("DTCalibValidation") << " Booking histos for SL: " << slId;
483 
484  // Compose the chamber name
485  stringstream wheel; wheel << slId.wheel();
486  stringstream station; station << slId.station();
487  stringstream sector; sector << slId.sector();
488  stringstream superLayer; superLayer << slId.superlayer();
489  // Define the step
490  stringstream Step; Step << step;
491 
492 
493  string slHistoName =
494  "_STEP" + Step.str() +
495  "_W" + wheel.str() +
496  "_St" + station.str() +
497  "_Sec" + sector.str() +
498  "_SL" + superLayer.str();
499 
500  theDbe->setCurrentFolder("DT/DTCalibValidation/Wheel" + wheel.str() +
501  "/Station" + station.str() +
502  "/Sector" + sector.str());
503  // Create the monitor elements
504  vector<MonitorElement *> histos;
505  // Note hte order matters
506  histos.push_back(theDbe->book1D("hResDist"+slHistoName,
507  "Residuals on the distance from wire (rec_hit - segm_extr) (cm)",
508  200, -0.4, 0.4));
509  histos.push_back(theDbe->book2D("hResDistVsDist"+slHistoName,
510  "Residuals on the distance (cm) from wire (rec_hit - segm_extr) vs distance (cm)",
511  100, 0, 2.5, 200, -0.4, 0.4));
512  if(detailedAnalysis){
513  histos.push_back(theDbe->book1D("hResPos"+slHistoName,
514  "Residuals on the position from wire (rec_hit - segm_extr) (cm)",
515  200, -0.4, 0.4));
516  histos.push_back(theDbe->book2D("hResPosVsPos"+slHistoName,
517  "Residuals on the position (cm) from wire (rec_hit - segm_extr) vs distance (cm)",
518  200, -2.5, 2.5, 200, -0.4, 0.4));
519  }
520 
521  histosPerSL[make_pair(slId, step)] = histos;
522 }
523 
524 
525 
526 // Fill a set of histograms for a given SL
528  float distance,
529  float residualOnDistance,
530  float position,
531  float residualOnPosition,
532  int step) {
533  // FIXME: optimization of the number of searches
534  vector<MonitorElement *> histos = histosPerSL[make_pair(slId,step)];
535  histos[0]->Fill(residualOnDistance);
536  histos[1]->Fill(distance, residualOnDistance);
537  if(detailedAnalysis){
538  histos[2]->Fill(residualOnPosition);
539  histos[3]->Fill(position, residualOnPosition);
540  }
541 }
542 
543 
544 // Local Variables:
545 // show-trailing-whitespace: t
546 // truncate-lines: t
547 // End:
type
Definition: HCALResponse.h:21
float wirePosition(int wireNumber) const
Returns the x position in the layer of a given wire number.
Definition: DTTopology.cc:86
dictionary parameters
Definition: Parameters.py:2
float recHitDistFromWire(const DTRecHit1DPair &hitPair, const DTLayer *layer)
const DTChamberRecSegment2D * phiSegment() const
The superPhi segment: 0 if no phi projection available.
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:52
T y() const
Definition: PV3DBase.h:63
void endJob()
Endjob.
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:67
void beginRun(const edm::Run &, const edm::EventSetup &)
BeginRun.
void bookHistos()
Definition: Histogram.h:33
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
const DTLayer * layer(DTLayerId id) const
Return a layer given its id.
Definition: DTGeometry.cc:110
int nevent
Definition: AMPTWrapper.h:74
virtual LocalVector localDirection() const
Local direction in Chamber frame.
void beginJob()
BeginJob.
void compute(const DTGeometry *dtGeom, const DTRecSegment4D &segment, const std::map< DTWireId, std::vector< type > > &recHitsPerWire, int step)
const DTTopology & specificTopology() const
Definition: DTLayer.cc:42
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::map< DTWireId, std::vector< DTRecHit1DPair > > map1DRecHitsPerWire(const DTRecHitCollection *dt1DRecHitPairs)
void fillHistos(DTSuperLayerId slId, float distance, float residualOnDistance, float position, float residualOnPosition, int step)
void analyze(const edm::Event &event, const edm::EventSetup &setup)
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
virtual LocalPoint localPosition() const
Local position in Chamber frame.
#define LogTrace(id)
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
virtual int dimension() const
Dimension (in parameter space)
virtual LocalPoint localPosition() const
Return the 3-dimensional local position.
Definition: DTRecHit1D.h:60
int wire() const
Return the wire number.
Definition: DTWireId.h:56
int superlayer() const
Return the superlayer number (deprecated method name)
const DTSLRecSegment2D * zSegment() const
The Z segment: 0 if not zed projection available.
T const * product() const
Definition: Handle.h:81
virtual ~DTCalibValidation()
Destructor.
const T & get() const
Definition: EventSetup.h:55
const DTChamber * chamber(DTChamberId id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:100
void bookHistos(DTSuperLayerId slId, int step)
int sector() const
Definition: DTChamberId.h:61
static int position[264][3]
Definition: ReadPGInfo.cc:509
DTCalibValidation(const edm::ParameterSet &pset)
Constructor.
int station() const
Return the station number.
Definition: DTChamberId.h:51
float recHitPosition(const DTRecHit1DPair &hitPair, const DTLayer *layer, const DTChamber *chamber, float segmPos, int sl)
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
T x() const
Definition: PV3DBase.h:62
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
Definition: Run.h:41
virtual LocalPoint localPosition() const
const type * findBestRecHit(const DTLayer *layer, DTWireId wireId, const std::vector< type > &recHits, const float simHitDist)
DTWireId wireId() const
Return the wireId.
Definition: DTRecHit1D.h:107