CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GenJet.cc
Go to the documentation of this file.
1 // GenJet.cc
2 // Fedor Ratnikov, UMd
3 // $Id: GenJet.cc,v 1.12 2008/05/26 11:22:12 arizzi Exp $
4 
5 #include <sstream>
6 
8 
9 //Own header file
11 
12 using namespace reco;
13 
14 GenJet::GenJet (const LorentzVector& fP4, const Point& fVertex,
15  const Specific& fSpecific)
16  : Jet (fP4, fVertex),
17  m_specific (fSpecific)
18 {}
19 
20 GenJet::GenJet (const LorentzVector& fP4, const Point& fVertex,
21  const Specific& fSpecific,
22  const Jet::Constituents& fConstituents)
23  : Jet (fP4, fVertex, fConstituents),
24  m_specific (fSpecific)
25 {}
26 
28  const Specific& fSpecific,
29  const Jet::Constituents& fConstituents)
30  : Jet (fP4, Point(0,0,0), fConstituents),
31  m_specific (fSpecific)
32 {}
33 
35 float GenJet::detectorEta (float fZVertex) const {
36  return Jet::detectorEta (fZVertex, eta());
37 }
38 
39 const GenParticle* GenJet::genParticle (const Candidate* fConstituent) {
40  const Candidate* base = fConstituent;
41  if (fConstituent->hasMasterClone ()) base = fConstituent->masterClone().get ();
42  const GenParticle* result = dynamic_cast<const GenParticle*> (base);
43  if (!result) throw cms::Exception("Invalid Constituent") << "GenJet constituent is not of the type GenParticle";
44  return result;
45 }
46 
47 const GenParticle* GenJet::getGenConstituent (unsigned fIndex) const {
48  // no direct access, have to iterate for now
49  int index (fIndex);
51  for (; --index >= 0 && daugh != end (); daugh++) {}
52  if (daugh != end ()) { // in range
53  const Candidate* constituent = &*daugh; // deref
54  return genParticle (constituent);
55  }
56  return 0;
57 }
58 
59 std::vector <const GenParticle*> GenJet::getGenConstituents () const {
60  std::vector <const GenParticle*> result;
61  for (unsigned i = 0; i < numberOfDaughters (); i++) result.push_back (getGenConstituent (i));
62  return result;
63 }
64 
65 GenJet* GenJet::clone () const {
66  return new GenJet (*this);
67 }
68 
69 bool GenJet::overlap( const Candidate & ) const {
70  return false;
71 }
72 
74  std::ostringstream out;
75  out << Jet::print () // generic jet info
76  << " GenJet specific:" << std::endl
77  << " em/had/invisible/aux energies: "
78  << emEnergy() << '/' << hadEnergy() << '/' << invisibleEnergy() << '/' << auxiliaryEnergy() << std::endl;
79  out << " MC particles:" << std::endl;
80  std::vector <const GenParticle*> mcparts = getGenConstituents ();
81  for (unsigned i = 0; i < mcparts.size (); i++) {
82  const GenParticle* mcpart = mcparts[i];
83  if (mcpart) {
84  out << " #" << i << " PDG code:" << mcpart->pdgId()
85  << ", p/pt/eta/phi: " << mcpart->p() << '/' << mcpart->pt() << '/' << mcpart->eta() << '/' << mcpart->phi() << std::endl;
86  }
87  else {
88  out << " #" << i << " No information about constituent" << std::endl;
89  }
90  }
91  return out.str ();
92 }
tuple base
Main Program
Definition: newFWLiteAna.py:92
int i
Definition: DBlmapReader.cc:9
float auxiliaryEnergy() const
Definition: GenJet.h:64
float hadEnergy() const
Definition: GenJet.h:60
static float detectorEta(float fZVertex, float fPhysicsEta)
static function to convert physics eta to detector eta
Definition: Jet.cc:326
virtual double p() const GCC11_FINAL
magnitude of momentum vector
float emEnergy() const
Definition: GenJet.h:58
virtual int pdgId() const GCC11_FINAL
PDG identifier.
Base class for all types of Jets.
Definition: Jet.h:21
std::vector< Constituent > Constituents
Definition: Jet.h:24
float detectorEta(float fZVertex) const
Detector Eta (use reference Z and jet kinematics only)
Definition: GenJet.cc:35
virtual float phi() const GCC11_FINAL
momentum azimuthal angle
virtual bool hasMasterClone() const =0
virtual size_t numberOfDaughters() const
number of daughters
virtual std::string print() const
Print object.
Definition: GenJet.cc:73
tuple result
Definition: query.py:137
Jets made from MC generator particles.
Definition: GenJet.h:25
unsigned int index
Definition: LeafCandidate.h:34
virtual std::vector< const GenParticle * > getGenConstituents() const
get all constituents
Definition: GenJet.cc:59
float invisibleEnergy() const
Definition: GenJet.h:62
virtual float eta() const GCC11_FINAL
momentum pseudorapidity
tuple out
Definition: dbtoconf.py:99
virtual const_iterator begin() const
first daughter const_iterator
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
virtual bool overlap(const Candidate &) const
Polymorphic overlap.
Definition: GenJet.cc:69
virtual const GenParticle * getGenConstituent(unsigned fIndex) const
get specific constituent
Definition: GenJet.cc:47
virtual const_iterator end() const
last daughter const_iterator
math::XYZPoint Point
point in the space
Definition: Candidate.h:45
static const GenParticle * genParticle(const reco::Candidate *fConstituent)
convert generic constituent to specific type
Definition: GenJet.cc:39
virtual GenJet * clone() const
Polymorphic clone.
Definition: GenJet.cc:65
virtual std::string print() const
Print object.
Definition: Jet.cc:425
virtual float pt() const GCC11_FINAL
transverse momentum
value_type const * get() const
Definition: RefToBase.h:212
virtual const CandidateBaseRef & masterClone() const =0