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 edm::EDConsumerBase

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
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

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
 
- 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::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

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

Author
S. Bolognesi and G. Cerminara - INFN Torino

Definition at line 30 of file DTSegment2DQuality.h.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 37 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().

37  {
38  // Get the debug parameter for verbose output
39  debug = pset.getUntrackedParameter<bool>("debug");
41  // the name of the simhit collection
42  simHitLabel = pset.getUntrackedParameter<InputTag>("simHitLabel");
43  // the name of the 2D rec hit collection
44  segment2DLabel = pset.getUntrackedParameter<InputTag>("segment2DLabel");
45 
46  //sigma resolution on position
47  sigmaResPos = pset.getParameter<double>("sigmaResPos");
48  //sigma resolution on angle
49  sigmaResAngle = pset.getParameter<double>("sigmaResAngle");
50 
51  // Create the root file
52  //theFile = new TFile(rootFileName.c_str(), "RECREATE");
53  //theFile->cd();
54 
55  if(debug)
56  cout << "[DTSegment2DQuality] Constructor called " << endl;
57 
58  // ----------------------
59  // get hold of back-end interface
60  dbe_ = 0;
62  if ( dbe_ ) {
63  if (debug) {
64  dbe_->setVerbose(1);
65  } else {
66  dbe_->setVerbose(0);
67  }
68  }
69  if ( dbe_ ) {
70  if ( debug ) dbe_->showDirStructure();
71  }
72 
73  h2DHitRPhi = new HRes2DHit ("RPhi",dbe_,true,true);
74  h2DHitRZ= new HRes2DHit ("RZ",dbe_,true,true);
75  h2DHitRZ_W0= new HRes2DHit ("RZ_W0",dbe_,true,true);
76  h2DHitRZ_W1= new HRes2DHit ("RZ_W1",dbe_,true,true);
77  h2DHitRZ_W2= new HRes2DHit ("RZ_W2",dbe_,true,true);
78 
79  h2DHitEff_RPhi= new HEff2DHit ("RPhi",dbe_);
80  h2DHitEff_RZ= new HEff2DHit ("RZ",dbe_);
81  h2DHitEff_RZ_W0= new HEff2DHit ("RZ_W0",dbe_);
82  h2DHitEff_RZ_W1= new HEff2DHit ("RZ_W1",dbe_);
83  h2DHitEff_RZ_W2= new HEff2DHit ("RZ_W2",dbe_);
84  if(debug)
85  cout << "[DTSegment2DQuality] hitsos created " << endl;
86 }
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:548
HEff2DHit * h2DHitEff_RZ_W1
HEff2DHit * h2DHitEff_RPhi
edm::InputTag segment2DLabel
tuple cout
Definition: gather_cfg.py:121
void showDirStructure(void) const
Definition: DQMStore.cc:2961
HEff2DHit * h2DHitEff_RZ_W2
DTSegment2DQuality::~DTSegment2DQuality ( )
virtual

Destructor.

Definition at line 89 of file DTSegment2DQuality.cc.

89  {
90 
91 }

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 120 of file DTSegment2DQuality.cc.

References funct::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().

