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