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
DTSegment2DQuality Class Reference

#include <DTSegment2DQuality.h>

Inheritance diagram for DTSegment2DQuality:
edm::EDAnalyzer

Public Member Functions

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

Private Attributes

DQMStoredbe_
 
bool debug
 
HEff2DHith2DHitEff_RPhi
 
HEff2DHith2DHitEff_RZ
 
HEff2DHith2DHitEff_RZ_W0
 
HEff2DHith2DHitEff_RZ_W1
 
HEff2DHith2DHitEff_RZ_W2
 
HRes2DHith2DHitRPhi
 
HRes2DHith2DHitRZ
 
HRes2DHith2DHitRZ_W0
 
HRes2DHith2DHitRZ_W1
 
HRes2DHith2DHitRZ_W2
 
std::string rootFileName
 
edm::InputTag segment2DLabel
 
double sigmaResAngle
 
double sigmaResPos
 
edm::InputTag simHitLabel
 

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)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Basic analyzer class which accesses 2D DTSegments and plot resolution comparing reconstructed and simulated quantities

Date:
2010/09/13 09:49:18
Revision:
1.4
Author
S. Bolognesi and G. Cerminara - INFN Torino

Definition at line 32 of file DTSegment2DQuality.h.

Constructor & Destructor Documentation

DTSegment2DQuality::DTSegment2DQuality ( const edm::ParameterSet pset)

Constructor.

Definition at line 39 of file DTSegment2DQuality.cc.

References gather_cfg::cout, dbe_, debug, DTHitQualityUtils::debug, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), cppFunctionSkipper::operator, DQMStore::setVerbose(), and DQMStore::showDirStructure().

39  {
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 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::InputTag simHitLabel
HEff2DHit * h2DHitEff_RZ_W0
void setVerbose(unsigned level)
Definition: DQMStore.cc:393
HEff2DHit * h2DHitEff_RZ_W1
HEff2DHit * h2DHitEff_RPhi
edm::InputTag segment2DLabel
tuple cout
Definition: gather_cfg.py:121
void showDirStructure(void) const
Definition: DQMStore.cc:2761
HEff2DHit * h2DHitEff_RZ_W2
DTSegment2DQuality::~DTSegment2DQuality ( )
virtual

Destructor.

Definition at line 91 of file DTSegment2DQuality.cc.

91  {
92 
93 }

Member Function Documentation

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

Perform the real analysis.

Implements edm::EDAnalyzer.

Definition at line 122 of file DTSegment2DQuality.cc.

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

122  {
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 }
virtual LocalError localPositionError() const
local position error in SL frame
float xx() const
Definition: LocalError.h:24
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:52
edm::InputTag simHitLabel
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:68
#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.
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:46
void Fill(float angleSimSegment, float angleRecSegment, float posSimSegment, float posRecSegment, float etaSimSegment, float phiSimSegment, float sigmaPos, float sigmaAngle)
Definition: Histograms.h:388
HEff2DHit * h2DHitEff_RZ_W0
bool isValid() const
Definition: HandleBase.h:76
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...
HEff2DHit * h2DHitEff_RZ_W1
HEff2DHit * h2DHitEff_RPhi
void Fill(float etaSimSegm, float phiSimSegm, float posSimSegm, float angleSimSegm, bool fillRecHit)
Definition: Histograms.h:529
const T & get() const
Definition: EventSetup.h:55
tuple simHits
Definition: trackerHits.py:16
T eta() const
Definition: PV3DBase.h:75
edm::InputTag segment2DLabel
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)
T x() const
Definition: PV3DBase.h:61
HEff2DHit * h2DHitEff_RZ_W2
void DTSegment2DQuality::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 95 of file DTSegment2DQuality.cc.

95  {
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 
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 }
void ComputeEfficiency()
Definition: Histograms.h:549
HEff2DHit * h2DHitEff_RZ_W0
HEff2DHit * h2DHitEff_RZ_W1
HEff2DHit * h2DHitEff_RPhi
HEff2DHit * h2DHitEff_RZ_W2

Member Data Documentation

DQMStore* DTSegment2DQuality::dbe_
private

Definition at line 76 of file DTSegment2DQuality.h.

bool DTSegment2DQuality::debug
private

Definition at line 54 of file DTSegment2DQuality.h.

HEff2DHit* DTSegment2DQuality::h2DHitEff_RPhi
private

Definition at line 71 of file DTSegment2DQuality.h.

HEff2DHit* DTSegment2DQuality::h2DHitEff_RZ
private

Definition at line 72 of file DTSegment2DQuality.h.

HEff2DHit* DTSegment2DQuality::h2DHitEff_RZ_W0
private

Definition at line 73 of file DTSegment2DQuality.h.

HEff2DHit* DTSegment2DQuality::h2DHitEff_RZ_W1
private

Definition at line 74 of file DTSegment2DQuality.h.

HEff2DHit* DTSegment2DQuality::h2DHitEff_RZ_W2
private

Definition at line 75 of file DTSegment2DQuality.h.

HRes2DHit* DTSegment2DQuality::h2DHitRPhi
private

Definition at line 65 of file DTSegment2DQuality.h.

HRes2DHit* DTSegment2DQuality::h2DHitRZ
private

Definition at line 66 of file DTSegment2DQuality.h.

HRes2DHit* DTSegment2DQuality::h2DHitRZ_W0
private

Definition at line 67 of file DTSegment2DQuality.h.

HRes2DHit* DTSegment2DQuality::h2DHitRZ_W1
private

Definition at line 68 of file DTSegment2DQuality.h.

HRes2DHit* DTSegment2DQuality::h2DHitRZ_W2
private

Definition at line 69 of file DTSegment2DQuality.h.

std::string DTSegment2DQuality::rootFileName
private

Definition at line 56 of file DTSegment2DQuality.h.

edm::InputTag DTSegment2DQuality::segment2DLabel
private

Definition at line 59 of file DTSegment2DQuality.h.

double DTSegment2DQuality::sigmaResAngle
private

Definition at line 63 of file DTSegment2DQuality.h.

double DTSegment2DQuality::sigmaResPos
private

Definition at line 61 of file DTSegment2DQuality.h.

edm::InputTag DTSegment2DQuality::simHitLabel
private

Definition at line 58 of file DTSegment2DQuality.h.