CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtFullLeptonicEvent.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("TtFullLeptonicEvent");
13 
14  log << "+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ \n";
15 
16  // get some information from the genEvent
17  log << " TtGenEvent says: ";
18  if( !this->genEvent()->isTtBar() ) log << "Not TtBar";
19  else if( this->genEvent()->isFullHadronic() ) log << "Fully Hadronic TtBar";
20  else if( this->genEvent()->isSemiLeptonic() ) log << "Semi-leptonic TtBar";
21  else if( this->genEvent()->isFullLeptonic() ) {
22  log << "Fully Leptonic TtBar, ";
23  switch( this->genEvent()->fullLeptonicChannel().first ) {
24  case WDecay::kElec : log << "Electron-"; break;
25  case WDecay::kMuon : log << "Muon-" ; break;
26  case WDecay::kTau : log << "Tau-" ; break;
27  default : log << "Unknown-" ; break;
28  }
29  switch( this->genEvent()->fullLeptonicChannel().second ) {
30  case WDecay::kElec : log << "Electron-"; break;
31  case WDecay::kMuon : log << "Muon-" ; break;
32  case WDecay::kTau : log << "Tau-" ; break;
33  default : log << "Unknown-" ; break;
34  }
35  log << "Channel";
36  }
37  log << "\n";
38 
39  // get number of available hypothesis classes
40  log << " Number of available event hypothesis classes: " << this->numberOfAvailableHypoClasses() << " \n";
41 
42  // create a legend for the jetLepComb
43  log << " - JetLepComb: ";
44  log << " b ";
45  log << " bbar ";
46  log << " e1(+) ";
47  log << " e2(-) ";
48  log << " mu1(+)";
49  log << " mu2(-)";
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 not (yet) applicable to TtFullLeptonicEvent --> skipping" ; continue;
60  case kWMassMaxSumPt : log << " WMassMaxSumPt not (yet) applicable to TtFullLeptonicEvent --> skipping" ; continue;
61  case kMaxSumPtWMass : log << " MaxSumPtWMass not (yet) applicable to TtFullLeptonicEvent --> skipping" ; continue;
62  case kGenMatch : log << " GenMatch" ; break;
63  case kMVADisc : log << " MVADisc not (yet) applicable to TtFullLeptonicEvent --> skipping" ; continue;
64  case kKinFit : log << " KinFit not (yet) applicable to TtFullLeptonicEvent --> skipping" ; continue;
65  case kKinSolution : log << " KinSolution" ; break;
66  case kWMassDeltaTopMass : log << " WMassDeltaTopMass not (yet) applicable to TtFullLeptonicEvent --> skipping"; continue;
67  case kHitFit : log << " HitFit not (yet) applicable to TtFullLeptonicEvent --> skipping" ; continue;
68  default : log << " Unknown TtEvent::HypoClassKey provided --> skipping" ; continue;
69  }
70  log << "-Hypothesis: \n";
71  unsigned nOfHyp = this->numberOfAvailableHypos(hypKey);
72  if(nOfHyp > 1) {
73  log << " * Number of available jet combinations: " << nOfHyp << "\n";
74  if(verbosity < 10)
75  log << " The following was found to be the best one:\n";
76  }
77  // if verbosity level is smaller than 10, never show more than the best jet combination
78  if(verbosity < 10)
79  nOfHyp = 1;
80  for(unsigned cmb=0; cmb<nOfHyp; cmb++) {
81  // check if hypothesis is valid
82  if( !this->isHypoValid(hypKey, cmb) )
83  log << " * Not valid! \n";
84  // get meta information for valid hypothesis
85  else {
86  // jetLepComb
87  log << " * JetLepComb:";
88  std::vector<int> jets = this->jetLeptonCombination(hypKey, cmb);
89  for(unsigned int iJet = 0; iJet < jets.size(); iJet++) {
90  log << " " << jets[iJet] << " ";
91  }
92  log << "\n";
93  // specialties for some hypotheses
94  switch(hypKey) {
95  case kGenMatch : log << " * Sum(DeltaR) : " << this->genMatchSumDR(cmb) << " \n"
96  << " * Sum(DeltaPt): " << this->genMatchSumPt(cmb) << " \n"; break;
97  case kKinSolution : log << " * Weight : " << this->solWeight(cmb) << " \n"
98  << " * isWrongCharge: " << this->isWrongCharge() << " \n"; break;
99  default : break;
100  }
101  // kinematic quantities of particles (if last digit of verbosity level > 1)
102  if(verbosity%10 >= 2) {
103  log << " * Candidates (pt; eta; phi; mass):\n";
104  if(verbosity%10 >= 3)
105  printParticle(log, "top pair", this->topPair(hypKey, cmb));
106  printParticle(log, "top ", this->top (hypKey, cmb));
107  printParticle(log, "W plus ", this->wPlus(hypKey, cmb));
108  if(verbosity%10 >= 3) {
109  printParticle(log, "b ", this->b (hypKey, cmb));
110  printParticle(log, "leptonBar ", this->leptonBar(hypKey, cmb));
111  printParticle(log, "neutrino ", this->neutrino (hypKey, cmb));
112  }
113  printParticle(log, "topBar ", this->topBar(hypKey, cmb));
114  printParticle(log, "W minus ", this->wMinus(hypKey, cmb));
115  if(verbosity%10 >= 3) {
116  printParticle(log, "bBar ", this->bBar (hypKey, cmb));
117  printParticle(log, "lepton ", this->lepton (hypKey, cmb));
118  printParticle(log, "neutrinoBar ", this->neutrinoBar(hypKey, cmb));
119  }
120  }
121  }
122  }
123  }
124 
125  log << "+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++";
126 }
HypoClassKey
supported classes of event hypotheses
Definition: TtEvent.h:31
const reco::Candidate * bBar(const std::string &key, const unsigned &cmb=0) const
get anti-b of the given hypothesis
const edm::RefProd< TtGenEvent > & genEvent() const
get TtGenEvent
Definition: TtEvent.h:53
const math::XYZTLorentzVector * topPair() const
get combined 4-vector of top and topBar from the TtGenEvent
Definition: TtEvent.h:111
unsigned int numberOfAvailableHypos(const std::string &key) const
return number of available hypotheses within a given hypothesis class
Definition: TtEvent.h:70
const reco::Candidate * top(const std::string &key, const unsigned &cmb=0) const
get top of the given hypothesis
unsigned int numberOfAvailableHypoClasses() const
return number of available hypothesis classes
Definition: TtEvent.h:68
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:26
const reco::Candidate * lepton(const std::string &key, const unsigned &cmb=0) const
get lepton of the given hypothesis
U second(std::pair< T, U > const &p)
double genMatchSumDR(const unsigned &cmb=0) const
return the sum dr of the generator match if available; -1 else
Definition: TtEvent.h:84
void print(const int verbosity=1) const
const reco::Candidate * neutrino(const std::string &key, const unsigned &cmb=0) const
get neutrino of the given hypothesis
vector< PseudoJet > jets
bool isWrongCharge() const
return if the kinematic solution of hypothesis &#39;cmb&#39; is right or wrong charge if available; -1 else ...
bool first
Definition: L1TdeRCT.cc:75
const reco::Candidate * neutrinoBar(const std::string &key, const unsigned &cmb=0) const
get anti-neutrino of the given hypothesis
const reco::Candidate * b(const std::string &key, const unsigned &cmb=0) const
get b of the given hypothesis
const reco::Candidate * topBar(const std::string &key, const unsigned &cmb=0) const
get anti-top of the given hypothesis
double solWeight(const unsigned &cmb=0) const
return the weight of the kinematic solution of hypothesis &#39;cmb&#39; if available; -1 else ...
const reco::Candidate * wPlus(const std::string &key, const unsigned &cmb=0) const
get Wplus of the given hypothesis
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
double genMatchSumPt(const unsigned &cmb=0) const
return the sum pt of the generator match if available; -1 else
Definition: TtEvent.h:82
const reco::Candidate * wMinus(const std::string &key, const unsigned &cmb=0) const
get Wminus of the given hypothesis
std::map< HypoClassKey, std::vector< HypoCombPair > > evtHyp_
Definition: TtEvent.h:153
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:78
const reco::Candidate * leptonBar(const std::string &key, const unsigned &cmb=0) const
get anti-lepton of the given hypothesis