CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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  inputDTRecSegment4DCollection_(parset.getParameter<edm::InputTag>("inputDTRecSegment4DCollection")),
17  inputCSCSegmentCollection_(parset.getParameter<edm::InputTag>("inputCSCSegmentCollection"))
18 {
19 }
20 
22 }
23 
25 {
28 
29  edm::LogVerbatim("TrackerMuonHitExtractor") <<"\nThere are "<< dtSegmentCollectionH_->size()<<" DT segments.";
30  unsigned int index_dt_segment = 0;
32  segment != dtSegmentCollectionH_->end(); ++segment , index_dt_segment++) {
33  LocalPoint segmentLocalPosition = segment->localPosition();
34  LocalVector segmentLocalDirection = segment->localDirection();
35  LocalError segmentLocalPositionError = segment->localPositionError();
36  LocalError segmentLocalDirectionError = segment->localDirectionError();
37  DetId geoid = segment->geographicalId();
38  DTChamberId dtdetid = DTChamberId(geoid);
39  int wheel = dtdetid.wheel();
40  int station = dtdetid.station();
41  int sector = dtdetid.sector();
42 
43  float segmentX = segmentLocalPosition.x();
44  float segmentY = segmentLocalPosition.y();
45  float segmentdXdZ = segmentLocalDirection.x()/segmentLocalDirection.z();
46  float segmentdYdZ = segmentLocalDirection.y()/segmentLocalDirection.z();
47  float segmentXerr = sqrt(segmentLocalPositionError.xx());
48  float segmentYerr = sqrt(segmentLocalPositionError.yy());
49  float segmentdXdZerr = sqrt(segmentLocalDirectionError.xx());
50  float segmentdYdZerr = sqrt(segmentLocalDirectionError.yy());
51 
52  edm::LogVerbatim("TrackerMuonHitExtractor")
53  <<"\nDT segment index :"<<index_dt_segment
54  <<"\nchamber Wh:"<<wheel<<",St:"<<station<<",Se:"<<sector
55  <<"\nLocal Position (X,Y)=("<<segmentX<<","<<segmentY<<") +/- ("<<segmentXerr<<","<<segmentYerr<<"), "
56  <<"Local Direction (dXdZ,dYdZ)=("<<segmentdXdZ<<","<<segmentdYdZ<<") +/- ("<<segmentdXdZerr<<","<<segmentdYdZerr<<")";
57  }
58 
59  edm::LogVerbatim("TrackerMuonHitExtractor") <<"\nThere are "<< cscSegmentCollectionH_->size()<<" CSC segments.";
60  unsigned int index_csc_segment = 0;
62  segment != cscSegmentCollectionH_->end(); ++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
86  <<"\nchamber Endcap:"<<endcap<<",St:"<<station<<",Ri:"<<ring<<",Ch:"<<chamber
87  <<"\nLocal Position (X,Y)=("<<segmentX<<","<<segmentY<<") +/- ("<<segmentXerr<<","<<segmentYerr<<"), "
88  <<"Local Direction (dXdZ,dYdZ)=("<<segmentdXdZ<<","<<segmentdYdZ<<") +/- ("<<segmentdXdZerr<<","<<segmentdYdZerr<<")";
89  }
90 
91 }
92 std::vector<const TrackingRecHit *>
94  std::vector<const TrackingRecHit *> ret;
95 
96  int wheel, station, sector;
97  int endcap, /*station, */ ring, chamber;
98 
99  edm::LogVerbatim("TrackerMuonHitExtractor") <<"Number of chambers: "<<mu.matches().size()
101  unsigned int index_chamber = 0;
102 
103  for(std::vector<reco::MuonChamberMatch>::const_iterator chamberMatch = mu.matches().begin();
104  chamberMatch != mu.matches().end(); ++chamberMatch, index_chamber++) {
105  std::stringstream chamberStr;
106  chamberStr <<"\nchamber index: "<<index_chamber;
107 
108  int subdet = chamberMatch->detector();
109  DetId did = chamberMatch->id;
110 
111  if (subdet == MuonSubdetId::DT) {
112  DTChamberId dtdetid = DTChamberId(did);
113  wheel = dtdetid.wheel();
114  station = dtdetid.station();
115  sector = dtdetid.sector();
116  chamberStr << ", DT chamber Wh:"<<wheel<<",St:"<<station<<",Se:"<<sector;
117  }
118  else if (subdet == MuonSubdetId::CSC) {
119  CSCDetId cscdetid = CSCDetId(did);
120  endcap = cscdetid.endcap();
121  station = cscdetid.station();
122  ring = cscdetid.ring();
123  chamber = cscdetid.chamber();
124  chamberStr << ", CSC chamber End:"<<endcap<<",St:"<<station<<",Ri:"<<ring<<",Ch:"<<chamber;
125  }
126 
127  chamberStr << ", Number of segments: "<<chamberMatch->segmentMatches.size();
128  edm::LogVerbatim("TrackerMuonHitExtractor") << chamberStr.str();
129 
130  unsigned int index_segment = 0;
131 
132  for(std::vector<reco::MuonSegmentMatch>::const_iterator segmentMatch = chamberMatch->segmentMatches.begin();
133  segmentMatch != chamberMatch->segmentMatches.end(); ++segmentMatch, index_segment++) {
134 
135  float segmentX = segmentMatch->x;
136  float segmentY = segmentMatch->y ;
137  float segmentdXdZ = segmentMatch->dXdZ;
138  float segmentdYdZ = segmentMatch->dYdZ;
139  float segmentXerr = segmentMatch->xErr;
140  float segmentYerr = segmentMatch->yErr;
141  float segmentdXdZerr = segmentMatch->dXdZErr;
142  float segmentdYdZerr = segmentMatch->dYdZErr;
143 
144  CSCSegmentRef segmentCSC = segmentMatch->cscSegmentRef;
145  DTRecSegment4DRef segmentDT = segmentMatch->dtSegmentRef;
146 
147  bool segment_arbitrated_Ok = (segmentMatch->isMask(reco::MuonSegmentMatch::BestInChamberByDR) &&
148  segmentMatch->isMask(reco::MuonSegmentMatch::BelongsToTrackByDR));
149 
150  std::string ARBITRATED(" ***Arbitrated Off*** ");
151  if (segment_arbitrated_Ok) ARBITRATED = " ***ARBITRATED OK*** ";
152 
153  if (subdet == MuonSubdetId::DT) {
154  edm::LogVerbatim("TrackerMuonHitExtractor")
155  <<"\n\t segment index: "<<index_segment << ARBITRATED
156  <<"\n\t Local Position (X,Y)=("<<segmentX<<","<<segmentY<<") +/- ("<<segmentXerr<<","<<segmentYerr<<"), "
157  <<"\n\t Local Direction (dXdZ,dYdZ)=("<<segmentdXdZ<<","<<segmentdYdZ<<") +/- ("<<segmentdXdZerr<<","<<segmentdYdZerr<<")";
158 
159  if (!segment_arbitrated_Ok) continue;
160 
161  if (segmentDT.get() != 0) {
162  const DTRecSegment4D* segment = segmentDT.get();
163 
164  edm::LogVerbatim("TrackerMuonHitExtractor")<<"\t ===> MATCHING with DT segment with index = "<<segmentDT.key();
165 
166  if(segment->hasPhi()) {
167  const DTChamberRecSegment2D* phiSeg = segment->phiSegment();
168  std::vector<const TrackingRecHit*> phiHits = phiSeg->recHits();
169  for(std::vector<const TrackingRecHit*>::const_iterator ihit = phiHits.begin();
170  ihit != phiHits.end(); ++ihit) {
171  ret.push_back(*ihit);
172  }
173  }
174 
175  if(segment->hasZed()) {
176  const DTSLRecSegment2D* zSeg = (*segment).zSegment();
177  std::vector<const TrackingRecHit*> zedHits = zSeg->recHits();
178  for(std::vector<const TrackingRecHit*>::const_iterator ihit = zedHits.begin();
179  ihit != zedHits.end(); ++ihit) {
180  ret.push_back(*ihit);
181  }
182  }
183  } else edm::LogWarning("TrackerMuonHitExtractor")<<"\n***WARNING: UNMATCHED DT segment ! \n";
184  } // if (subdet == MuonSubdetId::DT)
185 
186  else if (subdet == MuonSubdetId::CSC) {
187  edm::LogVerbatim("TrackerMuonHitExtractor")
188  <<"\n\t segment index: "<<index_segment << ARBITRATED
189  <<"\n\t Local Position (X,Y)=("<<segmentX<<","<<segmentY<<") +/- ("<<segmentXerr<<","<<segmentYerr<<"), "
190  <<"\n\t Local Direction (dXdZ,dYdZ)=("<<segmentdXdZ<<","<<segmentdYdZ<<") +/- ("<<segmentdXdZerr<<","<<segmentdYdZerr<<")";
191 
192  if (!segment_arbitrated_Ok) continue;
193 
194  if (segmentCSC.get() != 0) {
195  const CSCSegment* segment = segmentCSC.get();
196 
197  edm::LogVerbatim("TrackerMuonHitExtractor")<<"\t ===> MATCHING with CSC segment with index = "<<segmentCSC.key();
198 
199  std::vector<const TrackingRecHit*> hits = segment->recHits();
200  for(std::vector<const TrackingRecHit*>::const_iterator ihit = hits.begin();
201  ihit != hits.end(); ++ihit) {
202  ret.push_back(*ihit);
203  }
204  } else edm::LogWarning("TrackerMuonHitExtractor")<<"\n***WARNING: UNMATCHED CSC segment ! \n";
205  } // else if (subdet == MuonSubdetId::CSC)
206 
207  } // loop on vector<MuonSegmentMatch>
208  } // loop on vector<MuonChamberMatch>
209 
210  return ret;
211 }
212 
int chamber() const
Definition: CSCDetId.h:70
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:62
edm::Handle< CSCSegmentCollection > cscSegmentCollectionH_
TrackerMuonHitExtractor(const edm::ParameterSet &)
int endcap() const
Definition: CSCDetId.h:95
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:45
int iEvent
Definition: GenABIO.cc:243
static const int CSC
Definition: MuonSubdetId.h:15
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:46
T z() const
Definition: PV3DBase.h:63
const int mu
Definition: Constants.h:23
static const unsigned int BestInChamberByDR
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:31
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
int numberOfMatches(ArbitrationType type=SegmentAndTrackArbitration) const
get number of chambers with matched segments
Definition: Muon.cc:56
void init(const edm::Event &, const edm::EventSetup &)
edm::Handle< DTRecSegment4DCollection > dtSegmentCollectionH_
int ring() const
Definition: CSCDetId.h:77
bool hasZed() const
Does it have the Z projection?
Definition: DetId.h:20
key_type key() const
Accessor for product key.
Definition: Ref.h:266
std::vector< MuonChamberMatch > & matches()
get muon matching information
Definition: Muon.h:141
static const unsigned int BelongsToTrackByDR
int sector() const
Definition: DTChamberId.h:63
int station() const
Definition: CSCDetId.h:88
static const int DT
Definition: MuonSubdetId.h:14
int station() const
Return the station number.
Definition: DTChamberId.h:53
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:242
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:47
T x() const
Definition: PV3DBase.h:61