CMS 3D CMS Logo

TrackerMuonHitExtractor.cc
Go to the documentation of this file.
1 //
2 // modified & integrated by Giovanni Abbiendi
3 // from code by Arun Luthra: UserCode/luthra/MuonTrackSelector/src/MuonTrackSelector.cc
4 //
13 #include <sstream>
14 
16  edm::ConsumesCollector && ic) :
17  inputDTRecSegment4DToken_(ic.consumes<DTRecSegment4DCollection>(parset.getParameter<edm::InputTag>("inputDTRecSegment4DCollection"))),
18  inputCSCSegmentToken_(ic.consumes<CSCSegmentCollection>(parset.getParameter<edm::InputTag>("inputCSCSegmentCollection"))),
19  inputDTRecSegment4DCollection_(parset.getParameter<edm::InputTag>("inputDTRecSegment4DCollection")),
20  inputCSCSegmentCollection_(parset.getParameter<edm::InputTag>("inputCSCSegmentCollection"))
21 {
22 }
23 
25  inputDTRecSegment4DCollection_(parset.getParameter<edm::InputTag>("inputDTRecSegment4DCollection")),
26  inputCSCSegmentCollection_(parset.getParameter<edm::InputTag>("inputCSCSegmentCollection"))
27 {
28 }
29 
31 }
32 
34 {
37 
38  edm::LogVerbatim("TrackerMuonHitExtractor") <<"\nThere are "<< dtSegmentCollectionH_->size()<<" DT segments.";
39  unsigned int index_dt_segment = 0;
41  segment != dtSegmentCollectionH_->end(); ++segment , index_dt_segment++) {
42  LocalPoint segmentLocalPosition = segment->localPosition();
43  LocalVector segmentLocalDirection = segment->localDirection();
44  LocalError segmentLocalPositionError = segment->localPositionError();
45  LocalError segmentLocalDirectionError = segment->localDirectionError();
46  DetId geoid = segment->geographicalId();
47  DTChamberId dtdetid = DTChamberId(geoid);
48  int wheel = dtdetid.wheel();
49  int station = dtdetid.station();
50  int sector = dtdetid.sector();
51 
52  float segmentX = segmentLocalPosition.x();
53  float segmentY = segmentLocalPosition.y();
54  float segmentdXdZ = segmentLocalDirection.x()/segmentLocalDirection.z();
55  float segmentdYdZ = segmentLocalDirection.y()/segmentLocalDirection.z();
56  float segmentXerr = sqrt(segmentLocalPositionError.xx());
57  float segmentYerr = sqrt(segmentLocalPositionError.yy());
58  float segmentdXdZerr = sqrt(segmentLocalDirectionError.xx());
59  float segmentdYdZerr = sqrt(segmentLocalDirectionError.yy());
60 
61  edm::LogVerbatim("TrackerMuonHitExtractor")
62  <<"\nDT segment index :"<<index_dt_segment
63  <<"\nchamber Wh:"<<wheel<<",St:"<<station<<",Se:"<<sector
64  <<"\nLocal Position (X,Y)=("<<segmentX<<","<<segmentY<<") +/- ("<<segmentXerr<<","<<segmentYerr<<"), "
65  <<"Local Direction (dXdZ,dYdZ)=("<<segmentdXdZ<<","<<segmentdYdZ<<") +/- ("<<segmentdXdZerr<<","<<segmentdYdZerr<<")";
66  }
67 
68  edm::LogVerbatim("TrackerMuonHitExtractor") <<"\nThere are "<< cscSegmentCollectionH_->size()<<" CSC segments.";
69  unsigned int index_csc_segment = 0;
71  segment != cscSegmentCollectionH_->end(); ++segment , index_csc_segment++) {
72  LocalPoint segmentLocalPosition = segment->localPosition();
73  LocalVector segmentLocalDirection = segment->localDirection();
74  LocalError segmentLocalPositionError = segment->localPositionError();
75  LocalError segmentLocalDirectionError = segment->localDirectionError();
76 
77  DetId geoid = segment->geographicalId();
78  CSCDetId cscdetid = CSCDetId(geoid);
79  int endcap = cscdetid.endcap();
80  int station = cscdetid.station();
81  int ring = cscdetid.ring();
82  int chamber = cscdetid.chamber();
83 
84  float segmentX = segmentLocalPosition.x();
85  float segmentY = segmentLocalPosition.y();
86  float segmentdXdZ = segmentLocalDirection.x()/segmentLocalDirection.z();
87  float segmentdYdZ = segmentLocalDirection.y()/segmentLocalDirection.z();
88  float segmentXerr = sqrt(segmentLocalPositionError.xx());
89  float segmentYerr = sqrt(segmentLocalPositionError.yy());
90  float segmentdXdZerr = sqrt(segmentLocalDirectionError.xx());
91  float segmentdYdZerr = sqrt(segmentLocalDirectionError.yy());
92 
93  edm::LogVerbatim("TrackerMuonHitExtractor")
94  <<"\nCSC segment index :"<<index_csc_segment
95  <<"\nchamber Endcap:"<<endcap<<",St:"<<station<<",Ri:"<<ring<<",Ch:"<<chamber
96  <<"\nLocal Position (X,Y)=("<<segmentX<<","<<segmentY<<") +/- ("<<segmentXerr<<","<<segmentYerr<<"), "
97  <<"Local Direction (dXdZ,dYdZ)=("<<segmentdXdZ<<","<<segmentdYdZ<<") +/- ("<<segmentdXdZerr<<","<<segmentdYdZerr<<")";
98  }
99 
100 }
101 std::vector<const TrackingRecHit *>
103  std::vector<const TrackingRecHit *> ret;
104 
105  int wheel, station, sector;
106  int endcap, /*station, */ ring, chamber;
107 
108  edm::LogVerbatim("TrackerMuonHitExtractor") <<"Number of chambers: "<<mu.matches().size()
110  unsigned int index_chamber = 0;
111 
112  for(std::vector<reco::MuonChamberMatch>::const_iterator chamberMatch = mu.matches().begin();
113  chamberMatch != mu.matches().end(); ++chamberMatch, index_chamber++) {
114  std::stringstream chamberStr;
115  chamberStr <<"\nchamber index: "<<index_chamber;
116 
117  int subdet = chamberMatch->detector();
118  DetId did = chamberMatch->id;
119 
120  if (subdet == MuonSubdetId::DT) {
121  DTChamberId dtdetid = DTChamberId(did);
122  wheel = dtdetid.wheel();
123  station = dtdetid.station();
124  sector = dtdetid.sector();
125  chamberStr << ", DT chamber Wh:"<<wheel<<",St:"<<station<<",Se:"<<sector;
126  }
127  else if (subdet == MuonSubdetId::CSC) {
128  CSCDetId cscdetid = CSCDetId(did);
129  endcap = cscdetid.endcap();
130  station = cscdetid.station();
131  ring = cscdetid.ring();
132  chamber = cscdetid.chamber();
133  chamberStr << ", CSC chamber End:"<<endcap<<",St:"<<station<<",Ri:"<<ring<<",Ch:"<<chamber;
134  }
135 
136  chamberStr << ", Number of segments: "<<chamberMatch->segmentMatches.size();
137  edm::LogVerbatim("TrackerMuonHitExtractor") << chamberStr.str();
138 
139  unsigned int index_segment = 0;
140 
141  for(std::vector<reco::MuonSegmentMatch>::const_iterator segmentMatch = chamberMatch->segmentMatches.begin();
142  segmentMatch != chamberMatch->segmentMatches.end(); ++segmentMatch, index_segment++) {
143 
144  float segmentX = segmentMatch->x;
145  float segmentY = segmentMatch->y ;
146  float segmentdXdZ = segmentMatch->dXdZ;
147  float segmentdYdZ = segmentMatch->dYdZ;
148  float segmentXerr = segmentMatch->xErr;
149  float segmentYerr = segmentMatch->yErr;
150  float segmentdXdZerr = segmentMatch->dXdZErr;
151  float segmentdYdZerr = segmentMatch->dYdZErr;
152 
153  CSCSegmentRef segmentCSC = segmentMatch->cscSegmentRef;
154  DTRecSegment4DRef segmentDT = segmentMatch->dtSegmentRef;
155 
156  bool segment_arbitrated_Ok = (segmentMatch->isMask(reco::MuonSegmentMatch::BestInChamberByDR) &&
157  segmentMatch->isMask(reco::MuonSegmentMatch::BelongsToTrackByDR));
158 
159  std::string ARBITRATED(" ***Arbitrated Off*** ");
160  if (segment_arbitrated_Ok) ARBITRATED = " ***ARBITRATED OK*** ";
161 
162  if (subdet == MuonSubdetId::DT) {
163  edm::LogVerbatim("TrackerMuonHitExtractor")
164  <<"\n\t segment index: "<<index_segment << ARBITRATED
165  <<"\n\t Local Position (X,Y)=("<<segmentX<<","<<segmentY<<") +/- ("<<segmentXerr<<","<<segmentYerr<<"), "
166  <<"\n\t Local Direction (dXdZ,dYdZ)=("<<segmentdXdZ<<","<<segmentdYdZ<<") +/- ("<<segmentdXdZerr<<","<<segmentdYdZerr<<")";
167 
168  if (!segment_arbitrated_Ok) continue;
169 
170  if (segmentDT.get() != 0) {
171  const DTRecSegment4D* segment = segmentDT.get();
172 
173  edm::LogVerbatim("TrackerMuonHitExtractor")<<"\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();
179  ihit != phiHits.end(); ++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();
188  ihit != zedHits.end(); ++ihit) {
189  ret.push_back(*ihit);
190  }
191  }
192  } else edm::LogWarning("TrackerMuonHitExtractor")<<"\n***WARNING: UNMATCHED DT segment ! \n";
193  } // if (subdet == MuonSubdetId::DT)
194 
195  else if (subdet == MuonSubdetId::CSC) {
196  edm::LogVerbatim("TrackerMuonHitExtractor")
197  <<"\n\t segment index: "<<index_segment << ARBITRATED
198  <<"\n\t Local Position (X,Y)=("<<segmentX<<","<<segmentY<<") +/- ("<<segmentXerr<<","<<segmentYerr<<"), "
199  <<"\n\t Local Direction (dXdZ,dYdZ)=("<<segmentdXdZ<<","<<segmentdYdZ<<") +/- ("<<segmentdXdZerr<<","<<segmentdYdZerr<<")";
200 
201  if (!segment_arbitrated_Ok) continue;
202 
203  if (segmentCSC.get() != 0) {
204  const CSCSegment* segment = segmentCSC.get();
205 
206  edm::LogVerbatim("TrackerMuonHitExtractor")<<"\t ===> MATCHING with CSC segment with index = "<<segmentCSC.key();
207 
208  std::vector<const TrackingRecHit*> hits = segment->recHits();
209  for(std::vector<const TrackingRecHit*>::const_iterator ihit = hits.begin();
210  ihit != hits.end(); ++ihit) {
211  ret.push_back(*ihit);
212  }
213  } else edm::LogWarning("TrackerMuonHitExtractor")<<"\n***WARNING: UNMATCHED CSC segment ! \n";
214  } // else if (subdet == MuonSubdetId::CSC)
215 
216  } // loop on vector<MuonSegmentMatch>
217  } // loop on vector<MuonChamberMatch>
218 
219  return ret;
220 }
221 
int chamber() const
Definition: CSCDetId.h:68
edm::InputTag inputDTRecSegment4DCollection_
float xx() const
Definition: LocalError.h:24
std::vector< const TrackingRecHit * > getMuonHits(const reco::Muon &mu) const
const DTChamberRecSegment2D * phiSegment() const
The superPhi segment: 0 if no phi projection available.
T y() const
Definition: PV3DBase.h:63
key_type key() const
Accessor for product key.
Definition: Ref.h:264
edm::Handle< CSCSegmentCollection > cscSegmentCollectionH_
int endcap() const
Definition: CSCDetId.h:93
int iEvent
Definition: GenABIO.cc:230
static const int CSC
Definition: MuonSubdetId.h:13
float yy() const
Definition: LocalError.h:26
virtual std::vector< const TrackingRecHit * > recHits() const
Access to component RecHits (if any)
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
const int mu
Definition: Constants.h:22
static const unsigned int BestInChamberByDR
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:244
bool hasPhi() const
Does it have the Phi projection?
virtual std::vector< const TrackingRecHit * > recHits() const
Access to component RecHits (if any)
Definition: CSCSegment.cc:30
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:413
void init(const edm::Event &, const edm::EventSetup &)
edm::Handle< DTRecSegment4DCollection > dtSegmentCollectionH_
int ring() const
Definition: CSCDetId.h:75
bool hasZed() const
Does it have the Z projection?
Definition: DetId.h:18
std::vector< MuonChamberMatch > & matches()
get muon matching information
Definition: Muon.h:144
int numberOfMatches(ArbitrationType type=SegmentAndTrackArbitration) const
get number of chambers with matched segments
static const unsigned int BelongsToTrackByDR
TrackerMuonHitExtractor(const edm::ParameterSet &, edm::ConsumesCollector &&ic)
HLT enums.
int sector() const
Definition: DTChamberId.h:61
int station() const
Definition: CSCDetId.h:86
static const int DT
Definition: MuonSubdetId.h:12
int station() const
Return the station number.
Definition: DTChamberId.h:51
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
T x() const
Definition: PV3DBase.h:62