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