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