CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
DTSegmentSelector Class Reference

#include <DTSegmentSelector.h>

Public Member Functions

 DTSegmentSelector (edm::ParameterSet const &pset, edm::ConsumesCollector &iC)
 
bool operator() (DTRecSegment4D const &, edm::Event const &, edm::EventSetup const &)
 
 ~DTSegmentSelector ()
 

Private Member Functions

bool checkNoisySegment (edm::ESHandle< DTStatusFlag > const &, std::vector< DTRecHit1D > const &)
 

Private Attributes

bool checkNoisyChannels_
 
double maxAnglePhi_
 
double maxAngleZ_
 
double maxChi2_
 
int minHitsPhi_
 
int minHitsZ_
 
edm::InputTag muonTags_
 
edm::EDGetTokenT< reco::MuonCollectionmuonToken_
 

Detailed Description

Definition at line 24 of file DTSegmentSelector.h.

Constructor & Destructor Documentation

DTSegmentSelector::DTSegmentSelector ( edm::ParameterSet const &  pset,
edm::ConsumesCollector iC 
)

Definition at line 16 of file DTSegmentSelector.cc.

References edm::ConsumesCollector::consumes(), muonTags_, and muonToken_.

16  :
17  muonTags_(pset.getParameter<edm::InputTag>("Muons")),
18  checkNoisyChannels_(pset.getParameter<bool>("checkNoisyChannels")),
19  minHitsPhi_(pset.getParameter<int>("minHitsPhi")),
20  minHitsZ_(pset.getParameter<int>("minHitsZ")),
21  maxChi2_(pset.getParameter<double>("maxChi2")),
22  maxAnglePhi_(pset.getParameter<double>("maxAnglePhi")),
23  maxAngleZ_(pset.getParameter<double>("maxAngleZ"))
24 {
26 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
edm::EDGetTokenT< reco::MuonCollection > muonToken_
edm::InputTag muonTags_
DTSegmentSelector::~DTSegmentSelector ( )
inline

Definition at line 27 of file DTSegmentSelector.h.

References checkNoisySegment(), and operator()().

27 {}

Member Function Documentation

bool DTSegmentSelector::checkNoisySegment ( edm::ESHandle< DTStatusFlag > const &  statusMap,
std::vector< DTRecHit1D > const &  dtHits 
)
private

Definition at line 119 of file DTSegmentSelector.cc.

References DTStatusFlag::cellStatus(), and LogTrace.

Referenced by operator()(), and ~DTSegmentSelector().

119  {
120 
121  bool segmentNoisy = false;
122 
123  std::vector<DTRecHit1D>::const_iterator dtHit = dtHits.begin();
124  std::vector<DTRecHit1D>::const_iterator dtHits_end = dtHits.end();
125  for(; dtHit != dtHits_end; ++dtHit){
126  DTWireId wireId = dtHit->wireId();
127  // Check for noisy channels to skip them
128  bool isNoisy = false, isFEMasked = false, isTDCMasked = false, isTrigMask = false,
129  isDead = false, isNohv = false;
130  statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
131  if(isNoisy) {
132  LogTrace("Calibration") << "Wire: " << wireId << " is noisy, skipping!";
133  segmentNoisy = true; break;
134  }
135  }
136  return segmentNoisy;
137 }
#define LogTrace(id)
int cellStatus(int wheelId, int stationId, int sectorId, int slId, int layerId, int cellId, bool &noiseFlag, bool &feMask, bool &tdcMask, bool &trigMask, bool &deadFlag, bool &nohvFlag) const
get content
Definition: DTStatusFlag.h:100
bool DTSegmentSelector::operator() ( DTRecSegment4D const &  segment,
edm::Event const &  event,
edm::EventSetup const &  setup 
)

Definition at line 29 of file DTSegmentSelector.cc.

References checkNoisyChannels_, checkNoisySegment(), DTRecSegment4D::chi2(), DTRecSegment4D::degreesOfFreedom(), runTauDisplay::dr, PVValHelper::dx, PVValHelper::dy, edm::EventSetup::get(), DTRecSegment4D::hasPhi(), DTRecSegment4D::hasZed(), edm::InputTag::label(), DTRecSegment4D::localDirection(), DTRecSegment4D::localPosition(), electrons_cff::matched, maxAnglePhi_, maxAngleZ_, maxChi2_, minHitsPhi_, minHitsZ_, extraflags_cff::muons, muonTags_, muonToken_, DTRecSegment4D::phiSegment(), Geom::pi(), edm::Handle< T >::product(), TrackingRecHit::rawId(), mps_fire::result, DTRecSegment2D::specificRecHits(), mathSSE::sqrt(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), PV3DBase< T, PVType, FrameType >::z(), and DTRecSegment4D::zSegment().

Referenced by ~DTSegmentSelector().

29  {
30 
31  bool result = true;
32 
33  /* get the muon collection if one is specified
34  (not specifying a muon collection switches off muon matching */
35  if (!muonTags_.label().empty()) {
37  event.getByToken(muonToken_,muonH);
38  const std::vector<reco::Muon>* muons = muonH.product();
39  //std::cout << " Muon collection size: " << muons->size() << std::endl;
40  if (muons->empty()) return false;
41 
42  DTChamberId segId(segment.rawId());
43 
44  bool matched = false;
45  for (auto &imuon : *muons)
46  for (const auto &ch : imuon.matches()) {
47  DetId chId(ch.id.rawId());
48  if (chId!=segId) continue;
49  if ( imuon.pt()<15 || !imuon.isGlobalMuon()) continue;
50 
51  int nsegs=ch.segmentMatches.size();
52  if (!nsegs) continue;
53 
54  LocalPoint posHit = segment.localPosition();
55  float dx = (posHit.x() ? posHit.x()-ch.x : 0 );
56  float dy = (posHit.y() ? posHit.y()-ch.y : 0 );
57  float dr = sqrt(dx*dx+dy*dy);
58  if (dr<5) matched=true;
59  }
60 
61  if (!matched) result=false;
62  }
63 
64 
65 
67  if(checkNoisyChannels_) setup.get<DTStatusFlagRcd>().get(statusMap);
68 
69  // Get the Phi 2D segment
70  int nPhiHits = -1;
71  bool segmentNoisyPhi = false;
72  if( segment.hasPhi() ){
73  const DTChamberRecSegment2D* phiSeg = segment.phiSegment(); // phiSeg lives in the chamber RF
74  std::vector<DTRecHit1D> phiRecHits = phiSeg->specificRecHits();
75  nPhiHits = phiRecHits.size();
76  if(checkNoisyChannels_) segmentNoisyPhi = checkNoisySegment(statusMap,phiRecHits);
77  //} else result = false;
78  }
79  // Get the Theta 2D segment
80  int nZHits = -1;
81  bool segmentNoisyZ = false;
82  if( segment.hasZed() ){
83  const DTSLRecSegment2D* zSeg = segment.zSegment(); // zSeg lives in the SL RF
84  std::vector<DTRecHit1D> zRecHits = zSeg->specificRecHits();
85  nZHits = zRecHits.size();
86  if(checkNoisyChannels_) segmentNoisyZ = checkNoisySegment(statusMap,zRecHits);
87  }
88 
89  // Segment selection
90  // Discard segment if it has a noisy cell
91  if(segmentNoisyPhi || segmentNoisyZ)
92  result = false;
93 
94  // 2D-segment number of hits
95  if(segment.hasPhi() && nPhiHits < minHitsPhi_)
96  result = false;
97 
98  if(segment.hasZed() && nZHits < minHitsZ_)
99  result = false;
100 
101  // Segment chi2
102  double chiSquare = segment.chi2()/segment.degreesOfFreedom();
103  if(chiSquare > maxChi2_)
104  result = false;
105 
106  // Segment angle
107  LocalVector segment4DLocalDir = segment.localDirection();
108  double angleZ = fabs( atan(segment4DLocalDir.y()/segment4DLocalDir.z())*180./Geom::pi() );
109  if( angleZ > maxAngleZ_)
110  result = false;
111 
112  double anglePhi = fabs( atan(segment4DLocalDir.x()/segment4DLocalDir.z())*180./Geom::pi() );
113  if( anglePhi > maxAnglePhi_)
114  result = false;
115 
116  return result;
117 }
T y() const
Definition: PV3DBase.h:63
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:2
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
Definition: DetId.h:18
edm::EDGetTokenT< reco::MuonCollection > muonToken_
T const * product() const
Definition: Handle.h:74
bool checkNoisySegment(edm::ESHandle< DTStatusFlag > const &, std::vector< DTRecHit1D > const &)
std::string const & label() const
Definition: InputTag.h:36
constexpr double pi()
Definition: Pi.h:31
T x() const
Definition: PV3DBase.h:62
edm::InputTag muonTags_

Member Data Documentation

bool DTSegmentSelector::checkNoisyChannels_
private

Definition at line 35 of file DTSegmentSelector.h.

Referenced by operator()().

double DTSegmentSelector::maxAnglePhi_
private

Definition at line 39 of file DTSegmentSelector.h.

Referenced by operator()().

double DTSegmentSelector::maxAngleZ_
private

Definition at line 40 of file DTSegmentSelector.h.

Referenced by operator()().

double DTSegmentSelector::maxChi2_
private

Definition at line 38 of file DTSegmentSelector.h.

Referenced by operator()().

int DTSegmentSelector::minHitsPhi_
private

Definition at line 36 of file DTSegmentSelector.h.

Referenced by operator()().

int DTSegmentSelector::minHitsZ_
private

Definition at line 37 of file DTSegmentSelector.h.

Referenced by operator()().

edm::InputTag DTSegmentSelector::muonTags_
private

Definition at line 33 of file DTSegmentSelector.h.

Referenced by DTSegmentSelector(), and operator()().

edm::EDGetTokenT<reco::MuonCollection> DTSegmentSelector::muonToken_
private

Definition at line 34 of file DTSegmentSelector.h.

Referenced by DTSegmentSelector(), and operator()().