ATLAS Offline Software
Loading...
Searching...
No Matches
eFEXegAlgo.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//***************************************************************************
6// eFEXegAlgo - description
7// -------------------
8// begin : 24 02 2020
9// email : antonio.jacques.costa@cern.ch ulla.blumenschein@cern.ch tong.qiu@cern.ch
10// ***************************************************************************/
11#include <vector>
12#include <mutex>
13
17#include "L1CaloFEXSim/eTower.h"
18
19
20namespace LVL1 {
21
22 // default constructor for persistency
23 eFEXegAlgo::eFEXegAlgo(const std::string& type, const std::string& name, const IInterface* parent):
24 AthAlgTool(type, name, parent)
25 {
26 declareInterface<eFEXegAlgo>(this);
27 }
28
33
35
36 ATH_CHECK(m_eTowerContainerKey.initialize());
38
39 return StatusCode::SUCCESS;
40 }
41
42
43 StatusCode eFEXegAlgo::safetyTest() const {
44
46 if(!eTowerContainer.isValid()){
47 ATH_MSG_FATAL("Could not retrieve container " << m_eTowerContainerKey.key() );
48 return StatusCode::FAILURE;
49 }
50
51 return StatusCode::SUCCESS;
52 }
53
54 void eFEXegAlgo::setup(int inputTable[3][3], int efex_id, int fpga_id, int central_eta) {
55
56 std::copy(&inputTable[0][0], &inputTable[0][0] + 9, &m_eFEXegAlgoTowerID[0][0]);
57
58 m_efexid = efex_id;
59 m_fpgaid = fpga_id;
60 m_central_eta = central_eta;
61
62 setSeed();
63
64 }
65
67
69
71 et = tmpTower->getLayerTotalET(0) + tmpTower->getLayerTotalET(1) + tmpTower->getLayerTotalET(2) + tmpTower->getLayerTotalET(3);
72 }
73
81
88
96
97 void eFEXegAlgo::getReta(std::vector<unsigned int> & retavec) {
98
99 unsigned int coresum = 0; // 3x2 L2 sum : core
100 unsigned int totalsum = 0; // 7x3 L2 sum : total
101 unsigned int envsum = 0; // total - core : env
102
103 retavec.clear(); // clear the output vector before starting
104
105 // window limits
106 int iTotalStart = m_seedID-3;
107 int iTotalEnd = m_seedID+3;
108 int iCoreStart = m_seedID-1;
109 int iCoreEnd = m_seedID+1;
110 int phiStart = -999;
111 int phiEnd = -99;
112 if (m_seed_UnD) {
113 phiStart = 1;
114 phiEnd = 2;
115 } else {
116 phiStart = 0;
117 phiEnd = 1;
118 }
119
120 // 3x2 and 7x3 L2 sum
121 for (int i=iTotalStart; i<=iTotalEnd; ++i) { // eta
122 for(int j=0; j<=2; ++j) { // phi
123 if (i>=iCoreStart && i <= iCoreEnd && j>=phiStart && j<=phiEnd) {
124 unsigned int tmp_et; getWindowET(2,j,i,tmp_et);
125 coresum += tmp_et;
126 }
127
128 unsigned int tmptot_et; getWindowET(2,j,i,tmptot_et);
129 totalsum += tmptot_et;
130 }
131 }
132
133 // get environment
134 envsum = totalsum - coresum;
135
136 // Overflow handling
137 if (coresum > 0xffff) coresum = 0xffff;
138 if (envsum > 0xffff) envsum = 0xffff;
139
140 // Return results
141 retavec.push_back(coresum);
142 retavec.push_back(envsum);
143
144 }
145
146 void eFEXegAlgo::getRhad(std::vector<unsigned int> & rhadvec) {
147
148 unsigned int hadsum = 0; // 3x3 Towers Had
149 unsigned int emsum = 0; // (1x3 + 3x3 + 3x3 + 1x3) SCs EM
150
151 rhadvec.clear(); // clear the output vector before starting
152
153 int iCoreStart = m_seedID-1;
154 int iCoreEnd = m_seedID+1;
155
157
158 if(m_algoVersion==0) {
159 // 3x3 Towers Had ; 1x3 L0 + 1x3 L3 EM
160 for (int i=0; i<3; ++i) { // phi
161 for (int j=0; j<=2; ++j) { // eta
162 if (((m_efexid%3 == 0) && (m_fpgaid == 0) && (m_central_eta == 0) && (j == 0)) || ((m_efexid%3 == 2) && (m_fpgaid == 3) && (m_central_eta == 5) && (j == 2))) {
163 continue;
164 } else {
165 const eTower * tTower = eTowerContainer->findTower(m_eFEXegAlgoTowerID[i][j]);
166 hadsum += tTower->getLayerTotalET(4);
167 if (j==1) {
168 emsum += ( tTower->getLayerTotalET(0) + tTower->getLayerTotalET(3) );
169 }
170 }
171 }
172 }
173
174 // 3x3 SCs L1 and L2 sum
175 for (int i=iCoreStart; i<=iCoreEnd; ++i) { // eta
176 for(int j=0; j<=2; ++j) { // phi
177 unsigned int tmp_et_a, tmp_et_b;
178 getWindowET(1,j,i,tmp_et_a);
179 getWindowET(2,j,i,tmp_et_b);
180 emsum += ( tmp_et_a + tmp_et_b );
181 }
182 }
183 } else {
184
185 // see slide 4 of https://indico.cern.ch/event/1513502/contributions/6389265/attachments/3019474/5326755/Rate%20Reduction%20Methods%20Overview.pdf
186 // Update 29-Apr-2025: removing 4 corners of EM3 from sum
187
188 // EM cluster depends on UnD flag
189 int phi2 = (m_seed_UnD > 0 ? 2 : 0);
190
191 // 3x3 Towers Had + EM3 (excluding corners in EM3); 1x2 L0
192 for (int i = 0; i < 3; ++i) { // phi
193 for (int j = 0; j <= 2; ++j) { // eta
194 if (((m_efexid % 3 == 0) && (m_fpgaid == 0) && (m_central_eta == 0) && (j == 0)) ||
195 ((m_efexid % 3 == 2) && (m_fpgaid == 3) && (m_central_eta == 5) && (j == 2))) {
196 continue;
197 } else {
198 const eTower *tTower = eTowerContainer->findTower(m_eFEXegAlgoTowerID[i][j]);
199 hadsum += tTower->getLayerTotalET(4) + (((i==0&&j==0)||(i==2&&j==2)||(i==0&&j==2)||(i==2&&j==0)) ? 0 : tTower->getLayerTotalET(3));
200 // For PS add central tower + UnD phi neighbour
201 if (j == 1 && (i == 1 || i == phi2)) {
202 emsum += (tTower->getLayerTotalET(0));
203 }
204 }
205 }
206 }
207
208 // 3x2 SCs L1 and L2 sum, so phi range depends on UnD flag
209 int phiStart = 0;
210 int phiEnd = 1;
211 if (m_seed_UnD > 0) {
212 phiStart = 1;
213 phiEnd = 2;
214 }
215 for (int i = iCoreStart; i <= iCoreEnd; ++i) { // eta
216 for (int j = phiStart; j <= phiEnd; ++j) { // phi
217 unsigned int tmp_et_a, tmp_et_b;
218 getWindowET(1, j, i, tmp_et_a);
219 getWindowET(2, j, i, tmp_et_b);
220 emsum += (tmp_et_a + tmp_et_b);
221 }
222 }
223 }
224 // Overflow handling
225 if (emsum > 0xffff) emsum = 0xffff;
226 if (hadsum > 0xffff) hadsum = 0xffff;
227
228 // Return results
229 rhadvec.push_back(emsum);
230 rhadvec.push_back(hadsum);
231
232 }
233
234 void LVL1::eFEXegAlgo::getWstot(std::vector<unsigned int> & output){
235 unsigned int numer = 0;
236 unsigned int den = 0;
237
238 output.clear(); // clear the output vector before starting
239
240 int iStart = m_seedID - 2;
241 int iEnd = m_seedID + 2;
242
243 for (int i = iStart; i <= iEnd; ++i) { // eta
244 int diff = i - m_seedID;
245 unsigned int weight = diff*diff;
246 for (int j = 0; j <= 2; ++j) { // phi
247 unsigned int eT;
248 getWindowET(1, j, i, eT);
249 // NB need to be careful as numer and den are defined such that wstot=numer/den,
250 // but in the firmware (and therefore this bitwise code) we actually
251 // check that den/numer < Threshold
252 numer += eT*weight;
253 den += eT;
254 }
255 }
256
257 // Overflow handling
258 //if (den > 0xffff) den = 0xffff; - commented out so that denom can overflow, will then automatically pass all thresholds (see eFEXFPGA::SetIsoWP)
259 if (numer > 0xffff) numer = 0xffff;
260
261 // Return results
262 output.push_back(den);
263 output.push_back(numer);
264
265 }
266
269 void LVL1::eFEXegAlgo::getClusterCells(std::vector<unsigned int> &cellETs) {
270
271 int phiUpDownID = 0;
272 if (m_seed_UnD) phiUpDownID = 2;
273
274 // Initialise results vector
275 cellETs.resize(16,0);
276 // Fill vector with 2 PS cells, 6 L1, 6 L2, 2 L3
277 // Presampler
278 getWindowET(0, 1, 0, cellETs[0]);
279 getWindowET(0, phiUpDownID, 0, cellETs[1]);
280 // central phi Layer 1
281 getWindowET(1, 1, m_seedID, cellETs[2]);
282 getWindowET(1, 1, m_seedID - 1, cellETs[3]);
283 getWindowET(1, 1, m_seedID + 1, cellETs[4]);
284 // top/bottom phi Layer 1
285 getWindowET(1, phiUpDownID, m_seedID, cellETs[5]);
286 getWindowET(1, phiUpDownID, m_seedID - 1, cellETs[6]);
287 getWindowET(1, phiUpDownID, m_seedID + 1, cellETs[7]);
288 // central phi Layer 2
289 getWindowET(2, 1, m_seedID, cellETs[8]);
290 getWindowET(2, 1, m_seedID - 1, cellETs[9]);
291 getWindowET(2, 1, m_seedID + 1, cellETs[10]);
292 // top/bottom phi Layer 2
293 getWindowET(2, phiUpDownID, m_seedID, cellETs[11]);
294 getWindowET(2, phiUpDownID, m_seedID - 1, cellETs[12]);
295 getWindowET(2, phiUpDownID, m_seedID + 1, cellETs[13]);
296 // Layer 3
297 getWindowET(3, 1, 0, cellETs[14]);
298 getWindowET(3, phiUpDownID, 0, cellETs[15]);
299
300 return;
301
302 }
303
304 unsigned int LVL1::eFEXegAlgo::getET() {
305
307 std::vector<unsigned int> clusterCells;
308 getClusterCells(clusterCells);
309
310
312 unsigned int PS_ET = 0;
313
314 if(m_algoVersion==0) {
315 PS_ET = dmCorrection(clusterCells[0], 0)
316 + dmCorrection(clusterCells[1], 0);
317 } else {
318 // 2025 algoVersion only uses 1 PS scell, except at most extreme eta values
319 PS_ET = dmCorrection(clusterCells[0], 0);
320 if ( ((m_efexid%3) == 0 && m_fpgaid == 0) || ((m_efexid%3) == 2 && m_fpgaid == 3)) {
321 PS_ET += dmCorrection(clusterCells[1], 0);
322 }
323 }
324 unsigned int L1_ET = dmCorrection(clusterCells[2], 1)
325 + dmCorrection(clusterCells[3], 1)
326 + dmCorrection(clusterCells[4], 1)
327 + dmCorrection(clusterCells[5], 1)
328 + dmCorrection(clusterCells[6], 1)
329 + dmCorrection(clusterCells[7], 1);
330 unsigned int L2_ET = dmCorrection(clusterCells[8], 2)
331 + dmCorrection(clusterCells[9], 2)
332 + dmCorrection(clusterCells[10], 2)
333 + dmCorrection(clusterCells[11], 2)
334 + dmCorrection(clusterCells[12], 2)
335 + dmCorrection(clusterCells[13], 2);
336 unsigned int L3_ET = clusterCells[14] + clusterCells[15];
337
339 unsigned int totET = PS_ET + L1_ET + L2_ET + L3_ET;
340
341 // overflow handling
342 if (totET > 0xffff) totET = 0xffff;
343
344 return totET;
345
346 }
347
348
350 MsgStream& msg)
351 {
352 for (const auto& p : dmCorrections) {
353 if (p.first < 25 || p.first >= 50) continue;
354 m_corrections[0][p.first - 25] = p.second["EmPS"].data<int>();
355 m_corrections[1][p.first - 25] = p.second["EmFR"].data<int>();
356 m_corrections[2][p.first - 25] = p.second["EmMD"].data<int>();
357 if (msg.level() <= MSG::DEBUG) {
358 msg << MSG::DEBUG << "DM Correction for etaIdx=" << (p.first - 25) << " : [" << m_corrections[0][p.first - 25] << ","
359 << m_corrections[1][p.first - 25] << "," << m_corrections[2][p.first - 25] << "]"
360 << endmsg;
361 }
362 }
363 }
364
365 unsigned int LVL1::eFEXegAlgo::dmCorrection (unsigned int ET, unsigned int layer) {
367 if ( !m_dmCorr || layer > 2 ) return ET;
368
371 int efexEta = m_efexid%3;
372 int ieta = 0;
373 if (efexEta == 2) { // Rightmost eFEX
374 // m_central_eta has range 1-4 or 1-5
375 ieta = 8 + m_fpgaid*4 + m_central_eta - 1;
376 }
377 else if (efexEta == 1 && m_fpgaid > 1) { // central eFEX, eta > 0
378 // m_central_eta has range 1-4
379 ieta = (m_fpgaid-2)*4 + m_central_eta - 1;
380 }
381 else if (efexEta == 1) { // central eFEX, eta < 0
382 // m_central_eta had range 1-4
383 ieta = (1-m_fpgaid)*4 + (4-m_central_eta);
384 }
385 else { // Leftmost eFEX
386 // m_central_eta has range 0-4 or 1-4
387 ieta = 8 + 4*(3-m_fpgaid) + (4-m_central_eta);
388 }
389
390 static std::once_flag flag;
391 std::call_once(flag, [&]() {
392
393 if (!m_dmCorrectionsKey.empty()) {
394 // replace m_corrections values with values from database ... only try this once
395 SG::ReadCondHandle <CondAttrListCollection> dmCorrections{m_dmCorrectionsKey/*, ctx*/ };
396 if (dmCorrections.isValid()) {
397 if(dmCorrections->size()==0 && Gaudi::Hive::currentContext().eventID().time_stamp()>1672527600) { // not an error for data before 2023 (will include MC21 and MC23a)
398 ATH_MSG_ERROR("No dead material corrections found in conditions database for this event in folder " << m_dmCorrectionsKey.key());
399 throw std::runtime_error("No dead material corrections found in database for this event");
400 }
401 m_corrections = Corrections (*dmCorrections.cptr(), msg());
402 }
403 ATH_MSG_INFO("Loaded DM Corrections from database");
404 }
405 });
406
408 unsigned int factor = m_corrections.corr(layer, ieta);
409
414
415 unsigned int correction = ET;
416 for (int bit = 6; bit >= 0; bit--) {
417 correction /= 2;
418 if (factor & (1<<bit))
419 ET += correction;
420 }
422 return ET;
423 }
424
425
426 std::unique_ptr<eFEXegTOB> LVL1::eFEXegAlgo::geteFEXegTOB() {
427
428 std::unique_ptr<eFEXegTOB> out = std::make_unique<eFEXegTOB>();
429 out->setET(getET());
430
431 std::vector<unsigned int> temvector;
432 getWstot(temvector);
433 // For wstot, num and den seem switched around, but this is because the 'numerator' and 'denominator'
434 // mean different things at different points in the processing chain
435 // When the threshold comparison is done in the SetIsoWP function, we actually check Den/Num
436 out->setWstotNum(temvector[1]);
437 out->setWstotDen(temvector[0]);
438 getRhad(temvector);
439 out->setRhadEM(temvector[0]);
440 out->setRhadHad(temvector[1]);
441 getReta(temvector);
442 out->setRetaCore(temvector[0]);
443 out->setRetaEnv(temvector[1]);
444 out->setSeedUnD(m_seed_UnD);
445 out->setSeed(m_seedID);
446 return out;
447 }
448
449 void LVL1::eFEXegAlgo::getWindowET(int layer, int jPhi, int SCID, unsigned int & outET) {
450
452
453 if (SCID<0) { // left towers in eta
454 if ((m_efexid%3 == 0) && (m_fpgaid == 0) && (m_central_eta == 0)) {
455 outET = 0;
456 } else {
457 int etaID = 4+SCID;
458 const eTower * tmpTower = eTowerContainer->findTower(m_eFEXegAlgoTowerID[jPhi][0]);
459 if (layer==1 || layer==2) {
460 outET = tmpTower->getET(layer,etaID);
461 } else if (layer==0 || layer==3 || layer==4) {
462 outET = tmpTower->getLayerTotalET(layer);
463 }
464 }
465 } else if (SCID>=0 && SCID<4) { // central towers in eta
466 const eTower * tmpTower = eTowerContainer->findTower(m_eFEXegAlgoTowerID[jPhi][1]);
467 if (layer==1 || layer==2) {
468 outET = tmpTower->getET(layer,SCID);
469 } else if (layer==0 || layer==3 || layer==4) {
470 outET = tmpTower->getLayerTotalET(layer);
471 }
472 } else if (SCID>=4){ // right towers in eta
473 if ((m_efexid%3 == 2) && (m_fpgaid == 3) && (m_central_eta == 5)) {
474 outET = 0;
475 } else {
476 int etaID = SCID-4;
477 const eTower * tmpTower = eTowerContainer->findTower(m_eFEXegAlgoTowerID[jPhi][2]);
478 if (layer==1 || layer==2) {
479 outET = tmpTower->getET(layer,etaID);
480 } else if (layer==0 || layer==3 || layer==4) {
481 outET = tmpTower->getLayerTotalET(layer);
482 }
483 }
484 }
485
486 // overflow handling
487 if (outET > 0xffff) outET = 0xffff;
488
489 }
490
491
492// Utility function to calculate and return jet discriminant sums for specified location
493// Intended to allow xAOD TOBs to be decorated with this information
494 void eFEXegAlgo::getSums(unsigned int seed, bool UnD,
495 std::vector<unsigned int> & RetaSums,
496 std::vector<unsigned int> & RhadSums,
497 std::vector<unsigned int> & WstotSums)
498 {
499 // Set seed parameters to supplied values
500 m_seed_UnD = UnD;
501 m_seedID = seed;
502 m_hasSeed = true;
503
504 // Now just call the 3 discriminant calculation methods
505 getReta(RetaSums);
506 getRhad(RhadSums);
507 getWstot(WstotSums);
508
509 }
510
511
512// Find seed and UnD flag
514
515 m_hasSeed = false;
516 m_seed_UnD = false;
517 unsigned int tmpID = 999;
518 unsigned int maxET = 0;
519
520 for (int i=0; i<4 ; ++i) {
521 int iSeedL = i-1;
522 int iSeedR = i+1;
523
524 // eta ID of candidate seed
525 unsigned int cETUp;
526 getWindowET(2,2,i,cETUp);
527 unsigned int iSeedET;
528 getWindowET(2,1,i, iSeedET);
529 unsigned int cETDown;
530 getWindowET(2,0,i, cETDown);
531
532 // left of candidate seed
533 unsigned int lETUp;
534 getWindowET(2,2,iSeedL,lETUp);
535 unsigned int lET;
536 getWindowET(2,1,iSeedL,lET);
537 unsigned int lETDown;
538 getWindowET(2,0,iSeedL,lETDown);
539
540 // right of candidate seed
541 unsigned int rETUp;
542 getWindowET(2,2,iSeedR,rETUp);
543 unsigned int rET;
544 getWindowET(2,1,iSeedR,rET);
545 unsigned int rETDown;
546 getWindowET(2,0,iSeedR,rETDown);
547
548 // greater or equal than for left and down cells, greater than for right and up ones
549 if (iSeedET>=lET && iSeedET>rET
550 && iSeedET>=lETUp && iSeedET>cETUp && iSeedET>rETUp
551 && iSeedET>=lETDown && iSeedET>=cETDown && iSeedET>rETDown) {
552 if (iSeedET>=maxET) { // if two maxima exist and have the same ET, keep the one to the right
553 maxET = iSeedET;
554 tmpID = i;
555 }
556 }
557 }
558
559 if(tmpID!=999) {
560 m_seedID = tmpID;
561 m_hasSeed = true;
562 unsigned int tmp_et_up, tmp_et_down;
563 getWindowET(2,2,m_seedID,tmp_et_up);
564 getWindowET(2,0,m_seedID,tmp_et_down);
565 if (tmp_et_up >= tmp_et_down) {
566 m_seed_UnD = true; // go up if energy greater or equal to bottom
567 }
568 }
569 }
570
571} // namespace LVL1
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define endmsg
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
void diff(const Jet &rJet1, const Jet &rJet2, std::map< std::string, double > varDiff)
Difference between jets - Non-Class function required by trigger.
Definition Jet.cxx:631
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
MsgStream & msg() const
This class is a collection of AttributeLists where each one is associated with a channel number.
virtual ~eFEXegAlgo()
Destructor.
virtual unsigned int getET()
virtual void getReta(std::vector< unsigned int > &)
virtual void getRealPhi(float &phi)
int m_eFEXegAlgoTowerID[3][3]
Definition eFEXegAlgo.h:70
SG::ReadHandleKey< LVL1::eTowerContainer > m_eTowerContainerKey
Definition eFEXegAlgo.h:97
Gaudi::Property< int > m_algoVersion
Definition eFEXegAlgo.h:94
virtual void setup(int inputTable[3][3], int efex_id, int fpga_id, int central_eta)
virtual void getRealEta(float &eta)
virtual std::unique_ptr< eFEXegTOB > geteFEXegTOB()
virtual void getSums(unsigned int seed, bool UnD, std::vector< unsigned int > &RetaSums, std::vector< unsigned int > &RhadSums, std::vector< unsigned int > &WstotSums)
virtual StatusCode initialize()
standard Athena-Algorithm method
SG::ReadCondHandleKey< CondAttrListCollection > m_dmCorrectionsKey
Definition eFEXegAlgo.h:100
eFEXegAlgo(const std::string &type, const std::string &name, const IInterface *parent)
Constructors.
virtual unsigned int dmCorrection(unsigned int ET, unsigned int layer)
virtual void getCoreEMTowerET(unsigned int &et)
virtual void getWstot(std::vector< unsigned int > &)
unsigned int m_seedID
Definition eFEXegAlgo.h:69
virtual void getWindowET(int layer, int jPhi, int SCID, unsigned int &)
virtual void getClusterCells(std::vector< unsigned int > &cellETs)
Return cell ET values used in cluster.
Gaudi::Property< bool > m_dmCorr
Definition eFEXegAlgo.h:93
virtual void getRhad(std::vector< unsigned int > &)
virtual StatusCode safetyTest() const
virtual void getCoreHADTowerET(unsigned int &et)
const LVL1::eTower * findTower(int towerID) const
fast find method given identifier.
The eTower class is an interface object for eFEX trigger algorithms The purposes are twofold:
Definition eTower.h:38
int getLayerTotalET(unsigned int layer) const
Get total ET sum of all cells in a given layer in MeV.
Definition eTower.cxx:249
float phi() const
Definition eTower.h:65
int getPosNeg() const
Definition eTower.h:121
int getET(unsigned int layer, int cell=0) const
Get ET of a specified cell in MeV.
Definition eTower.cxx:172
float eta() const
Definition eTower.h:64
const_pointer_type cptr()
eFexTowerBuilder creates xAOD::eFexTowerContainer from supercells (LATOME) and triggerTowers (TREX) i...
Extra patterns decribing particle interation process.