120  {
121  //theFile->cd();
122 
123  // Get the DT Geometry
124  ESHandle<DTGeometry> dtGeom;
125  eventSetup.get<MuonGeometryRecord>().get(dtGeom);
126 
127  // Get the SimHit collection from the event
129  event.getByLabel(simHitLabel, simHits); //FIXME: second string to be removed
130 
131  //Map simHits by sl
132  map<DTSuperLayerId, PSimHitContainer > simHitsPerSl;
133  for(PSimHitContainer::const_iterator simHit = simHits->begin();
134  simHit != simHits->end(); simHit++){
135  // Create the id of the sl (the simHits in the DT known their wireId)
136  DTSuperLayerId slId = ((DTWireId(simHit->detUnitId())).layerId()).superlayerId();
137  // Fill the map
138  simHitsPerSl[slId].push_back(*simHit);
139  }
140 
141  // Get the 2D rechits from the event
143  event.getByLabel(segment2DLabel, segment2Ds);
144 
145  if(!segment2Ds.isValid()) {
146  if(debug) cout << "[DTSegment2DQuality]**Warning: no 2DSegments with label: " << segment2DLabel
147  << " in this event, skipping!" << endl;
148  return;
149  }
150 
151  // Loop over all superlayers containing a segment
152  DTRecSegment2DCollection::id_iterator slId;
153  for (slId = segment2Ds->id_begin();
154  slId != segment2Ds->id_end();
155  ++slId){
156 
157  //------------------------- simHits ---------------------------//
158  //Get simHits of each superlayer
159  PSimHitContainer simHits = simHitsPerSl[(*slId)];
160 
161  // Map simhits per wire
162  map<DTWireId, PSimHitContainer > simHitsPerWire = DTHitQualityUtils::mapSimHitsPerWire(simHits);
163  map<DTWireId, const PSimHit*> muSimHitPerWire = DTHitQualityUtils::mapMuSimHitsPerWire(simHitsPerWire);
164  int nMuSimHit = muSimHitPerWire.size();
165  if(nMuSimHit == 0 || nMuSimHit == 1) {
166  if(debug && nMuSimHit == 1)
167  cout << "[DTSegment2DQuality] Only " << nMuSimHit << " mu SimHit in this SL, skipping!" << endl;
168  continue; // If no or only one mu SimHit is found skip this SL
169  }
170  if(debug)
171  cout << "=== SL " << (*slId) << " has " << nMuSimHit << " SimHits" << endl;
172 
173  //Find outer and inner mu SimHit to build a segment
174  pair<const PSimHit*, const PSimHit*> inAndOutSimHit = DTHitQualityUtils::findMuSimSegment(muSimHitPerWire);
175  //Check that outermost and innermost SimHit are not the same
176  if(inAndOutSimHit.first ==inAndOutSimHit.second ) {
177  cout << "[DTHitQualityUtils]***Warning: outermost and innermost SimHit are the same!" << endl;
178  continue;
179  }
180 
181  //Find direction and position of the sim Segment in SL RF
182  pair<LocalVector, LocalPoint> dirAndPosSimSegm = DTHitQualityUtils::findMuSimSegmentDirAndPos(inAndOutSimHit,
183  (*slId),&(*dtGeom));
184 
185  LocalVector simSegmLocalDir = dirAndPosSimSegm.first;
186  LocalPoint simSegmLocalPos = dirAndPosSimSegm.second;
187  if(debug)
188  cout<<" Simulated segment: local direction "<<simSegmLocalDir<<endl
189  <<" local position "<<simSegmLocalPos<<endl;
190  const DTSuperLayer* superLayer = dtGeom->superLayer(*slId);
191  GlobalPoint simSegmGlobalPos = superLayer->toGlobal(simSegmLocalPos);
192 
193  //Atan(x/z) angle and x position in SL RF
194  float angleSimSeg = DTHitQualityUtils::findSegmentAlphaAndBeta(simSegmLocalDir).first;
195  float posSimSeg = simSegmLocalPos.x();
196  //Position (in eta,phi coordinates) in the global RF
197  float etaSimSeg = simSegmGlobalPos.eta();
198  float phiSimSeg = simSegmGlobalPos.phi();
199 
200 
201  //---------------------------- recHits --------------------------//
202  // Get the range of rechit for the corresponding slId
203  bool recHitFound = false;
204  DTRecSegment2DCollection::range range = segment2Ds->get(*slId);
205  int nsegm = distance(range.first, range.second);
206  if(debug)
207  cout << " Sl: " << *slId << " has " << nsegm
208  << " 2D segments" << endl;
209 
210  if (nsegm!=0) {
211  // Find the best RecHit: look for the 2D RecHit with the angle closest
212  // to that of segment made of SimHits.
213  // RecHits must have delta alpha and delta position within 5 sigma of
214  // the residual distribution (we are looking for residuals of segments
215  // usefull to the track fit) for efficency purpose
216  const DTRecSegment2D* bestRecHit = 0;
217  bool bestRecHitFound = false;
218  double deltaAlpha = 99999;
219 
220  // Loop over the recHits of this slId
221  for (DTRecSegment2DCollection::const_iterator segment2D = range.first;
222  segment2D!=range.second;
223  ++segment2D){
224  // Check the dimension
225  if((*segment2D).dimension() != 2) {
226  if(debug) cout << "[DTSegment2DQuality]***Error: This is not 2D segment!!!" << endl;
227  abort();
228  }
229  // Segment Local Direction and position (in SL RF)
230  LocalVector recSegDirection = (*segment2D).localDirection();
231  LocalPoint recSegPosition = (*segment2D).localPosition();
232 
233  float recSegAlpha = DTHitQualityUtils::findSegmentAlphaAndBeta(recSegDirection).first;
234  if(debug)
235  cout << " RecSegment direction: " << recSegDirection << endl
236  << " position : " << recSegPosition << endl
237  << " alpha : " << recSegAlpha << endl;
238 
239  if(fabs(recSegAlpha - angleSimSeg) < deltaAlpha) {
240  deltaAlpha = fabs(recSegAlpha - angleSimSeg);
241  bestRecHit = &(*segment2D);
242  bestRecHitFound = true;
243  }
244  } // End of Loop over all 2D RecHits
245 
246  if(bestRecHitFound) {
247  // Best rechit direction and position in SL RF
248  LocalPoint bestRecHitLocalPos = bestRecHit->localPosition();
249  LocalVector bestRecHitLocalDir = bestRecHit->localDirection();
250 
251  LocalError bestRecHitLocalPosErr = bestRecHit->localPositionError();
252  LocalError bestRecHitLocalDirErr = bestRecHit->localDirectionError();
253 
254  float angleBestRHit = DTHitQualityUtils::findSegmentAlphaAndBeta(bestRecHitLocalDir).first;
255 
256  if(fabs(angleBestRHit - angleSimSeg) < 5*sigmaResAngle &&
257  fabs(bestRecHitLocalPos.x() - posSimSeg) < 5*sigmaResPos) {
258  recHitFound = true;
259  }
260 
261  // Fill Residual histos
262  HRes2DHit *hRes = 0;
263  if((*slId).superlayer() == 1 || (*slId).superlayer() == 3) { //RPhi SL
264  hRes = h2DHitRPhi;
265  } else if((*slId).superlayer() == 2) { //RZ SL
266  h2DHitRZ->Fill(angleSimSeg,
267  angleBestRHit,
268  posSimSeg,
269  bestRecHitLocalPos.x(),
270  etaSimSeg,
271  phiSimSeg,
272  sqrt(bestRecHitLocalPosErr.xx()),
273  sqrt(bestRecHitLocalDirErr.xx())
274  );
275  if(abs((*slId).wheel()) == 0)
276  hRes = h2DHitRZ_W0;
277  else if(abs((*slId).wheel()) == 1)
278  hRes = h2DHitRZ_W1;
279  else if(abs((*slId).wheel()) == 2)
280  hRes = h2DHitRZ_W2;
281  }
282  hRes->Fill(angleSimSeg,
283  angleBestRHit,
284  posSimSeg,
285  bestRecHitLocalPos.x(),
286  etaSimSeg,
287  phiSimSeg,
288  sqrt(bestRecHitLocalPosErr.xx()),
289  sqrt(bestRecHitLocalDirErr.xx())
290  );
291  }
292  } //end of if(nsegm!=0)
293 
294  // Fill Efficiency plot
295  HEff2DHit *hEff = 0;
296  if((*slId).superlayer() == 1 || (*slId).superlayer() == 3) { //RPhi SL
297  hEff = h2DHitEff_RPhi;
298  } else if((*slId).superlayer() == 2) { //RZ SL
299  h2DHitEff_RZ->Fill(etaSimSeg, phiSimSeg, posSimSeg, angleSimSeg, recHitFound);
300  if(abs((*slId).wheel()) == 0)
301  hEff = h2DHitEff_RZ_W0;
302  else if(abs((*slId).wheel()) == 1)
303  hEff = h2DHitEff_RZ_W1;
304  else if(abs((*slId).wheel()) == 2)
305  hEff = h2DHitEff_RZ_W2;
306  }
307  hEff->Fill(etaSimSeg, phiSimSeg, posSimSeg, angleSimSeg, recHitFound);
308  } // End of loop over superlayers
309 }
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:50
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:69
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) ...
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:388
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
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:76
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:62
HEff2DHit * h2DHitEff_RZ_W2
void DTSegment2DQuality::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 93 of file DTSegment2DQuality.cc.

