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

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 ()
 
void endLuminosityBlock (edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &c)
 
virtual ~DTSegment2DSLPhiQuality ()
 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
 
bool doall
 
HEff2DHith2DHitEff_SuperPhi
 
HRes2DHith2DHitSuperPhi
 
bool local
 
std::string rootFileName
 
edm::InputTag segment4DLabel
 
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 reconstructed with both SL Phi and plot resolution comparing reconstructed and simulated quantities

Author
S. Bolognesi and G. Cerminara - INFN Torino

Definition at line 30 of file DTSegment2DSLPhiQuality.h.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 40 of file DTSegment2DSLPhiQuality.cc.

References dbe_, debug, DTHitQualityUtils::debug, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), and cppFunctionSkipper::operator.

40  {
41  // Get the debug parameter for verbose output
42  debug = pset.getUntrackedParameter<bool>("debug");
44 
45  // the name of the simhit collection
46  simHitLabel = pset.getUntrackedParameter<InputTag>("simHitLabel");
47  // the name of the 4D rec hit collection
48  segment4DLabel = pset.getUntrackedParameter<InputTag>("segment4DLabel");
49 
50  //sigma resolution on position
51  sigmaResPos = pset.getParameter<double>("sigmaResPos");
52  //sigma resolution on angle
53  sigmaResAngle = pset.getParameter<double>("sigmaResAngle");
54  doall = pset.getUntrackedParameter<bool>("doall", false);
55  local = pset.getUntrackedParameter<bool>("local", false);
56 
57  // Create the root file
58  //theFile = new TFile(rootFileName.c_str(), "RECREATE");
59  //theFile->cd();
60  // ----------------------
61  // get hold of back-end interface
62  dbe_ = 0;
63  dbe_ = Service<DQMStore>().operator->();
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  // Book the histos
76  h2DHitSuperPhi = new HRes2DHit ("SuperPhi",dbe_,doall,local);
77  if(doall) h2DHitEff_SuperPhi = new HEff2DHit ("SuperPhi",dbe_);
78 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
void setVerbose(unsigned level)
Definition: DQMStore.cc:549
void showDirStructure(void) const
Definition: DQMStore.cc:2962
DTSegment2DSLPhiQuality::~DTSegment2DSLPhiQuality ( )
virtual

Destructor.

Definition at line 81 of file DTSegment2DSLPhiQuality.cc.

81  {
82 }

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 104 of file DTSegment2DSLPhiQuality.cc.

References gather_cfg::cout, debug, PV3DBase< T, PVType, FrameType >::eta(), 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().

104  {
105  //theFile->cd();
106 
107  // Get the DT Geometry
108  ESHandle<DTGeometry> dtGeom;
109  eventSetup.get<MuonGeometryRecord>().get(dtGeom);
110 
111  // Get the SimHit collection from the event
113  event.getByLabel(simHitLabel, simHits); //FIXME: second string to be removed
114 
115  //Map simHits by chamber
116  map<DTChamberId, PSimHitContainer > simHitsPerCh;
117  for(PSimHitContainer::const_iterator simHit = simHits->begin();
118  simHit != simHits->end(); simHit++){
119  // Create the id of the chamber (the simHits in the DT known their wireId)
120  DTChamberId chamberId = (((DTWireId(simHit->detUnitId())).layerId()).superlayerId()).chamberId();
121  // Fill the map
122  simHitsPerCh[chamberId].push_back(*simHit);
123  }
124 
125  // Get the 4D rechits from the event
127  event.getByLabel(segment4DLabel, segment4Ds);
128 
129  if(!segment4Ds.isValid()) {
130  if(debug) cout << "[DTSegment2DSLPhiQuality]**Warning: no 4D Segments with label: " << segment4DLabel
131  << " in this event, skipping!" << endl;
132  return;
133  }
134 
135  // Loop over all chambers containing a segment
136  DTRecSegment4DCollection::id_iterator chamberId;
137  for (chamberId = segment4Ds->id_begin();
138  chamberId != segment4Ds->id_end();
139  ++chamberId){
140 
141  //------------------------- simHits ---------------------------//
142  //Get simHits of each chamber
143  PSimHitContainer simHits = simHitsPerCh[(*chamberId)];
144 
145  // Map simhits per wire
146  map<DTWireId, PSimHitContainer > simHitsPerWire = DTHitQualityUtils::mapSimHitsPerWire(simHits);
147  map<DTWireId, const PSimHit*> muSimHitPerWire = DTHitQualityUtils::mapMuSimHitsPerWire(simHitsPerWire);
148  int nMuSimHit = muSimHitPerWire.size();
149  if(nMuSimHit == 0 || nMuSimHit == 1) {
150  if(debug && nMuSimHit == 1)
151  cout << "[DTSegment2DSLPhiQuality] Only " << nMuSimHit << " mu SimHit in this chamber, skipping!" << endl;
152  continue; // If no or only one mu SimHit is found skip this chamber
153  }
154  if(debug)
155  cout << "=== Chamber " << (*chamberId) << " has " << nMuSimHit << " SimHits" << endl;
156 
157  //Find outer and inner mu SimHit to build a segment
158  pair<const PSimHit*, const PSimHit*> inAndOutSimHit = DTHitQualityUtils::findMuSimSegment(muSimHitPerWire);
159 
160  //Find direction and position of the sim Segment in Chamber RF
161  pair<LocalVector, LocalPoint> dirAndPosSimSegm = DTHitQualityUtils::findMuSimSegmentDirAndPos(inAndOutSimHit,
162  (*chamberId),&(*dtGeom));
163 
164  LocalVector simSegmLocalDir = dirAndPosSimSegm.first;
165  LocalPoint simSegmLocalPos = dirAndPosSimSegm.second;
166  const DTChamber* chamber = dtGeom->chamber(*chamberId);
167  GlobalPoint simSegmGlobalPos = chamber->toGlobal(simSegmLocalPos);
168 
169  //Atan(x/z) angle and x position in SL RF
170  float angleSimSeg = DTHitQualityUtils::findSegmentAlphaAndBeta(simSegmLocalDir).first;
171  float posSimSeg = simSegmLocalPos.x();
172  //Position (in eta,phi coordinates) in lobal RF
173  float etaSimSeg = simSegmGlobalPos.eta();
174  float phiSimSeg = simSegmGlobalPos.phi();
175 
176  if(debug)
177  cout<<" Simulated segment: local direction "<<simSegmLocalDir<<endl
178  <<" local position "<<simSegmLocalPos<<endl
179  <<" angle "<<angleSimSeg<<endl;
180 
181  //---------------------------- recHits --------------------------//
182  // Get the range of rechit for the corresponding chamberId
183  bool recHitFound = false;
184  DTRecSegment4DCollection::range range = segment4Ds->get(*chamberId);
185  int nsegm = distance(range.first, range.second);
186  if(debug)
187  cout << " Chamber: " << *chamberId << " has " << nsegm
188  << " 4D segments" << endl;
189 
190  if (nsegm!=0) {
191  // Find the best RecHit: look for the 4D RecHit with the phi angle closest
192  // to that of segment made of SimHits.
193  // RecHits must have delta alpha and delta position within 5 sigma of
194  // the residual distribution (we are looking for residuals of segments
195  // usefull to the track fit) for efficency purpose
196  const DTRecSegment2D* bestRecHit = 0;
197  bool bestRecHitFound = false;
198  double deltaAlpha = 99999;
199 
200  // Loop over the recHits of this chamberId
201  for (DTRecSegment4DCollection::const_iterator segment4D = range.first;
202  segment4D!=range.second;
203  ++segment4D){
204  // Check the dimension
205  if((*segment4D).dimension() != 4) {
206  if(debug) cout << "[DTSegment2DSLPhiQuality]***Error: This is not 4D segment!!!" << endl;
207  continue;
208  }
209 
210  //Get 2D superPhi segments from 4D segments
211  const DTChamberRecSegment2D* phiSegment2D = (*segment4D).phiSegment();
212  if((*phiSegment2D).dimension() != 2) {
213  if(debug) cout << "[DTSegment2DQuality]***Error: This is not 2D segment!!!" << endl;
214  abort();
215  }
216 
217  // Segment Local Direction and position (in Chamber RF)
218  LocalVector recSegDirection = (*phiSegment2D).localDirection();
219 
220  float recSegAlpha = DTHitQualityUtils::findSegmentAlphaAndBeta(recSegDirection).first;
221  if(debug)
222  cout << " RecSegment direction: " << recSegDirection << endl
223  << " position : " << (*phiSegment2D).localPosition() << endl
224  << " alpha : " << recSegAlpha << endl;
225 
226  if(fabs(recSegAlpha - angleSimSeg) < deltaAlpha) {
227  deltaAlpha = fabs(recSegAlpha - angleSimSeg);
228  bestRecHit = &(*phiSegment2D);
229  bestRecHitFound = true;
230  }
231  } // End of Loop over all 4D RecHits of this chambers
232 
233  if(bestRecHitFound) {
234  // Best rechit direction and position in Chamber RF
235  LocalPoint bestRecHitLocalPos = bestRecHit->localPosition();
236  LocalVector bestRecHitLocalDir = bestRecHit->localDirection();
237 
238  LocalError bestRecHitLocalPosErr = bestRecHit->localPositionError();
239  LocalError bestRecHitLocalDirErr = bestRecHit->localDirectionError();
240 
241  float angleBestRHit = DTHitQualityUtils::findSegmentAlphaAndBeta(bestRecHitLocalDir).first;
242  if(fabs(angleBestRHit - angleSimSeg) < 5*sigmaResAngle &&
243  fabs(bestRecHitLocalPos.x() - posSimSeg) < 5*sigmaResPos) {
244  recHitFound = true;
245  }
246 
247  // Fill Residual histos
248  h2DHitSuperPhi->Fill(angleSimSeg,
249  angleBestRHit,
250  posSimSeg,
251  bestRecHitLocalPos.x(),
252  etaSimSeg,
253  phiSimSeg,
254  sqrt(bestRecHitLocalPosErr.xx()),
255  sqrt(bestRecHitLocalDirErr.xx())
256  );
257  }
258  } //end of if(nsegm!=0)
259 
260  // Fill Efficiency plot
261  if(doall) {h2DHitEff_SuperPhi->Fill(etaSimSeg,
262  phiSimSeg,
263  posSimSeg,
264  angleSimSeg,
265  recHitFound);}
266  } // End of loop over chambers
267 }
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
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
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...
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
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
void DTSegment2DSLPhiQuality::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 91 of file DTSegment2DSLPhiQuality.cc.

91  {
92  // Write the histos to file
93  //theFile->cd();
94  //h2DHitSuperPhi->Write();
95 
97  //h2DHitEff_SuperPhi->Write();
98 
99  //if ( rootFileName.size() != 0 && dbe_ ) dbe_->save(rootFileName);
100  //theFile->Close();
101 }
void ComputeEfficiency()
Definition: Histograms.h:549
void DTSegment2DSLPhiQuality::endLuminosityBlock ( edm::LuminosityBlock const &  lumiSeg,
edm::EventSetup const &  c 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 84 of file DTSegment2DSLPhiQuality.cc.

85  {
86 }

Member Data Documentation

DQMStore* DTSegment2DSLPhiQuality::dbe_
private

Definition at line 68 of file DTSegment2DSLPhiQuality.h.

bool DTSegment2DSLPhiQuality::debug
private

Definition at line 55 of file DTSegment2DSLPhiQuality.h.

bool DTSegment2DSLPhiQuality::doall
private

Definition at line 69 of file DTSegment2DSLPhiQuality.h.

HEff2DHit* DTSegment2DSLPhiQuality::h2DHitEff_SuperPhi
private

Definition at line 67 of file DTSegment2DSLPhiQuality.h.

HRes2DHit* DTSegment2DSLPhiQuality::h2DHitSuperPhi
private

Definition at line 66 of file DTSegment2DSLPhiQuality.h.

bool DTSegment2DSLPhiQuality::local
private

Definition at line 70 of file DTSegment2DSLPhiQuality.h.

std::string DTSegment2DSLPhiQuality::rootFileName
private

Definition at line 57 of file DTSegment2DSLPhiQuality.h.

edm::InputTag DTSegment2DSLPhiQuality::segment4DLabel
private

Definition at line 60 of file DTSegment2DSLPhiQuality.h.

double DTSegment2DSLPhiQuality::sigmaResAngle
private

Definition at line 64 of file DTSegment2DSLPhiQuality.h.

double DTSegment2DSLPhiQuality::sigmaResPos
private

Definition at line 62 of file DTSegment2DSLPhiQuality.h.

edm::InputTag DTSegment2DSLPhiQuality::simHitLabel
private

Definition at line 59 of file DTSegment2DSLPhiQuality.h.