CMS 3D CMS Logo

L3MuonTrajectoryBuilder.cc
Go to the documentation of this file.
1 
24 
27 
29 
32 
36 
38 
41 
44 
49 
50 
51 //----------------
52 // Constructors --
53 //----------------
56  edm::ConsumesCollector& iC) : GlobalTrajectoryBuilderBase(par, service, iC) {
58  theTkCollName = par.getParameter<edm::InputTag>("tkTrajLabel");
59  theBeamSpotInputTag = par.getParameter<edm::InputTag>("tkTrajBeamSpot");
60  theMaxChi2 = par.getParameter<double>("tkTrajMaxChi2");
61  theDXYBeamSpot = par.getParameter<double>("tkTrajMaxDXYBeamSpot");
62  theUseVertex = par.getParameter<bool>("tkTrajUseVertex");
63  theVertexCollInputTag = par.getParameter<edm::InputTag>("tkTrajVertex");
65 }
66 
67 //--------------
68 // Destructor --
69 //--------------
72 }
73 
77  desc.add("MuonTrackingRegionBuilder",descTRB);
78 }
79 
80 //
81 // Get information from event
82 //
84  const std::string category = "Muon|RecoMuon|L3MuonTrajectoryBuilder|setEvent";
85 
87 
88  // get tracker TrackCollection from Event
89  event.getByToken(theTrackToken,allTrackerTracks);
90  LogDebug(category)
91  << "Found " << allTrackerTracks->size()
92  << " tracker Tracks with label "<< theTkCollName;
93 
94  if( theUseVertex ) {
95  // PV
97  if ( pvHandle.isValid() ) {
98  vtx = pvHandle->front();
99  }
100  else {
101  edm::LogInfo(category) << "No Primary Vertex available from EventSetup \n";
102  }
103  }
104  else {
105  // BS
106  event.getByLabel(theBeamSpotInputTag, beamSpotHandle);
107  if( beamSpotHandle.isValid() ) {
109  }
110  else {
111  edm::LogInfo(category) << "No beam spot available from EventSetup \n";
112  }
113  }
114 }
115 
116 
117 //
118 // reconstruct trajectories
119 //
121 
122  const std::string category = "Muon|RecoMuon|L3MuonTrajectoryBuilder|trajectories";
123 
124  // cut on muons with low momenta
125  if ( (staCandIn).second->pt() < thePtCut || (staCandIn).second->innerMomentum().Rho() < thePtCut || (staCandIn).second->innerMomentum().R() < 2.5 ) return CandidateContainer();
126 
127  // convert the STA track into a Trajectory if Trajectory not already present
128  TrackCand staCand(staCandIn);
129 
130  std::vector<TrackCand> trackerTracks;
131 
132  std::vector<TrackCand> regionalTkTracks = makeTkCandCollection(staCand);
133  LogDebug(category) << "Found " << regionalTkTracks.size() << " tracks within region of interest";
134 
135  // match tracker tracks to muon track
136  trackerTracks = trackMatcher()->match(staCand, regionalTkTracks);
137 
138  LogDebug(category) << "Found " << trackerTracks.size() << " matching tracker tracks within region of interest";
139  if ( trackerTracks.empty() ) return CandidateContainer();
140 
141  // build a combined tracker-muon MuonCandidate
142  // turn tkMatchedTracks into MuonCandidates
143  LogDebug(category) << "turn tkMatchedTracks into MuonCandidates";
144  CandidateContainer tkTrajs;
145  for (std::vector<TrackCand>::const_iterator tkt = trackerTracks.begin(); tkt != trackerTracks.end(); tkt++) {
146  if ((*tkt).first != nullptr && (*tkt).first->isValid()) {
147  MuonCandidate* muonCand = new MuonCandidate( nullptr ,staCand.second,(*tkt).second, new Trajectory(*(*tkt).first));
148  tkTrajs.push_back(muonCand);
149  } else {
150  MuonCandidate* muonCand = new MuonCandidate( nullptr ,staCand.second,(*tkt).second, nullptr);
151  tkTrajs.push_back(muonCand);
152  }
153  }
154 
155  if ( tkTrajs.empty() ) {
156  LogDebug(category) << "tkTrajs empty";
157  return CandidateContainer();
158  }
159 
160  CandidateContainer result = build(staCand, tkTrajs);
161  LogDebug(category) << "Found "<< result.size() << " L3Muons from one L2Cand";
162 
163  // free memory
164  if ( staCandIn.first == nullptr) delete staCand.first;
165 
166  for( CandidateContainer::const_iterator it = tkTrajs.begin(); it != tkTrajs.end(); ++it) {
167  if ( (*it)->trajectory() ) delete (*it)->trajectory();
168  if ( (*it)->trackerTrajectory() ) delete (*it)->trackerTrajectory();
169  if ( *it ) delete (*it);
170  }
171  tkTrajs.clear();
172 
173  for ( std::vector<TrackCand>::const_iterator is = regionalTkTracks.begin(); is != regionalTkTracks.end(); ++is) {
174  delete (*is).first;
175  }
176 
177  return result;
178 }
179 
180 
181 //
182 // make a TrackCand collection using tracker Track, Trajectory information
183 //
184 std::vector<L3MuonTrajectoryBuilder::TrackCand> L3MuonTrajectoryBuilder::makeTkCandCollection(const TrackCand& staCand) {
185  const std::string category = "Muon|RecoMuon|L3MuonTrajectoryBuilder|makeTkCandCollection";
186  std::vector<TrackCand> tkCandColl;
187  std::vector<TrackCand> tkTrackCands;
188 
189 // for (auto&& tkTrack: allTrackerTracks){
190 // auto tkCand = TrackCand((Trajectory*)(0),tkTrack);
191  for ( unsigned int position = 0; position != allTrackerTracks->size(); ++position ) {
193  TrackCand tkCand = TrackCand((Trajectory*)nullptr,tkTrackRef);
194  tkCandColl.push_back(tkCand);
195  }
196 
197  //Loop over TrackCand collection made from allTrackerTracks in previous step
198  for (auto&& tkCand: tkCandColl){
199  auto& tk = tkCand.second;
200  bool canUseL3MTS = false;
201  // check the seedRef is non-null first; and then
202  if (tk->seedRef().isNonnull()){
203  auto a = dynamic_cast<const L3MuonTrajectorySeed*>(tk->seedRef().get());
204  canUseL3MTS = a != nullptr;
205  }
206  if (canUseL3MTS){
208  // May still need provenance here, so using trackref:
209  reco::TrackRef staTrack = l3seedRef->l2Track();
210  if( staTrack == (staCand.second) ) {
211  // Apply filters (dxy, chi2 cut)
212  double tk_vtx;
213  if( theUseVertex ) tk_vtx = tk->dxy(vtx.position());
214  else tk_vtx = tk->dxy(beamSpot.position());
215  if( fabs(tk_vtx) > theDXYBeamSpot || tk->normalizedChi2() > theMaxChi2 ) continue;
216  tkTrackCands.push_back(tkCand);
217  }
218  }
219  else{
220  // We will try to match all tracker tracks with the muon:
221  double tk_vtx;
222  if( theUseVertex ) tk_vtx = tk->dxy(vtx.position());
223  else tk_vtx = tk->dxy(beamSpot.position());
224  if( fabs(tk_vtx) > theDXYBeamSpot || tk->normalizedChi2() > theMaxChi2 ) continue;
225  tkTrackCands.push_back(tkCand);
226  }
227  }
228 
229  return tkTrackCands;
230 }
231 
#define LogDebug(id)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
MuonCandidate::CandidateContainer CandidateContainer
GlobalMuonTrackMatcher * trackMatcher() const
MuonTrajectoryBuilder::CandidateContainer trajectories(const TrackCand &) override
return a container reconstructed muons starting from a given track
std::pair< const Trajectory *, reco::TrackRef > TrackCand
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
static void fillDescriptions(edm::ParameterSetDescription &descriptions)
Add default values for fillDescriptions.
const Point & position() const
position
Definition: Vertex.h:109
U second(std::pair< T, U > const &p)
edm::Handle< reco::TrackCollection > allTrackerTracks
edm::Handle< reco::VertexCollection > pvHandle
void setEvent(const edm::Event &) override
Pass the Event to the algo at each event.
void setEvent(const edm::Event &) override
pass the Event to the algo at each event
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::vector< MuonCandidate * > CandidateContainer
Definition: MuonCandidate.h:20
double match(const TrackCand &sta, const TrackCand &track, int matchOption=0, int surfaceOption=1) const
bool isValid() const
Definition: HandleBase.h:74
edm::EDGetTokenT< reco::TrackCollection > theTrackToken
MuonTrajectoryBuilder::CandidateContainer build(const TrackCand &, MuonTrajectoryBuilder::CandidateContainer &) const
build combined trajectory from sta Track and tracker RecHits
~L3MuonTrajectoryBuilder() override
Destructor.
TrajectoryCleaner * theTrajectoryCleaner
edm::Handle< reco::BeamSpot > beamSpotHandle
double a
Definition: hdecay.h:121
static int position[264][3]
Definition: ReadPGInfo.cc:509
const Point & position() const
position
Definition: BeamSpot.h:62
std::vector< TrackCand > makeTkCandCollection(const TrackCand &) override
Make a TrackCand collection using tracker Track, Trajectory information.
L3MuonTrajectoryBuilder(const edm::ParameterSet &, const MuonServiceProxy *, edm::ConsumesCollector &)
Constructor with Parameter Set and MuonServiceProxy.
static void fillDescriptionsHLT(edm::ParameterSetDescription &descriptions)
Definition: event.py:1