CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
DTSegment2DSLPhiQuality Class Reference

#include <DTSegment2DSLPhiQuality.h>

Inheritance diagram for DTSegment2DSLPhiQuality:
edm::EDAnalyzer

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &eventSetup)
 Perform the real analysis. More...
 
 DTSegment2DSLPhiQuality (const edm::ParameterSet &pset)
 Constructor. More...
 
void endJob ()
 
virtual ~DTSegment2DSLPhiQuality ()
 Destructor. More...
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Attributes

bool debug
 
HEff2DHith2DHitEff_SuperPhi
 
HRes2DHith2DHitSuperPhi
 
std::string rootFileName
 
std::string segment4DLabel
 
double sigmaResAngle
 
double sigmaResPos
 
std::string simHitLabel
 
TFile * theFile
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Basic analyzer class which accesses 2D DTSegments reconstructed with both SL Phi and plot resolution comparing reconstructed and simulated quantities

Date:
2007/06/08 15:50:31
Revision:
1.2
Author
S. Bolognesi and G. Cerminara - INFN Torino

Definition at line 28 of file DTSegment2DSLPhiQuality.h.

Constructor & Destructor Documentation

DTSegment2DSLPhiQuality::DTSegment2DSLPhiQuality ( const edm::ParameterSet pset)

Constructor.

Definition at line 38 of file DTSegment2DSLPhiQuality.cc.

References debug, DTHitQualityUtils::debug, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), dtT0WireCalibration_cfg::rootFileName, and interactiveExample::theFile.

38  {
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 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
DTSegment2DSLPhiQuality::~DTSegment2DSLPhiQuality ( )
virtual

Destructor.

Definition at line 64 of file DTSegment2DSLPhiQuality.cc.

64  {
65 }

Member Function Documentation

void DTSegment2DSLPhiQuality::analyze ( const edm::Event event,
const edm::EventSetup eventSetup 
)
virtual

Perform the real analysis.

Implements edm::EDAnalyzer.

Definition at line 80 of file DTSegment2DSLPhiQuality.cc.

References gather_cfg::cout, debug, PV3DBase< T, PVType, FrameType >::eta(), DTHitQualityUtils::findMuSimSegment(), DTHitQualityUtils::findMuSimSegmentDirAndPos(), DTHitQualityUtils::findSegmentAlphaAndBeta(), edm::EventSetup::get(), DTRecSegment2D::localDirection(), DTRecSegment2D::localDirectionError(), DTRecSegment2D::localPosition(), DTRecSegment2D::localPositionError(), DTHitQualityUtils::mapMuSimHitsPerWire(), DTHitQualityUtils::mapSimHitsPerWire(), PV3DBase< T, PVType, FrameType >::phi(), trackerHits::simHits, mathSSE::sqrt(), interactiveExample::theFile, GeomDet::toGlobal(), PV3DBase< T, PVType, FrameType >::x(), and LocalError::xx().

80  {
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 }
virtual LocalError localPositionError() const
local position error in SL frame
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
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.
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
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
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)
T x() const
Definition: PV3DBase.h:56
void DTSegment2DSLPhiQuality::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 67 of file DTSegment2DSLPhiQuality.cc.

References interactiveExample::theFile.

67  {
68  // Write the histos to file
69  theFile->cd();
70 
72 
75 
76  theFile->Close();
77 }
void ComputeEfficiency()
Definition: Histograms.h:523
void Write()
Definition: Histograms.h:591
void Write()
Definition: Histograms.h:390

Member Data Documentation

bool DTSegment2DSLPhiQuality::debug
private

Definition at line 50 of file DTSegment2DSLPhiQuality.h.

HEff2DHit* DTSegment2DSLPhiQuality::h2DHitEff_SuperPhi
private

Definition at line 62 of file DTSegment2DSLPhiQuality.h.

HRes2DHit* DTSegment2DSLPhiQuality::h2DHitSuperPhi
private

Definition at line 61 of file DTSegment2DSLPhiQuality.h.

std::string DTSegment2DSLPhiQuality::rootFileName
private

Definition at line 52 of file DTSegment2DSLPhiQuality.h.

std::string DTSegment2DSLPhiQuality::segment4DLabel
private

Definition at line 55 of file DTSegment2DSLPhiQuality.h.

double DTSegment2DSLPhiQuality::sigmaResAngle
private

Definition at line 59 of file DTSegment2DSLPhiQuality.h.

double DTSegment2DSLPhiQuality::sigmaResPos
private

Definition at line 57 of file DTSegment2DSLPhiQuality.h.

std::string DTSegment2DSLPhiQuality::simHitLabel
private

Definition at line 54 of file DTSegment2DSLPhiQuality.h.

TFile* DTSegment2DSLPhiQuality::theFile
private

Definition at line 48 of file DTSegment2DSLPhiQuality.h.