93  {
94  // Write the histos to file
95  //theFile->cd();
96 
97  //h2DHitRPhi->Write();
98  //h2DHitRZ->Write();
99  //h2DHitRZ_W0->Write();
100  //h2DHitRZ_W1->Write();
101  //h2DHitRZ_W2->Write();
102 
108 
109  //h2DHitEff_RPhi->Write();
110  //h2DHitEff_RZ->Write();
111  //h2DHitEff_RZ_W0->Write();
112  //h2DHitEff_RZ_W1->Write();
113  //h2DHitEff_RZ_W2->Write();
114  //if ( rootFileName.size() != 0 && dbe_ ) dbe_->save(rootFileName);
115 
116  //theFile->Close();
117 }
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 74 of file DTSegment2DQuality.h.

bool DTSegment2DQuality::debug
private

Definition at line 52 of file DTSegment2DQuality.h.

HEff2DHit* DTSegment2DQuality::h2DHitEff_RPhi
private

Definition at line 69 of file DTSegment2DQuality.h.

HEff2DHit* DTSegment2DQuality::h2DHitEff_RZ
private

Definition at line 70 of file DTSegment2DQuality.h.

HEff2DHit* DTSegment2DQuality::h2DHitEff_RZ_W0
private

Definition at line 71 of file DTSegment2DQuality.h.

