CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/DataFormats/JetReco/src/PFJet.cc

Go to the documentation of this file.
00001 // PFJet.cc
00002 // Fedor Ratnikov UMd
00003 // $Id: PFJet.cc,v 1.17 2010/03/10 21:52:18 pandolf Exp $
00004 #include <sstream>
00005 #include <typeinfo>
00006 
00007 #include "FWCore/Utilities/interface/Exception.h"
00008 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00009 #include "DataFormats/ParticleFlowReco/interface/PFBlock.h"
00010 
00011 //Own header file
00012 #include "DataFormats/JetReco/interface/PFJet.h"
00013 
00014 using namespace reco;
00015 
00016 PFJet::PFJet (const LorentzVector& fP4, const Point& fVertex, 
00017                   const Specific& fSpecific, 
00018                   const Jet::Constituents& fConstituents)
00019   : Jet (fP4, fVertex, fConstituents),
00020     m_specific (fSpecific)
00021 {}
00022 
00023 PFJet::PFJet (const LorentzVector& fP4, const Point& fVertex, 
00024               const Specific& fSpecific)
00025   : Jet (fP4, fVertex),
00026     m_specific (fSpecific)
00027 {}
00028 
00029 PFJet::PFJet (const LorentzVector& fP4, 
00030                   const Specific& fSpecific, 
00031                   const Jet::Constituents& fConstituents)
00032   : Jet (fP4, Point(0,0,0), fConstituents),
00033     m_specific (fSpecific)
00034 {}
00035 
00036 reco::PFCandidatePtr PFJet::getPFConstituent (unsigned fIndex) const {
00037 
00038   Constituent dau = daughterPtr (fIndex);
00039   if ( dau.isNonnull() && dau.isAvailable() ) {
00040     const PFCandidate* pfCandidate = dynamic_cast <const PFCandidate*> (dau.get());
00041     if (pfCandidate) {
00042       return edm::Ptr<PFCandidate> (dau.id(), pfCandidate, dau.key() );
00043     }
00044     else {
00045       throw cms::Exception("Invalid Constituent") << "PFJet constituent is not of PFCandidate type";
00046     }
00047    }
00048    else {
00049      return PFCandidatePtr();
00050    }
00051 }
00052 
00053 std::vector <reco::PFCandidatePtr> PFJet::getPFConstituents () const {
00054   std::vector <PFCandidatePtr> result;
00055   for (unsigned i = 0;  i <  numberOfDaughters (); i++) result.push_back (getPFConstituent(i));
00056   return result;
00057 }
00058 
00059 
00060 reco::TrackRefVector PFJet::getTrackRefs() const {
00061   // result will contain chargedMultiplicity() elements
00062   reco::TrackRefVector result;
00063   result.reserve( chargedMultiplicity() );
00064   for (unsigned i = 0;  i <  numberOfDaughters (); i++) {
00065     const reco::PFCandidatePtr pfcand = getPFConstituent (i);
00066     reco::TrackRef trackref = pfcand->trackRef();
00067     if( trackref.isNonnull() ) {
00068       result.push_back( trackref );
00069     }
00070   }
00071 
00072   return result;
00073 }
00074 
00075 
00076 PFJet* PFJet::clone () const {
00077   return new PFJet (*this);
00078 }
00079 
00080 bool PFJet::overlap( const Candidate & ) const {
00081   return false;
00082 }
00083 
00084 std::string PFJet::print () const {
00085   std::ostringstream out;
00086   out << Jet::print () // generic jet info
00087       << "    PFJet specific:" << std::endl
00088       << "      charged hadron energy/multiplicity: " << chargedHadronEnergy () << '/' << chargedHadronMultiplicity () << std::endl
00089       << "      neutral hadron energy/multiplicity: " << neutralHadronEnergy () << '/' << neutralHadronMultiplicity () << std::endl
00090       << "      photon energy/multiplicity: " << photonEnergy () << '/' << photonMultiplicity () << std::endl
00091       << "      electron energy/multiplicity: " << electronEnergy () << '/' << electronMultiplicity () << std::endl
00092       << "      muon energy/multiplicity: " << muonEnergy () << '/' << muonMultiplicity () << std::endl
00093       << "      HF Hadron energy/multiplicity: " << HFHadronEnergy () << '/' << HFHadronMultiplicity () << std::endl
00094       << "      HF EM particle energy/multiplicity: " << HFEMEnergy () << '/' << HFEMMultiplicity () << std::endl
00095       << "      charged/neutral hadrons energy: " << chargedHadronEnergy () << '/' << neutralHadronEnergy () << std::endl
00096       << "      charged/neutral em energy: " << chargedEmEnergy () << '/' << neutralEmEnergy () << std::endl
00097       << "      charged muon energy: " << chargedMuEnergy () << '/' << std::endl
00098       << "      charged/neutral multiplicity: " << chargedMultiplicity () << '/' << neutralMultiplicity () << std::endl;
00099   out << "      PFCandidate constituents:" << std::endl;
00100   std::vector <PFCandidatePtr> constituents = getPFConstituents ();
00101   for (unsigned i = 0; i < constituents.size (); ++i) {
00102     if (constituents[i].get()) {
00103       out << "      #" << i << " " << *(constituents[i]) << std::endl;
00104     }
00105     else {
00106       out << "      #" << i << " PFCandidate is not available in the event"  << std::endl;
00107     }
00108   }
00109   return out.str ();
00110 }
00111 
00112 std::ostream& reco::operator<<(std::ostream& out, const reco::PFJet& jet) {
00113 
00114   if(out) {
00115     out<<"PFJet "
00116        <<"(pt, eta, phi) = "<<jet.pt()<<","<<jet.eta()<<","<<jet.phi()
00117        <<"  (Rch,Rnh,Rgamma,Re,Rmu,RHFHad,RHFEM) = "
00118        <<jet.chargedHadronEnergyFraction()<<","
00119        <<jet.neutralHadronEnergyFraction()<<","
00120        <<jet.photonEnergyFraction()<<","
00121        <<jet.electronEnergyFraction()<<","
00122        <<jet.muonEnergyFraction()<<","
00123        <<jet.HFHadronEnergyFraction()<<","
00124        <<jet.HFEMEnergyFraction();
00125   }
00126   return out;
00127 }