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 
17 
22 
24 
25  public:
26 
27  typedef boost::shared_ptr<fastjet::ClusterSequence> ClusterSequencePtr;
28  typedef boost::shared_ptr<fastjet::GhostedAreaSpec> ActiveAreaSpecPtr;
29  typedef boost::shared_ptr<fastjet::RangeDefinition> RangeDefPtr;
30  typedef boost::shared_ptr<fastjet::JetDefinition> JetDefPtr;
31 
33  virtual ~PileUpSubtractor(){;}
34 
35  virtual void setDefinition(JetDefPtr const & jetDef);
36  virtual void reset(std::vector<edm::Ptr<reco::Candidate> >& input,
37  std::vector<fastjet::PseudoJet>& towers,
38  std::vector<fastjet::PseudoJet>& output);
39  virtual void setupGeometryMap(edm::Event& iEvent,const edm::EventSetup& iSetup);
40  virtual void calculatePedestal(std::vector<fastjet::PseudoJet> const & coll);
41  virtual void subtractPedestal(std::vector<fastjet::PseudoJet> & coll);
42  virtual void calculateOrphanInput(std::vector<fastjet::PseudoJet> & orphanInput);
43  virtual void offsetCorrectJets();
44  virtual double getMeanAtTower(const reco::CandidatePtr & in) const;
45  virtual double getSigmaAtTower(const reco::CandidatePtr & in) const;
46  virtual double getPileUpAtTower(const reco::CandidatePtr & in) const;
47  virtual double getPileUpEnergy(int ijet) const {return jetOffset_[ijet];}
48  virtual double getCone(double cone, double eta, double phi, double& et, double& pu);
49  int getN(const reco::CandidatePtr & in) const;
50  int getNwithJets(const reco::CandidatePtr & in) const;
51 
52  int ieta(const reco::CandidatePtr & in) const;
53  int iphi(const reco::CandidatePtr & in) const;
54 
55  protected:
56 
57  // From jet producer
58  JetDefPtr fjJetDefinition_; // fastjet jet definition
59  ClusterSequencePtr fjClusterSeq_; // fastjet cluster sequence
60  std::vector<edm::Ptr<reco::Candidate> >* inputs_; // input candidates
61  std::vector<fastjet::PseudoJet>* fjInputs_; // fastjet inputs
62  std::vector<fastjet::PseudoJet>* fjJets_; // fastjet jets
63  std::vector<fastjet::PseudoJet> fjOriginalInputs_; // to back-up unsubtracted fastjet inputs
64 
65  // PU subtraction parameters
66  bool reRunAlgo_;
69  double jetPtMin_;
70  double puPtMin_;
71 
72  double ghostEtaMax;
74  double ghostArea;
75 
76  double nSigmaPU_; // number of sigma for pileup
77  double radiusPU_; // pileup radius
78  ActiveAreaSpecPtr fjActiveArea_; // fastjet active area definition
79  CaloGeometry const * geo_; // geometry
80  int ietamax_; // maximum eta in geometry
81  int ietamin_; // minimum eta in geometry
82  std::vector<HcalDetId> allgeomid_; // all det ids in the geometry
83  std::map<int,int> geomtowers_; // map of geometry towers to det id
84  std::map<int,int> ntowersWithJets_; // number of towers with jets
85  std::map<int,double> esigma_; // energy sigma
86  std::map<int,double> emean_; // energy mean
87 
88  std::vector<double> jetOffset_;
89 
90 };
91 
93 namespace edm {class ParameterSet; class EventSetup; class ConsumesCollector;}
95 
96 #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:45
int ieta(const reco::CandidatePtr &in) const
virtual void offsetCorrectJets()
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()