CMS 3D CMS Logo

GlobalMuonTrajectoryBuilder.cc
Go to the documentation of this file.
1 
25 
26 //---------------
27 // C++ Headers --
28 //---------------
29 
30 #include <iostream>
31 #include <iomanip>
32 #include <algorithm>
33 
34 //-------------------------------
35 // Collaborating Class Headers --
36 //-------------------------------
37 
40 
42 
44 
48 
49 using namespace std;
50 using namespace edm;
51 
52 //----------------
53 // Constructors --
54 //----------------
55 
60 
61 {
62  theTkTrackLabel = par.getParameter<edm::InputTag>("TrackerCollectionLabel");
64  thePrimaryVtxLabel = par.getParameter<edm::InputTag>("VertexCollectionLabel");
66  selectHighPurity_ = par.getParameter<bool>("selectHighPurity");
67 }
68 
69 //--------------
70 // Destructor --
71 //--------------
72 
74 
75 //
76 // get information from event
77 //
79  const std::string category = "Muon|RecoMuon|GlobalMuonTrajectoryBuilder|setEvent";
80 
82 
83  // get tracker TrackCollection from Event
84  event.getByToken(allTrackerTracksToken, allTrackerTracks);
85  LogDebug(category) << " Found " << allTrackerTracks->size() << " tracker Tracks with label " << theTkTrackLabel;
86 
87  // get primary vertex from Event
88  event.getByToken(primaryVertexToken, vertexCollection);
89  LogDebug(category) << " Found " << vertexCollection->size() << " tracker Tracks with label " << thePrimaryVtxLabel;
90 }
91 
92 //
93 // reconstruct trajectories
94 //
96  const std::string category = "Muon|RecoMuon|GlobalMuonTrajectoryBuilder|trajectories";
97 
98  // cut on muons with low momenta
99  LogTrace(category) << " STA pt " << staCandIn.second->pt() << " rho " << staCandIn.second->innerMomentum().Rho()
100  << " R " << staCandIn.second->innerMomentum().R() << " theCut " << thePtCut;
101 
102  // convert the STA track into a Trajectory if Trajectory not already present
103  TrackCand staCand(staCandIn);
104 
105  vector<TrackCand> regionalTkTracks = makeTkCandCollection(staCand);
106  LogTrace(category) << " Found " << regionalTkTracks.size() << " tracks within region of interest";
107 
108  // match tracker tracks to muon track
109  vector<TrackCand> trackerTracks = trackMatcher()->match(staCand, regionalTkTracks);
110  LogTrace(category) << " Found " << trackerTracks.size() << " matching tracker tracks within region of interest";
111 
112  if (trackerTracks.empty()) {
113  if (staCandIn.first == nullptr)
114  delete staCand.first;
115 
116  return CandidateContainer();
117  }
118 
119  // build a combined tracker-muon MuonCandidate
120  //
121  // turn tkMatchedTracks into MuonCandidates
122  //
123  LogTrace(category) << " Turn tkMatchedTracks into MuonCandidates";
124  CandidateContainer tkTrajs;
125  for (vector<TrackCand>::const_iterator tkt = trackerTracks.begin(); tkt != trackerTracks.end(); tkt++) {
126  tkTrajs.push_back(std::make_unique<MuonCandidate>(nullptr, staCand.second, (*tkt).second, nullptr));
127  }
128 
129  if (tkTrajs.empty()) {
130  LogTrace(category) << " tkTrajs empty";
131  if (staCandIn.first == nullptr)
132  delete staCand.first;
133 
134  return CandidateContainer();
135  }
136 
137  CandidateContainer result = build(staCand, tkTrajs);
138  LogTrace(category) << " Found " << result.size() << " GLBMuons from one STACand";
139 
140  // free memory
141  if (staCandIn.first == nullptr)
142  delete staCand.first;
143 
144  return result;
145 }
146 
147 //
148 // make a TrackCand collection using tracker Track, Trajectory information
149 //
150 vector<GlobalMuonTrajectoryBuilder::TrackCand> GlobalMuonTrajectoryBuilder::makeTkCandCollection(
151  const TrackCand& staCand) {
152  const std::string category = "Muon|RecoMuon|GlobalMuonTrajectoryBuilder|makeTkCandCollection";
153 
154  vector<TrackCand> tkCandColl;
155 
156  vector<TrackCand> tkTrackCands;
157 
158  for (unsigned int position = 0; position != allTrackerTracks->size(); ++position) {
160  TrackCand tkCand = TrackCand((Trajectory*)nullptr, tkTrackRef);
161  if (selectHighPurity_ && !tkTrackRef->quality(reco::TrackBase::highPurity)) {
162  const reco::VertexCollection* recoVertices = vertexCollection.product();
163  if (!(*recoVertices)[0].isFake())
164  continue;
165  }
166  tkTrackCands.push_back(tkCand);
167  }
168 
169  tkCandColl = chooseRegionalTrackerTracks(staCand, tkTrackCands);
170 
171  return tkCandColl;
172 }
MuonTrajectoryBuilder::CandidateContainer build(const TrackCand &, MuonTrajectoryBuilder::CandidateContainer &) const
build combined trajectory from sta Track and tracker RecHits
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
MuonCandidate::CandidateContainer CandidateContainer
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::pair< const Trajectory *, reco::TrackRef > TrackCand
EDGetTokenT< ProductType > mayConsume(edm::InputTag const &tag)
MuonTrajectoryBuilder::CandidateContainer trajectories(const TrackCand &) override
reconstruct trajectories from standalone and tracker only Tracks
void setEvent(const edm::Event &) override
pass the Event to the algo at each event
T const * product() const
Definition: Handle.h:70
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
std::vector< TrackCand > makeTkCandCollection(const TrackCand &) override
make a TrackCand collection using tracker Track, Trajectory information
std::vector< TrackCand > chooseRegionalTrackerTracks(const TrackCand &, const std::vector< TrackCand > &)
choose tracker tracks within region of interest
edm::EDGetTokenT< reco::TrackCollection > allTrackerTracksToken
#define LogTrace(id)
edm::EDGetTokenT< reco::VertexCollection > primaryVertexToken
~GlobalMuonTrajectoryBuilder() override
destructor
void setEvent(const edm::Event &) override
pass the Event to the algo at each event
double match(const TrackCand &sta, const TrackCand &track, int matchOption=0, int surfaceOption=1) const
GlobalMuonTrajectoryBuilder(const edm::ParameterSet &, const MuonServiceProxy *, edm::ConsumesCollector &)
constructor with Parameter Set and MuonServiceProxy
std::vector< std::unique_ptr< MuonCandidate > > CandidateContainer
Definition: MuonCandidate.h:18
edm::Handle< reco::TrackCollection > allTrackerTracks
edm::Handle< reco::VertexCollection > vertexCollection
HLT enums.
static int position[264][3]
Definition: ReadPGInfo.cc:289
GlobalMuonTrackMatcher * trackMatcher() const
Definition: event.py:1
#define LogDebug(id)