CMS 3D CMS Logo

TtDecayChannelSelector.cc
Go to the documentation of this file.
6 
7 // static const string for status check in
8 // TtDecayChannelSelector::search functions
9 static const std::string kGenParticles = "genParticles";
10 
11 // maximal number of possible leptonic decay
12 // channels
13 static const unsigned int kDecayChannels = 3;
14 
16  : invert_(cfg.getParameter<bool>("invert")),
17  allowElectron_(false),
18  allowMuon_(false),
19  allow1Prong_(false),
20  allow3Prong_(false) {
21  // tau decays are not restricted if this PSet does not exist at all
22  restrictTauDecays_ = cfg.existsAs<edm::ParameterSet>("restrictTauDecays");
23  // determine allowed tau decays
24  if (restrictTauDecays_) {
25  edm::ParameterSet allowedTauDecays = cfg.getParameter<edm::ParameterSet>("restrictTauDecays");
26  // tau decays are not restricted if none of the following parameters exists
27  restrictTauDecays_ = (allowedTauDecays.existsAs<bool>("electron") || allowedTauDecays.existsAs<bool>("muon") ||
28  allowedTauDecays.existsAs<bool>("oneProng") || allowedTauDecays.existsAs<bool>("threeProng"));
29  // specify the different possible restrictions of the tau decay channels
31  (allowedTauDecays.existsAs<bool>("electron") ? allowedTauDecays.getParameter<bool>("electron") : false);
32  allowMuon_ = (allowedTauDecays.existsAs<bool>("muon") ? allowedTauDecays.getParameter<bool>("muon") : false);
33  allow1Prong_ =
34  (allowedTauDecays.existsAs<bool>("oneProng") ? allowedTauDecays.getParameter<bool>("oneProng") : false);
35  allow3Prong_ =
36  (allowedTauDecays.existsAs<bool>("threeProng") ? allowedTauDecays.getParameter<bool>("threeProng") : false);
37  }
38  // allowed top decays PSet
39  edm::ParameterSet allowedTopDecays = cfg.getParameter<edm::ParameterSet>("allowedTopDecays");
40 
41  // fill decayBranchA_
42  edm::ParameterSet decayBranchA = allowedTopDecays.getParameter<edm::ParameterSet>("decayBranchA");
43  decayBranchA_.push_back(decayBranchA.getParameter<bool>("electron"));
44  decayBranchA_.push_back(decayBranchA.getParameter<bool>("muon"));
45  decayBranchA_.push_back(decayBranchA.getParameter<bool>("tau"));
46 
47  // fill decay branchB_
48  edm::ParameterSet decayBranchB = allowedTopDecays.getParameter<edm::ParameterSet>("decayBranchB");
49  decayBranchB_.push_back(decayBranchB.getParameter<bool>("electron"));
50  decayBranchB_.push_back(decayBranchB.getParameter<bool>("muon"));
51  decayBranchB_.push_back(decayBranchB.getParameter<bool>("tau"));
52 
53  // fill allowedDecays_
54  for (unsigned int d = 0; d < kDecayChannels; ++d) {
56  }
57 }
58 
59 bool TtDecayChannelSelector::operator()(const reco::GenParticleCollection& parts, std::string_view inputType) const {
60  bool verbose =
61  false; // set this to true for debugging and add TtDecayChannelSelector category to the MessageLogger in your cfg file
62  unsigned int iLep = 0;
63  unsigned int iTop = 0, iBeauty = 0, iElec = 0, iMuon = 0, iTau = 0;
64  for (reco::GenParticleCollection::const_iterator top = parts.begin(); top != parts.end(); ++top) {
65  if (search(top, TopDecayID::tID, inputType)) {
66  ++iTop;
67  for (reco::GenParticle::const_iterator td = top->begin(); td != top->end(); ++td) {
68  if (search(td, TopDecayID::bID, inputType)) {
69  ++iBeauty;
70  }
71  if (search(td, TopDecayID::WID, inputType)) {
72  for (reco::GenParticle::const_iterator wd = td->begin(); wd != td->end(); ++wd) {
73  if (std::abs(wd->pdgId()) == TopDecayID::elecID) {
74  ++iElec;
75  }
76  if (std::abs(wd->pdgId()) == TopDecayID::muonID) {
77  ++iMuon;
78  }
79  if (std::abs(wd->pdgId()) == TopDecayID::tauID) {
80  if (restrictTauDecays_) {
81  // count as iTau if it is leptonic, one-prong
82  // or three-prong and ignore increasing iLep
83  // though else
84  if (tauDecay(*wd)) {
85  ++iTau;
86  } else {
87  ++iLep;
88  }
89  } else {
90  ++iTau;
91  }
92  }
93  }
94  }
95  }
96  }
97  }
98  if (verbose) {
99  edm::LogVerbatim log("TtDecayChannelSelector");
100  log << "----------------------"
101  << "\n"
102  << " iTop : " << iTop << "\n"
103  << " iBeauty : " << iBeauty << "\n"
104  << " iElec : " << iElec << "\n"
105  << " iMuon : " << iMuon << "\n"
106  << " iTau : " << iTau + iLep;
107  if (restrictTauDecays_ && (iTau + iLep) > 0) {
108  log << " (" << iTau << ")\n";
109  } else {
110  log << "\n";
111  }
112  log << "- - - - - - - - - - - "
113  << "\n";
114  }
115  iLep += iElec + iMuon + iTau;
116 
117  bool accept = false;
118  unsigned int channel = decayChannel();
119  if ((iTop == 2) && (iBeauty == 2)) {
120  if (channel == iLep) {
121  if (channel == 0) {
122  // no lepton: accept without restriction we already
123  // know that the number of leptons is correct
124  accept = true;
125  }
126  if (channel == 1) {
127  // one lepton: check that this one is allowed
128  accept = (iElec && allowedDecays_[Elec]) || (iMuon && allowedDecays_[Muon]) || (iTau && allowedDecays_[Tau]);
129  }
130  if (channel == 2) {
131  if (checkSum(allowedDecays_) == channel) {
132  // no redundancy
133  accept = (allowedDecays_[Elec] == (int)iElec) && (allowedDecays_[Muon] == (int)iMuon) &&
134  (allowedDecays_[Tau] == (int)iTau);
135  } else {
136  // reject events with wrong tau decays
137  if (iElec + iMuon + iTau != channel) {
138  accept = false;
139  } else {
140  if ((iElec == 2) || (iMuon == 2) || (iTau == 2)) {
141  // same lepton twice: check that this is allowed.
142  accept = (allowedDecays_[Elec] == (int)iElec) || (allowedDecays_[Muon] == (int)iMuon) ||
143  (allowedDecays_[Tau] == (int)iTau);
144  } else {
145  // two different leptons: look if there is a possible combination
146  accept = (((iElec && decayBranchA_[Elec]) &&
147  ((iMuon && decayBranchB_[Muon]) || (iTau && decayBranchB_[Tau]))) ||
148  ((iMuon && decayBranchA_[Muon]) &&
149  ((iElec && decayBranchB_[Elec]) || (iTau && decayBranchB_[Tau]))) ||
150  ((iTau && decayBranchA_[Tau]) &&
151  ((iElec && decayBranchB_[Elec]) || (iMuon && decayBranchB_[Muon]))));
152  }
153  }
154  }
155  }
156  }
157  accept = ((!invert_ && accept) || (!(!invert_) && !accept));
158  } else {
159  edm::LogWarning("NoVtbDecay") << "Decay is not via Vtb";
160  }
161  if (verbose)
162  edm::LogVerbatim("TtDecayChannelSelector") << " accept : " << accept;
163  return accept;
164 }
165 
166 bool TtDecayChannelSelector::search(reco::GenParticleCollection::const_iterator& part,
167  int pdgId,
168  std::string_view inputType) const {
169  if (inputType == kGenParticles) {
170  return (std::abs(part->pdgId()) == pdgId && part->status() == TopDecayID::unfrag) ? true : false;
171  } else {
172  return (std::abs(part->pdgId()) == pdgId) ? true : false;
173  }
174 }
175 
177  int pdgId,
178  std::string_view inputType) const {
179  if (inputType == kGenParticles) {
180  return (std::abs(part->pdgId()) == pdgId && part->status() == TopDecayID::unfrag) ? true : false;
181  } else {
182  return (std::abs(part->pdgId()) == pdgId) ? true : false;
183  }
184 }
185 
187  // if stable, return 1 or 0
188  if (part.status() == 1) {
189  return (part.charge() != 0);
190  }
191  // if unstable, call recursively on daughters
192  int prong = 0;
193  for (reco::Candidate::const_iterator daughter = part.begin(); daughter != part.end(); ++daughter) {
194  prong += countProngs(*daughter);
195  }
196  return prong;
197 }
198 
200  bool electronTau = false;
201  bool muonTau = false;
202  unsigned int nch = 0;
203  // loop on tau decays, check for an elec
204  // or muon and count charged particles
205  for (reco::Candidate::const_iterator daughter = tau.begin(); daughter != tau.end(); ++daughter) {
206  // if the tau daughter is again a tau, this means that the particle has
207  // still to be propagated; in that case, return the result of the same
208  // method applied on the daughter of the current particle
209  if (daughter->pdgId() == tau.pdgId()) {
210  return tauDecay(*daughter);
211  }
212  // check for electron from tau decay
213  electronTau |= (std::abs(daughter->pdgId()) == TopDecayID::elecID);
214  // check for muon from tau decay
215  muonTau |= (std::abs(daughter->pdgId()) == TopDecayID::muonID);
216  // count charged particles
217  nch += countProngs(*daughter);
218  }
219  return ((allowElectron_ && electronTau) || (allowMuon_ && muonTau) ||
220  (allow1Prong_ && !electronTau && !muonTau && nch == 1) ||
221  (allow3Prong_ && !electronTau && !muonTau && nch == 3));
222 }
Log< level::Info, true > LogVerbatim
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
static const int bID
Definition: TopGenEvent.h:13
bool verbose
static const int unfrag
Definition: TopGenEvent.h:11
Decay decayBranchA_
top decay branch 1
unsigned int checkSum(const Decay &vec) const
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:172
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31
static const int tID
Definition: TopGenEvent.h:12
static const int tauID
Definition: TopGenEvent.h:20
static const unsigned int kDecayChannels
bool tauDecay(const reco::Candidate &) const
check tau decay to be leptonic, 1-prong or 3-prong
unsigned int countProngs(const reco::Candidate &part) const
count the number of charged particles for tau decays
Definition: Muon.py:1
bool restrictTauDecays_
restrict tau decays
const_iterator end() const
last daughter const_iterator
Definition: Candidate.h:145
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
TtDecayChannelSelector(const edm::ParameterSet &)
std contructor
d
Definition: ztail.py:151
unsigned int decayChannel() const
return decay channel to select for from configuration
Decay decayBranchB_
top decay branch 2
bool allowElectron_
allow tau decays into electron
bool allow3Prong_
allow 2-prong tau decays
Definition: Tau.py:1
bool allowMuon_
allow tau decays into muon
static const int muonID
Definition: TopGenEvent.h:19
part
Definition: HCALResponse.h:20
allowedTopDecays
allow given lepton in corresponding decay branch for a given decay-channel selection; all leptons to ...
bool operator()(const reco::GenParticleCollection &parts, std::string_view inputType) const
operator for decay channel selection
static const int WID
Definition: TopGenEvent.h:17
static const std::string kGenParticles
static const int elecID
Definition: TopGenEvent.h:18
bool allow1Prong_
allow 1-prong tau decays
Log< level::Warning, false > LogWarning
bool search(reco::GenParticleCollection::const_iterator &part, int pdgId, std::string_view inputType) const
search for particle with pdgId in given listing (for top)
bool invert_
invert selection
const_iterator begin() const
first daughter const_iterator
Definition: Candidate.h:143