CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EgammaTowerIsolation.h
Go to the documentation of this file.
1 #ifndef EgammaTowerIsolation_h
2 #define EgammaTowerIsolation_h
3 
4 //*****************************************************************************
5 // File: EgammaTowerIsolation.h
6 // ----------------------------------------------------------------------------
7 // OrigAuth: Matthias Mozer
8 // Institute: IIHE-VUB
9 // Adding feature to exclude towers used by H/E
10 //
11 // 11/11/12 Hack by VI to make it 100 times faster
12 //=============================================================================
13 //*****************************************************************************
14 
15 #include<vector>
16 
17 
18 //CMSSW includes
23 
24 
25 #include <cmath>
26 #include <algorithm>
27 #include <cstdint>
28 
30 
31 
32 
33 /*
34  for each set of cuts it will compute Et for all, depth1 and depth2 twice:
35  one between inner and outer and once inside outer vetoid the tower to excude
36 
37  */
38 template<unsigned int NC>
40  public:
41 
42  struct Sum {
43  Sum(): he{0},h2{0},heBC{0},h2BC{0}
44  {}
45  float he[NC];
46  float h2[NC];
47  float heBC[NC];
48  float h2BC[NC];
49  };
50 
51  // number of cuts
52  constexpr static unsigned int NCuts = NC;
53 
54  //constructors
55 
57  EgammaTowerIsolationNew (float extRadius[NC],
58  float intRadius[NC],
59  CaloTowerCollection const & towers) ;
60 
61 
63 
64  void compute(bool et, Sum&sum, reco::Candidate const & cand, CaloTowerDetId const * first, CaloTowerDetId const * last) const {
65  reco::SuperCluster const & sc = *cand.get<reco::SuperClusterRef>().get();
66  return compute(et,sum,sc,first,last);
67  }
68  void compute(bool et, Sum &sum, reco::SuperCluster const & sc, CaloTowerDetId const * first, CaloTowerDetId const * last) const;
69 
70  void setRadius(float const extRadius[NC],float const intRadius[NC]) {
71  for (std::size_t i=0; i!=NCuts; ++i) {
72  extRadius2_[i]=extRadius[i]*extRadius[i];
73  intRadius2_[i]=intRadius[i]*intRadius[i];
74  }
75  maxEta = *std::max_element(extRadius,extRadius+NC);
76  }
77 
78 public:
79 
80  float extRadius2_[NCuts] ;
81  float intRadius2_[NCuts] ;
82 
83  float maxEta;
84  //SOA
85  const uint32_t nt;
86  float * eta;
87  float * phi;
88  float * he;
89  float * h2;
90  float * st;
91  uint32_t * id;
92  uint32_t * mem=nullptr;
93  void initSoa() {
94  mem = new uint32_t[nt*6];
95  eta = (float*)(mem); phi = eta+nt; he = phi+nt; h2 = he+nt; st = h2+nt;
96  id = (uint32_t*)(st) + nt;
97  }
98 
99 
100 };
101 
102 namespace etiStat {
103 
104  struct Count {
105  uint32_t create=0;
106  uint32_t comp=0;
107  uint32_t span=0;
108  static Count count;
109  ~Count();
110  };
111 
112 }
113 
114 
115 template<unsigned int NC>
116 inline
118  float intRadius[NC],
119  CaloTowerCollection const & towers) :
120  maxEta(*std::max_element(extRadius,extRadius+NC)),
121  nt(towers.size()) {
122  if (nt==0) return;
123  initSoa();
124 
126 
127 
128  for (std::size_t i=0; i!=NCuts; ++i) {
129  extRadius2_[i]=extRadius[i]*extRadius[i];
130  intRadius2_[i]=intRadius[i]*intRadius[i];
131  }
132 
133  // sort in eta (kd-tree anoverkill,does not vectorize...)
134  uint32_t index[nt];
135  float e[nt];
136  for (std::size_t k=0; k!=nt; ++k) {
137  e[k]=towers[k].eta();
138  index[k]=k;
139  std::push_heap(index,index+k+1,[&e](uint32_t i, uint32_t j){ return e[i]<e[j];});
140  }
141  std::sort_heap(index,index+nt,[&e](uint32_t i, uint32_t j){ return e[i]<e[j];});
142 
143 
144  for (std::size_t i=0;i!=nt; ++i) {
145  auto j = index[i];
146  eta[i]=towers[j].eta();
147  phi[i]=towers[j].phi();
148  id[i]=towers[j].id();
149  st[i] = 1.f/std::cosh(eta[i]); // std::sin(towers[j].theta()); //;
150  he[i] = towers[j].hadEnergy();
151  h2[i] = towers[j].hadEnergyHeOuterLayer();
152  }
153 }
154 
155 
156 template<unsigned int NC>
157 inline
158 void
160  if (nt==0) return;
161 
163 
164  float candEta = sc.eta();
165  float candPhi = sc.phi();
166 
167 
168  auto lb = std::lower_bound(eta,eta+nt,candEta-maxEta);
169  auto ub = std::upper_bound(lb,eta+nt,candEta+maxEta);
170  uint32_t il = lb-eta;
171  uint32_t iu = std::min(nt,uint32_t(ub-eta+1));
172 
173  etiStat::Count::count.span += (iu-il);
174 
175  bool ok[iu-il];
176  for (std::size_t i=il;i!=iu; ++i)
177  ok[i-il] = (std::find(first,last,id[i])==last);
178 
179 
180 
181  // should be restricted in eta....
182  for (std::size_t i=il;i!=iu; ++i) {
183  float dr2 = reco::deltaR2(candEta,candPhi,eta[i], phi[i]);
184  float tt = et ? st[i] : 1.f;
185  for (std::size_t j=0; j!=NCuts; ++j) {
186  if (dr2<extRadius2_[j]) {
187  if (dr2>=intRadius2_[j]) {
188  sum.he[j] +=he[i]*tt;
189  sum.h2[j] +=h2[i]*tt;
190  }
191  if(ok[i-il]) {
192  sum.heBC[j] +=he[i]*tt;
193  sum.h2BC[j] +=h2[i]*tt;
194  }
195  }
196  }
197  }
198 }
199 
201 public:
202 
204 
205  //constructors
206  EgammaTowerIsolation (float extRadiusI,
207  float intRadiusI,
208  float etLow,
209  signed int depth,
210  const CaloTowerCollection* towers );
211 
212  double getTowerEtSum (const reco::Candidate* cand, const std::vector<CaloTowerDetId> * detIdToExclude=0 ) const{
213  reco::SuperCluster const & sc = *cand->get<reco::SuperClusterRef>().get();
214  return getSum(true,sc,detIdToExclude);
215  }
216  double getTowerESum (const reco::Candidate* cand, const std::vector<CaloTowerDetId> * detIdToExclude=0) const{
217  reco::SuperCluster const & sc = *cand->get<reco::SuperClusterRef>().get();
218  return getSum(false,sc,detIdToExclude);
219  }
220  double getTowerEtSum (reco::SuperCluster const * sc, const std::vector<CaloTowerDetId> * detIdToExclude=0 ) const{
221  return getSum(true,*sc,detIdToExclude);
222  }
223  double getTowerESum (reco::SuperCluster const * sc, const std::vector<CaloTowerDetId> * detIdToExclude=0) const{
224  return getSum(false,*sc,detIdToExclude);
225  }
226  private:
227  double getSum (bool et, reco::SuperCluster const & sc, const std::vector<CaloTowerDetId> * detIdToExclude) const;
228 
229 
230 private:
233  static uint32_t id15;
234  signed int depth_;
235  float extRadius;
236  float intRadius;
237 };
238 
239 
240 
241 #endif
int i
Definition: DBlmapReader.cc:9
double getSum(bool et, reco::SuperCluster const &sc, const std::vector< CaloTowerDetId > *detIdToExclude) const
double maxEta
static constexpr unsigned int NCuts
T eta() const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
double getTowerESum(reco::SuperCluster const *sc, const std::vector< CaloTowerDetId > *detIdToExclude=0) const
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:161
static EgammaTowerIsolationNew< 1 > * newAlgo
void setRadius(float const extRadius[NC], float const intRadius[NC])
int j
Definition: DBlmapReader.cc:9
static const CaloTowerCollection * oldTowers
bool first
Definition: L1TdeRCT.cc:79
double getTowerESum(const reco::Candidate *cand, const std::vector< CaloTowerDetId > *detIdToExclude=0) const
int nt
Definition: AMPTWrapper.h:32
int k[5][pyjets_maxn]
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:58
double getTowerEtSum(const reco::Candidate *cand, const std::vector< CaloTowerDetId > *detIdToExclude=0) const
T get() const
get a component
Definition: Candidate.h:219
EgammaTowerIsolation(float extRadiusI, float intRadiusI, float etLow, signed int depth, const CaloTowerCollection *towers)
double getTowerEtSum(reco::SuperCluster const *sc, const std::vector< CaloTowerDetId > *detIdToExclude=0) const
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:164
tuple size
Write out results.
#define constexpr
void compute(bool et, Sum &sum, reco::Candidate const &cand, CaloTowerDetId const *first, CaloTowerDetId const *last) const
Definition: DDAxes.h:10