CMS 3D CMS Logo

GenJet.cc
Go to the documentation of this file.
1 // GenJet.cc
2 // Fedor Ratnikov, UMd
3 
4 #include <sstream>
5 
7 
8 //Own header file
10 
11 using namespace reco;
12 
13 GenJet::GenJet(const LorentzVector& fP4, const Point& fVertex, const Specific& fSpecific)
14  : Jet(fP4, fVertex), m_specific(fSpecific) {}
15 
17  const Point& fVertex,
18  const Specific& fSpecific,
19  const Jet::Constituents& fConstituents)
20  : Jet(fP4, fVertex, fConstituents), m_specific(fSpecific) {}
21 
22 GenJet::GenJet(const LorentzVector& fP4, const Specific& fSpecific, const Jet::Constituents& fConstituents)
23  : Jet(fP4, Point(0, 0, 0), fConstituents), m_specific(fSpecific) {}
24 
26 float GenJet::detectorEta(float fZVertex) const { return Jet::detectorEta(fZVertex, eta()); }
27 
28 const GenParticle* GenJet::genParticle(const Candidate* fConstituent) {
29  const Candidate* base = fConstituent;
30  if (fConstituent->hasMasterClone())
31  base = fConstituent->masterClone().get();
32  const GenParticle* result = dynamic_cast<const GenParticle*>(base);
33  if (!result)
34  throw cms::Exception("Invalid Constituent") << "GenJet constituent is not of the type GenParticle";
35  return result;
36 }
37 
38 const GenParticle* GenJet::getGenConstituent(unsigned fIndex) const {
39  // no direct access, have to iterate for now
40  int index(fIndex);
42  for (; --index >= 0 && daugh != end(); daugh++) {
43  }
44  if (daugh != end()) { // in range
45  const Candidate* constituent = &*daugh; // deref
46  return genParticle(constituent);
47  }
48  return nullptr;
49 }
50 
51 std::vector<const GenParticle*> GenJet::getGenConstituents() const {
52  std::vector<const GenParticle*> result;
53  for (unsigned i = 0; i < numberOfDaughters(); i++)
54  result.push_back(getGenConstituent(i));
55  return result;
56 }
57 
58 GenJet* GenJet::clone() const { return new GenJet(*this); }
59 
60 bool GenJet::overlap(const Candidate&) const { return false; }
61 
63  std::ostringstream out;
64  out << Jet::print() // generic jet info
65  << " GenJet specific:" << std::endl
66  << " em/had/invisible/aux energies: " << emEnergy() << '/' << hadEnergy() << '/' << invisibleEnergy() << '/'
67  << auxiliaryEnergy() << std::endl;
68  out << " MC particles:" << std::endl;
69  std::vector<const GenParticle*> mcparts = getGenConstituents();
70  for (unsigned i = 0; i < mcparts.size(); i++) {
71  const GenParticle* mcpart = mcparts[i];
72  if (mcpart) {
73  out << " #" << i << " PDG code:" << mcpart->pdgId() << ", p/pt/eta/phi: " << mcpart->p() << '/'
74  << mcpart->pt() << '/' << mcpart->eta() << '/' << mcpart->phi() << std::endl;
75  } else {
76  out << " #" << i << " No information about constituent" << std::endl;
77  }
78  }
79  return out.str();
80 }
float detectorEta(float fZVertex) const
Detector Eta (use reference Z and jet kinematics only)
Definition: GenJet.cc:26
static float detectorEta(float fZVertex, float fPhysicsEta)
static function to convert physics eta to detector eta
float hadEnergy() const
Definition: GenJet.h:88
float auxiliaryEnergy() const
Definition: GenJet.h:92
double pt() const final
transverse momentum
virtual const GenParticle * getGenConstituent(unsigned fIndex) const
get specific constituent
Definition: GenJet.cc:38
std::string print() const override
Print object.
Definition: GenJet.cc:62
base
Main Program
Definition: newFWLiteAna.py:92
std::vector< Constituent > Constituents
Definition: Jet.h:23
int pdgId() const final
PDG identifier.
virtual bool hasMasterClone() const =0
double p() const final
magnitude of momentum vector
Definition: Jet.py:1
virtual std::string print() const
Print object.
const_iterator end() const
last daughter const_iterator
Definition: Candidate.h:145
float invisibleEnergy() const
Definition: GenJet.h:90
Jets made from MC generator particles.
Definition: GenJet.h:23
unsigned int index
Definition: LeafCandidate.h:31
size_t numberOfDaughters() const override
number of daughters
unsigned int index
index type
Definition: Candidate.h:50
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
fixed size matrix
Structure Point Contains parameters of Gaussian fits to DMRs.
virtual std::vector< const GenParticle * > getGenConstituents() const
get all constituents
Definition: GenJet.cc:51
GenJet * clone() const override
Polymorphic clone.
Definition: GenJet.cc:58
bool overlap(const Candidate &) const override
Polymorphic overlap.
Definition: GenJet.cc:60
static const GenParticle * genParticle(const reco::Candidate *fConstituent)
convert generic constituent to specific type
Definition: GenJet.cc:28
float emEnergy() const
Definition: GenJet.h:86
double phi() const final
momentum azimuthal angle
value_type const * get() const
Definition: RefToBase.h:209
const_iterator begin() const
first daughter const_iterator
Definition: Candidate.h:143
virtual const CandidateBaseRef & masterClone() const =0
double eta() const final
momentum pseudorapidity