CMS 3D CMS Logo

ZtoMMEventSelector.cc
Go to the documentation of this file.
7 #include "TLorentzVector.h"
8 
11 
12 using namespace std;
13 using namespace edm;
14 
16  : verbose_(ps.getUntrackedParameter<bool>("verbose", false)),
17  muonTag_(ps.getUntrackedParameter<edm::InputTag>("muonInputTag", edm::InputTag("muons"))),
18  bsTag_(ps.getUntrackedParameter<edm::InputTag>("offlineBeamSpot", edm::InputTag("offlineBeamSpot"))),
19  muonToken_(consumes<reco::MuonCollection>(muonTag_)),
20  bsToken_(consumes<reco::BeamSpot>(bsTag_)),
21  maxEta_(ps.getUntrackedParameter<double>("maxEta", 2.1)),
22  minPt_(ps.getUntrackedParameter<double>("minPt", 5)),
23  maxNormChi2_(ps.getUntrackedParameter<double>("maxNormChi2", 10)),
24  maxD0_(ps.getUntrackedParameter<double>("maxD0", 0.02)),
25  maxDz_(ps.getUntrackedParameter<double>("maxDz", 20.)),
26  minPixelHits_(ps.getUntrackedParameter<uint32_t>("minPixelHits", 1)),
27  minStripHits_(ps.getUntrackedParameter<uint32_t>("minStripHits", 8)),
28  minChambers_(ps.getUntrackedParameter<uint32_t>("minChambers", 2)),
29  minMatches_(ps.getUntrackedParameter<uint32_t>("minMatches", 2)),
30  minMatchedStations_(ps.getUntrackedParameter<double>("minMatchedStations", 2)),
31  maxIso_(ps.getUntrackedParameter<double>("maxIso", 0.3)),
32  minPtHighest_(ps.getUntrackedParameter<double>("minPtHighest", 24)),
33  minInvMass_(ps.getUntrackedParameter<double>("minInvMass", 60)),
34  maxInvMass_(ps.getUntrackedParameter<double>("maxInvMass", 120)) {}
36  // Read Muon Collection
38  iEvent.getByToken(muonToken_, muonColl);
39 
40  // and the beamspot
42  iEvent.getByToken(bsToken_, beamSpot);
43 
44  std::vector<TLorentzVector> list;
45  if (muonColl.isValid()) {
46  for (auto const& mu : *muonColl) {
47  if (!mu.isGlobalMuon())
48  continue;
49  if (!mu.isPFMuon())
50  continue;
51  if (std::fabs(mu.eta()) >= maxEta_)
52  continue;
53  if (mu.pt() < minPt_)
54  continue;
55 
56  reco::TrackRef gtk = mu.globalTrack();
57  double chi2 = gtk->chi2();
58  double ndof = gtk->ndof();
59  double chbyndof = (ndof > 0) ? chi2 / ndof : 0;
60  if (chbyndof >= maxNormChi2_)
61  continue;
62 
63  if (beamSpot.isValid()) {
64  reco::TrackRef tk = mu.innerTrack();
65  double abstrkd0 = std::abs(tk->dxy(beamSpot->position()));
66  if (abstrkd0 >= maxD0_)
67  continue;
68  double abstrkdz = std::abs(tk->dz(beamSpot->position()));
69  if (abstrkdz >= maxDz_)
70  continue;
71  } else {
72  edm::LogError("ZtoMMEventSelector") << "Error >> Failed to get BeamSpot for label: " << bsTag_;
73  }
74 
75  const reco::HitPattern& hitp = gtk->hitPattern();
77  continue;
79  continue;
80 
81  // Hits/section in the muon chamber
82  if (mu.numberOfChambers() < minChambers_)
83  continue;
84  if (mu.numberOfMatches() < minMatches_)
85  continue;
86  if (mu.numberOfMatchedStations() < minMatchedStations_)
87  continue;
89  continue;
90 
91  // PF Isolation
92  const reco::MuonPFIsolation& pfIso04 = mu.pfIsolationR04();
93  double absiso = pfIso04.sumChargedParticlePt +
94  std::max(0.0, pfIso04.sumNeutralHadronEt + pfIso04.sumPhotonEt - 0.5 * pfIso04.sumPUPt);
95  if (absiso / mu.pt() > maxIso_)
96  continue;
97 
98  TLorentzVector lv;
99  lv.SetPtEtaPhiE(mu.pt(), mu.eta(), mu.phi(), mu.energy());
100  list.push_back(lv);
101  }
102  } else {
103  edm::LogError("ZtoMMEventSelector") << "Error >> Failed to get MuonCollection for label: " << muonTag_;
104  return false;
105  }
106 
107  if (list.size() < 2)
108  return false;
109  if (list[0].Pt() < minPtHighest_)
110  return false;
111  TLorentzVector zv = list[0] + list[1];
112  double mass = zv.M();
113  if (mass < minInvMass_ || mass > maxInvMass_)
114  return false;
115 
116  return true;
117 }
118 
119 // Define this as a plug-in
int numberOfValidPixelHits() const
Definition: HitPattern.h:831
const edm::InputTag muonTag_
const double maxNormChi2_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
float *__restrict__ zv
const float maxD0_
Definition: Constants.h:72
Log< level::Error, false > LogError
float sumPhotonEt
sum pt of PF photons
const double minPtHighest_
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
int numberOfValidStripHits() const
Definition: HitPattern.h:843
float sumNeutralHadronEt
sum pt of neutral hadrons
float sumChargedParticlePt
sum-pt of charged Particles(inludes e/mu)
float sumPUPt
sum pt of charged Particles not from PV (for Pu corrections)
int iEvent
Definition: GenABIO.cc:224
const edm::EDGetTokenT< reco::MuonCollection > muonToken_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isGoodMuon(const reco::Muon &muon, SelectionType type, reco::Muon::ArbitrationType arbitrationType=reco::Muon::SegmentAndTrackArbitration)
main GoodMuon wrapper call
ZtoMMEventSelector(const edm::ParameterSet &)
const edm::InputTag bsTag_
bool isValid() const
Definition: HandleBase.h:70
const edm::EDGetTokenT< reco::BeamSpot > bsToken_
bool filter(edm::Event &, edm::EventSetup const &) override
fixed size matrix
HLT enums.
const double maxInvMass_