CMS 3D CMS Logo

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

#include <TrackerMuonHitExtractor.h>

Public Member Functions

std::vector< const TrackingRecHit * > getMuonHits (const reco::Muon &mu) const
 
void init (const edm::Event &)
 
 TrackerMuonHitExtractor (const edm::ParameterSet &, edm::ConsumesCollector &&ic)
 
 ~TrackerMuonHitExtractor ()=default
 

Private Attributes

const edm::EDGetTokenT< CSCSegmentCollectioninputCSCSegmentToken_
 
const edm::EDGetTokenT< DTRecSegment4DCollectioninputDTRecSegment4DToken_
 

Detailed Description

Definition at line 18 of file TrackerMuonHitExtractor.h.

Constructor & Destructor Documentation

◆ TrackerMuonHitExtractor()

TrackerMuonHitExtractor::TrackerMuonHitExtractor ( const edm::ParameterSet parset,
edm::ConsumesCollector &&  ic 
)
explicit

Definition at line 16 of file TrackerMuonHitExtractor.cc.

18  ic.consumes<DTRecSegment4DCollection>(parset.getParameter<edm::InputTag>("inputDTRecSegment4DCollection"))),
20  ic.consumes<CSCSegmentCollection>(parset.getParameter<edm::InputTag>("inputCSCSegmentCollection"))) {}
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const edm::EDGetTokenT< DTRecSegment4DCollection > inputDTRecSegment4DToken_
const edm::EDGetTokenT< CSCSegmentCollection > inputCSCSegmentToken_

◆ ~TrackerMuonHitExtractor()

TrackerMuonHitExtractor::~TrackerMuonHitExtractor ( )
default

Member Function Documentation

◆ getMuonHits()

std::vector< const TrackingRecHit * > TrackerMuonHitExtractor::getMuonHits ( const reco::Muon mu) const

Definition at line 92 of file TrackerMuonHitExtractor.cc.

References reco::MuonSegmentMatch::BelongsToTrackByDR, reco::MuonSegmentMatch::BestInChamberByDR, relativeConstraints::chamber, CSCDetId::chamber(), MuonSubdetId::CSC, MuonSubdetId::DT, CSCDetId::endcap(), makeMuonMisalignmentScenario::endcap, edm::Ref< C, T, F >::get(), DTRecSegment4D::hasPhi(), DTRecSegment4D::hasZed(), hfClusterShapes_cfi::hits, edm::Ref< C, T, F >::key(), amptDefaultParameters_cff::mu, reco::Muon::NoArbitration, DTRecSegment4D::phiSegment(), CSCSegment::recHits(), DTRecSegment2D::recHits(), runTheMatrix::ret, CSCDetId::ring(), relativeConstraints::ring, hgcalTBTopologyTester_cfi::sector, DTChamberId::sector(), reco::Muon::SegmentAndTrackArbitration, DTChamberId::station(), relativeConstraints::station, CSCDetId::station(), AlCaHLTBitMon_QueryRunRegistry::string, DTChamberId::wheel(), and makeMuonMisalignmentScenario::wheel.

Referenced by MuonToTrackingParticleAssociatorByHitsImpl::associateMuons().