HEff2DHit* DTSegment2DQuality::h2DHitEff_RZ_W1
private

Definition at line 72 of file DTSegment2DQuality.h.

HEff2DHit* DTSegment2DQuality::h2DHitEff_RZ_W2
private

Definition at line 73 of file DTSegment2DQuality.h.

HRes2DHit* DTSegment2DQuality::h2DHitRPhi
private

Definition at line 63 of file DTSegment2DQuality.h.

HRes2DHit* DTSegment2DQuality::h2DHitRZ
private

Definition at line 64 of file DTSegment2DQuality.h.

HRes2DHit* DTSegment2DQuality::h2DHitRZ_W0
private

Definition at line 65 of file DTSegment2DQuality.h.

HRes2DHit* DTSegment2DQuality::h2DHitRZ_W1
private

Definition at line 66 of file DTSegment2DQuality.h.

HRes2DHit* DTSegment2DQuality::h2DHitRZ_W2
private

Definition at line 67 of file DTSegment2DQuality.h.

std::string DTSegment2DQuality::rootFileName
private

Definition at line 54 of file DTSegment2DQuality.h.

edm::InputTag DTSegment2DQuality::segment2DLabel
private

Definition at line 57 of file DTSegment2DQuality.h.

double DTSegment2DQuality::sigmaResAngle
private

Definition at line 61 of file DTSegment2DQuality.h.

double DTSegment2DQuality::sigmaResPos
private

Definition at line 59 of file DTSegment2DQuality.h.

edm::InputTag DTSegment2DQuality::simHitLabel
private

Definition at line 56 of file DTSegment2DQuality.h.