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_
 
edm::ESGetToken< DTStatusFlag, DTStatusFlagRcdtheStatusMapToken_
 

Detailed Description

Definition at line 26 of file DTSegmentSelector.h.

Constructor & Destructor Documentation

◆ DTSegmentSelector()

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

Definition at line 16 of file DTSegmentSelector.cc.

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

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")) {
26 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
edm::ESGetToken< DTStatusFlag, DTStatusFlagRcd > theStatusMapToken_
edm::EDGetTokenT< reco::MuonCollection > muonToken_
edm::InputTag muonTags_

◆ ~DTSegmentSelector()

DTSegmentSelector::~DTSegmentSelector ( )
inline

Definition at line 29 of file DTSegmentSelector.h.

29 {}

Member Function Documentation

◆ checkNoisySegment()

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

Definition at line 126 of file DTSegmentSelector.cc.

References DTStatusFlag::cellStatus(), and LogTrace.

Referenced by operator()().

127  {
128  bool segmentNoisy = false;
129 
130  std::vector<DTRecHit1D>::const_iterator dtHit = dtHits.begin();
131  std::vector<DTRecHit1D>::const_iterator dtHits_end = dtHits.end();
132  for (; dtHit != dtHits_end; ++dtHit) {
133  DTWireId wireId = dtHit->wireId();
134  // Check for noisy channels to skip them
135  bool isNoisy = false, isFEMasked = false, isTDCMasked = false, isTrigMask = false, isDead = false, isNohv = false;
136  statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
137  if (isNoisy) {
138  LogTrace("Calibration") << "Wire: " << wireId << " is noisy, skipping!";
139  segmentNoisy = true;
140  break;
141  }
142  }
143  return segmentNoisy;
144 }
#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:90

◆ operator()()

bool DTSegmentSelector::operator() ( DTRecSegment4D const &  segment,
edm::Event const &  event,
edm::EventSetup const &  setup 
)

Definition at line 28 of file DTSegmentSelector.cc.

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

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())
41  return false;
42 
43  DTChamberId segId(segment.rawId());
44 
45  bool matched = false;
46  for (auto& imuon : *muons)
47  for (const auto& ch : imuon.matches()) {
48  DetId chId(ch.id.rawId());
49  if (chId != segId)
50  continue;
51  if (imuon.pt() < 15 || !imuon.isGlobalMuon())
52  continue;
53 
54  int nsegs = ch.segmentMatches.size();
55  if (!nsegs)
56  continue;
57 
58  LocalPoint posHit = segment.localPosition();
59  float dx = (posHit.x() ? posHit.x() - ch.x : 0);
60  float dy = (posHit.y() ? posHit.y() - ch.y : 0);
61  float dr = sqrt(dx * dx + dy * dy);
62  if (dr < 5)
63  matched = true;
64  }
65 
66  if (!matched)
67  result = false;
68  }
69 
72  statusMap = setup.getHandle(theStatusMapToken_);
73 
74  // Get the Phi 2D segment
75  int nPhiHits = -1;
76  bool segmentNoisyPhi = false;
77  if (segment.hasPhi()) {
78  const DTChamberRecSegment2D* phiSeg = segment.phiSegment(); // phiSeg lives in the chamber RF
79  std::vector<DTRecHit1D> phiRecHits = phiSeg->specificRecHits();
80  nPhiHits = phiRecHits.size();
82  segmentNoisyPhi = checkNoisySegment(statusMap, phiRecHits);
83  //} else result = false;
84  }
85  // Get the Theta 2D segment
86  int nZHits = -1;
87  bool segmentNoisyZ = false;
88  if (segment.hasZed()) {
89  const DTSLRecSegment2D* zSeg = segment.zSegment(); // zSeg lives in the SL RF
90  std::vector<DTRecHit1D> zRecHits = zSeg->specificRecHits();
91  nZHits = zRecHits.size();
93  segmentNoisyZ = checkNoisySegment(statusMap, zRecHits);
94  }
95 
96  // Segment selection
97  // Discard segment if it has a noisy cell
98  if (segmentNoisyPhi || segmentNoisyZ)
99  result = false;
100 
101  // 2D-segment number of hits
102  if (segment.hasPhi() && nPhiHits < minHitsPhi_)
103  result = false;
104 
105  if (segment.hasZed() && nZHits < minHitsZ_)
106  result = false;
107 
108  // Segment chi2
109  double chiSquare = segment.chi2() / segment.degreesOfFreedom();
110  if (chiSquare > maxChi2_)
111  result = false;
112 
113  // Segment angle
114  LocalVector segment4DLocalDir = segment.localDirection();
115  double angleZ = fabs(atan(segment4DLocalDir.y() / segment4DLocalDir.z()) * 180. / Geom::pi());
116  if (angleZ > maxAngleZ_)
117  result = false;
118 
119  double anglePhi = fabs(atan(segment4DLocalDir.x() / segment4DLocalDir.z()) * 180. / Geom::pi());
120  if (anglePhi > maxAnglePhi_)
121  result = false;
122 
123  return result;
124 }
T z() const
Definition: PV3DBase.h:61
T const * product() const
Definition: Handle.h:70
std::string const & label() const
Definition: InputTag.h:36
muons
the two sets of parameters below are mutually exclusive, depending if RECO or ALCARECO is used the us...
Definition: DiMuonV_cfg.py:212
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
T sqrt(T t)
Definition: SSEVec.h:19
edm::ESGetToken< DTStatusFlag, DTStatusFlagRcd > theStatusMapToken_
Definition: DetId.h:17
edm::EDGetTokenT< reco::MuonCollection > muonToken_
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
bool checkNoisySegment(edm::ESHandle< DTStatusFlag > const &, std::vector< DTRecHit1D > const &)
constexpr double pi()
Definition: Pi.h:31
edm::InputTag muonTags_

Member Data Documentation

◆ checkNoisyChannels_

bool DTSegmentSelector::checkNoisyChannels_
private

Definition at line 37 of file DTSegmentSelector.h.

Referenced by operator()().

◆ maxAnglePhi_

double DTSegmentSelector::maxAnglePhi_
private

Definition at line 41 of file DTSegmentSelector.h.

Referenced by operator()().

◆ maxAngleZ_

double DTSegmentSelector::maxAngleZ_
private

Definition at line 42 of file DTSegmentSelector.h.

Referenced by operator()().

◆ maxChi2_

double DTSegmentSelector::maxChi2_
private

Definition at line 40 of file DTSegmentSelector.h.

Referenced by operator()().

◆ minHitsPhi_

int DTSegmentSelector::minHitsPhi_
private

Definition at line 38 of file DTSegmentSelector.h.

Referenced by operator()().

◆ minHitsZ_

int DTSegmentSelector::minHitsZ_
private

Definition at line 39 of file DTSegmentSelector.h.

Referenced by operator()().

◆ muonTags_

edm::InputTag DTSegmentSelector::muonTags_
private

Definition at line 35 of file DTSegmentSelector.h.

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

◆ muonToken_

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

Definition at line 36 of file DTSegmentSelector.h.

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

◆ theStatusMapToken_

edm::ESGetToken<DTStatusFlag, DTStatusFlagRcd> DTSegmentSelector::theStatusMapToken_
private

Definition at line 43 of file DTSegmentSelector.h.

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