92  {
93  std::vector<const TrackingRecHit *> ret;
94 
95  int wheel, station, sector;
96  int endcap, /*station, */ ring, chamber;
97 
98  edm::LogVerbatim("TrackerMuonHitExtractor")
99  << "Number of chambers: " << mu.matches().size()
100  << ", arbitrated: " << mu.numberOfMatches(reco::Muon::SegmentAndTrackArbitration)
101  << ", tot (no arbitration): " << mu.numberOfMatches(reco::Muon::NoArbitration);
102  unsigned int index_chamber = 0;
103  int n_segments_noArb = 0;
104 
105  for (std::vector<reco::MuonChamberMatch>::const_iterator chamberMatch = mu.matches().begin();
106  chamberMatch != mu.matches().end();
107  ++chamberMatch, index_chamber++) {
108  std::stringstream chamberStr;
109  chamberStr << "\nchamber index: " << index_chamber;
110 
111  int subdet = chamberMatch->detector();
112  DetId did = chamberMatch->id;
113 
114  if (subdet == MuonSubdetId::DT) {
115  DTChamberId dtdetid = DTChamberId(did);
116  wheel = dtdetid.wheel();
117  station = dtdetid.station();
118  sector = dtdetid.sector();
119  chamberStr << ", DT chamber Wh:" << wheel << ",St:" << station << ",Se:" << sector;
120  } else if (subdet == MuonSubdetId::CSC) {
121  CSCDetId cscdetid = CSCDetId(did);
122  endcap = cscdetid.endcap();
123  station = cscdetid.station();
124  ring = cscdetid.ring();
125  chamber = cscdetid.chamber();
126  chamberStr << ", CSC chamber End:" << endcap << ",St:" << station << ",Ri:" << ring << ",Ch:" << chamber;
127  }
128 
129  chamberStr << ", Number of segments: " << chamberMatch->segmentMatches.size();
130  edm::LogVerbatim("TrackerMuonHitExtractor") << chamberStr.str();
131  n_segments_noArb = n_segments_noArb + chamberMatch->segmentMatches.size();
132 
133  unsigned int index_segment = 0;
134 
135  for (std::vector<reco::MuonSegmentMatch>::const_iterator segmentMatch = chamberMatch->segmentMatches.begin();
136  segmentMatch != chamberMatch->segmentMatches.end();
137  ++segmentMatch, index_segment++) {
138  float segmentX = segmentMatch->x;
139  float segmentY = segmentMatch->y;
140  float segmentdXdZ = segmentMatch->dXdZ;
141  float segmentdYdZ = segmentMatch->dYdZ;
142  float segmentXerr = segmentMatch->xErr;
143  float segmentYerr = segmentMatch->yErr;
144  float segmentdXdZerr = segmentMatch->dXdZErr;
145  float segmentdYdZerr = segmentMatch->dYdZErr;
146 
147  CSCSegmentRef segmentCSC = segmentMatch->cscSegmentRef;
148  DTRecSegment4DRef segmentDT = segmentMatch->dtSegmentRef;
149 
150  bool segment_arbitrated_Ok = (segmentMatch->isMask(reco::MuonSegmentMatch::BestInChamberByDR) &&
151  segmentMatch->isMask(reco::MuonSegmentMatch::BelongsToTrackByDR));
152 
153  std::string ARBITRATED(" ***Arbitrated Off*** ");
154  if (segment_arbitrated_Ok)
155  ARBITRATED = " ***ARBITRATED OK*** ";
156 
157  if (subdet == MuonSubdetId::DT) {
158  edm::LogVerbatim("TrackerMuonHitExtractor")
159  << "\n\t segment index: " << index_segment << ARBITRATED << "\n\t Local Position (X,Y)=(" << segmentX
160  << "," << segmentY << ") +/- (" << segmentXerr << "," << segmentYerr << "), "
161  << "\n\t Local Direction (dXdZ,dYdZ)=(" << segmentdXdZ << "," << segmentdYdZ << ") +/- (" << segmentdXdZerr
162  << "," << segmentdYdZerr << ")";
163 
164  // with the following line the DT segments failing standard arbitration
165  // are skipped
166  if (!segment_arbitrated_Ok)
167  continue;
168 
169  if (segmentDT.get() != nullptr) {
170  const DTRecSegment4D *segment = segmentDT.get();
171 
172  edm::LogVerbatim("TrackerMuonHitExtractor")
173  << "\t ===> MATCHING with DT segment with index = " << segmentDT.key();
174 
175  if (segment->hasPhi()) {
176  const DTChamberRecSegment2D *phiSeg = segment->phiSegment();
177  std::vector<const TrackingRecHit *> phiHits = phiSeg->recHits();
178  for (std::vector<const TrackingRecHit *>::const_iterator ihit = phiHits.begin(); ihit != phiHits.end();
179  ++ihit) {
180  ret.push_back(*ihit);
181  }
182  }
183 
184  if (segment->hasZed()) {
185  const DTSLRecSegment2D *zSeg = (*segment).zSegment();
186  std::vector<const TrackingRecHit *> zedHits = zSeg->recHits();
187  for (std::vector<const TrackingRecHit *>::const_iterator ihit = zedHits.begin(); ihit != zedHits.end();
188  ++ihit) {
189  ret.push_back(*ihit);
190  }
191  }
192  } else
193  edm::LogWarning("TrackerMuonHitExtractor") << "\n***WARNING: UNMATCHED DT segment ! \n";
194  } // if (subdet == MuonSubdetId::DT)
195 
196  else if (subdet == MuonSubdetId::CSC) {
197  edm::LogVerbatim("TrackerMuonHitExtractor")
198  << "\n\t segment index: " << index_segment << ARBITRATED << "\n\t Local Position (X,Y)=(" << segmentX
199  << "," << segmentY << ") +/- (" << segmentXerr << "," << segmentYerr << "), "
200  << "\n\t Local Direction (dXdZ,dYdZ)=(" << segmentdXdZ << "," << segmentdYdZ << ") +/- (" << segmentdXdZerr
201  << "," << segmentdYdZerr << ")";
202 
203  // with the following line the CSC segments failing standard arbitration
204  // are skipped
205  if (!segment_arbitrated_Ok)
206  continue;
207 
208  if (segmentCSC.get() != nullptr) {
209  const CSCSegment *segment = segmentCSC.get();
210 
211  edm::LogVerbatim("TrackerMuonHitExtractor")
212  << "\t ===> MATCHING with CSC segment with index = " << segmentCSC.key();
213 
214  std::vector<const TrackingRecHit *> hits = segment->recHits();
215  for (std::vector<const TrackingRecHit *>::const_iterator ihit = hits.begin(); ihit != hits.end(); ++ihit) {
216  ret.push_back(*ihit);
217  }
218  } else
219  edm::LogWarning("TrackerMuonHitExtractor") << "\n***WARNING: UNMATCHED CSC segment ! \n";
220  } // else if (subdet == MuonSubdetId::CSC)
221 
222  } // loop on vector<MuonSegmentMatch>
223  } // loop on vector<MuonChamberMatch>
224 
225  edm::LogVerbatim("TrackerMuonHitExtractor") << "\n N. matched Segments before arbitration = " << n_segments_noArb;
226 
227  return ret;
228 }
Log< level::Info, true > LogVerbatim
int station() const
Return the station number.
Definition: DTChamberId.h:42
std::vector< const TrackingRecHit * > recHits() const override
Access to component RecHits (if any)
Definition: CSCSegment.cc:32
bool hasPhi() const
Does it have the Phi projection?
ret
prodAgent to be discontinued
key_type key() const
Accessor for product key.
Definition: Ref.h:250
int chamber() const
Definition: CSCDetId.h:62
std::vector< const TrackingRecHit * > recHits() const override
Access to component RecHits (if any)
static const unsigned int BestInChamberByDR
Definition: DetId.h:17
int station() const
Definition: CSCDetId.h:79
const DTChamberRecSegment2D * phiSegment() const
The superPhi segment: 0 if no phi projection available.
int endcap() const
Definition: CSCDetId.h:85
static const unsigned int BelongsToTrackByDR
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
int sector() const
Definition: DTChamberId.h:49
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
static constexpr int DT
Definition: MuonSubdetId.h:11
int ring() const
Definition: CSCDetId.h:68
Log< level::Warning, false > LogWarning
bool hasZed() const
Does it have the Z projection?
static constexpr int CSC
Definition: MuonSubdetId.h:12

