CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTRecHitQuality.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 
8 #include "DTRecHitQuality.h"
10 
14 
20 
25 
26 
27 #include <iostream>
28 #include <map>
29 
30 using namespace std;
31 using namespace edm;
32 
33 
34 
35 
36 // Constructor
38  // Get the debug parameter for verbose output
39  debug = pset.getUntrackedParameter<bool>("debug");
40  // the name of the simhit collection
41  simHitLabel = pset.getUntrackedParameter<InputTag>("simHitLabel");
42  // the name of the 1D rec hit collection
43  recHitLabel = pset.getUntrackedParameter<InputTag>("recHitLabel");
44  // the name of the 2D rec hit collection
45  segment2DLabel = pset.getUntrackedParameter<InputTag>("segment2DLabel");
46  // the name of the 4D rec hit collection
47  segment4DLabel = pset.getUntrackedParameter<InputTag>("segment4DLabel");
48 
49  // Switches for analysis at various steps
50  doStep1 = pset.getUntrackedParameter<bool>("doStep1", false);
51  doStep2 = pset.getUntrackedParameter<bool>("doStep2", false);
52  doStep3 = pset.getUntrackedParameter<bool>("doStep3", false);
53  doall = pset.getUntrackedParameter<bool>("doall", false);
54  local = pset.getUntrackedParameter<bool>("local", true);
55  // if(doall) doStep1
56  // Create the root file
57  //theFile = new TFile(rootFileName.c_str(), "RECREATE");
58  //theFile->cd();
59 
60 
61  // ----------------------
62  // get hold of back-end interface
63  dbe_ = 0;
64  dbe_ = Service<DQMStore>().operator->();
65  /*if ( dbe_ ) {
66  if (debug) {
67  dbe_->setVerbose(1);
68  } else {
69  dbe_->setVerbose(0);
70  }
71  }*/
72  dbe_->setVerbose(0);
73  /*if ( dbe_ ) {
74  if ( debug ) dbe_->showDirStructure();
75  }*/
76  if(doall && doStep1){
77  hRes_S1RPhi= new HRes1DHit("S1RPhi",dbe_,true,local); // RecHits, 1. step, RPhi
78  hRes_S1RPhi_W0= new HRes1DHit("S1RPhi_W0",dbe_,true,local); // RecHits, 1. step, RZ, wheel 0
79  hRes_S1RPhi_W1= new HRes1DHit("S1RPhi_W1",dbe_,true,local); // RecHits, 1. step, RZ, wheel +-1
80  hRes_S1RPhi_W2= new HRes1DHit("S1RPhi_W2",dbe_,true,local); // RecHits, 1. step, RZ, wheel +-2
81  hRes_S1RZ= new HRes1DHit("S1RZ",dbe_,true,local); // RecHits, 1. step, RZ
82  hRes_S1RZ_W0= new HRes1DHit("S1RZ_W0",dbe_,true,local); // RecHits, 1. step, RZ, wheel 0
83  hRes_S1RZ_W1= new HRes1DHit("S1RZ_W1",dbe_,true,local); // RecHits, 1. step, RZ, wheel +-1
84  hRes_S1RZ_W2= new HRes1DHit("S1RZ_W2",dbe_,true,local); // RecHits, 1. step, RZ, wheel +-2
85  hEff_S1RPhi= new HEff1DHit("S1RPhi",dbe_); // RecHits, 1. step, RPhi
86  hEff_S1RZ= new HEff1DHit("S1RZ",dbe_); // RecHits, 1. step, RZ
87  hEff_S1RZ_W0= new HEff1DHit("S1RZ_W0",dbe_); // RecHits, 1. step, RZ, wheel 0
88  hEff_S1RZ_W1= new HEff1DHit("S1RZ_W1",dbe_); // RecHits, 1. step, RZ, wheel +-1
89  hEff_S1RZ_W2= new HEff1DHit("S1RZ_W2",dbe_); // RecHits, 1. step, RZ, wheel +-2
90  }
91  if(doall && doStep2){
92  hRes_S2RPhi= new HRes1DHit("S2RPhi",dbe_,true,local); // RecHits, 2. step, RPhi
93  hRes_S2RPhi_W0= new HRes1DHit("S2RPhi_W0",dbe_,true,local); // RecHits, 2. step, RPhi, wheel 0
94  hRes_S2RPhi_W1= new HRes1DHit("S2RPhi_W1",dbe_,true,local); // RecHits, 2. step, RPhi, wheel +-1
95  hRes_S2RPhi_W2= new HRes1DHit("S2RPhi_W2",dbe_,true,local); // RecHits, 2. step, RPhi, wheel +-2
96  hRes_S2RZ= new HRes1DHit("S2RZ",dbe_,true,local); // RecHits, 2. step, RZ
97  hRes_S2RZ_W0= new HRes1DHit("S2RZ_W0",dbe_,true,local); // RecHits, 2. step, RZ, wheel 0
98  hRes_S2RZ_W1= new HRes1DHit("S2RZ_W1",dbe_,true,local); // RecHits, 2. step, RZ, wheel +-1
99  hRes_S2RZ_W2= new HRes1DHit("S2RZ_W2",dbe_,true,local); // RecHits, 2. step, RZ, wheel +-2
100  hEff_S2RPhi= new HEff1DHit("S2RPhi",dbe_); // RecHits, 2. step, RPhi
101  hEff_S2RZ_W0= new HEff1DHit("S2RZ_W0",dbe_); // RecHits, 2. step, RZ, wheel 0
102  hEff_S2RZ_W1= new HEff1DHit("S2RZ_W1",dbe_); // RecHits, 2. step, RZ, wheel +-1
103  hEff_S2RZ_W2= new HEff1DHit("S2RZ_W2",dbe_); // RecHits, 2. step, RZ, wheel +-2
104  hEff_S2RZ= new HEff1DHit("S2RZ",dbe_); // RecHits, 2. step, RZ
105  }
106  if(doStep3){
107  hRes_S3RPhi= new HRes1DHit("S3RPhi",dbe_,doall,local); // RecHits, 3. step, RPhi
108  hRes_S3RPhi_W0= new HRes1DHit("S3RPhi_W0",dbe_,doall,local); // RecHits, 3. step, RPhi, wheel 0
109  hRes_S3RPhi_W1= new HRes1DHit("S3RPhi_W1",dbe_,doall,local); // RecHits, 3. step, RPhi, wheel +-1
110  hRes_S3RPhi_W2= new HRes1DHit("S3RPhi_W2",dbe_,doall,local); // RecHits, 3. step, RPhi, wheel +-2
111  hRes_S3RZ= new HRes1DHit("S3RZ",dbe_,doall,local); // RecHits, 3. step, RZ
112  hRes_S3RZ_W0= new HRes1DHit("S3RZ_W0",dbe_,doall,local); // RecHits, 3. step, RZ, wheel 0
113  hRes_S3RZ_W1= new HRes1DHit("S3RZ_W1",dbe_,doall,local); // RecHits, 3. step, RZ, wheel +-1
114  hRes_S3RZ_W2= new HRes1DHit("S3RZ_W2",dbe_,doall,local); // RecHits, 3. step, RZ, wheel +-2
115 
116  if (local) {
117  // Plots with finer granularity, not to be included in DQM
118  TString name1="RPhi_W";
119  TString name2="RZ_W";
120  for (long w=0;w<=2;++w) {
121  for (long s=1;s<=4;++s){
122  hRes_S3RPhiWS[w][s-1] = new HRes1DHit(("S3"+name1+w+"_St"+s).Data(),dbe_,doall,local);
123  hEff_S1RPhiWS[w][s-1] = new HEff1DHit(("S1"+name1+w+"_St"+s).Data(),dbe_);
124  hEff_S3RPhiWS[w][s-1] = new HEff1DHit(("S3"+name1+w+"_St"+s).Data(),dbe_);
125  if (s!=4) {
126  hRes_S3RZWS[w][s-1] = new HRes1DHit(("S3"+name2+w+"_St"+s).Data(),dbe_,doall,local);
127  hEff_S1RZWS[w][s-1] = new HEff1DHit(("S1"+name2+w+"_St"+s).Data(),dbe_);
128  hEff_S3RZWS[w][s-1] = new HEff1DHit(("S3"+name2+w+"_St"+s).Data(),dbe_);
129  }
130  }
131  }
132  }
133 
134 
135  if(doall){
136  hEff_S3RPhi= new HEff1DHit("S3RPhi",dbe_); // RecHits, 3. step, RPhi
137  hEff_S3RZ= new HEff1DHit("S3RZ",dbe_); // RecHits, 3. step, RZ
138  hEff_S3RZ_W0= new HEff1DHit("S3RZ_W0",dbe_); // RecHits, 3. step, RZ, wheel 0
139  hEff_S3RZ_W1= new HEff1DHit("S3RZ_W1",dbe_); // RecHits, 3. step, RZ, wheel +-1
140  hEff_S3RZ_W2= new HEff1DHit("S3RZ_W2",dbe_); // RecHits, 3. step, RZ, wheel +-2
141  }
142  }
143 }
144 
145 
146 // Destructor
148 }
149 
150 
152  edm::EventSetup const& c){
153 
154 }
155 
157  // Write the histos to file
158  if(doall){
159  if(doStep1){
160  hEff_S1RPhi->ComputeEfficiency();
161  hEff_S1RZ->ComputeEfficiency();
162  hEff_S1RZ_W0->ComputeEfficiency();
163  hEff_S1RZ_W1->ComputeEfficiency();
164  hEff_S1RZ_W2->ComputeEfficiency();
165  }
166  if(doStep2){
167  hEff_S2RPhi->ComputeEfficiency();
168  hEff_S2RZ->ComputeEfficiency();
169  hEff_S2RZ_W0->ComputeEfficiency();
170  hEff_S2RZ_W1->ComputeEfficiency();
171  hEff_S2RZ_W2->ComputeEfficiency();
172  }
173  if(doStep3){
174  hEff_S3RPhi->ComputeEfficiency();
175  hEff_S3RZ->ComputeEfficiency();
176  hEff_S3RZ_W0->ComputeEfficiency();
177  hEff_S3RZ_W1->ComputeEfficiency();
178  hEff_S3RZ_W2->ComputeEfficiency();
179  }
180  }
181 }
182 
183 // The real analysis
184  void DTRecHitQuality::analyze(const Event & event, const EventSetup& eventSetup){
185  if(debug)
186  cout << "--- [DTRecHitQuality] Analysing Event: #Run: " << event.id().run()
187  << " #Event: " << event.id().event() << endl;
188  //theFile->cd();
189  // Get the DT Geometry
190  ESHandle<DTGeometry> dtGeom;
191  eventSetup.get<MuonGeometryRecord>().get(dtGeom);
192 
193  // Get the SimHit collection from the event
195  event.getByLabel(simHitLabel, simHits);
196 
197  // Map simhits per wire
198  map<DTWireId, PSimHitContainer > simHitsPerWire =
200 
201 
202 
203  //=======================================================================================
204  // RecHit analysis at Step 1
205  if(doStep1 && doall) {
206  if(debug)
207  cout << " -- DTRecHit S1: begin analysis:" << endl;
208  // Get the rechit collection from the event
209  Handle<DTRecHitCollection> dtRecHits;
210  event.getByLabel(recHitLabel, dtRecHits);
211 
212  if(!dtRecHits.isValid()) {
213  if(debug) cout << "[DTRecHitQuality]**Warning: no 1DRechits with label: " << recHitLabel << " in this event, skipping!" << endl;
214  return;
215  }
216 
217  // Map rechits per wire
218  map<DTWireId,vector<DTRecHit1DPair> > recHitsPerWire =
219  map1DRecHitsPerWire(dtRecHits.product());
220 
221  compute(dtGeom.product(), simHitsPerWire, recHitsPerWire, 1);
222  }
223 
224 
225  //=======================================================================================
226  // RecHit analysis at Step 2
227  if(doStep2 && doall) {
228  if(debug)
229  cout << " -- DTRecHit S2: begin analysis:" << endl;
230 
231  // Get the 2D rechits from the event
233  event.getByLabel(segment2DLabel, segment2Ds);
234 
235  if(!segment2Ds.isValid()) {
236  if(debug) cout << "[DTRecHitQuality]**Warning: no 2DSegments with label: " << segment2DLabel
237  << " in this event, skipping!" << endl;
238 
239  }
240  else{
241  // Map rechits per wire
242  map<DTWireId,vector<DTRecHit1D> > recHitsPerWire =
243  map1DRecHitsPerWire(segment2Ds.product());
244 
245  compute(dtGeom.product(), simHitsPerWire, recHitsPerWire, 2);
246  }
247  }
248 
249  //=======================================================================================
250  // RecHit analysis at Step 3
251  if(doStep3) {
252  if(debug)
253  cout << " -- DTRecHit S3: begin analysis:" << endl;
254 
255  // Get the 4D rechits from the event
257  event.getByLabel(segment4DLabel, segment4Ds);
258 
259  if(!segment4Ds.isValid()) {
260  if(debug) cout << "[DTRecHitQuality]**Warning: no 4D Segments with label: " << segment4DLabel
261  << " in this event, skipping!" << endl;
262  return;
263  }
264 
265  // Map rechits per wire
266  map<DTWireId,vector<DTRecHit1D> > recHitsPerWire =
267  map1DRecHitsPerWire(segment4Ds.product());
268 
269  compute(dtGeom.product(), simHitsPerWire, recHitsPerWire, 3);
270  }
271 
272  }
273 
274 
275 
276 // Return a map between DTRecHit1DPair and wireId
277 map<DTWireId, vector<DTRecHit1DPair> >
279  map<DTWireId, vector<DTRecHit1DPair> > ret;
280 
281  for(DTRecHitCollection::const_iterator rechit = dt1DRecHitPairs->begin();
282  rechit != dt1DRecHitPairs->end(); rechit++) {
283  ret[(*rechit).wireId()].push_back(*rechit);
284  }
285 
286  return ret;
287 }
288 
289 
290 // Return a map between DTRecHit1D at S2 and wireId
291 map<DTWireId, vector<DTRecHit1D> >
293  map<DTWireId, vector<DTRecHit1D> > ret;
294 
295  // Loop over all 2D segments
296  for(DTRecSegment2DCollection::const_iterator segment = segment2Ds->begin();
297  segment != segment2Ds->end();
298  segment++) {
299  vector<DTRecHit1D> component1DHits= (*segment).specificRecHits();
300  // Loop over all component 1D hits
301  for(vector<DTRecHit1D>::const_iterator hit = component1DHits.begin();
302  hit != component1DHits.end();
303  hit++) {
304  ret[(*hit).wireId()].push_back(*hit);
305  }
306  }
307  return ret;
308 }
309 
310 
311 
312 // Return a map between DTRecHit1D at S3 and wireId
313 map<DTWireId, std::vector<DTRecHit1D> >
315  map<DTWireId, vector<DTRecHit1D> > ret;
316  // Loop over all 4D segments
317  for(DTRecSegment4DCollection::const_iterator segment = segment4Ds->begin();
318  segment != segment4Ds->end();
319  segment++) {
320  // Get component 2D segments
321  vector<const TrackingRecHit*> segment2Ds = (*segment).recHits();
322  // Loop over 2D segments:
323  for(vector<const TrackingRecHit*>::const_iterator segment2D = segment2Ds.begin();
324  segment2D != segment2Ds.end();
325  segment2D++) {
326  // Get 1D component rechits
327  vector<const TrackingRecHit*> hits = (*segment2D)->recHits();
328  // Loop over them
329  for(vector<const TrackingRecHit*>::const_iterator hit = hits.begin();
330  hit != hits.end(); hit++) {
331  const DTRecHit1D* hit1D = dynamic_cast<const DTRecHit1D*>(*hit);
332  ret[hit1D->wireId()].push_back(*hit1D);
333  }
334  }
335  }
336 
337  return ret;
338 }
339 
340 // Compute SimHit distance from wire (cm)
342  DTWireId wireId,
343  const PSimHit& hit) {
344  float xwire = layer->specificTopology().wirePosition(wireId.wire());
345  LocalPoint entryP = hit.entryPoint();
346  LocalPoint exitP = hit.exitPoint();
347  float xEntry = entryP.x()-xwire;
348  float xExit = exitP.x()-xwire;
349 
350  return fabs(xEntry - (entryP.z()*(xExit-xEntry))/(exitP.z()-entryP.z()));//FIXME: check...
351 }
352 
353 // Compute SimHit impact angle (in direction perp to wire)
355  DTWireId wireId,
356  const PSimHit& hit) {
357  LocalPoint entryP = hit.entryPoint();
358  LocalPoint exitP = hit.exitPoint();
359  float theta=(exitP.x()-entryP.x())/(exitP.z()-entryP.z());
360  return atan(theta);
361 }
362 
363 // Compute SimHit distance from FrontEnd
365  DTWireId wireId,
366  const PSimHit& hit) {
367  LocalPoint entryP = hit.entryPoint();
368  LocalPoint exitP = hit.exitPoint();
369  float wireLenght=layer->specificTopology().cellLenght();
370  return (entryP.y()+exitP.y())/2.+wireLenght;
371 }
372 
373 
374 // Find the RecHit closest to the muon SimHit
375 template <typename type>
376 const type*
378  DTWireId wireId,
379  const vector<type>& recHits,
380  const float simHitDist) {
381  float res = 99999;
382  const type* theBestRecHit = 0;
383  // Loop over RecHits within the cell
384  for(typename vector<type>::const_iterator recHit = recHits.begin();
385  recHit != recHits.end();
386  recHit++) {
387  float distTmp = recHitDistFromWire(*recHit, layer);
388  if(fabs(distTmp-simHitDist) < res) {
389  res = fabs(distTmp-simHitDist);
390  theBestRecHit = &(*recHit);
391  }
392  } // End of loop over RecHits within the cell
393 
394  return theBestRecHit;
395 }
396 
397 
398 // Compute the distance from wire (cm) of a hits in a DTRecHit1DPair
399 float
401  // Compute the rechit distance from wire
402  return fabs(hitPair.localPosition(DTEnums::Left).x() -
403  hitPair.localPosition(DTEnums::Right).x())/2.;
404 }
405 
406 
407 
408 // Compute the distance from wire (cm) of a hits in a DTRecHit1D
409 float
411  return fabs(recHit.localPosition().x() - layer->specificTopology().wirePosition(recHit.wireId().wire()));
412 }
413 
414 
415 template <typename type>
417  const std::map<DTWireId, std::vector<PSimHit> >& simHitsPerWire,
418  const std::map<DTWireId, std::vector<type> >& recHitsPerWire,
419  int step) {
420  // Loop over cells with a muon SimHit
421  for(map<DTWireId, vector<PSimHit> >::const_iterator wireAndSHits = simHitsPerWire.begin();
422  wireAndSHits != simHitsPerWire.end();
423  wireAndSHits++) {
424  DTWireId wireId = (*wireAndSHits).first;
425  vector<PSimHit> simHitsInCell = (*wireAndSHits).second;
426 
427  // Get the layer
428  const DTLayer* layer = dtGeom->layer(wireId);
429 
430  // Look for a mu hit in the cell
431  const PSimHit* muSimHit = DTHitQualityUtils::findMuSimHit(simHitsInCell);
432  if (muSimHit==0) {
433  if (debug)
434  cout << " No mu SimHit in channel: " << wireId << ", skipping! " << endl;
435  continue; // Skip this cell
436  }
437 
438  // Find the distance of the simhit from the wire
439  float simHitWireDist = simHitDistFromWire(layer, wireId, *muSimHit);
440  // Skip simhits out of the cell
441  if(simHitWireDist>2.1) {
442  if(debug)
443  cout << " [DTRecHitQuality]###Warning: The mu SimHit in out of the cell, skipping!" << endl;
444  continue; // Skip this cell
445  }
446  GlobalPoint simHitGlobalPos = layer->toGlobal(muSimHit->localPosition());
447 
448  // find SH impact angle
449  float simHitTheta = simHitImpactAngle(layer, wireId, *muSimHit);
450 
451  // find SH distance from FE
452  float simHitFEDist = simHitDistFromFE(layer, wireId, *muSimHit);
453 
454  bool recHitReconstructed = false;
455 
456  // Look for RecHits in the same cell
457  if(recHitsPerWire.find(wireId) == recHitsPerWire.end()) {
458  // No RecHit found in this cell
459  if(debug)
460  cout << " No RecHit found at Step: " << step << " in cell: " << wireId << endl;
461  } else {
462  recHitReconstructed = true;
463  // vector<type> recHits = (*wireAndRecHits).second;
464  vector<type> recHits = recHitsPerWire.at(wireId);
465  if(debug)
466  cout << " " << recHits.size() << " RecHits, Step " << step << " in channel: " << wireId << endl;
467 
468  // Find the best RecHit
469  const type* theBestRecHit = findBestRecHit(layer, wireId, recHits, simHitWireDist);
470 
471 
472  float recHitWireDist = recHitDistFromWire(*theBestRecHit, layer);
473  if(debug)
474  cout << " SimHit distance from wire: " << simHitWireDist << endl
475  << " SimHit distance from FE: " << simHitFEDist << endl
476  << " SimHit distance angle " << simHitTheta << endl
477  << " RecHit distance from wire: " << recHitWireDist << endl;
478  float recHitErr = recHitPositionError(*theBestRecHit);
479  HRes1DHit *hRes = 0;
480  HRes1DHit *hResTot = 0;
481 
482  // Fill residuals and pulls
483  // Select the histo to be filled
484  if(step == 1) {
485  // Step 1
486  if(wireId.superLayer() != 2) {
487  hResTot = hRes_S1RPhi;
488  if(wireId.wheel() == 0)
489  hRes = hRes_S1RPhi_W0;
490  if(abs(wireId.wheel()) == 1)
491  hRes = hRes_S1RPhi_W1;
492  if(abs(wireId.wheel()) == 2)
493  hRes = hRes_S1RPhi_W2;
494  } else {
495  hResTot = hRes_S1RZ;
496  if(wireId.wheel() == 0)
497  hRes = hRes_S1RZ_W0;
498  if(abs(wireId.wheel()) == 1)
499  hRes = hRes_S1RZ_W1;
500  if(abs(wireId.wheel()) == 2)
501  hRes = hRes_S1RZ_W2;
502  }
503 
504  } else if(step == 2) {
505  // Step 2
506  if(wireId.superlayer() != 2) {
507  hRes = hRes_S2RPhi;
508  if(wireId.wheel() == 0)
509  hRes = hRes_S2RPhi_W0;
510  if(abs(wireId.wheel()) == 1)
511  hRes = hRes_S2RPhi_W1;
512  if(abs(wireId.wheel()) == 2)
513  hRes = hRes_S2RPhi_W2;
514  } else {
515  hResTot = hRes_S2RZ;
516  if(wireId.wheel() == 0)
517  hRes = hRes_S2RZ_W0;
518  if(abs(wireId.wheel()) == 1)
519  hRes = hRes_S2RZ_W1;
520  if(abs(wireId.wheel()) == 2)
521  hRes = hRes_S2RZ_W2;
522  }
523 
524  } else if(step == 3) {
525  // Step 3
526  if(wireId.superlayer() != 2) {
527  hResTot = hRes_S3RPhi;
528  if(wireId.wheel() == 0)
529  hRes = hRes_S3RPhi_W0;
530  if(abs(wireId.wheel()) == 1)
531  hRes = hRes_S3RPhi_W1;
532  if(abs(wireId.wheel()) == 2)
533  hRes = hRes_S3RPhi_W2;
534  if (local) hRes_S3RPhiWS[abs(wireId.wheel())][wireId.station()-1]->Fill(simHitWireDist, simHitTheta, simHitFEDist, recHitWireDist, simHitGlobalPos.eta(),simHitGlobalPos.phi(),recHitErr,wireId.station());
535 
536  } else {
537  hResTot = hRes_S3RZ;
538  if(wireId.wheel() == 0)
539  hRes = hRes_S3RZ_W0;
540  if(abs(wireId.wheel()) == 1)
541  hRes = hRes_S3RZ_W1;
542  if(abs(wireId.wheel()) == 2)
543  hRes = hRes_S3RZ_W2;
544 
545  if (local) hRes_S3RZWS[abs(wireId.wheel())][wireId.station()-1]->Fill(simHitWireDist, simHitTheta, simHitFEDist, recHitWireDist, simHitGlobalPos.eta(),simHitGlobalPos.phi(),recHitErr,wireId.station());
546  }
547  }
548  // Fill
549  hRes->Fill(simHitWireDist, simHitTheta, simHitFEDist, recHitWireDist, simHitGlobalPos.eta(),
550  simHitGlobalPos.phi(),recHitErr,wireId.station());
551  if(hResTot != 0)
552  hResTot->Fill(simHitWireDist, simHitTheta, simHitFEDist, recHitWireDist, simHitGlobalPos.eta(),
553  simHitGlobalPos.phi(),recHitErr,wireId.station());
554  }
555 
556  // Fill Efficiencies
557  if(doall){
558  HEff1DHit *hEff = 0;
559  HEff1DHit *hEffTot = 0;
560  if(step == 1) {
561  // Step 1
562  if(wireId.superlayer() != 2) {
563  hEff = hEff_S1RPhi;
564  if (local) hEff_S1RPhiWS[abs(wireId.wheel())][wireId.station()-1]->Fill(simHitWireDist, simHitGlobalPos.eta(), simHitGlobalPos.phi(), recHitReconstructed);
565  } else {
566  hEffTot = hEff_S1RZ;
567  if(wireId.wheel() == 0)
568  hEff = hEff_S1RZ_W0;
569  if(abs(wireId.wheel()) == 1)
570  hEff = hEff_S1RZ_W1;
571  if(abs(wireId.wheel()) == 2)
572  hEff = hEff_S1RZ_W2;
573  if (local) hEff_S1RZWS[abs(wireId.wheel())][wireId.station()-1]->Fill(simHitWireDist, simHitGlobalPos.eta(), simHitGlobalPos.phi(), recHitReconstructed);
574  }
575 
576  } else if(step == 2) {
577  // Step 2
578  if(wireId.superlayer() != 2) {
579  hEff = hEff_S2RPhi;
580  } else {
581  hEffTot = hEff_S2RZ;
582  if(wireId.wheel() == 0)
583  hEff = hEff_S2RZ_W0;
584  if(abs(wireId.wheel()) == 1)
585  hEff = hEff_S2RZ_W1;
586  if(abs(wireId.wheel()) == 2)
587  hEff = hEff_S2RZ_W2;
588  }
589 
590  } else if(step == 3) {
591  // Step 3
592  if(wireId.superlayer() != 2) {
593  hEff = hEff_S3RPhi;
594  if (local) hEff_S3RPhiWS[abs(wireId.wheel())][wireId.station()-1]->Fill(simHitWireDist, simHitGlobalPos.eta(), simHitGlobalPos.phi(), recHitReconstructed);
595  } else {
596  hEffTot = hEff_S3RZ;
597  if(wireId.wheel() == 0)
598  hEff = hEff_S3RZ_W0;
599  if(abs(wireId.wheel()) == 1)
600  hEff = hEff_S3RZ_W1;
601  if(abs(wireId.wheel()) == 2)
602  hEff = hEff_S3RZ_W2;
603  if (local) hEff_S3RZWS[abs(wireId.wheel())][wireId.station()-1]->Fill(simHitWireDist, simHitGlobalPos.eta(), simHitGlobalPos.phi(), recHitReconstructed);
604  }
605 
606  }
607  // Fill
608  hEff->Fill(simHitWireDist, simHitGlobalPos.eta(), simHitGlobalPos.phi(), recHitReconstructed);
609  if(hEffTot != 0)
610  hEffTot->Fill(simHitWireDist, simHitGlobalPos.eta(), simHitGlobalPos.phi(), recHitReconstructed);
611  }
612  }
613 }
614 
615 // Return the error on the measured (cm) coordinate
617  return sqrt(recHit.localPositionError(DTEnums::Left).xx());
618 }
619 
620 // Return the error on the measured (cm) coordinate
622  return sqrt(recHit.localPositionError().xx());
623 }
624 
type
Definition: HCALResponse.h:21
T getUntrackedParameter(std::string const &, T const &) const
float wirePosition(int wireNumber) const
Returns the x position in the layer of a given wire number.
Definition: DTTopology.cc:86
void Fill(float distSimHit, float etaSimHit, float phiSimHit, bool fillRecHit)
Definition: Histograms.h:217
float xx() const
Definition: LocalError.h:24
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:47
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
void Fill(float distSimHit, float thetaSimHit, float distFESimHit, float distRecHit, float etaSimHit, float phiSimHit, float errRecHit, int station)
Definition: Histograms.h:95
Geom::Theta< T > theta() const
T y() const
Definition: PV3DBase.h:63
static const PSimHit * findMuSimHit(const edm::PSimHitContainer &hits)
Select the SimHit from a muon in a vector of SimHits.
void compute(const DTGeometry *dtGeom, const std::map< DTWireId, std::vector< PSimHit > > &simHitsPerWire, const std::map< DTWireId, std::vector< type > > &recHitsPerWire, int step)
const type * findBestRecHit(const DTLayer *layer, DTWireId wireId, const std::vector< type > &recHits, const float simHitDist)
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup)
Perform the real analysis.
virtual ~DTRecHitQuality()
Destructor.
float recHitPositionError(const DTRecHit1DPair &recHit)
virtual LocalError localPositionError() const
const DTLayer * layer(DTLayerId id) const
Return a layer given its id.
Definition: DTGeometry.cc:110
const DTTopology & specificTopology() const
Definition: DTLayer.cc:42
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:38
float simHitDistFromWire(const DTLayer *layer, DTWireId wireId, const PSimHit &hit)
DTRecHitQuality(const edm::ParameterSet &pset)
Constructor.
Local3DPoint localPosition() const
Definition: PSimHit.h:44
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
float recHitDistFromWire(const DTRecHit1DPair &hitPair, const DTLayer *layer)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int superLayer() const
Return the superlayer number.
std::map< DTWireId, std::vector< DTRecHit1DPair > > map1DRecHitsPerWire(const DTRecHitCollection *dt1DRecHitPairs)
A set of histograms of residuals and pulls for 1D RecHits.
Definition: Histograms.h:26
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
bool isValid() const
Definition: HandleBase.h:76
virtual LocalError localPositionError() const
Return the 3-dimensional error on the local position.
Definition: DTRecHit1D.h:66
float simHitImpactAngle(const DTLayer *layer, DTWireId wireId, const PSimHit &hit)
DQMStore * dbe_
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)
#define debug
Definition: HDRShower.cc:19
float simHitDistFromFE(const DTLayer *layer, DTWireId wireId, const PSimHit &hit)
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
tuple simHits
Definition: trackerHits.py:16
T const * product() const
Definition: Handle.h:81
T eta() const
Definition: PV3DBase.h:76
tuple cout
Definition: gather_cfg.py:121
T w() const
static std::map< DTWireId, edm::PSimHitContainer > mapSimHitsPerWire(const edm::PSimHitContainer &simhits)
int station() const
Return the station number.
Definition: DTChamberId.h:51
float cellLenght() const
Definition: DTTopology.h:73
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
T x() const
Definition: PV3DBase.h:62
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:35
void endLuminosityBlock(edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &c)
virtual LocalPoint localPosition() const
DTWireId wireId() const
Return the wireId.
Definition: DTRecHit1D.h:107