CMS 3D CMS Logo

ZtoMMEventSelector.cc
Go to the documentation of this file.
1 // user includes
16 
17 // ROOT includes
18 #include "TLorentzVector.h"
19 
21 public:
22  explicit ZtoMMEventSelector(const edm::ParameterSet&);
23  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
24  bool filter(edm::Event&, edm::EventSetup const&) override;
25 
26 private:
27  bool verbose_;
32 
33  const double maxEta_;
34  const double minPt_;
35  const double maxNormChi2_;
36  const double maxD0_;
37  const double maxDz_;
38  const int minPixelHits_;
39  const int minStripHits_;
40  const int minChambers_;
41  const int minMatches_;
43  const double maxIso_;
44  const double minPtHighest_;
45  const double minInvMass_;
46  const double maxInvMass_;
47 };
48 
49 using namespace std;
50 using namespace edm;
51 
54  desc.addUntracked<bool>("verbose", false);
55  desc.addUntracked<edm::InputTag>("muonInputTag", edm::InputTag("muons"));
56  desc.addUntracked<edm::InputTag>("offlineBeamSpot", edm::InputTag("offlineBeamSpot"));
57  desc.addUntracked<double>("maxEta", 2.4);
58  desc.addUntracked<double>("minPt", 5);
59  desc.addUntracked<double>("maxNormChi2", 1000);
60  desc.addUntracked<double>("maxD0", 0.02);
61  desc.addUntracked<double>("maxDz", 20.);
62  desc.addUntracked<uint32_t>("minPixelHits", 1);
63  desc.addUntracked<uint32_t>("minStripHits", 8);
64  desc.addUntracked<uint32_t>("minChambers", 2);
65  desc.addUntracked<uint32_t>("minMatches", 2);
66  desc.addUntracked<double>("minMatchedStations", 2);
67  desc.addUntracked<double>("maxIso", 0.3);
68  desc.addUntracked<double>("minPtHighest", 24);
69  desc.addUntracked<double>("minInvMass", 75);
70  desc.addUntracked<double>("maxInvMass", 105);
71  descriptions.addWithDefaultLabel(desc);
72 }
73 
75  : verbose_(ps.getUntrackedParameter<bool>("verbose", false)),
76  muonTag_(ps.getUntrackedParameter<edm::InputTag>("muonInputTag", edm::InputTag("muons"))),
77  bsTag_(ps.getUntrackedParameter<edm::InputTag>("offlineBeamSpot", edm::InputTag("offlineBeamSpot"))),
78  muonToken_(consumes<reco::MuonCollection>(muonTag_)),
79  bsToken_(consumes<reco::BeamSpot>(bsTag_)),
80  maxEta_(ps.getUntrackedParameter<double>("maxEta", 2.4)),
81  minPt_(ps.getUntrackedParameter<double>("minPt", 5)),
82  maxNormChi2_(ps.getUntrackedParameter<double>("maxNormChi2", 1000)),
83  maxD0_(ps.getUntrackedParameter<double>("maxD0", 0.02)),
84  maxDz_(ps.getUntrackedParameter<double>("maxDz", 20.)),
85  minPixelHits_(ps.getUntrackedParameter<uint32_t>("minPixelHits", 1)),
86  minStripHits_(ps.getUntrackedParameter<uint32_t>("minStripHits", 8)),
87  minChambers_(ps.getUntrackedParameter<uint32_t>("minChambers", 2)),
88  minMatches_(ps.getUntrackedParameter<uint32_t>("minMatches", 2)),
89  minMatchedStations_(ps.getUntrackedParameter<double>("minMatchedStations", 2)),
90  maxIso_(ps.getUntrackedParameter<double>("maxIso", 0.3)),
91  minPtHighest_(ps.getUntrackedParameter<double>("minPtHighest", 24)),
92  minInvMass_(ps.getUntrackedParameter<double>("minInvMass", 75)),
93  maxInvMass_(ps.getUntrackedParameter<double>("maxInvMass", 105)) {}
94 
96  // Read Muon Collection
98  iEvent.getByToken(muonToken_, muonColl);
99 
100  // and the beamspot
102  iEvent.getByToken(bsToken_, beamSpot);
103 
104  std::vector<TLorentzVector> list;
105  if (muonColl.isValid()) {
106  for (auto const& mu : *muonColl) {
107  if (!mu.isGlobalMuon())
108  continue;
109  if (!mu.isPFMuon())
110  continue;
111  if (std::fabs(mu.eta()) >= maxEta_)
112  continue;
113  if (mu.pt() < minPt_)
114  continue;
115 
116  reco::TrackRef gtk = mu.globalTrack();
117  double chi2 = gtk->chi2();
118  double ndof = gtk->ndof();
119  double chbyndof = (ndof > 0) ? chi2 / ndof : 0;
120  if (chbyndof >= maxNormChi2_)
121  continue;
122 
123  if (beamSpot.isValid()) {
124  reco::TrackRef tk = mu.innerTrack();
125  double abstrkd0 = std::abs(tk->dxy(beamSpot->position()));
126  if (abstrkd0 >= maxD0_)
127  continue;
128  double abstrkdz = std::abs(tk->dz(beamSpot->position()));
129  if (abstrkdz >= maxDz_)
130  continue;
131  } else {
132  edm::LogError("ZtoMMEventSelector") << "Error >> Failed to get BeamSpot for label: " << bsTag_;
133  }
134 
135  const reco::HitPattern& hitp = gtk->hitPattern();
137  continue;
139  continue;
140 
141  // Hits/section in the muon chamber
142  if (mu.numberOfChambers() < minChambers_)
143  continue;
144  if (mu.numberOfMatches() < minMatches_)
145  continue;
146  if (mu.numberOfMatchedStations() < minMatchedStations_)
147  continue;
149  continue;
150 
151  // PF Isolation
152  const reco::MuonPFIsolation& pfIso04 = mu.pfIsolationR04();
153  double absiso = pfIso04.sumChargedParticlePt +
154  std::max(0.0, pfIso04.sumNeutralHadronEt + pfIso04.sumPhotonEt - 0.5 * pfIso04.sumPUPt);
155  if (absiso / mu.pt() > maxIso_)
156  continue;
157 
158  TLorentzVector lv;
159  lv.SetPtEtaPhiE(mu.pt(), mu.eta(), mu.phi(), mu.energy());
160  list.push_back(lv);
161  }
162  } else {
163  edm::LogError("ZtoMMEventSelector") << "Error >> Failed to get MuonCollection for label: " << muonTag_;
164  return false;
165  }
166 
167  if (list.size() < 2)
168  return false;
169  if (list[0].Pt() < minPtHighest_)
170  return false;
171  TLorentzVector zv = list[0] + list[1];
172  double mass = zv.M();
173  if (mass < minInvMass_ || mass > maxInvMass_)
174  return false;
175 
176  return true;
177 }
178 
179 // Define this as a plug-in
int numberOfValidPixelHits() const
Definition: HitPattern.h:831
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
const edm::InputTag muonTag_
const float maxD0_
Definition: Constants.h:82
Log< level::Error, false > LogError
const double minPtHighest_
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
int numberOfValidStripHits() const
Definition: HitPattern.h:843
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:224
const edm::EDGetTokenT< reco::MuonCollection > muonToken_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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.