CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtSemiLeptonicEvent.cc
Go to the documentation of this file.
4 
5 // print info via MessageLogger
6 void
8 {
9  if(verbosity%10 <= 0)
10  return;
11 
12  edm::LogInfo log("TtSemiLeptonicEvent");
13 
14  log << "++++++++++++++++++++++++++++++++++++++++++++++++++ \n";
15 
16  // get some information from the genEvent (if available)
17  if( !genEvt_ ) log << " TtGenEvent not available! \n";
18  else {
19  log << " TtGenEvent says: ";
20  if( !this->genEvent()->isTtBar() ) log << "Not TtBar";
21  else if( this->genEvent()->isFullHadronic() ) log << "Fully Hadronic TtBar";
22  else if( this->genEvent()->isFullLeptonic() ) log << "Fully Leptonic TtBar";
23  else if( this->genEvent()->isSemiLeptonic() ) {
24  log << "Semi-leptonic TtBar, ";
25  switch( this->genEvent()->semiLeptonicChannel() ) {
26  case WDecay::kElec : log << "Electron"; break;
27  case WDecay::kMuon : log << "Muon" ; break;
28  case WDecay::kTau : log << "Tau" ; break;
29  default : log << "Unknown" ; break;
30  }
31  log << " Channel";
32  }
33  log << "\n";
34  }
35 
36  // get number of available hypothesis classes
37  log << " Number of available event hypothesis classes: " << this->numberOfAvailableHypoClasses() << " \n";
38 
39  // create a legend for the jetLepComb
40  log << " - JetLepComb: ";
41  for(unsigned idx = 0; idx < 5; idx++) {
42  switch(idx) {
43  case TtSemiLepEvtPartons::LightQ : log << "LightP "; break;
44  case TtSemiLepEvtPartons::LightQBar : log << "LightQ "; break;
45  case TtSemiLepEvtPartons::HadB : log << " HadB "; break;
46  case TtSemiLepEvtPartons::LepB : log << " LepB "; break;
47  case TtSemiLepEvtPartons::Lepton : log << "Lepton "; break;
48  }
49  }
50  log << "\n";
51 
52  // get details from the hypotheses
53  typedef std::map<HypoClassKey, std::vector<HypoCombPair> >::const_iterator EventHypo;
54  for(EventHypo hyp = evtHyp_.begin(); hyp != evtHyp_.end(); ++hyp) {
55  HypoClassKey hypKey = (*hyp).first;
56  // header for each hypothesis
57  log << "-------------------------------------------------- \n";
58  switch(hypKey) {
59  case kGeom : log << " Geom" ; break;
60  case kWMassDeltaTopMass : log << " WMassDeltaTopMass"; break;
61  case kWMassMaxSumPt : log << " WMassMaxSumPt" ; break;
62  case kMaxSumPtWMass : log << " MaxSumPtWMass" ; break;
63  case kGenMatch : log << " GenMatch" ; break;
64  case kMVADisc : log << " MVADisc" ; break;
65  case kKinFit : log << " KinFit" ; break;
66  default : log << " Unknown";
67  }
68  log << "-Hypothesis: \n";
69  log << " * Number of real neutrino solutions: " << this->numberOfRealNeutrinoSolutions(hypKey) << "\n";
70  unsigned nOfHyp = this->numberOfAvailableHypos(hypKey);
71  if(nOfHyp > 1) {
72  log << " * Number of available jet combinations: " << nOfHyp << "\n";
73  if(verbosity < 10)
74  log << " The following was found to be the best one:\n";
75  }
76  // if verbosity level is smaller than 10, never show more than the best jet combination
77  if(verbosity < 10)
78  nOfHyp = 1;
79  for(unsigned cmb=0; cmb<nOfHyp; cmb++) {
80  // check if hypothesis is valid
81  if( !this->isHypoValid(hypKey, cmb) )
82  log << " * Not valid! \n";
83  // get meta information for valid hypothesis
84  else {
85  // jetLepComb
86  log << " * JetLepComb:";
87  std::vector<int> jets = this->jetLeptonCombination(hypKey, cmb);
88  for(unsigned int iJet = 0; iJet < jets.size(); iJet++) {
89  log << " " << jets[iJet] << " ";
90  }
91  log << "\n";
92  // specialties for some hypotheses
93  switch(hypKey) {
94  case kGenMatch : log << " * Sum(DeltaR) : " << this->genMatchSumDR(cmb) << " \n"
95  << " * Sum(DeltaPt): " << this->genMatchSumPt(cmb) << " \n"; break;
96  case kMVADisc : log << " * Method : " << this->mvaMethod() << " \n"
97  << " * Discrim.: " << this->mvaDisc(cmb) << " \n"; break;
98  case kKinFit : log << " * Chi^2 : " << this->fitChi2(cmb) << " \n"
99  << " * Prob(Chi^2): " << this->fitProb(cmb) << " \n"; break;
100  default : break;
101  }
102  // kinematic quantities of particles (if last digit of verbosity level > 1)
103  if(verbosity%10 >= 2) {
104  log << " * Candidates (pt; eta; phi; mass):\n";
105  printParticle(log, "hadronic top", this->hadronicDecayTop(hypKey, cmb));
106  printParticle(log, "hadronic W ", this->hadronicDecayW (hypKey, cmb));
107  if(verbosity%10 >= 3) {
108  printParticle(log, "hadronic b ", this->hadronicDecayB (hypKey, cmb));
109  printParticle(log, "hadronic p ", this->hadronicDecayQuark (hypKey, cmb));
110  printParticle(log, "hadronic q ", this->hadronicDecayQuarkBar(hypKey, cmb));
111  }
112  printParticle(log, "leptonic top", this->leptonicDecayTop(hypKey, cmb));
113  printParticle(log, "leptonic W ", this->leptonicDecayW (hypKey, cmb));
114  if(verbosity%10 >= 3) {
115  printParticle(log, "leptonic b ", this->leptonicDecayB(hypKey, cmb));
116  printParticle(log, "lepton ", this->singleLepton (hypKey, cmb));
117  printParticle(log, "neutrino ", this->singleNeutrino(hypKey, cmb));
118  }
119  }
120  }
121  }
122  }
123 
124  log << "++++++++++++++++++++++++++++++++++++++++++++++++++";
125 }
const reco::GenParticle * leptonicDecayTop() const
get leptonic top of the TtGenEvent
HypoClassKey
supported classes of event hypotheses
Definition: TtEvent.h:31
const edm::RefProd< TtGenEvent > & genEvent() const
get TtGenEvent
Definition: TtEvent.h:53
unsigned int numberOfAvailableHypos(const std::string &key) const
return number of available hypotheses within a given hypothesis class
Definition: TtEvent.h:70
const reco::GenParticle * leptonicDecayB() const
get leptonic b of the TtGenEvent
unsigned int numberOfAvailableHypoClasses() const
return number of available hypothesis classes
Definition: TtEvent.h:68
const reco::GenParticle * hadronicDecayB() const
get hadronic b of the TtGenEvent
void printParticle(edm::LogInfo &log, const char *name, const reco::Candidate *cand) const
print pt, eta, phi, mass of a given candidate into an existing LogInfo
Definition: TtEvent.cc:50
double fitChi2(const unsigned &cmb=0) const
return the chi2 of the kinemtaic fit of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:86
const reco::GenParticle * singleLepton() const
get lepton top of the TtGenEvent
double genMatchSumDR(const unsigned &cmb=0) const
return the sum dr of the generator match if available; -1 else
Definition: TtEvent.h:80
const reco::GenParticle * leptonicDecayW() const
get leptonic W of the TtGenEvent
const reco::GenParticle * hadronicDecayW() const
get hadronic W of the TtGenEvent
void print(const int verbosity=1) const
double mvaDisc(const unsigned &cmb=0) const
return the mva discriminant value of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:84
const reco::GenParticle * hadronicDecayQuarkBar() const
get hadronic light quark of the TtGenEvent
const int verbosity
Log< T >::type log(const T &t)
Definition: Log.h:22
edm::RefProd< TtGenEvent > genEvt_
reference to TtGenEvent (has to be kept in the event!)
Definition: TtEvent.h:121
double fitProb(const unsigned &cmb=0) const
return the fit probability of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:88
bool isHypoValid(const std::string &key, const unsigned &cmb=0) const
check if hypothesis &#39;cmb&#39; within the hypothesis class was valid; if not it lead to an empty Composite...
Definition: TtEvent.h:64
const reco::GenParticle * hadronicDecayTop() const
get hadronic top of the TtGenEvent
double genMatchSumPt(const unsigned &cmb=0) const
return the sum pt of the generator match if available; -1 else
Definition: TtEvent.h:78
std::string mvaMethod() const
return the label of the mva method in use for the jet parton association (if kMVADisc is not availabl...
Definition: TtEvent.h:82
const reco::GenParticle * hadronicDecayQuark() const
get hadronic light quark of the TtGenEvent
std::map< HypoClassKey, std::vector< HypoCombPair > > evtHyp_
Definition: TtEvent.h:124
const int numberOfRealNeutrinoSolutions(const HypoClassKey &key) const
get number of real neutrino solutions for a given hypo class
std::vector< int > jetLeptonCombination(const std::string &key, const unsigned &cmb=0) const
return the vector of jet lepton combinatorics for a given hypothesis and class
Definition: TtEvent.h:74
const reco::GenParticle * singleNeutrino() const
get neutrino of the TtGenEvent