CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTSegment2DSLPhiQuality.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * $Date: 2007/10/25 11:58:37 $
5  * $Revision: 1.1 $
6  * \author S. Bolognesi and G. Cerminara - INFN Torino
7  */
8 
11 
16 
20 
25 
26 #include "Histograms.h"
27 
28 #include "TFile.h"
29 
30 #include <iostream>
31 #include <map>
32 
33 
34 using namespace std;
35 using namespace edm;
36 
37 // Constructor
39  // Get the debug parameter for verbose output
40  debug = pset.getUntrackedParameter<bool>("debug");
42 
43  rootFileName = pset.getUntrackedParameter<string>("rootFileName");
44  // the name of the simhit collection
45  simHitLabel = pset.getUntrackedParameter<string>("simHitLabel", "SimG4Object");
46  // the name of the 4D rec hit collection
47  segment4DLabel = pset.getUntrackedParameter<string>("segment4DLabel");
48 
49  //sigma resolution on position
50  sigmaResPos = pset.getParameter<double>("sigmaResPos");
51  //sigma resolution on angle
52  sigmaResAngle = pset.getParameter<double>("sigmaResAngle");
53 
54  // Create the root file
55  theFile = new TFile(rootFileName.c_str(), "RECREATE");
56  theFile->cd();
57 
58  // Book the histos
59  h2DHitSuperPhi = new HRes2DHit ("SuperPhi");
60  h2DHitEff_SuperPhi = new HEff2DHit ("SuperPhi");
61 }
62 
63 // Destructor
65 }
66 
68  // Write the histos to file
69  theFile->cd();
70 
71  h2DHitSuperPhi->Write();
72 
73  h2DHitEff_SuperPhi->ComputeEfficiency();
74  h2DHitEff_SuperPhi->Write();
75 
76  theFile->Close();
77 }
78 
79 // The real analysis
80 void DTSegment2DSLPhiQuality::analyze(const Event & event, const EventSetup& eventSetup){
81  theFile->cd();
82 
83  // Get the DT Geometry
84  ESHandle<DTGeometry> dtGeom;
85  eventSetup.get<MuonGeometryRecord>().get(dtGeom);
86 
87  // Get the SimHit collection from the event
89  event.getByLabel(simHitLabel, "MuonDTHits", simHits); //FIXME: second string to be removed
90 
91  //Map simHits by chamber
92  map<DTChamberId, PSimHitContainer > simHitsPerCh;
93  for(PSimHitContainer::const_iterator simHit = simHits->begin();
94  simHit != simHits->end(); simHit++){
95  // Create the id of the chamber (the simHits in the DT known their wireId)
96  DTChamberId chamberId = (((DTWireId(simHit->detUnitId())).layerId()).superlayerId()).chamberId();
97  // Fill the map
98  simHitsPerCh[chamberId].push_back(*simHit);
99  }
100 
101  // Get the 4D rechits from the event
103  event.getByLabel(segment4DLabel, segment4Ds);
104 
105  // Loop over all chambers containing a segment
107  for (chamberId = segment4Ds->id_begin();
108  chamberId != segment4Ds->id_end();
109  ++chamberId){
110 
111  //------------------------- simHits ---------------------------//
112  //Get simHits of each chamber
113  PSimHitContainer simHits = simHitsPerCh[(*chamberId)];
114 
115  // Map simhits per wire
116  map<DTWireId, PSimHitContainer > simHitsPerWire = DTHitQualityUtils::mapSimHitsPerWire(simHits);
117  map<DTWireId, const PSimHit*> muSimHitPerWire = DTHitQualityUtils::mapMuSimHitsPerWire(simHitsPerWire);
118  int nMuSimHit = muSimHitPerWire.size();
119  if(nMuSimHit == 0 || nMuSimHit == 1) {
120  if(debug && nMuSimHit == 1)
121  cout << "[DTSegment2DSLPhiQuality] Only " << nMuSimHit << " mu SimHit in this chamber, skipping!" << endl;
122  continue; // If no or only one mu SimHit is found skip this chamber
123  }
124  if(debug)
125  cout << "=== Chamber " << (*chamberId) << " has " << nMuSimHit << " SimHits" << endl;
126 
127  //Find outer and inner mu SimHit to build a segment
128  pair<const PSimHit*, const PSimHit*> inAndOutSimHit = DTHitQualityUtils::findMuSimSegment(muSimHitPerWire);
129 
130  //Find direction and position of the sim Segment in Chamber RF
131  pair<LocalVector, LocalPoint> dirAndPosSimSegm = DTHitQualityUtils::findMuSimSegmentDirAndPos(inAndOutSimHit,
132  (*chamberId),&(*dtGeom));
133 
134  LocalVector simSegmLocalDir = dirAndPosSimSegm.first;
135  LocalPoint simSegmLocalPos = dirAndPosSimSegm.second;
136  const DTChamber* chamber = dtGeom->chamber(*chamberId);
137  GlobalPoint simSegmGlobalPos = chamber->toGlobal(simSegmLocalPos);
138 
139  //Atan(x/z) angle and x position in SL RF
140  float angleSimSeg = DTHitQualityUtils::findSegmentAlphaAndBeta(simSegmLocalDir).first;
141  float posSimSeg = simSegmLocalPos.x();
142  //Position (in eta,phi coordinates) in lobal RF
143  float etaSimSeg = simSegmGlobalPos.eta();
144  float phiSimSeg = simSegmGlobalPos.phi();
145 
146  if(debug)
147  cout<<" Simulated segment: local direction "<<simSegmLocalDir<<endl
148  <<" local position "<<simSegmLocalPos<<endl
149  <<" angle "<<angleSimSeg<<endl;
150 
151  //---------------------------- recHits --------------------------//
152  // Get the range of rechit for the corresponding chamberId
153  bool recHitFound = false;
154  DTRecSegment4DCollection::range range = segment4Ds->get(*chamberId);
155  int nsegm = distance(range.first, range.second);
156  if(debug)
157  cout << " Chamber: " << *chamberId << " has " << nsegm
158  << " 4D segments" << endl;
159 
160  if (nsegm!=0) {
161  // Find the best RecHit: look for the 4D RecHit with the phi angle closest
162  // to that of segment made of SimHits.
163  // RecHits must have delta alpha and delta position within 5 sigma of
164  // the residual distribution (we are looking for residuals of segments
165  // usefull to the track fit) for efficency purpose
166  const DTRecSegment2D* bestRecHit = 0;
167  bool bestRecHitFound = false;
168  double deltaAlpha = 99999;
169 
170  // Loop over the recHits of this chamberId
171  for (DTRecSegment4DCollection::const_iterator segment4D = range.first;
172  segment4D!=range.second;
173  ++segment4D){
174  // Check the dimension
175  if((*segment4D).dimension() != 4) {
176  if(debug) cout << "[DTSegment2DSLPhiQuality]***Error: This is not 4D segment!!!" << endl;
177  continue;
178  }
179 
180  //Get 2D superPhi segments from 4D segments
181  const DTChamberRecSegment2D* phiSegment2D = (*segment4D).phiSegment();
182  if((*phiSegment2D).dimension() != 2) {
183  if(debug) cout << "[DTSegment2DQuality]***Error: This is not 2D segment!!!" << endl;
184  abort();
185  }
186 
187  // Segment Local Direction and position (in Chamber RF)
188  LocalVector recSegDirection = (*phiSegment2D).localDirection();
189 
190  float recSegAlpha = DTHitQualityUtils::findSegmentAlphaAndBeta(recSegDirection).first;
191  if(debug)
192  cout << " RecSegment direction: " << recSegDirection << endl
193  << " position : " << (*phiSegment2D).localPosition() << endl
194  << " alpha : " << recSegAlpha << endl;
195 
196  if(fabs(recSegAlpha - angleSimSeg) < deltaAlpha) {
197  deltaAlpha = fabs(recSegAlpha - angleSimSeg);
198  bestRecHit = &(*phiSegment2D);
199  bestRecHitFound = true;
200  }
201  } // End of Loop over all 4D RecHits of this chambers
202 
203  if(bestRecHitFound) {
204  // Best rechit direction and position in Chamber RF
205  LocalPoint bestRecHitLocalPos = bestRecHit->localPosition();
206  LocalVector bestRecHitLocalDir = bestRecHit->localDirection();
207 
208  LocalError bestRecHitLocalPosErr = bestRecHit->localPositionError();
209  LocalError bestRecHitLocalDirErr = bestRecHit->localDirectionError();
210 
211  float angleBestRHit = DTHitQualityUtils::findSegmentAlphaAndBeta(bestRecHitLocalDir).first;
212  if(fabs(angleBestRHit - angleSimSeg) < 5*sigmaResAngle &&
213  fabs(bestRecHitLocalPos.x() - posSimSeg) < 5*sigmaResPos) {
214  recHitFound = true;
215  }
216 
217  // Fill Residual histos
218  h2DHitSuperPhi->Fill(angleSimSeg,
219  angleBestRHit,
220  posSimSeg,
221  bestRecHitLocalPos.x(),
222  etaSimSeg,
223  phiSimSeg,
224  sqrt(bestRecHitLocalPosErr.xx()),
225  sqrt(bestRecHitLocalDirErr.xx())
226  );
227  }
228  } //end of if(nsegm!=0)
229 
230  // Fill Efficiency plot
231  h2DHitEff_SuperPhi->Fill(etaSimSeg,
232  phiSimSeg,
233  posSimSeg,
234  angleSimSeg,
235  recHitFound);
236  } // End of loop over chambers
237 }
238 
239 
240 
241 
242 
243 
T getParameter(std::string const &) const
virtual LocalError localPositionError() const
local position error in SL frame
DTSegment2DSLPhiQuality(const edm::ParameterSet &pset)
Constructor.
T getUntrackedParameter(std::string const &, T const &) const
float xx() const
Definition: LocalError.h:19
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:52
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup)
Perform the real analysis.
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
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.
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
T sqrt(T t)
Definition: SSEVec.h:28
tuple pset
Definition: CrabTask.py:85
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...
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
virtual ~DTSegment2DSLPhiQuality()
Destructor.