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 #include <atomic>
29 
31 
32 
33 
34 /*
35  for each set of cuts it will compute Et for all, depth1 and depth2 twice:
36  one between inner and outer and once inside outer vetoid the tower to excude
37 
38  */
39 template<unsigned int NC>
41  public:
42 
43  struct Sum {
44  Sum(): he{0},h2{0},heBC{0},h2BC{0}
45  {}
46  float he[NC];
47  float h2[NC];
48  float heBC[NC];
49  float h2BC[NC];
50  };
51 
52  // number of cuts
53  constexpr static unsigned int NCuts = NC;
54 
55  //constructors
56 
58  EgammaTowerIsolationNew (float extRadius[NC],
59  float intRadius[NC],
60  CaloTowerCollection const & towers) ;
61 
62 
64 
65  void compute(bool et, Sum&sum, reco::Candidate const & cand, CaloTowerDetId const * first, CaloTowerDetId const * last) const {
66  reco::SuperCluster const & sc = *cand.get<reco::SuperClusterRef>().get();
67  return compute(et,sum,sc,first,last);
68  }
69  void compute(bool et, Sum &sum, reco::SuperCluster const & sc, CaloTowerDetId const * first, CaloTowerDetId const * last) const;
70 
71  void setRadius(float const extRadius[NC],float const intRadius[NC]) {
72  for (std::size_t i=0; i!=NCuts; ++i) {
73  extRadius2_[i]=extRadius[i]*extRadius[i];
74  intRadius2_[i]=intRadius[i]*intRadius[i];
75  }
76  maxEta = *std::max_element(extRadius,extRadius+NC);
77  }
78 
79 public:
80 
81  float extRadius2_[NCuts] ;
82  float intRadius2_[NCuts] ;
83 
84  float maxEta;
85  //SOA
86  const uint32_t nt;
87  float * eta;
88  float * phi;
89  float * he;
90  float * h2;
91  float * st;
92  uint32_t * id;
93  uint32_t * mem=nullptr;
94  void initSoa() {
95  mem = new uint32_t[nt*6];
96  eta = (float*)(mem); phi = eta+nt; he = phi+nt; h2 = he+nt; st = h2+nt;
97  id = (uint32_t*)(st) + nt;
98  }
99 
100 
101 };
102 
103 /*#define ETISTATDEBUG*/
104 #ifdef ETISTATDEBUG
105 namespace etiStat {
106 
107  struct Count {
108  std::atomic<uint32_t> create=0;
109  std::atomic<uint32_t> comp=0;
110  std::atomic<uint32_t> span=0;
111  static Count count;
112  ~Count();
113  };
114 
115 }
116 #endif
117 
118 template<unsigned int NC>
119 inline
121  float intRadius[NC],
122  CaloTowerCollection const & towers) :
123  maxEta(*std::max_element(extRadius,extRadius+NC)),
124  nt(towers.size()) {
125  if (nt==0) return;
126  initSoa();
127 
128 #ifdef ETISTATDEBUG
129  etiStat::Count::count.create++;
130 #endif
131 
132  for (std::size_t i=0; i!=NCuts; ++i) {
133  extRadius2_[i]=extRadius[i]*extRadius[i];
134  intRadius2_[i]=intRadius[i]*intRadius[i];
135  }
136 
137  // sort in eta (kd-tree anoverkill,does not vectorize...)
138  uint32_t index[nt];
139  #ifdef __clang__
140  std::vector<float> e(nt);
141  #else
142  float e[nt];
143  #endif
144  for (std::size_t k=0; k!=nt; ++k) {
145  e[k]=towers[k].eta();
146  index[k]=k;
147  std::push_heap(index,index+k+1,[&e](uint32_t i, uint32_t j){ return e[i]<e[j];});
148  }
149  std::sort_heap(index,index+nt,[&e](uint32_t i, uint32_t j){ return e[i]<e[j];});
150 
151 
152  for (std::size_t i=0;i!=nt; ++i) {
153  auto j = index[i];
154  eta[i]=towers[j].eta();
155  phi[i]=towers[j].phi();
156  id[i]=towers[j].id();
157  st[i] = 1.f/std::cosh(eta[i]); // std::sin(towers[j].theta()); //;
158  he[i] = towers[j].hadEnergy();
159  h2[i] = towers[j].hadEnergyHeOuterLayer();
160  }
161 }
162 
163 
164 template<unsigned int NC>
165 inline
166 void
168  if (nt==0) return;
169 
170 #ifdef ETISTATDEBUG
171  etiStat::Count::count.comp++;
172 #endif
173 
174  float candEta = sc.eta();
175  float candPhi = sc.phi();
176 
177 
178  auto lb = std::lower_bound(eta,eta+nt,candEta-maxEta);
179  auto ub = std::upper_bound(lb,eta+nt,candEta+maxEta);
180  uint32_t il = lb-eta;
181  uint32_t iu = std::min(nt,uint32_t(ub-eta+1));
182 
183 #ifdef ETISTATDEBUG
184  etiStat::Count::count.span += (iu-il);
185 #endif
186 
187  bool ok[iu-il];
188  for (std::size_t i=il;i!=iu; ++i)
189  ok[i-il] = (std::find(first,last,id[i])==last);
190 
191 
192 
193  // should be restricted in eta....
194  for (std::size_t i=il;i!=iu; ++i) {
195  float dr2 = reco::deltaR2(candEta,candPhi,eta[i], phi[i]);
196  float tt = et ? st[i] : 1.f;
197  for (std::size_t j=0; j!=NCuts; ++j) {
198  if (dr2<extRadius2_[j]) {
199  if (dr2>=intRadius2_[j]) {
200  sum.he[j] +=he[i]*tt;
201  sum.h2[j] +=h2[i]*tt;
202  }
203  if(ok[i-il]) {
204  sum.heBC[j] +=he[i]*tt;
205  sum.h2BC[j] +=h2[i]*tt;
206  }
207  }
208  }
209  }
210 }
211 
213 public:
214 
216 
217  //constructors
218  EgammaTowerIsolation (float extRadiusI,
219  float intRadiusI,
220  float etLow,
221  signed int depth,
222  const CaloTowerCollection* towers );
223 
224  double getTowerEtSum (const reco::Candidate* cand, const std::vector<CaloTowerDetId> * detIdToExclude=0 ) const{
225  reco::SuperCluster const & sc = *cand->get<reco::SuperClusterRef>().get();
226  return getSum(true,sc,detIdToExclude);
227  }
228  double getTowerESum (const reco::Candidate* cand, const std::vector<CaloTowerDetId> * detIdToExclude=0) const{
229  reco::SuperCluster const & sc = *cand->get<reco::SuperClusterRef>().get();
230  return getSum(false,sc,detIdToExclude);
231  }
232  double getTowerEtSum (reco::SuperCluster const * sc, const std::vector<CaloTowerDetId> * detIdToExclude=0 ) const{
233  return getSum(true,*sc,detIdToExclude);
234  }
235  double getTowerESum (reco::SuperCluster const * sc, const std::vector<CaloTowerDetId> * detIdToExclude=0) const{
236  return getSum(false,*sc,detIdToExclude);
237  }
238  private:
239  double getSum (bool et, reco::SuperCluster const & sc, const std::vector<CaloTowerDetId> * detIdToExclude) const;
240 
241 
242 private:
243  signed int depth_;
244  float extRadius;
245  float intRadius;
246 };
247 
248 
249 
250 #endif
int i
Definition: DBlmapReader.cc:9
double getSum(bool et, reco::SuperCluster const &sc, const std::vector< CaloTowerDetId > *detIdToExclude) const
double maxEta
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
#define constexpr
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:163
void setRadius(float const extRadius[NC], float const intRadius[NC])
int j
Definition: DBlmapReader.cc:9
static unsigned int NCuts
T min(T a, T b)
Definition: MathUtil.h:58
double getTowerESum(const reco::Candidate *cand, const std::vector< CaloTowerDetId > *detIdToExclude=0) const
int nt
Definition: AMPTWrapper.h:32
Geom::Phi< T > phi() const
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
double getTowerEtSum(const reco::Candidate *cand, const std::vector< CaloTowerDetId > *detIdToExclude=0) const
T get() const
get a component
Definition: Candidate.h:217
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:166
SurfaceDeformation * create(int type, const std::vector< double > &params)
tuple size
Write out results.
void compute(bool et, Sum &sum, reco::Candidate const &cand, CaloTowerDetId const *first, CaloTowerDetId const *last) const