ATLAS Offline Software
Loading...
Searching...
No Matches
GTowerMappingDataCondAlg.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3*/
4
6#include "TMath.h"
7#include "TVector2.h"
8#include "CaloDetDescr/CaloDetDescrElement.h"
9#include "GTowerHelpers.h"
10
11namespace LVL1
12{
13 GTowerMappingDataCondAlg::GTowerMappingDataCondAlg(const std::string &name, ISvcLocator *pSvcLocator)
15 {
16 }
17
19
21 {
23 ATH_CHECK(detStore()->retrieve(m_gtowerID));
24 return StatusCode::SUCCESS;
25 }
26
29 const CaloSuperCellDetDescrManager *scDetMgr) const
30 {
31 // Define GTower geometry
32 std::array<float, 5> etaBounds{3.1, 3.5, 4.0, 4.45, 4.9};
33 constexpr std::size_t nEtaTowers = etaBounds.size() - 1;
34 constexpr std::size_t nPhiTowers = 16;
35 float dPhi = 2 * TMath::Pi() / nPhiTowers;
36 for (std::size_t iEta = 0; iEta < etaBounds.size() - 1; ++iEta)
37 {
38 float eta = (etaBounds[iEta] + etaBounds[iEta + 1]) / 2;
39 float dEta = etaBounds[iEta + 1] - etaBounds[iEta];
40
41 for (int sign : {1, -1})
42 {
43 for (std::size_t iPhi = 0; iPhi < nPhiTowers; ++iPhi)
44 {
45 std::size_t index = data.size();
46 float phi = dPhi * (iPhi + 0.5);
47 if (phi > TMath::Pi())
48 phi -= 2 * TMath::Pi();
49 JGTowerHelper helper(eta * sign, dEta, phi, dPhi);
50 helper.SetSampling(2);
51 for (std::size_t scIdx = 0; scIdx < m_scid->calo_cell_hash_max(); ++scIdx)
52 {
53 Identifier scid = m_scid->cell_id(scIdx);
54 if (m_scid->is_tile(scid) && m_scid->sampling(scid) != 2)
55 continue;
56 const CaloDetDescrElement *dde = scDetMgr->get_element(scid);
57 if (!dde)
58 {
59 ATH_MSG_ERROR("Failed to load CaloDetDescrElement");
60 return StatusCode::FAILURE;
61 }
62 if (std::abs(dde->eta_raw()) < 3.2)
63 continue;
64 if (helper.inBox(dde->eta_raw(), dde->phi_raw()) && sign * m_scid->pos_neg(scid) > 0)
65 helper.SetSCIndices(scIdx);
66 }
67 data.push_back(std::move(helper));
68 if (iEta > 0)
69 {
70 if (sign > 0)
71 helper.setPreviousEtaIndex(index - 2 * nPhiTowers);
72 else
73 helper.setNextEtaIndex(index - 2 * nPhiTowers);
74 }
75 if (iEta < nEtaTowers - 1)
76 {
77 if (sign > 0)
78 helper.setNextEtaIndex(index - 2 * nPhiTowers);
79 else
80 helper.setPreviousEtaIndex(index - 2 * nPhiTowers);
81 }
82 helper.setPreviousPhiIndex(iPhi == 0 ? index + nPhiTowers - 1 : index - 1);
83 helper.setNextPhiIndex(iPhi == nPhiTowers - 1 ? index - (nPhiTowers - 1) : index + 1);
84 }
85 }
86 }
87 return StatusCode::SUCCESS;
88 }
89
91 {
92 for (JGTowerHelper &helper : data)
93 helper.setArea(gFEX::towerArea(helper.Eta()));
94 return StatusCode::SUCCESS;
95 }
96} // namespace LVL1
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
char data[hepevt_bytes_allocation_ATLAS]
Definition HepEvt.cxx:11
int sign(int a)
const ServiceHandle< StoreGateSvc > & detStore() const
This class groups all DetDescr information related to a CaloCell.
const CaloDetDescrElement * get_element(const Identifier &cellId) const
get element by its identifier
virtual StatusCode buildForwardMapping(JGTowerMappingData &data, const CaloSuperCellDetDescrManager *mgr) const override
GTowerMappingDataCondAlg(const std::string &name, ISvcLocator *pSvcLocator)
virtual StatusCode loadTowerAreas(JGTowerMappingData &data) const override
virtual StatusCode initialize() override
JGTowerMappingDataCondAlgBase(const std::string &name, ISvcLocator *pSvcLocator)
float towerArea(float eta)
Get the GCaloTower areas from their eta bins.
eFexTowerBuilder creates xAOD::eFexTowerContainer from supercells (LATOME) and triggerTowers (TREX) i...
std::vector< JGTowerHelper > JGTowerMappingData
Definition index.py:1