ATLAS Offline Software
Loading...
Searching...
No Matches
jTowerBuilder.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5
8#include "jTowerBuilder.h"
10
11// TOWER IS A COLLECTION OF SUPER CELLS
12// IT SHOULD HAVE A UNIQUE ID
13// IT SHOULD BE ABLE TO RETURN LIST OF SUPER CELLS BELONGING TO IT
14
15// THIS IS A CLASS DESIGNED TO BUILD AN JTOWER USING THE JTOWER CLASS AND THEN PRINT THE RELEVANT INFORMATION TO THE SCREEN USING FUNCTION CALLS FROM THE JTOWER CLASS
16
17namespace LVL1 {
18
19jTowerBuilder::jTowerBuilder(const std::string& type,const std::string& name,const IInterface* parent): AthAlgTool(type,name,parent) {
20 declareInterface<IjTowerBuilder>(this);
21}
22
24{
25 ATH_CHECK( m_BDToolKey.initialize() );
26
27 return StatusCode::SUCCESS;
28}
29
30
31
32void jTowerBuilder::init(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) const {
33
34 execute(jTowerContainerRaw);
35
36 jTowerContainerRaw->clearContainerMap();
37 jTowerContainerRaw->fillContainerMap();
38
39}
40
41
43
44}
45
46
47void jTowerBuilder::execute(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) const {
48 BuildAllTowers(jTowerContainerRaw);
49}
50
51 // TOWER IDs FOR CLARITY
52 // Left Barrel IETower = 100000 + X
53 // Right Barrel IETower = 200000 + X
54 // Left Transition ID Tower = 300000 + X;
55 // Right Transition ID Tower = 400000 + X;
56 // Left Endcap ID Tower = 500000 + X
57 // Right Endcap ID Tower = 600000 + X
58 // Left Hadronic Endcap ID Tower = 11100000 + X --> These are just Layer 5 of Endcap Towers. They will never be generated as standalone jTowers.
59 // Right Haronic Endcap ID Tower = 22200000 + X --> These are just Layer 5 of Endcap Towers. They will never be generated as standalone jTowers.
60
61void jTowerBuilder::BuildEMBjTowers(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) const {
62 // Regions 0 only. Region 1 is 'transition region'.
63 for (int ieta = 0; ieta < 14; ++ieta) { // loop over 14 eta steps (ignoring last step as it is transition region)
64 float centre_eta = (0.1*ieta) + (0.05) ;
65 for (int iphi = 0; iphi < 64; ++iphi) { // loop over 64 phi steps
66 int key_eta = ieta;
67 float centre_phi = (m_TT_Size_phi*iphi) + (m_TT_Size_phi/2);
68 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 100000, -1, -1*centre_eta, centre_phi);
69 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 200000, 1, centre_eta, centre_phi);
70 }
71 }
72
73}
74
75void jTowerBuilder::BuildTRANSjTowers(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) const {
76 int TRANS_MODIFIER = 14;
77 int tmpVal = TRANS_MODIFIER;
78
79 for (int ieta = tmpVal; ieta < tmpVal + 1; ieta++) { // loop over eta steps
80 float centre_eta = (0.1*ieta) + (0.05);
81 for (int iphi = 0; iphi < 64; ++iphi) { // loop over 64 phi steps
82 int key_eta = ieta;
83 float centre_phi = (m_TT_Size_phi*iphi) + (m_TT_Size_phi/2);
84 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 300000, -1,-1*centre_eta, centre_phi);
85 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 400000, 1, centre_eta, centre_phi);
86 }
87 }
88
89}
90
91void jTowerBuilder::BuildEMEjTowers(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) const {
92 // Region 1
93 int EME_MODIFIER = 15;
94 int tmpVal = EME_MODIFIER;
95
96 for (int ieta = tmpVal; ieta < tmpVal + 3; ++ieta) { // loop over eta steps
97 float centre_eta =(0.1*ieta) + (0.05) ;
98 for (int iphi = 0; iphi < 64; ++iphi) { // loop over 64 phi steps
99 int key_eta = ieta;
100 float centre_phi = (m_TT_Size_phi*iphi) + (m_TT_Size_phi/2);
101 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 500000, -1, -1*centre_eta, centre_phi);
102 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 600000, 1, centre_eta, centre_phi);
103 }
104 EME_MODIFIER++;
105 }
106
107 // Region 2
108 tmpVal = EME_MODIFIER;
109 for (int ieta = tmpVal; ieta < tmpVal + 2; ++ieta) { // loop over eta steps
110 float centre_eta = (0.1*ieta) + (0.05);
111 for (int iphi = 0; iphi < 64; ++iphi) { // loop over 64 phi steps
112 int key_eta = ieta;
113 float centre_phi = (m_TT_Size_phi*iphi) + (m_TT_Size_phi/2);
114 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 500000, -1,-1*centre_eta, centre_phi);
115 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 600000, 1, centre_eta, centre_phi);
116 }
117 EME_MODIFIER++;
118 }
119
120 // Region 3
121 tmpVal = EME_MODIFIER;
122 for (int ieta = tmpVal; ieta < tmpVal + 4; ++ieta) { // loop over eta steps
123 float centre_eta= (0.1*ieta) + (0.05) ;
124 for (int iphi = 0; iphi < 64; ++iphi) { // loop over 64 phi steps
125 int key_eta = ieta;
126 float centre_phi = (m_TT_Size_phi*iphi) + (m_TT_Size_phi/2);
127 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 500000, -1,-1*centre_eta, centre_phi);
128 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 600000, 1, centre_eta, centre_phi);
129 }
130 EME_MODIFIER++;
131 }
132
133 // Region 4
134 tmpVal = EME_MODIFIER;
135 for (int ieta = tmpVal; ieta < tmpVal + 1; ++ieta) { // loop over eta steps
136 float centre_eta = (0.1*ieta) + (0.05);
137 for (int iphi = 0; iphi < 64; ++iphi) { // loop over 64 phi steps
138 int key_eta = ieta;
139 //float centre_phi =(TT_Size*iphi) + (0.5*TT_Size) ;
140 float centre_phi = (m_TT_Size_phi*iphi) + (m_TT_Size_phi/2);
141 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 500000, -1, -1*centre_eta, centre_phi);
142 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 600000, 1, centre_eta, centre_phi);
143 }
144 EME_MODIFIER++;
145 }
146
147
148}
149
150 // EMIE = Electromagnetic Inner ECAL - i.e. the forward ECAL region at high eta
151void jTowerBuilder::BuildEMIEjTowers(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) const {
152 int EMIE_MODIFIER = 25;
153 int tmpVal = EMIE_MODIFIER;
154 int cellCountEta = 0;
155
156 for (int ieta = tmpVal; ieta < tmpVal + 3; ++ieta) { // loop over eta steps (there are 3 here, 2.5-2.7, 2.7-2.9, 2.9-3.1)
157 cellCountEta++;
158 float centre_eta =(0.1*ieta) + (0.1*cellCountEta) ;
159 for (int iphi = 0; iphi < 32; ++iphi) { // loop over 32 phi steps
160 int key_eta = ieta;
161 float centre_phi = (2*m_TT_Size_phi*iphi) + m_TT_Size_phi;
162 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, /*7*/500000, -1, -1*centre_eta, centre_phi);
163 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, /*8*/600000, 1, centre_eta, centre_phi);
164 }
165 EMIE_MODIFIER++;
166 }
167
168 tmpVal = EMIE_MODIFIER;
169 for (int ieta = tmpVal; ieta < tmpVal + 1; ++ieta) { // loop over eta steps (there are 1 here, 3.1-3.2)
170 float centre_eta = (0.1*ieta + 0.3) + (0.05);
171 for (int iphi = 0; iphi < 32; ++iphi) { // loop over 32 phi steps
172 int key_eta = ieta;
173 float centre_phi = (2*m_TT_Size_phi*iphi) + m_TT_Size_phi;
174 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, /*7*/500000, -1, -1*centre_eta, centre_phi);
175 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, /*8*/600000, 1, centre_eta, centre_phi);
176 }
177 EMIE_MODIFIER++;
178 }
179
180}
181
182void jTowerBuilder::BuildFCALjTowers(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) const {
183 int FCAL_MODIFIER = 29; // there's 0.1 overlap with EMIE here in eta, but in counting we pretend it's the next one along.
184 int tmpVal = FCAL_MODIFIER;
185
186 //These jTowers split between all of the layers as well (FCAL0,1,2) but we treat them as though they are a single flat layer of 24 supercells and also pretend that they do not overlap - when they definitely do...
187 //This means that from a tower numbering perspective we start with FCAL0 and work upwards (in numbers?), but real ordering in eta is different and this has to be built into the jTower internal properties
188 //Right now we are unfortunately using hard-coded eta and phi positions to do this, but maybe these should be drawn from a database someday
189
190 // THIS REGION DEFINES 1 jTOWER PER SUPERCELL! AS SUCH THE jTOWER AND SUPERCELL POSITIONS IN ETA-PHI SPACE WILL MATCH EXACTY
191 // (That's good right? Supposed to make life easier?)
192
193 // 21/01/21 IN THE MC:
194 // FCAL 0 Region [NOT SPECIFIED IN MC] has 12 supercells in 3.2 < eta < 4.88, and 16 supercells in phi. Supercells are 0.14 x 0.4. posneg +-2
195 // FCAL 1 Region [NOT SPECIFIED IN MC] has 8 supercells in 3.2 < eta < 4.48, and 16 supercells in phi. Supercells are 0.16 x 0.4. posneg +-2
196 // FCAL 2 Region [NOT SPECIFIED IN MC] has 4 supercells in 3.2 < eta < 4.48, and 16 supercells in phi. Supercells are 0.32 x 0.4. posneg +-2
197
198 //FCAL0
199 float eta_width = 1.4;
200 int cellCountEta = 0;
201 int FCAL0_INITIAL = FCAL_MODIFIER;
202 std::vector<int> TT_etapos{31,33,34,36,37,39,40,42,43,45,46,48}; // Eta position of each supercell, need to be change for the real coords. Future MR
203 for (int ieta = tmpVal; ieta < tmpVal + 12; ++ieta) { // loop over eta steps (there are 12 here in varying positions for FCAL0)
204 int key_eta = ieta - FCAL0_INITIAL;
205 float centre_eta = (TT_etapos[cellCountEta]+eta_width/2)/10.0;
206 cellCountEta++;
207
208 for (int iphi = 0; iphi < 16; ++iphi) { // loop over 16 phi steps
209 float centre_phi = (2*m_TT_Size_phi_FCAL*iphi) + m_TT_Size_phi_FCAL;
210 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 700000, -1, -1*centre_eta, centre_phi, 0);
211 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 800000, 1, centre_eta, centre_phi, 0);
212 }
213 FCAL_MODIFIER++;
214 }
215
216 //FCAL1
217 eta_width = 1.6;
218 cellCountEta = 0;
219 tmpVal = FCAL_MODIFIER;
220 TT_etapos = {31,33,35,37,39,41,43,44};// Eta position of each supercell, need to be change for the real coords. Future MR
221 int FCAL1_INITIAL = FCAL_MODIFIER;
222 for (int ieta = tmpVal; ieta < tmpVal + 8; ++ieta) { // loop over eta steps (there are 8 here in varying positions for FCAL1)
223 int key_eta = ieta - FCAL1_INITIAL;
224 float centre_eta = (TT_etapos[cellCountEta]+eta_width/2)/10.0;
225 cellCountEta++;
226 for (int iphi = 0; iphi < 16; ++iphi) { // loop over 16 phi steps
227 float centre_phi = (2*m_TT_Size_phi_FCAL*iphi) + m_TT_Size_phi_FCAL;
228 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 900000, -1, -1*centre_eta, centre_phi, 1);
229 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 1000000, 1, centre_eta, centre_phi, 1);
230 }
231 FCAL_MODIFIER++;
232 }
233
234
235 //FCAL2
236 eta_width = 3.2;
237 cellCountEta = 0;
238 tmpVal = FCAL_MODIFIER;
239 TT_etapos = {31,34,37,41};// Eta position of each supercell, need to be change for the real coords. Future MR
240 int FCAL2_INITIAL = FCAL_MODIFIER;
241 for (int ieta = tmpVal; ieta < tmpVal + 4; ++ieta) { // loop over eta steps (there are 4 here in varying positions for FCAL2)
242 int key_eta = ieta - FCAL2_INITIAL;
243 float centre_eta = (TT_etapos[cellCountEta]+eta_width/2)/10.0;
244 cellCountEta++;
245 for (int iphi = 0; iphi < 16; ++iphi) { // loop over 16 phi steps
246 float centre_phi = (2*m_TT_Size_phi_FCAL*iphi) + m_TT_Size_phi_FCAL;
247 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 1100000, -1, -1*centre_eta, centre_phi, 2);
248 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 1200000, 1, centre_eta, centre_phi, 2);
249 }
250 FCAL_MODIFIER++;
251 }
252
253
254}
255
256
257
258 void jTowerBuilder::BuildHECjTowers(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) const {
259 // Region 0
260 int HEC_MODIFIER = 29;
261 int tmpVal = HEC_MODIFIER;
262 for (int ieta = tmpVal; ieta < tmpVal + 10; ++ieta){ // loop over eta steps
263 for (int iphi = 0; iphi < 64; ++iphi){ // loop over 64 phi steps
264 int key_eta = ieta;
265 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 11100000, -1, ieta, iphi);
266 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 22200000, 1, ieta, iphi);
267 }
268 HEC_MODIFIER++;
269 }
270
271 // Region 1
272 tmpVal = HEC_MODIFIER;
273 for (int ieta = tmpVal; ieta < tmpVal + 4; ++ieta){ // loop over eta steps
274 for (int iphi = 0; iphi < 32; ++iphi){ // loop over 64 phi steps
275 int key_eta = ieta;
276 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 11100000, -1, ieta, iphi);
277 BuildSingleTower(jTowerContainerRaw, ieta, iphi, key_eta, 22200000, 1, ieta, iphi);
278 }
279 HEC_MODIFIER++;
280 }
281
282}
283//=================================================================================================================================================================
284
285void jTowerBuilder::BuildSingleTower(std::unique_ptr<jTowerContainer> & jTowerContainerRaw,float eta, float phi, int key_eta, float keybase, int posneg, float centre_eta, float centre_phi, int fcal_layer) const {
286 int towerID = keybase + phi + (64 * key_eta);
287 jTowerContainerRaw->push_back(eta, phi, towerID, posneg, centre_eta, centre_phi, fcal_layer);
288
289}
290
291
292
293StatusCode jTowerBuilder::AssignPileupAndNoiseValues(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) const{
294
296 if (!myDBTool.isValid()){
297 ATH_MSG_ERROR("Not able to read " << m_BDToolKey );
298 return StatusCode::FAILURE;
299 }
300
301 for(LVL1::jTower* jtower : *jTowerContainerRaw) {
302
303 auto [CutJetEM, CutJetHad, CutMetEM, CutMetHad] = myDBTool->get_NoiseCuts( jtower->OnlineID() );
304 auto [PileUpWeightEM, PileUpWeightHad, InverseWeightEM, InverseWeightHad] = myDBTool->get_PileUpValues( jtower->OnlineID() );
305
306 //Simulation used MeV not counts, those
307 int LSBscale_EM = 25; // cf LATOME
308 int LSBscale_HAD = 25; // cf LATOME
309
310 //TREX uses another conversion factor
311 if(std::abs(jtower->centreEta()) < 1.5){
312 LSBscale_HAD = 500;// cf TREX
313 }
314
315 //Since the COOL DB for FCAL individual towers are sharing the same OnlideID ( to save space)
316 //but in reality they are different towers. we need to set the parameters to 0
317 if(jtower->OfflineID() >= FEXAlgoSpaceDefs::jFEX_FCAL2_start){
318 PileUpWeightEM = 0;
319 InverseWeightEM = 0;
320 }
321 else if(jtower->OfflineID() >= FEXAlgoSpaceDefs::jFEX_FCAL1_start){
322 PileUpWeightHad = 0;
323 InverseWeightHad = 0;
324 }
325
326 jtower->setTTowerArea(PileUpWeightEM,0);
327 jtower->setTTowerArea(PileUpWeightHad,1);
328
329 jtower->setTTowerAreaInv(InverseWeightEM,0);
330 jtower->setTTowerAreaInv(InverseWeightHad,1);
331
332 jtower->setNoiseForMet(CutMetEM*LSBscale_EM,0);
333 jtower->setNoiseForMet(CutMetHad*LSBscale_HAD,1);
334
335 jtower->setNoiseForJet(CutJetEM*LSBscale_EM,0);
336 jtower->setNoiseForJet(CutJetHad*LSBscale_HAD,1);
337
338 }
339 return StatusCode::SUCCESS;
340}
341
342
343void jTowerBuilder::BuildAllTowers(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) const {
344 BuildEMBjTowers(jTowerContainerRaw);
345 BuildTRANSjTowers(jTowerContainerRaw);
346 BuildEMEjTowers(jTowerContainerRaw);
347 BuildEMIEjTowers(jTowerContainerRaw);
348 BuildFCALjTowers(jTowerContainerRaw);
349}
350
351} // end of LVL1 namespace
352
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)
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
static constexpr int jFEX_FCAL1_start
static constexpr int jFEX_FCAL2_start
void BuildTRANSjTowers(std::unique_ptr< jTowerContainer > &jTowerContainerRaw) const
void BuildAllTowers(std::unique_ptr< jTowerContainer > &jTowerContainerRaw) const
void BuildEMIEjTowers(std::unique_ptr< jTowerContainer > &jTowerContainerRaw) const
void BuildHECjTowers(std::unique_ptr< jTowerContainer > &jTowerContainerRaw) const
SG::ReadCondHandleKey< jFEXDBCondData > m_BDToolKey
jTowerBuilder(const std::string &type, const std::string &name, const IInterface *parent)
static constexpr float m_TT_Size_phi_FCAL
virtual void execute(std::unique_ptr< jTowerContainer > &jTowerContainerRaw) const override
virtual void init(std::unique_ptr< jTowerContainer > &jTowerContainerRaw) const override
void BuildFCALjTowers(std::unique_ptr< jTowerContainer > &jTowerContainerRaw) const
virtual StatusCode initialize() override
void BuildEMBjTowers(std::unique_ptr< jTowerContainer > &jTowerContainerRaw) const
static constexpr float m_TT_Size_phi
virtual StatusCode AssignPileupAndNoiseValues(std::unique_ptr< jTowerContainer > &jTowerContainerRaw) const override
void BuildSingleTower(std::unique_ptr< jTowerContainer > &jTowerContainerRaw, float eta, float phi, int key_eta, float keybase, int posneg, float centre_eta=0.0, float centre_phi=0.0, int fcal_layer=-1) const
virtual void reset() const override
void BuildEMEjTowers(std::unique_ptr< jTowerContainer > &jTowerContainerRaw) const
The jTower class is an interface object for jFEX trigger algorithms The purposes are twofold:
Definition jTower.h:36
eFexTowerBuilder creates xAOD::eFexTowerContainer from supercells (LATOME) and triggerTowers (TREX) i...