CMS 3D CMS Logo

HLTMuonL2SelectorForL3IO.cc
Go to the documentation of this file.
1 
15 
18  : l2Src_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("l2Src"))),
19  l3OISrc_(consumes<reco::RecoChargedCandidateCollection>(iConfig.getParameter<edm::InputTag>("l3OISrc"))),
20  l3linkToken_(consumes<reco::MuonTrackLinksCollection>(iConfig.getParameter<edm::InputTag>("InputLinks"))),
21  applyL3Filters_(iConfig.getParameter<bool>("applyL3Filters")),
22  max_NormalizedChi2_(iConfig.getParameter<double>("MaxNormalizedChi2")),
23  max_PtDifference_(iConfig.getParameter<double>("MaxPtDifference")),
24  min_Nhits_(iConfig.getParameter<int>("MinNhits")),
25  min_NmuonHits_(iConfig.getParameter<int>("MinNmuonHits")) {
26  LogTrace("Muon|RecoMuon|HLTMuonL2SelectorForL3IO") << "constructor called";
27  produces<reco::TrackCollection>();
28 }
29 
32 
35  const std::string metname = "Muon|RecoMuon|HLTMuonL2SelectorForL3IO";
36 
37  // IN:
39  iEvent.getByToken(l2Src_, l2muonH);
40 
42  iEvent.getByToken(l3OISrc_, l3muonH);
43 
44  // Read Links collection:
46  iEvent.getByToken(l3linkToken_, links);
47 
48  // OUT:
49  std::unique_ptr<reco::TrackCollection> result(new reco::TrackCollection());
50 
51  for (unsigned int il2 = 0; il2 != l2muonH->size(); ++il2) {
52  reco::TrackRef l2muRef(l2muonH, il2);
53  bool re_do_this_L2 = true;
54 
55  for (unsigned int il3 = 0; il3 != l3muonH->size(); ++il3) {
56  reco::RecoChargedCandidateRef cand(l3muonH, il3);
57  reco::TrackRef tk = cand->track();
58 
59  bool useThisLink = false;
60  for (unsigned int l(0); l < links->size() && !useThisLink; ++l) {
61  const reco::MuonTrackLinks* link = &links->at(l);
62 
63  // Check if the L3 link matches the L3 candidate
64  const reco::Track& globalTrack = *link->globalTrack();
65  float dR2 = deltaR2(tk->eta(), tk->phi(), globalTrack.eta(), globalTrack.phi());
66  if (dR2 < 0.02 * 0.02 and std::abs(tk->pt() - globalTrack.pt()) < 0.001 * tk->pt()) {
67  useThisLink = true;
68  }
69 
70  if (!useThisLink)
71  continue;
72 
73  // Check whether the stand-alone track matches a L2, if not, we will re-use this L2
74  const reco::TrackRef staTrack = link->standAloneTrack();
75  if (l2muRef == staTrack)
76  re_do_this_L2 = false;
77 
78  // Check the quality of the reconstructed L3, if poor quality, we will re-use this L2
79  if (staTrack == l2muRef && applyL3Filters_) {
80  re_do_this_L2 = true;
81  const reco::Track& globalTrack = *link->globalTrack();
82  if (globalTrack.numberOfValidHits() < min_Nhits_)
83  continue; // cut on number of hits
84  if (globalTrack.normalizedChi2() > max_NormalizedChi2_)
85  continue; //normalizedChi2 cut
86  if (globalTrack.hitPattern().numberOfValidMuonHits() < min_NmuonHits_)
87  continue; //min muon hits cut
88  if (std::abs(globalTrack.pt() - l2muRef->pt()) > max_PtDifference_ * globalTrack.pt())
89  continue; // pt difference
90  re_do_this_L2 = false;
91  }
92  }
93  }
94  if (re_do_this_L2)
95  result->push_back(*l2muRef); // used the L2 if no L3 if matched or if the matched L3 has poor quality cuts.
96  }
97  iEvent.put(std::move(result));
98 }
99 
102  desc.add<edm::InputTag>("l2Src", edm::InputTag("hltL2Muons", "UpdatedAtVtx"));
103  desc.add<edm::InputTag>("l3OISrc", edm::InputTag("hltNewOIL3MuonCandidates"));
104  desc.add<edm::InputTag>("InputLinks", edm::InputTag("hltNewOIL3MuonsLinksCombination"));
105  desc.add<bool>("applyL3Filters", true);
106  desc.add<int>("MinNhits", 1);
107  desc.add<double>("MaxNormalizedChi2", 20.0);
108  desc.add<int>("MinNmuonHits", 0);
109  desc.add<double>("MaxPtDifference", 999.0); //relative difference
110  descriptions.add("HLTMuonL2SelectorForL3IO", desc);
111 }
const edm::EDGetTokenT< reco::RecoChargedCandidateCollection > l3OISrc_
HLTMuonL2SelectorForL3IO(const edm::ParameterSet &)
constructor with config
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:572
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
const std::string metname
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:614
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
default values
std::vector< MuonTrackLinks > MuonTrackLinksCollection
collection of MuonTrackLinks
Definition: MuonFwd.h:22
const edm::EDGetTokenT< reco::MuonTrackLinksCollection > l3linkToken_
int iEvent
Definition: GenABIO.cc:224
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:617
double pt() const
track transverse momentum
Definition: TrackBase.h:602
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:740
ParameterDescriptionBase * add(U const &iLabel, T const &value)
#define LogTrace(id)
std::vector< RecoChargedCandidate > RecoChargedCandidateCollection
collectin of RecoChargedCandidate objects
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:483
const edm::EDGetTokenT< reco::TrackCollection > l2Src_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
void produce(edm::Event &, const edm::EventSetup &) override
select muons
~HLTMuonL2SelectorForL3IO() override
destructor
int numberOfValidMuonHits() const
Definition: HitPattern.h:793
def move(src, dest)
Definition: eostools.py:511