CMS 3D CMS Logo

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