CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTSegment2DQuality.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * $Date: 2009/11/04 17:22:32 $
5  * $Revision: 1.8 $
6  * \author S. Bolognesi and G. Cerminara - INFN Torino
7  */
8 
9 #include "DTSegment2DQuality.h"
11 
16 
20 
25 
26 #include "Histograms.h"
27 
28 #include "TFile.h"
29 
32 #include <iostream>
33 #include <map>
34 
35 using namespace std;
36 using namespace edm;
37 
38 // Constructor
40  // Get the debug parameter for verbose output
41  debug = pset.getUntrackedParameter<bool>("debug");
43  // the name of the simhit collection
44  simHitLabel = pset.getUntrackedParameter<InputTag>("simHitLabel");
45  // the name of the 2D rec hit collection
46  segment2DLabel = pset.getUntrackedParameter<InputTag>("segment2DLabel");
47 
48  //sigma resolution on position
49  sigmaResPos = pset.getParameter<double>("sigmaResPos");
50  //sigma resolution on angle
51  sigmaResAngle = pset.getParameter<double>("sigmaResAngle");
52 
53  // Create the root file
54  //theFile = new TFile(rootFileName.c_str(), "RECREATE");
55  //theFile->cd();
56 
57  if(debug)
58  cout << "[DTSegment2DQuality] Constructor called " << endl;
59 
60  // ----------------------
61  // get hold of back-end interface
62  dbe_ = 0;
64  if ( dbe_ ) {
65  if (debug) {
66  dbe_->setVerbose(1);
67  } else {
68  dbe_->setVerbose(0);
69  }
70  }
71  if ( dbe_ ) {
72  if ( debug ) dbe_->showDirStructure();
73  }
74 
75  h2DHitRPhi = new HRes2DHit ("RPhi",dbe_,true,true);
76  h2DHitRZ= new HRes2DHit ("RZ",dbe_,true,true);
77  h2DHitRZ_W0= new HRes2DHit ("RZ_W0",dbe_,true,true);
78  h2DHitRZ_W1= new HRes2DHit ("RZ_W1",dbe_,true,true);
79  h2DHitRZ_W2= new HRes2DHit ("RZ_W2",dbe_,true,true);
80 
81  h2DHitEff_RPhi= new HEff2DHit ("RPhi",dbe_);
82  h2DHitEff_RZ= new HEff2DHit ("RZ",dbe_);
83  h2DHitEff_RZ_W0= new HEff2DHit ("RZ_W0",dbe_);
84  h2DHitEff_RZ_W1= new HEff2DHit ("RZ_W1",dbe_);
85  h2DHitEff_RZ_W2= new HEff2DHit ("RZ_W2",dbe_);
86  if(debug)
87  cout << "[DTSegment2DQuality] hitsos created " << endl;
88 }
89 
90 // Destructor
92 
93 }
94 
96  // Write the histos to file
97  //theFile->cd();
98 
99  //h2DHitRPhi->Write();
100  //h2DHitRZ->Write();
101  //h2DHitRZ_W0->Write();
102  //h2DHitRZ_W1->Write();
103  //h2DHitRZ_W2->Write();
104 
105  h2DHitEff_RPhi->ComputeEfficiency();
106  h2DHitEff_RZ->ComputeEfficiency();
107  h2DHitEff_RZ_W0->ComputeEfficiency();
108  h2DHitEff_RZ_W1->ComputeEfficiency();
109  h2DHitEff_RZ_W2->ComputeEfficiency();
110 
111  //h2DHitEff_RPhi->Write();
112  //h2DHitEff_RZ->Write();
113  //h2DHitEff_RZ_W0->Write();
114  //h2DHitEff_RZ_W1->Write();
115  //h2DHitEff_RZ_W2->Write();
116  //if ( rootFileName.size() != 0 && dbe_ ) dbe_->save(rootFileName);
117 
118  //theFile->Close();
119 }
120 
121 // The real analysis
122 void DTSegment2DQuality::analyze(const Event & event, const EventSetup& eventSetup){
123  //theFile->cd();
124 
125  // Get the DT Geometry
126  ESHandle<DTGeometry> dtGeom;
127  eventSetup.get<MuonGeometryRecord>().get(dtGeom);
128 
129  // Get the SimHit collection from the event
131  event.getByLabel(simHitLabel, simHits); //FIXME: second string to be removed
132 
133  //Map simHits by sl
134  map<DTSuperLayerId, PSimHitContainer > simHitsPerSl;
135  for(PSimHitContainer::const_iterator simHit = simHits->begin();
136  simHit != simHits->end(); simHit++){
137  // Create the id of the sl (the simHits in the DT known their wireId)
138  DTSuperLayerId slId = ((DTWireId(simHit->detUnitId())).layerId()).superlayerId();
139  // Fill the map
140  simHitsPerSl[slId].push_back(*simHit);
141  }
142 
143  // Get the 2D rechits from the event
145  event.getByLabel(segment2DLabel, segment2Ds);
146 
147  if(!segment2Ds.isValid()) {
148  if(debug) cout << "[DTSegment2DQuality]**Warning: no 2DSegments with label: " << segment2DLabel
149  << " in this event, skipping!" << endl;
150  return;
151  }
152 
153  // Loop over all superlayers containing a segment
155  for (slId = segment2Ds->id_begin();
156  slId != segment2Ds->id_end();
157  ++slId){
158 
159  //------------------------- simHits ---------------------------//
160  //Get simHits of each superlayer
161  PSimHitContainer simHits = simHitsPerSl[(*slId)];
162 
163  // Map simhits per wire
164  map<DTWireId, PSimHitContainer > simHitsPerWire = DTHitQualityUtils::mapSimHitsPerWire(simHits);
165  map<DTWireId, const PSimHit*> muSimHitPerWire = DTHitQualityUtils::mapMuSimHitsPerWire(simHitsPerWire);
166  int nMuSimHit = muSimHitPerWire.size();
167  if(nMuSimHit == 0 || nMuSimHit == 1) {
168  if(debug && nMuSimHit == 1)
169  cout << "[DTSegment2DQuality] Only " << nMuSimHit << " mu SimHit in this SL, skipping!" << endl;
170  continue; // If no or only one mu SimHit is found skip this SL
171  }
172  if(debug)
173  cout << "=== SL " << (*slId) << " has " << nMuSimHit << " SimHits" << endl;
174 
175  //Find outer and inner mu SimHit to build a segment
176  pair<const PSimHit*, const PSimHit*> inAndOutSimHit = DTHitQualityUtils::findMuSimSegment(muSimHitPerWire);
177  //Check that outermost and innermost SimHit are not the same
178  if(inAndOutSimHit.first ==inAndOutSimHit.second ) {
179  cout << "[DTHitQualityUtils]***Warning: outermost and innermost SimHit are the same!" << endl;
180  continue;
181  }
182 
183  //Find direction and position of the sim Segment in SL RF
184  pair<LocalVector, LocalPoint> dirAndPosSimSegm = DTHitQualityUtils::findMuSimSegmentDirAndPos(inAndOutSimHit,
185  (*slId),&(*dtGeom));
186 
187  LocalVector simSegmLocalDir = dirAndPosSimSegm.first;
188  LocalPoint simSegmLocalPos = dirAndPosSimSegm.second;
189  if(debug)
190  cout<<" Simulated segment: local direction "<<simSegmLocalDir<<endl
191  <<" local position "<<simSegmLocalPos<<endl;
192  const DTSuperLayer* superLayer = dtGeom->superLayer(*slId);
193  GlobalPoint simSegmGlobalPos = superLayer->toGlobal(simSegmLocalPos);
194 
195  //Atan(x/z) angle and x position in SL RF
196  float angleSimSeg = DTHitQualityUtils::findSegmentAlphaAndBeta(simSegmLocalDir).first;
197  float posSimSeg = simSegmLocalPos.x();
198  //Position (in eta,phi coordinates) in the global RF
199  float etaSimSeg = simSegmGlobalPos.eta();
200  float phiSimSeg = simSegmGlobalPos.phi();
201 
202 
203  //---------------------------- recHits --------------------------//
204  // Get the range of rechit for the corresponding slId
205  bool recHitFound = false;
206  DTRecSegment2DCollection::range range = segment2Ds->get(*slId);
207  int nsegm = distance(range.first, range.second);
208  if(debug)
209  cout << " Sl: " << *slId << " has " << nsegm
210  << " 2D segments" << endl;
211 
212  if (nsegm!=0) {
213  // Find the best RecHit: look for the 2D RecHit with the angle closest
214  // to that of segment made of SimHits.
215  // RecHits must have delta alpha and delta position within 5 sigma of
216  // the residual distribution (we are looking for residuals of segments
217  // usefull to the track fit) for efficency purpose
218  const DTRecSegment2D* bestRecHit = 0;
219  bool bestRecHitFound = false;
220  double deltaAlpha = 99999;
221 
222  // Loop over the recHits of this slId
223  for (DTRecSegment2DCollection::const_iterator segment2D = range.first;
224  segment2D!=range.second;
225  ++segment2D){
226  // Check the dimension
227  if((*segment2D).dimension() != 2) {
228  if(debug) cout << "[DTSegment2DQuality]***Error: This is not 2D segment!!!" << endl;
229  abort();
230  }
231  // Segment Local Direction and position (in SL RF)
232  LocalVector recSegDirection = (*segment2D).localDirection();
233  LocalPoint recSegPosition = (*segment2D).localPosition();
234 
235  float recSegAlpha = DTHitQualityUtils::findSegmentAlphaAndBeta(recSegDirection).first;
236  if(debug)
237  cout << " RecSegment direction: " << recSegDirection << endl
238  << " position : " << recSegPosition << endl
239  << " alpha : " << recSegAlpha << endl;
240 
241  if(fabs(recSegAlpha - angleSimSeg) < deltaAlpha) {
242  deltaAlpha = fabs(recSegAlpha - angleSimSeg);
243  bestRecHit = &(*segment2D);
244  bestRecHitFound = true;
245  }
246  } // End of Loop over all 2D RecHits
247 
248  if(bestRecHitFound) {
249  // Best rechit direction and position in SL RF
250  LocalPoint bestRecHitLocalPos = bestRecHit->localPosition();
251  LocalVector bestRecHitLocalDir = bestRecHit->localDirection();
252 
253  LocalError bestRecHitLocalPosErr = bestRecHit->localPositionError();
254  LocalError bestRecHitLocalDirErr = bestRecHit->localDirectionError();
255 
256  float angleBestRHit = DTHitQualityUtils::findSegmentAlphaAndBeta(bestRecHitLocalDir).first;
257 
258  if(fabs(angleBestRHit - angleSimSeg) < 5*sigmaResAngle &&
259  fabs(bestRecHitLocalPos.x() - posSimSeg) < 5*sigmaResPos) {
260  recHitFound = true;
261  }
262 
263  // Fill Residual histos
264  HRes2DHit *hRes = 0;
265  if((*slId).superlayer() == 1 || (*slId).superlayer() == 3) { //RPhi SL
266  hRes = h2DHitRPhi;
267  } else if((*slId).superlayer() == 2) { //RZ SL
268  h2DHitRZ->Fill(angleSimSeg,
269  angleBestRHit,
270  posSimSeg,
271  bestRecHitLocalPos.x(),
272  etaSimSeg,
273  phiSimSeg,
274  sqrt(bestRecHitLocalPosErr.xx()),
275  sqrt(bestRecHitLocalDirErr.xx())
276  );
277  if(abs((*slId).wheel()) == 0)
278  hRes = h2DHitRZ_W0;
279  else if(abs((*slId).wheel()) == 1)
280  hRes = h2DHitRZ_W1;
281  else if(abs((*slId).wheel()) == 2)
282  hRes = h2DHitRZ_W2;
283  }
284  hRes->Fill(angleSimSeg,
285  angleBestRHit,
286  posSimSeg,
287  bestRecHitLocalPos.x(),
288  etaSimSeg,
289  phiSimSeg,
290  sqrt(bestRecHitLocalPosErr.xx()),
291  sqrt(bestRecHitLocalDirErr.xx())
292  );
293  }
294  } //end of if(nsegm!=0)
295 
296  // Fill Efficiency plot
297  HEff2DHit *hEff = 0;
298  if((*slId).superlayer() == 1 || (*slId).superlayer() == 3) { //RPhi SL
299  hEff = h2DHitEff_RPhi;
300  } else if((*slId).superlayer() == 2) { //RZ SL
301  h2DHitEff_RZ->Fill(etaSimSeg, phiSimSeg, posSimSeg, angleSimSeg, recHitFound);
302  if(abs((*slId).wheel()) == 0)
303  hEff = h2DHitEff_RZ_W0;
304  else if(abs((*slId).wheel()) == 1)
305  hEff = h2DHitEff_RZ_W1;
306  else if(abs((*slId).wheel()) == 2)
307  hEff = h2DHitEff_RZ_W2;
308  }
309  hEff->Fill(etaSimSeg, phiSimSeg, posSimSeg, angleSimSeg, recHitFound);
310  } // End of loop over superlayers
311 }
312 
313 
314 
315 
316 
317 
T getParameter(std::string const &) const
virtual LocalError localPositionError() const
local position error in SL frame
T getUntrackedParameter(std::string const &, T const &) const
virtual ~DTSegment2DQuality()
Destructor.
float xx() const
Definition: LocalError.h:24
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:52
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
#define abs(x)
Definition: mlp_lapack.h:159
identifier iterator
Definition: RangeMap.h:138
virtual LocalError localDirectionError() const
the local direction error (xx,xy,yy) in SL frame: only xx is not 0.
DTSegment2DQuality(const edm::ParameterSet &pset)
Constructor.
static std::map< DTWireId, const PSimHit * > mapMuSimHitsPerWire(const std::map< DTWireId, edm::PSimHitContainer > &simHitWireMap)
Create a map between the Mu SimHits and corresponding MuBarWireId ;.
static std::pair< const PSimHit *, const PSimHit * > findMuSimSegment(const std::map< DTWireId, const PSimHit * > &mapWireAndMuSimHit)
Find Innermost and outermost SimHit from Mu in a SL (they identify a simulated segment) ...
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:45
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup)
Perform the real analysis.
T sqrt(T t)
Definition: SSEVec.h:48
void Fill(float angleSimSegment, float angleRecSegment, float posSimSegment, float posRecSegment, float etaSimSegment, float phiSimSegment, float sigmaPos, float sigmaAngle)
Definition: Histograms.h:391
void setVerbose(unsigned level)
Definition: DQMStore.cc:398
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
static std::pair< double, double > findSegmentAlphaAndBeta(const LocalVector &direction)
Find the angles from a segment direction:
DQMStore * dbe_
virtual LocalPoint localPosition() const
local position in SL frame
static std::pair< LocalVector, LocalPoint > findMuSimSegmentDirAndPos(const std::pair< const PSimHit *, const PSimHit * > &inAndOutSimHit, const DetId detId, const DTGeometry *muonGeom)
Find direction and position of a segment (in local RF) from outer and inner mu SimHit in the RF of ob...
void Fill(float etaSimSegm, float phiSimSegm, float posSimSegm, float angleSimSegm, bool fillRecHit)
Definition: Histograms.h:532
const T & get() const
Definition: EventSetup.h:55
tuple simHits
Definition: trackerHits.py:16
T eta() const
Definition: PV3DBase.h:76
virtual LocalVector localDirection() const
the local direction in SL frame
tuple cout
Definition: gather_cfg.py:121
std::vector< PSimHit > PSimHitContainer
static std::map< DTWireId, edm::PSimHitContainer > mapSimHitsPerWire(const edm::PSimHitContainer &simhits)
void showDirStructure(void) const
Definition: DQMStore.cc:2766
#define debug
Definition: MEtoEDMFormat.h:34
T x() const
Definition: PV3DBase.h:62