CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtEvent.cc
Go to the documentation of this file.
3 #include <cstring>
4 
5 // find corresponding hypotheses based on JetLepComb
6 int
7 TtEvent::correspondingHypo(const HypoClassKey& key1, const unsigned& hyp1, const HypoClassKey& key2) const
8 {
9  for(unsigned hyp2 = 0; hyp2 < this->numberOfAvailableHypos(key2); ++hyp2) {
10  if( this->jetLeptonCombination(key1, hyp1) == this->jetLeptonCombination(key2, hyp2) )
11  return hyp2;
12  }
13  return -1; // if no corresponding hypothesis was found
14 }
15 
16 // return the corresponding enum value from a string
18 TtEvent::hypoClassKeyFromString(const std::string& label) const
19 {
20  static HypoClassKeyStringToEnum hypoClassKeyStringToEnumMap[] = {
21  { "kGeom", kGeom },
22  { "kWMassMaxSumPt", kWMassMaxSumPt },
23  { "kMaxSumPtWMass", kMaxSumPtWMass },
24  { "kGenMatch", kGenMatch },
25  { "kMVADisc", kMVADisc },
26  { "kKinFit", kKinFit },
27  { "kKinSolution", kKinSolution },
28  { "kWMassDeltaTopMass", kWMassDeltaTopMass},
29  { "kHitFit", kHitFit },
30  { 0, (HypoClassKey)-1 }
31  };
32 
33  bool found = false;
35  for(int i = 0; hypoClassKeyStringToEnumMap[i].label && (!found); ++i){
36  if(!strcmp(label.c_str(), hypoClassKeyStringToEnumMap[i].label)){
37  found = true;
38  value = hypoClassKeyStringToEnumMap[i].value;
39  }
40  }
41 
42  // in case of unrecognized selection type
43  if(!found){
44  throw cms::Exception("TtEventError") << label << " is not a recognized HypoClassKey";
45  }
46  return value;
47 }
48 
49 // print pt, eta, phi, mass of a given candidate into an existing LogInfo
50 void
51 TtEvent::printParticle(edm::LogInfo &log, const char* name, const reco::Candidate* cand) const
52 {
53  if(!cand) {
54  log << std::setw(15) << name << ": not available!\n";
55  return;
56  }
57  log << std::setprecision(3) << setiosflags(std::ios::fixed | std::ios::showpoint);
58  log << std::setw(15) << name << ": "
59  << std::setw( 7) << cand->pt() << "; "
60  << std::setw( 7) << cand->eta() << "; "
61  << std::setw( 7) << cand->phi() << "; "
62  << resetiosflags(std::ios::fixed | std::ios::showpoint) << setiosflags(std::ios::scientific)
63  << std::setw(10) << cand->mass() << "\n";
64  log << resetiosflags(std::ios::scientific);
65 }
int i
Definition: DBlmapReader.cc:9
HypoClassKey
supported classes of event hypotheses
Definition: TtEvent.h:31
unsigned int numberOfAvailableHypos(const std::string &key) const
return number of available hypotheses within a given hypothesis class
Definition: TtEvent.h:70
virtual double pt() const =0
transverse momentum
virtual double mass() const =0
mass
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:51
HypoClassKey hypoClassKeyFromString(const std::string &label) const
return the corresponding enum value from a string
Definition: TtEvent.cc:18
a lightweight map for selection type string label and enum value
Definition: TtEvent.h:37
int correspondingHypo(const std::string &key1, const unsigned &hyp1, const std::string &key2) const
return the hypothesis in hypothesis class &#39;key2&#39;, which corresponds to hypothesis &#39;hyp1&#39; in hypothesi...
Definition: TtEvent.h:102
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
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity