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 {
27  LogTrace("Muon|RecoMuon|HLTMuonL2SelectorForL3IO")<<"constructor called";
28  produces<reco::TrackCollection>();
29 }
30 
33 
36  const std::string metname = "Muon|RecoMuon|HLTMuonL2SelectorForL3IO";
37 
38  // IN:
40  iEvent.getByToken(l2Src_, l2muonH);
41 
43  iEvent.getByToken(l3OISrc_,l3muonH);
44 
45  // Read Links collection:
47  iEvent.getByToken(l3linkToken_, links);
48 
49  // OUT:
50  std::unique_ptr<reco::TrackCollection > result(new reco::TrackCollection());
51 
52  for (unsigned int il2=0; il2 != l2muonH->size(); ++il2){
53  reco::TrackRef l2muRef(l2muonH, il2);
54  bool re_do_this_L2=true;
55 
56  for (unsigned int il3=0; il3 != l3muonH->size(); ++il3){
58  reco::TrackRef tk = cand->track();
59 
60  bool useThisLink = false;
61  for (unsigned int l(0); l<links->size() && !useThisLink; ++l){
62  const reco::MuonTrackLinks* link = &links->at(l);
63 
64  // Check if the L3 link matches the L3 candidate
65  const reco::Track& globalTrack = *link->globalTrack();
66  float dR2 = deltaR2(tk->eta(),tk->phi(),globalTrack.eta(),globalTrack.phi());
67  if (dR2 < 0.02*0.02 and std::abs(tk->pt() - globalTrack.pt()) < 0.001*tk->pt()) {
68  useThisLink=true;
69  }
70 
71  if (!useThisLink) 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) re_do_this_L2=false;
76 
77  // Check the quality of the reconstructed L3, if poor quality, we will re-use this L2
78  if (staTrack==l2muRef && applyL3Filters_) {
79  re_do_this_L2 = true;
80  const reco::Track& globalTrack = *link->globalTrack();
81  if (globalTrack.numberOfValidHits()< min_Nhits_) continue; // cut on number of hits
82  if (globalTrack.normalizedChi2() > max_NormalizedChi2_ ) continue; //normalizedChi2 cut
83  if (globalTrack.hitPattern().numberOfValidMuonHits() < min_NmuonHits_ ) continue; //min muon hits cut
84  if (std::abs(globalTrack.pt()-l2muRef->pt()) > max_PtDifference_*globalTrack.pt()) continue; // pt difference
85  re_do_this_L2 = false;
86  }
87  }
88  }
89  if (re_do_this_L2) result->push_back(*l2muRef); // used the L2 if no L3 if matched or if the matched L3 has poor quality cuts.
90  }
91  iEvent.put(std::move(result));
92 }
93 
96  desc.add<edm::InputTag>("l2Src",edm::InputTag("hltL2Muons","UpdatedAtVtx"));
97  desc.add<edm::InputTag>("l3OISrc",edm::InputTag("hltNewOIL3MuonCandidates"));
98  desc.add<edm::InputTag>("InputLinks",edm::InputTag("hltNewOIL3MuonsLinksCombination"));
99  desc.add<bool>("applyL3Filters",true);
100  desc.add<int>("MinNhits",1);
101  desc.add<double>("MaxNormalizedChi2",20.0);
102  desc.add<int>("MinNmuonHits",0);
103  desc.add<double>("MaxPtDifference",999.0); //relative difference
104  descriptions.add("HLTMuonL2SelectorForL3IO",desc);
105 }
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:122
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:556
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
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:640
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:230
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:646
double pt() const
track transverse momentum
Definition: TrackBase.h:616
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:815
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:445
const edm::EDGetTokenT< reco::TrackCollection > l2Src_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
fixed size matrix
HLT enums.
int numberOfValidMuonHits() const
Definition: HitPattern.h:833
virtual void produce(edm::Event &, const edm::EventSetup &)
select muons
def move(src, dest)
Definition: eostools.py:510
virtual ~HLTMuonL2SelectorForL3IO()
destructor