CMS 3D CMS Logo

PileUpSubtractor.h
Go to the documentation of this file.
1 #ifndef __PUSubtractor__
2 #define __PUSubtractor__
3 
4 #include <vector>
5 #include "boost/shared_ptr.hpp"
6 #include "fastjet/PseudoJet.hh"
7 #include "fastjet/ClusterSequence.hh"
8 #include "fastjet/ClusterSequenceArea.hh"
9 #include "fastjet/GhostedAreaSpec.hh"
10 
15 
20 
22 
23  public:
24 
25  typedef boost::shared_ptr<fastjet::ClusterSequence> ClusterSequencePtr;
26  typedef boost::shared_ptr<fastjet::GhostedAreaSpec> ActiveAreaSpecPtr;
27  typedef boost::shared_ptr<fastjet::RangeDefinition> RangeDefPtr;
28  typedef boost::shared_ptr<fastjet::JetDefinition> JetDefPtr;
29 
31  virtual ~PileUpSubtractor(){;}
32 
33 virtual void setDefinition(JetDefPtr const & jetDef);
34 virtual void reset(std::vector<edm::Ptr<reco::Candidate> >& input,
35  std::vector<fastjet::PseudoJet>& towers,
36  std::vector<fastjet::PseudoJet>& output);
37 virtual void setupGeometryMap(edm::Event& iEvent,const edm::EventSetup& iSetup);
38 virtual void calculatePedestal(std::vector<fastjet::PseudoJet> const & coll);
39 virtual void subtractPedestal(std::vector<fastjet::PseudoJet> & coll);
40 virtual void calculateOrphanInput(std::vector<fastjet::PseudoJet> & orphanInput);
41 virtual void offsetCorrectJets();
42 virtual double getMeanAtTower(const reco::CandidatePtr & in) const;
43 virtual double getSigmaAtTower(const reco::CandidatePtr & in) const;
44 virtual double getPileUpAtTower(const reco::CandidatePtr & in) const;
45 virtual double getPileUpEnergy(int ijet) const {return jetOffset_[ijet];}
46  virtual double getCone(double cone, double eta, double phi, double& et, double& pu);
47  int getN(const reco::CandidatePtr & in) const;
48  int getNwithJets(const reco::CandidatePtr & in) const;
49 
50  int ieta(const reco::CandidatePtr & in) const;
51  int iphi(const reco::CandidatePtr & in) const;
52 
53  protected:
54 
55  // From jet producer
56  JetDefPtr fjJetDefinition_; // fastjet jet definition
57  ClusterSequencePtr fjClusterSeq_; // fastjet cluster sequence
58  std::vector<edm::Ptr<reco::Candidate> >* inputs_; // input candidates
59  std::vector<fastjet::PseudoJet>* fjInputs_; // fastjet inputs
60  std::vector<fastjet::PseudoJet>* fjJets_; // fastjet jets
61  std::vector<fastjet::PseudoJet> fjOriginalInputs_; // to back-up unsubtracted fastjet inputs
62 
63  // PU subtraction parameters
64  bool reRunAlgo_;
67  double jetPtMin_;
68  double puPtMin_;
69 
70  double nSigmaPU_; // number of sigma for pileup
71  double radiusPU_; // pileup radius
72  ActiveAreaSpecPtr fjActiveArea_; // fastjet active area definition
73  RangeDefPtr fjRangeDef_; // range definition
74 
75  CaloGeometry const * geo_; // geometry
76  int ietamax_; // maximum eta in geometry
77  int ietamin_; // minimum eta in geometry
78  std::vector<HcalDetId> allgeomid_; // all det ids in the geometry
79  std::map<int,int> geomtowers_; // map of geometry towers to det id
80  std::map<int,int> ntowersWithJets_; // number of towers with jets
81  std::map<int,double> esigma_; // energy sigma
82  std::map<int,double> emean_; // energy mean
83 
84  std::vector<double> jetOffset_;
85 
86 };
87 
89 namespace edm {class ParameterSet; class EventSetup; class ConsumesCollector;}
91 
92 #endif
virtual void calculateOrphanInput(std::vector< fastjet::PseudoJet > &orphanInput)
std::vector< double > jetOffset_
virtual double getMeanAtTower(const reco::CandidatePtr &in) const
boost::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
int getNwithJets(const reco::CandidatePtr &in) const
std::vector< fastjet::PseudoJet > * fjJets_
std::map< int, double > esigma_
std::map< int, int > geomtowers_
std::vector< fastjet::PseudoJet > fjOriginalInputs_
virtual double getPileUpEnergy(int ijet) const
virtual double getPileUpAtTower(const reco::CandidatePtr &in) const
static std::string const input
Definition: EdmProvDump.cc:44
int ieta(const reco::CandidatePtr &in) const
virtual void offsetCorrectJets()
RangeDefPtr fjRangeDef_
std::vector< fastjet::PseudoJet > * fjInputs_
std::map< int, int > ntowersWithJets_
int iEvent
Definition: GenABIO.cc:230
boost::shared_ptr< fastjet::GhostedAreaSpec > ActiveAreaSpecPtr
virtual void subtractPedestal(std::vector< fastjet::PseudoJet > &coll)
edmplugin::PluginFactory< PileUpSubtractor *(const edm::ParameterSet &, edm::ConsumesCollector &&)> PileUpSubtractorFactory
int iphi(const reco::CandidatePtr &in) const
JetCorrectorParametersCollection coll
Definition: classes.h:10
ActiveAreaSpecPtr fjActiveArea_
virtual double getSigmaAtTower(const reco::CandidatePtr &in) const
ClusterSequencePtr fjClusterSeq_
boost::shared_ptr< fastjet::RangeDefinition > RangeDefPtr
virtual void setDefinition(JetDefPtr const &jetDef)
et
define resolution functions of each parameter
HLT enums.
CaloGeometry const * geo_
std::map< int, double > emean_
PileUpSubtractor(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
virtual void calculatePedestal(std::vector< fastjet::PseudoJet > const &coll)
virtual void setupGeometryMap(edm::Event &iEvent, const edm::EventSetup &iSetup)
virtual double getCone(double cone, double eta, double phi, double &et, double &pu)
virtual void reset(std::vector< edm::Ptr< reco::Candidate > > &input, std::vector< fastjet::PseudoJet > &towers, std::vector< fastjet::PseudoJet > &output)
std::vector< HcalDetId > allgeomid_
JetDefPtr fjJetDefinition_
std::vector< edm::Ptr< reco::Candidate > > * inputs_
int getN(const reco::CandidatePtr &in) const
boost::shared_ptr< fastjet::JetDefinition > JetDefPtr
virtual ~PileUpSubtractor()