Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:31

0001 #ifndef RecoJets_JetProducers_interface_JetSpecific_h
0002 #define RecoJets_JetProducers_interface_JetSpecific_h
0003 
0004 #include "DataFormats/Candidate/interface/Candidate.h"
0005 #include "DataFormats/JetReco/interface/CaloJet.h"
0006 #include "DataFormats/JetReco/interface/PFJet.h"
0007 #include "DataFormats/JetReco/interface/GenJet.h"
0008 #include "DataFormats/JetReco/interface/TrackJet.h"
0009 #include "DataFormats/JetReco/interface/PFClusterJet.h"
0010 #include "DataFormats/JetReco/interface/BasicJet.h"
0011 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
0012 #include "Geometry/CaloTopology/interface/HcalTopology.h"
0013 
0014 class CaloSubdetectorGeometry;
0015 class CaloGeometry;
0016 
0017 namespace reco {
0018 
0019   //______________________________________________________________________________
0020   // Helper methods to write out specific types
0021 
0022   /// Make CaloJet specifics. Assumes PseudoJet is made from CaloTowerCandidates
0023   bool makeSpecific(std::vector<reco::CandidatePtr> const& towers,
0024                     const CaloSubdetectorGeometry* towerGeometry,
0025                     reco::CaloJet::Specific* caloJetSpecific,
0026                     const HcalTopology& topology);
0027 
0028   void writeSpecific(reco::CaloJet& jet,
0029                      reco::Particle::LorentzVector const& p4,
0030                      reco::Particle::Point const& point,
0031                      std::vector<reco::CandidatePtr> const& constituents,
0032                      CaloGeometry const& geometry,
0033                      HcalTopology const& topology);
0034 
0035   /// Make PFlowJet specifics. Assumes PseudoJet is made from ParticleFlowCandidates
0036   bool makeSpecific(std::vector<reco::CandidatePtr> const& particles,
0037                     reco::PFJet::Specific* pfJetSpecific,
0038                     edm::ValueMap<float> const* weights = nullptr);
0039 
0040   void writeSpecific(reco::PFJet& jet,
0041                      reco::Particle::LorentzVector const& p4,
0042                      reco::Particle::Point const& point,
0043                      std::vector<reco::CandidatePtr> const& constituents,
0044                      edm::ValueMap<float> const* weights = nullptr);
0045 
0046   /// Make GenJet specifics. Assumes PseudoJet is made from HepMCCandidate
0047   bool makeSpecific(std::vector<reco::CandidatePtr> const& mcparticles, reco::GenJet::Specific* genJetSpecific);
0048 
0049   void writeSpecific(reco::GenJet& jet,
0050                      reco::Particle::LorentzVector const& p4,
0051                      reco::Particle::Point const& point,
0052                      std::vector<reco::CandidatePtr> const& constituents);
0053 
0054   /// Make TrackJet. Assumes constituents point to tracks, through RecoChargedCandidates.
0055   void writeSpecific(reco::TrackJet& jet,
0056                      reco::Particle::LorentzVector const& p4,
0057                      reco::Particle::Point const& point,
0058                      std::vector<reco::CandidatePtr> const& constituents);
0059 
0060   /// Make PFClusterJet. Assumes PseudoJet is made from PFCluster
0061   void writeSpecific(reco::PFClusterJet& jet,
0062                      reco::Particle::LorentzVector const& p4,
0063                      reco::Particle::Point const& point,
0064                      std::vector<reco::CandidatePtr> const& constituents);
0065 
0066   /// Make BasicJet. Assumes nothing about the jet.
0067   void writeSpecific(reco::BasicJet& jet,
0068                      reco::Particle::LorentzVector const& p4,
0069                      reco::Particle::Point const& point,
0070                      std::vector<reco::CandidatePtr> const& constituents);
0071 
0072   /// converts eta to the corresponding HCAL subdetector.
0073   HcalSubdetector hcalSubdetector(int iEta, const HcalTopology& topology);
0074 
0075 }  // namespace reco
0076 
0077 #endif