ATLAS Offline Software
EgammaReSamp2Fex.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 // ********************************************************************
6 //
7 // NAME: EgammaReSamp2Fex.cxx
8 // PACKAGE: Trigger/TrigAlgorithms/TrigT2CaloEgamma
9 //
10 // AUTHOR: M.P. Casado, D.O. Damazio
11 //
12 //
13 // ********************************************************************
14 
15 #include "EgammaReSamp2Fex.h"
16 #include "CaloGeoHelpers/CaloSampling.h"
21 
22 EgammaReSamp2Fex::EgammaReSamp2Fex(const std::string& type, const std::string& name,
23  const IInterface* parent) :
25 {
26 }
27 
29  const CaloDetDescrElement*& caloDDE,
30  const EventContext& context) const
31 {
32  ATH_MSG_DEBUG("in execute(TrigEMCluster&)");
33 
34  // Region Selector, sampling 2
35  int sampling = 2;
36 
38  ATH_CHECK( m_dataSvc->loadCollections(context, roi, TTEM, sampling, sel) );
39 
40  double energyEta = 0.;
41  double energyPhi = 0.;
42 
43  // add these variables to take care of phi wrap-around
44  double energyNegPhi = 0.; // SRA
45  double energyNegPhiConv = 0.; // SRA
46  double energyPosPhi = 0.; // SRA
47 
48  // 1. Find seed cell (highest Et in ROI .. layer 2)
49  // 2. Find Et weighted eta, phi in 3*7 cell (layer 2) (photon + e id)
50  // 3. Find Et in cells of sizes 3*3, 3*7, 7*7 (layer 2 + strips)
51  // (3*7 for photon + e id)
52  // 4. Find cluster width in 3*5 cell, layer 2 (photon id, needs
53  // parabolic parametrisation)
54  // 5. Find strip energies and eta (2*5 window)
55  // 6. Find frac73 (photon id), (E1-E2)/(E1+E2) (e + photon id)
56 
57  double seedEnergy = 0.;
58  double seedPhi = 999.;
59  double seedEta = 999.;
60  double hotPhi = 999.;
61  double hotEta = 999.;
62  int ncells = 0;
63  // LVL1 positions
64  float etaL1 = rtrigEmCluster.eta();
65  float phiL1 = rtrigEmCluster.phi();
66 
67  const LArEM_ID* emID = m_larMgr->getEM_ID();
68  const LArCell* seedCell = nullptr;
69  const LArCell* hotCell = nullptr;
70  for (const LArCell* larcell : sel) {
71  if (larcell->energy() > seedEnergy) { // Hottest cell seach
72  float deta = std::abs(etaL1 - larcell->eta());
73  if (deta < m_maxHotCellDeta) { // Eta check is faster. Do it First
74  float dphi = std::abs(phiL1 - larcell->phi());
75  dphi = std::abs(M_PI - dphi);
76  dphi = std::abs(M_PI - dphi);
77  if (dphi < m_maxHotCellDphi) {
78  seedEnergy = larcell->energy();
79  seedCell = larcell;
80  } // End of dphi check
81  } // End of deta check
82  } // End of if energy
83  ncells++;
84  }
85  if (seedCell != nullptr) {
86  seedEta = seedCell->eta();
87  seedPhi = seedCell->phi();
88  // For the S-shape correction, we store the caloDDE of the hottest cell
89  caloDDE = (seedCell->caloDDE());
90  hotCell = seedCell;
91  hotEta = hotCell->eta();
92  hotPhi = hotCell->phi();
93  }
94  else {
95  return StatusCode::SUCCESS;
96  }
97 
98  std::map<const LArCell*, float> windows;
99  for (const LArCell* larcell : sel) {
100  float deta = std::abs(seedEta - larcell->eta());
101  if (deta < 0.025 + 0.002) { // Eta check is faster.
102  // Do it First 0.025 is cell, plus a little loose
103  float dphi = std::abs(seedPhi - larcell->phi());
104  dphi = std::abs(M_PI - dphi);
105  dphi = std::abs(M_PI - dphi);
106  if (dphi < 0.025 + 0.002) { // the same here (around 2*pi/64/4, but ok)
107  windows.emplace(larcell, 0.0);
108  } // End of dphi check
109  } // End of deta check
110  ncells++;
111  }
112 
113  for (const LArCell* larcell : sel) {
114  double etaCell = larcell->eta();
115  double phiCell = larcell->phi();
116  double energyCell = larcell->et();
117 
118  // Find position of current cell w.r.t. seed
119  float deta = std::abs(etaCell - seedEta);
120  if (deta > 0.05 + 0.002) continue; // find the energy weighted position with larger cluster
121  for (auto& [cell, energy] : windows) {
122  float deta1 = std::abs(etaCell - cell->eta());
123  float dphi = std::abs(phiCell - cell->phi());
124  if (dphi > M_PI) dphi = 2. * M_PI - dphi; // wrap (pi - 2pi)->(-pi - 0)
125 
126  // Electromagnetic measurements, done by layer
127  // layer 2: Et weighted eta, phi, summed Et in 3*7,
128  // cluster width in 3*5 cell
129  // emEnergy in 3*7 cell is also calculated, all samplings
130 
131  double cellSizeEta;
132  double cellSizePhi;
133  if (emID->is_em_barrel(larcell->ID())) {
134  cellSizeEta = etaSizeLArEMSamp2(etaCell, Calorimeter::EMBAR);
135  cellSizePhi = phiSizeLArEMSamp2(etaCell, Calorimeter::EMBAR);
136  }
137  else {
138  cellSizeEta = etaSizeLArEMSamp2(etaCell, Calorimeter::EMEND);
139  cellSizePhi = phiSizeLArEMSamp2(etaCell, Calorimeter::EMEND);
140  }
141 
142  // 3. do the 3*5 stuff
143  int nCellsEta = 3;
144  int nCellsPhi = 5;
145  if (deta1 <= 0.5 * double(nCellsEta - 1) * cellSizeEta + 0.01 &&
146  dphi <= 0.5 * double(nCellsPhi - 1) * cellSizePhi + 0.01) {
147  energy += energyCell;
148  } // End of do the 3*5 stuff
149  } // End of windows map loop
150  } // end of the loop
151 
152  float max = -999;
153  for (auto const& [cell, energy] : windows) {
154  if (energy > max) {
155  max = energy;
156  seedCell = cell;
157  }
158  }
159 
160  if (seedCell != nullptr) {
161  seedEta = seedCell->eta();
162  seedPhi = seedCell->phi();
163  // For the S-shape correction, we store the caloDDE of the hottest cell
164  caloDDE = (seedCell->caloDDE());
165  }
166  else {
167  return StatusCode::SUCCESS;
168  }
169 
170  bool cluster_in_barrel = true;
171  if (seedCell && seedCell->caloDDE()) cluster_in_barrel = seedCell->caloDDE()->is_lar_em_barrel();
172 
173  // Build desired quantities using the values previously obtained.
174 
175  double energy37Lay2 = 0.; // sum Et in 3*7 cell layer 2
176  double energy37Lay2NegPhi = 0.; // SRA for phi wrap-around
177  double energy37Lay2PosPhi = 0.; // SRA for phi wrap-around
178  double energy77Lay2 = 0.; // sum Et in 7*7 cell layer 2
179  double energy35Lay2 = 0.; // sum Et in 3*5 cell layer 2
180  double energy33Lay2 = 0.; // sum Et in 3*3 cell layer 2
181  double weightEta = 0.; // temporary variable for cluster width
182  double weightEta2 = 0.; // ditto
183  double clusterWidth35 = 0.; // cluster width in eta, 3*5 cell, layer 2
184  int nCellsEta = 0; // size of cell array in eta
185  int nCellsPhi = 0; // size of cell array in phi
186  double deta = 0.; // eta difference current cell - seed
187  double dphi = 0.; // phi difference current cell - seed
188 
189  double totalEnergy = 0;
191 
192  for (const LArCell* larcell : sel) {
193  double etaCell = larcell->eta();
194  double phiCell = larcell->phi();
195  double energyCell = larcell->energy();
196 
197  // Find position of current cell w.r.t. seed
198  deta = std::abs(etaCell - seedEta);
199  dphi = std::abs(phiCell - seedPhi);
200  if (dphi > M_PI) dphi = 2. * M_PI - dphi; // wrap (pi - 2pi)->(-pi - 0)
201 
202  double detaH = std::abs(etaCell - hotEta); // eta difference current cell - seed
203  double dphiH = std::abs(phiCell - hotPhi); // phi difference current cell - seed
204  if (dphiH > M_PI) dphiH = 2. * M_PI - dphiH; // wrap (pi - 2pi)->(-pi - 0)
205 
206  // Electromagnetic measurements, done by layer
207  // layer 2: Et weighted eta, phi, summed Et in 3*7,
208  // cluster width in 3*5 cell
209  // emEnergy in 3*7 cell is also calculated, all samplings
210 
211  double cellSizeEta;
212  double cellSizePhi;
213  if (emID->is_em_barrel(larcell->ID())) {
214  cellSizeEta = etaSizeLArEMSamp2(etaCell, Calorimeter::EMBAR);
215  cellSizePhi = phiSizeLArEMSamp2(etaCell, Calorimeter::EMBAR);
216  }
217  else {
218  cellSizeEta = etaSizeLArEMSamp2(etaCell, Calorimeter::EMEND);
219  cellSizePhi = phiSizeLArEMSamp2(etaCell, Calorimeter::EMEND);
220  }
221 
222  // hot cell stuff first
223  // 1. do the 7*7 stuff
224  nCellsEta = 7;
225  nCellsPhi = 7;
226  if (detaH <= 0.5 * double(nCellsEta - 1) * cellSizeEta + 0.01 &&
227  dphiH <= 0.5 * double(nCellsPhi - 1) * cellSizePhi + 0.01) {
228  // 2. do the 3*7 stuff
229  nCellsEta = 3;
230  nCellsPhi = 7;
231  if (detaH <= 0.5 * double(nCellsEta - 1) * cellSizeEta + 0.01 &&
232  dphiH <= 0.5 * double(nCellsPhi - 1) * cellSizePhi + 0.01) {
233  energy37Lay2 += energyCell;
234  } // End of do the 3*7 stuff
235  energy77Lay2 += energyCell;
236  } // End of do the 7*7 stuff
237 
238  // 1. do the 7*7 stuff
239  nCellsEta = 7;
240  nCellsPhi = 7;
241  if (deta <= 0.5 * double(nCellsEta - 1) * cellSizeEta + 0.005 &&
242  dphi <= 0.5 * double(nCellsPhi - 1) * cellSizePhi + 0.005) {
243  // 2. do the 3*7 stuff
244  nCellsEta = 3;
245  nCellsPhi = 7;
246  if (deta <= 0.5 * double(nCellsEta - 1) * cellSizeEta + 0.002 &&
247  dphi <= 0.5 * double(nCellsPhi - 1) * cellSizePhi + 0.002) {
248  energyEta += energyCell * etaCell;
249  // energy37Lay2 += energyCell;
250 
251  if (cluster_in_barrel) {
252  samp = larcell->caloDDE()->getSampling();
253  }
254 
255  if (phiCell > 0.) { // SRA phi wrap-around
256  energyPosPhi += double(energyCell) * phiCell;
257  energy37Lay2PosPhi += energyCell;
258  }
259  else {
260  energyNegPhi += double(energyCell) * phiCell;
261  energyNegPhiConv += double(energyCell) * (phiCell + 2.0 * M_PI);
262  energy37Lay2NegPhi += energyCell;
263  }
264  // 3. do the 3*5 stuff
265  nCellsEta = 3;
266  nCellsPhi = 5;
267  if (deta <= 0.5 * double(nCellsEta - 1) * cellSizeEta + 0.01 &&
268  dphi <= 0.5 * double(nCellsPhi - 1) * cellSizePhi + 0.01) {
269  weightEta += energyCell * etaCell;
270  weightEta2 += energyCell * etaCell * etaCell;
271  energy35Lay2 += energyCell;
272  // 3a. do the 3*3 stuff
273  // nCellsEta = 3;
274  nCellsPhi = 3;
275  if (deta <= 0.5 * double(nCellsEta - 1) * cellSizeEta + 0.01 &&
276  dphi <= 0.5 * double(nCellsPhi - 1) * cellSizePhi + 0.01) {
277  energy33Lay2 += energyCell;
278  } // End of do the 3*3 stuff
279  } // End of do the 3*5 stuff
280 
281  } // End of do the 3x7 stuff
282  // 4. do the 5*5 stuff
283  nCellsEta = 5;
284  nCellsPhi = 5;
285  if ((!cluster_in_barrel) && deta <= 0.5 * double(nCellsEta - 1) * cellSizeEta + 0.005 &&
286  dphi <= 0.5 * double(nCellsPhi - 1) * cellSizePhi + 0.005) {
287  samp = larcell->caloDDE()->getSampling();
288  } // End of do the 5*5 stuff
289  } // End of do the 7*7 stuff
290 
291  } // end of loop over sampling 2
292 
293  // Normalise energy weighted angles and calculate values
294 
295  if ((energy37Lay2PosPhi + energy37Lay2NegPhi) > 0.) { // dont divide by zero
296  // SRA phi wrap-around
297  double AvgNegPhi = 0.;
298  double AvgPosPhi = 0.;
299 
300  energyEta /= (energy37Lay2PosPhi + energy37Lay2NegPhi);
301 
302  if (energy37Lay2NegPhi > 0.) {
303  AvgNegPhi = energyNegPhi / energy37Lay2NegPhi;
304  }
305  else {
306  AvgNegPhi = -999.0;
307  }
308 
309  if (energy37Lay2PosPhi > 0.) {
310  AvgPosPhi = energyPosPhi / energy37Lay2PosPhi;
311  }
312  else {
313  AvgPosPhi = -999.0;
314  }
315 
316  if (AvgPosPhi == -999.0) {
317  if (AvgNegPhi != -999.0) {
318  energyPhi = AvgNegPhi;
319  }
320  }
321 
322  if (AvgNegPhi == -999.0) {
323  if (AvgPosPhi != -999.0) {
324  energyPhi = AvgPosPhi;
325  }
326  }
327 
328  if (AvgNegPhi != -999.0 && AvgPosPhi != -999.0) {
329  if ((AvgNegPhi > (-M_PI / 2.0)) && (AvgPosPhi < (M_PI / 2.0))) {
330  energyPhi = (energyNegPhi + energyPosPhi) / (energy37Lay2NegPhi + energy37Lay2PosPhi);
331  }
332  else {
333  if ((AvgNegPhi < (-M_PI / 2.0)) && (AvgPosPhi > (M_PI / 2.0))) {
334  energyPhi = (energyNegPhiConv + energyPosPhi) / (energy37Lay2NegPhi + energy37Lay2PosPhi);
335  if (energyPhi > M_PI) {
336  energyPhi = energyPhi - 2 * M_PI;
337  }
338  }
339  }
340  }
341  }
342  else {
343  energyEta = 99.;
344  energyPhi = 0.;
345  }
346 
347  for (const LArCell* larcell : sel) {
348  double etaCell = larcell->eta();
349  double phiCell = larcell->phi();
350  double energyCell = larcell->energy();
351 
352  // Find position of current cell w.r.t. seed
353  deta = std::abs(etaCell - energyEta);
354  dphi = std::abs(phiCell - energyPhi);
355  if (dphi > M_PI) dphi = 2. * M_PI - dphi; // wrap (pi - 2pi)->(-pi - 0)
356 
357  double cellSizeEta;
358  double cellSizePhi;
359  if (emID->is_em_barrel(larcell->ID())) {
360  cellSizeEta = etaSizeLArEMSamp2(etaCell, Calorimeter::EMBAR);
361  cellSizePhi = phiSizeLArEMSamp2(etaCell, Calorimeter::EMBAR);
362  }
363  else {
364  cellSizeEta = etaSizeLArEMSamp2(etaCell, Calorimeter::EMEND);
365  cellSizePhi = phiSizeLArEMSamp2(etaCell, Calorimeter::EMEND);
366  }
367 
368  nCellsEta = 3;
369  nCellsPhi = 7;
370  if ((deta <= 0.5 * double(nCellsEta - 1) * cellSizeEta + 0.01) &&
371  (dphi <= 0.5 * double(nCellsPhi - 1) * cellSizePhi + 0.01)) {
372 
373  if (cluster_in_barrel) {
374  totalEnergy += energyCell;
375  // samp = CaloSampling::getSampling(*larcell);
376  samp = larcell->caloDDE()->getSampling();
377  rtrigEmCluster.setEnergy(samp, rtrigEmCluster.energy(samp) + energyCell);
378  rtrigEmCluster.setRawEnergy(samp, rtrigEmCluster.rawEnergy(samp) + energyCell);
379  }
380  } // End of do the 3x7 stuff
381 
382  nCellsEta = 5;
383  nCellsPhi = 5;
384  if ((!cluster_in_barrel) && (deta <= 0.5 * double(nCellsEta - 1) * cellSizeEta + 0.01) &&
385  (dphi <= 0.5 * double(nCellsPhi - 1) * cellSizePhi + 0.01)) {
386  totalEnergy += energyCell;
387  // samp = CaloSampling::getSampling(*larcell);
388  samp = larcell->caloDDE()->getSampling();
389  rtrigEmCluster.setEnergy(samp, rtrigEmCluster.energy(samp) + energyCell);
390  rtrigEmCluster.setRawEnergy(samp, rtrigEmCluster.rawEnergy(samp) + energyCell);
391  } // End of do the 5*5 stuff
392 
393  } // end of loop over sampling 2
394 
395  // calculate cluster width
396 
397  if (energy35Lay2 > 0.) {
398  const double inv_energy35Lay2 = 1. / energy35Lay2;
399  clusterWidth35 = (weightEta2 * inv_energy35Lay2) -
400  (weightEta * inv_energy35Lay2) * (weightEta * inv_energy35Lay2);
401  clusterWidth35 > 0. ? clusterWidth35 = sqrt(clusterWidth35) : clusterWidth35 = 99.;
402  }
403  else {
404  clusterWidth35 = 99.;
405  }
406 
407  // Update cluster Variables
408 
409  // update stored variables
410  rtrigEmCluster.setE233(energy33Lay2);
411  rtrigEmCluster.setE237(energy37Lay2);
412  rtrigEmCluster.setE277(energy77Lay2);
413  rtrigEmCluster.setWeta2(clusterWidth35);
414  rtrigEmCluster.setRawEnergy(rtrigEmCluster.rawEnergy() + totalEnergy);
415  rtrigEmCluster.setEta(energyEta);
416  rtrigEmCluster.setPhi(energyPhi);
417  rtrigEmCluster.setRawEta(energyEta);
418  rtrigEmCluster.setRawPhi(energyPhi);
419  rtrigEmCluster.setNCells(ncells);
420 
421  return StatusCode::SUCCESS;
422 }
xAOD::TrigCaloCluster_v1::setRawPhi
void setRawPhi(float)
set Raw Phi (no calibration)
xAOD::TrigEMCluster_v1::setEta
void setEta(float)
set Eta (calibrated)
LArEM_ID.h
CaloCell::phi
virtual double phi() const override final
get phi (through CaloDetDescrElement)
Definition: CaloCell.h:359
xAOD::TrigEMCluster_v1::eta
float eta() const
get Eta (calibrated)
ReadCellNoiseFromCool.cell
cell
Definition: ReadCellNoiseFromCool.py:53
max
#define max(a, b)
Definition: cfImp.cxx:41
IReAlgToolCalo::m_larMgr
const CaloIdManager * m_larMgr
Calorimeter Id Manager for calorimeter part determination (Barrel versus EndCap)
Definition: IReAlgToolCalo.h:71
xAOD::TrigEMCluster_v1::setE237
void setE237(float)
set Energy in a 3x7 cluster (no calibration) around hottest cell
EgammaReSamp2Fex::phiSizeLArEMSamp2
double phiSizeLArEMSamp2(const double eta, const int calo) const
Definition: EgammaReSamp2Fex.h:69
xAOD::TrigEMCluster_v1::setEnergy
void setEnergy(float energy)
set Energy (calibrated)
CaloDetDescrElement
This class groups all DetDescr information related to a CaloCell. Provides a generic interface for al...
Definition: Calorimeter/CaloDetDescr/CaloDetDescr/CaloDetDescrElement.h:66
TrigEMCluster.h
CaloIdManager::getEM_ID
const LArEM_ID * getEM_ID(void) const
Definition: CaloIdManager.cxx:80
xAOD::TrigEMCluster_v1::setE233
void setE233(float)
set Energy in a 3x3 cluster (no calibration) around hottest cell
EgammaReSamp2Fex::EgammaReSamp2Fex
EgammaReSamp2Fex(const std::string &type, const std::string &name, const IInterface *parent)
Definition: EgammaReSamp2Fex.cxx:22
M_PI
#define M_PI
Definition: ActiveFraction.h:11
xAOD::TrigCaloCluster_v1::setRawEta
void setRawEta(float)
set Raw Eta (no calibration)
EgammaReSamp2Fex.h
EgammaReSamp2Fex::etaSizeLArEMSamp2
double etaSizeLArEMSamp2(const double eta, const int calo) const
Definition: EgammaReSamp2Fex.h:51
xAOD::TrigEMCluster_v1::setE277
void setE277(float)
set Energy in a 7x7 cluster (no calibration) around hottest cell
xAOD::TrigCaloCluster_v1::rawEnergy
float rawEnergy() const
get Raw Energy (no calibration)
Calorimeter::EMBAR
@ EMBAR
Definition: Calo_Def.h:21
xAOD::TrigEMCluster_v1::setPhi
void setPhi(float)
set Phi (calibrated)
IReAlgToolCalo::m_dataSvc
ServiceHandle< ITrigCaloDataAccessSvc > m_dataSvc
Definition: IReAlgToolCalo.h:84
xAOD::TrigEMCluster_v1::phi
float phi() const
get Phi (calibrated)
ParticleGun_FastCalo_ChargeFlip_Config.energy
energy
Definition: ParticleGun_FastCalo_ChargeFlip_Config.py:78
TTEM
@ TTEM
Definition: RegSelEnums.h:28
CaloSampling::CaloSample
CaloSample
Definition: Calorimeter/CaloGeoHelpers/CaloGeoHelpers/CaloSampling.h:22
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
IRoiDescriptor
Describes the API of the Region of Ineterest geometry.
Definition: IRoiDescriptor.h:23
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
CaloCell::caloDDE
const CaloDetDescrElement * caloDDE() const
get pointer to CaloDetDescrElement (data member)
Definition: CaloCell.h:305
sel
sel
Definition: SUSYToolsTester.cxx:92
EgammaReSamp2Fex::m_maxHotCellDphi
Gaudi::Property< float > m_maxHotCellDphi
Definition: EgammaReSamp2Fex.h:48
test_pyathena.parent
parent
Definition: test_pyathena.py:15
ATH_CHECK
#define ATH_CHECK
Definition: AthCheckMacros.h:40
xAOD::TrigEMCluster_v1::setWeta2
void setWeta2(float)
set cluster width (based on a 3x5 cluster - 2nd layer)
xAOD::double
double
Definition: CompositeParticle_v1.cxx:159
EgammaReSamp2Fex::m_maxHotCellDeta
Gaudi::Property< float > m_maxHotCellDeta
Definition: EgammaReSamp2Fex.h:47
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:195
LArCell
Data object for LAr calorimeter readout cell.
Definition: LArCell.h:53
LArTT_Selector< LArCellCont >
xAOD::TrigCaloCluster_v1::setNCells
void setNCells(int)
set number of cells used from RoI
IReAlgToolCalo
Base Class for Tools used for Egamma and Tau Feature Extraction Algorithms.
Definition: IReAlgToolCalo.h:37
CaloDetDescrElement::is_lar_em_barrel
bool is_lar_em_barrel() const
cell belongs to EM barrel
Definition: CaloDetDescrElement.cxx:98
LArEM_Base_ID::is_em_barrel
bool is_em_barrel(const Identifier id) const
test if the id belongs to the EM barrel
Calorimeter::EMEND
@ EMEND
Definition: Calo_Def.h:22
python.CaloScaleNoiseConfig.type
type
Definition: CaloScaleNoiseConfig.py:78
xAOD::TrigCaloCluster_v1::setRawEnergy
void setRawEnergy(float)
set Raw Energy (no calibration)
IRoiDescriptor.h
Calo_Def.h
LArEM_ID
Helper class for LArEM offline identifiers.
Definition: LArEM_ID.h:118
EgammaReSamp2Fex::execute
virtual StatusCode execute(xAOD::TrigEMCluster &rtrigEmCluster, const IRoiDescriptor &roi, const CaloDetDescrElement *&caloDDE, const EventContext &context) const override
execute feature extraction for the EM Calorimeter second layer
Definition: EgammaReSamp2Fex.cxx:28
xAOD::TrigEMCluster_v1
Description of a trigger EM cluster.
Definition: TrigEMCluster_v1.h:28
xAOD::TrigEMCluster_v1::energy
float energy() const
get Energy (calibrated)
CaloCell::eta
virtual double eta() const override final
get eta (through CaloDetDescrElement)
Definition: CaloCell.h:366
windows
int windows(float distance, float eta_pivot, int thr, int sector)
Definition: windows.cxx:14