◆ init()

void TrackerMuonHitExtractor::init ( const edm::Event iEvent)

Definition at line 22 of file TrackerMuonHitExtractor.cc.

References relativeConstraints::chamber, CSCDetId::chamber(), CSCDetId::endcap(), makeMuonMisalignmentScenario::endcap, iEvent, inputCSCSegmentToken_, inputDTRecSegment4DToken_, relativeConstraints::ring, CSCDetId::ring(), hgcalTBTopologyTester_cfi::sector, DTChamberId::sector(), mathSSE::sqrt(), DTChamberId::station(), relativeConstraints::station, CSCDetId::station(), DTChamberId::wheel(), makeMuonMisalignmentScenario::wheel, PV3DBase< T, PVType, FrameType >::x(), LocalError::xx(), PV3DBase< T, PVType, FrameType >::y(), LocalError::yy(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by MuonToTrackingParticleAssociatorEDProducer::produce().

22  {
23  const edm::Handle<DTRecSegment4DCollection> &dtSegmentCollectionH = iEvent.getHandle(inputDTRecSegment4DToken_);
24  const edm::Handle<CSCSegmentCollection> &cscSegmentCollectionH = iEvent.getHandle(inputCSCSegmentToken_);
25 
26  edm::LogVerbatim("TrackerMuonHitExtractor") << "\nThere are " << dtSegmentCollectionH->size() << " DT segments.";
27  unsigned int index_dt_segment = 0;
28  for (DTRecSegment4DCollection::const_iterator segment = dtSegmentCollectionH->begin();
29  segment != dtSegmentCollectionH->end();
30  ++segment, index_dt_segment++) {
31  LocalPoint segmentLocalPosition = segment->localPosition();
32  LocalVector segmentLocalDirection = segment->localDirection();
33  LocalError segmentLocalPositionError = segment->localPositionError();
34  LocalError segmentLocalDirectionError = segment->localDirectionError();
35  DetId geoid = segment->geographicalId();
36  DTChamberId dtdetid = DTChamberId(geoid);
37  int wheel = dtdetid.wheel();
38  int station = dtdetid.station();
39  int sector = dtdetid.sector();
40 
41  float segmentX = segmentLocalPosition.x();
42  float segmentY = segmentLocalPosition.y();
43  float segmentdXdZ = segmentLocalDirection.x() / segmentLocalDirection.z();
44  float segmentdYdZ = segmentLocalDirection.y() / segmentLocalDirection.z();
45  float segmentXerr = sqrt(segmentLocalPositionError.xx());
46  float segmentYerr = sqrt(segmentLocalPositionError.yy());
47  float segmentdXdZerr = sqrt(segmentLocalDirectionError.xx());
48  float segmentdYdZerr = sqrt(segmentLocalDirectionError.yy());
49 
50  edm::LogVerbatim("TrackerMuonHitExtractor")
51  << "\nDT segment index :" << index_dt_segment << "\nchamber Wh:" << wheel << ",St:" << station
52  << ",Se:" << sector << "\nLocal Position (X,Y)=(" << segmentX << "," << segmentY << ") +/- (" << segmentXerr
53  << "," << segmentYerr << "), "
54  << "Local Direction (dXdZ,dYdZ)=(" << segmentdXdZ << "," << segmentdYdZ << ") +/- (" << segmentdXdZerr << ","
55  << segmentdYdZerr << ")";
56  }
57 
58  edm::LogVerbatim("TrackerMuonHitExtractor") << "\nThere are " << cscSegmentCollectionH->size() << " CSC segments.";
59  unsigned int index_csc_segment = 0;
60  for (CSCSegmentCollection::const_iterator segment = cscSegmentCollectionH->begin();
61  segment != cscSegmentCollectionH->end();
62  ++segment, index_csc_segment++) {
63  LocalPoint segmentLocalPosition = segment->localPosition();
64  LocalVector segmentLocalDirection = segment->localDirection();
65  LocalError segmentLocalPositionError = segment->localPositionError();
66  LocalError segmentLocalDirectionError = segment->localDirectionError();
67 
68  DetId geoid = segment->geographicalId();
69  CSCDetId cscdetid = CSCDetId(geoid);
70  int endcap = cscdetid.endcap();
71  int station = cscdetid.station();
72  int ring = cscdetid.ring();
73  int chamber = cscdetid.chamber();
74 
75  float segmentX = segmentLocalPosition.x();
76  float segmentY = segmentLocalPosition.y();
77  float segmentdXdZ = segmentLocalDirection.x() / segmentLocalDirection.z();
78  float segmentdYdZ = segmentLocalDirection.y() / segmentLocalDirection.z();
79  float segmentXerr = sqrt(segmentLocalPositionError.xx());
80  float segmentYerr = sqrt(segmentLocalPositionError.yy());
81  float segmentdXdZerr = sqrt(segmentLocalDirectionError.xx());
82  float segmentdYdZerr = sqrt(segmentLocalDirectionError.yy());
83 
84  edm::LogVerbatim("TrackerMuonHitExtractor")
85  << "\nCSC segment index :" << index_csc_segment << "\nchamber Endcap:" << endcap << ",St:" << station
86  << ",Ri:" << ring << ",Ch:" << chamber << "\nLocal Position (X,Y)=(" << segmentX << "," << segmentY << ") +/- ("
87  << segmentXerr << "," << segmentYerr << "), "
88  << "Local Direction (dXdZ,dYdZ)=(" << segmentdXdZ << "," << segmentdYdZ << ") +/- (" << segmentdXdZerr << ","
89  << segmentdYdZerr << ")";
90  }
91 }
Log< level::Info, true > LogVerbatim
int station() const
Return the station number.
Definition: DTChamberId.h:42
const edm::EDGetTokenT< DTRecSegment4DCollection > inputDTRecSegment4DToken_
T z() const
Definition: PV3DBase.h:61
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
float yy() const
Definition: LocalError.h:24
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:19
int chamber() const
Definition: CSCDetId.h:62
const edm::EDGetTokenT< CSCSegmentCollection > inputCSCSegmentToken_
Definition: DetId.h:17
int station() const
Definition: CSCDetId.h:79
int endcap() const
Definition: CSCDetId.h:85
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
int sector() const
Definition: DTChamberId.h:49
int ring() const
Definition: CSCDetId.h:68
float xx() const
Definition: LocalError.h:22

Member Data Documentation

◆ inputCSCSegmentToken_

const edm::EDGetTokenT<CSCSegmentCollection> TrackerMuonHitExtractor::inputCSCSegmentToken_
private

Definition at line 28 of file TrackerMuonHitExtractor.h.

Referenced by init().

◆ inputDTRecSegment4DToken_

const edm::EDGetTokenT<DTRecSegment4DCollection> TrackerMuonHitExtractor::inputDTRecSegment4DToken_
private

Definition at line 27 of file TrackerMuonHitExtractor.h.

Referenced by init().