CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HiFJGridEmptyAreaCalculator.h
Go to the documentation of this file.
1 #ifndef HiJetBackground_HiFJGridEmptyAreaCalculator_h
2 #define HiJetBackground_HiFJGridEmptyAreaCalculator_h
3 
4 // system include files
5 #include <memory>
6 #include <sstream>
7 #include <string>
8 #include <vector>
9 
10 // user include files
13 
16 
19 
23 
25 public:
28 
29  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
30 
31 private:
32  void beginStream(edm::StreamID) override;
33  void produce(edm::Event&, const edm::EventSetup&) override;
34  void endStream() override;
35 
37  //\{
38  //----------------------------------------------------------------
39 
42 
43 private:
45  void setupGrid(double eta_min, double eta_max);
46  void setupGridJet(const reco::Jet* jet);
47 
49  int tileIndexJet(const reco::PFCandidate* pfCand);
50  int tileIndexEta(const reco::PFCandidate* pfCand);
51  int tileIndexEtaJet(const reco::PFCandidate* pfCand);
52  int tileIndexPhi(const reco::PFCandidate* pfCand);
53 
55  int numJetGridCells(std::vector<std::pair<int, int>>& indices);
56 
60  void calculateGridRho(const edm::Event& iEvent, const edm::EventSetup& iSetup);
61 
63  const double twopi_ = 2 * M_PI;
64 
66 
67  //parameters of grid covering the full acceptance
68  double ymin_;
69  double ymax_;
70  double dy_;
71  double dphi_;
72  double tileArea_;
73 
74  //parameters of grid around jets
75  double dyJet_;
76  double yminJet_;
77  double ymaxJet_;
79 
80  //leave bands at boundaries
81  double etaminJet_;
82  double etamaxJet_;
83 
84  //for the grid calculation covering the full acceptance
85  int ny_;
86  int nphi_;
87  int ntotal_;
88 
89  //for the grid calculation around each jet
91  int nyJet_;
92 
94  double gridWidth_;
95  double band_;
96  int hiBinCut_;
99 
100  std::vector<double> rhoVsEta_;
101  std::vector<double> meanRhoVsEta_;
102  std::vector<double> etaMaxGrid_;
103  std::vector<double> etaMinGrid_;
104 
105  int n_tiles() { return ntotal_; }
106 
113 
115 };
116 
117 #endif
edm::EDGetTokenT< std::vector< double > > mapRhoToken_
edm::EDGetTokenT< int > centralityBinToken_
int numJetGridCells(std::vector< std::pair< int, int >> &indices)
number of grid cells that overlap with jet constituents filling in the in between area ...
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double ymin_
internal parameters for grid
edm::EDGetTokenT< edm::View< reco::Jet > > jetsToken_
input tokens
int tileIndexEta(const reco::PFCandidate *pfCand)
edm::EDGetTokenT< reco::PFCandidateCollection > pfCandsToken_
Base class for all types of Jets.
Definition: Jet.h:20
void setupGrid(double eta_min, double eta_max)
configure the grid
edm::EDGetTokenT< std::vector< double > > mapRhoMToken_
int tileIndexEtaJet(const reco::PFCandidate *pfCand)
void setupGridJet(const reco::Jet *jet)
edm::EDGetTokenT< std::vector< double > > mapEtaToken_
int tileIndexPhi(const reco::PFCandidate *pfCand)
int iEvent
Definition: GenABIO.cc:224
void calculateAreaFractionOfJets(const edm::Event &iEvent, const edm::EventSetup &iSetup)
const double twopi_
information about the grid
void beginStream(edm::StreamID) override
#define M_PI
void calculateGridRho(const edm::Event &iEvent, const edm::EventSetup &iSetup)
HiFJGridEmptyAreaCalculator(const edm::ParameterSet &)
void produce(edm::Event &, const edm::EventSetup &) override
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
int tileIndexJet(const reco::PFCandidate *pfCand)
retrieve the grid cell index for a given PseudoJet
list indices
Definition: dqmdumpme.py:50