ATLAS Offline Software
Loading...
Searching...
No Matches
LVL1::gFEXSim Class Reference

The gFEXSim class defines the structure of the gFEX Its purpose is: More...

#include <gFEXSim.h>

Inheritance diagram for LVL1::gFEXSim:
Collaboration diagram for LVL1::gFEXSim:

Public Member Functions

 gFEXSim (const std::string &type, const std::string &name, const IInterface *parent)
 Constructors.
virtual ~gFEXSim ()
 Destructor.
virtual void reset () override
virtual void execute () override
virtual StatusCode initialize () override
virtual StatusCode executegFEXSim (const gTowersIDs &tmp, gFEXOutputCollection *gFEXOutputs) override
virtual std::vector< uint32_t > getgRhoTOBs () const override
virtual std::vector< uint32_t > getgBlockTOBs () const override
virtual std::vector< uint32_t > getgJetTOBs () const override
virtual std::vector< int32_t > getgScalarEJwojTOBs () const override
virtual std::vector< uint32_t > getgMETComponentsJwojTOBs () const override
virtual std::vector< uint32_t > getgMHTComponentsJwojTOBs () const override
virtual std::vector< uint32_t > getgMSTComponentsJwojTOBs () const override
virtual std::vector< uint32_t > getgMETComponentsNoiseCutTOBs () const override
virtual std::vector< uint32_t > getgMETComponentsRmsTOBs () const override
virtual std::vector< uint32_t > getgScalarENoiseCutTOBs () const override
virtual std::vector< uint32_t > getgScalarERmsTOBs () const override
ServiceHandle< StoreGateSvc > & evtStore ()
 The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
const ServiceHandle< StoreGateSvc > & detStore () const
 The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
virtual StatusCode sysInitialize () override
 Perform system initialization for an algorithm.
virtual StatusCode sysStart () override
 Handle START transition.
virtual std::vector< Gaudi::DataHandle * > inputHandles () const override
 Return this algorithm's input handles.
virtual std::vector< Gaudi::DataHandle * > outputHandles () const override
 Return this algorithm's output handles.
Gaudi::Details::PropertyBase & declareProperty (Gaudi::Property< T, V, H > &t)
void updateVHKA (Gaudi::Details::PropertyBase &)
MsgStream & msg () const
bool msgLvl (const MSG::Level lvl) const

Static Public Member Functions

static const InterfaceID & interfaceID ()

Protected Member Functions

void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce (T &h)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

gTowersIDs m_gTowersIDs
 Internal data.
CaloCellContainer m_sCellsCollection
std::vector< uint32_t > m_gRhoTobWords
std::vector< uint32_t > m_gBlockTobWords
std::vector< uint32_t > m_gJetTobWords
std::vector< int32_t > m_gScalarEJwojTobWords
std::vector< uint32_t > m_gMETComponentsJwojTobWords
std::vector< uint32_t > m_gMHTComponentsJwojTobWords
std::vector< uint32_t > m_gMSTComponentsJwojTobWords
std::vector< uint32_t > m_gMETComponentsNoiseCutTobWords
std::vector< uint32_t > m_gMETComponentsRmsTobWords
std::vector< uint32_t > m_gScalarENoiseCutTobWords
std::vector< uint32_t > m_gScalarERmsTobWords
ToolHandle< IgFEXFPGAm_gFEXFPGA_Tool {this, "gFEXFPGATool", "LVL1::gFEXFPGA", "Tool that simulates the FPGA hardware"}
ToolHandle< IgFEXJetAlgom_gFEXJetAlgoTool {this, "gFEXJetAlgoTool", "LVL1::gFEXJetAlgo", "Tool that runs the gFEX jet algorithm"}
ToolHandle< IgFEXJwoJAlgom_gFEXJwoJAlgoTool {this, "gFEXJwoJAlgoTool", "LVL1::gFEXJwoJAlgo", "Tool that runs the gFEX Jets without Jets algorithm"}
ToolHandle< IgFEXaltMetAlgom_gFEXaltMetAlgoTool {this, "gFEXaltMetAlgoTool", "LVL1::gFEXaltMetAlgo", "Tool that runs the gFEX noise cut and rho+RMS algorithms for MET"}
SG::ReadHandleKey< TrigConf::L1Menum_l1MenuKey {this, "L1TriggerMenu", "DetectorStore+L1TriggerMenu","Name of the L1Menu object to read configuration from"}
SG::WriteHandleKey< xAOD::gFexTowerContainerm_gTowersWriteKey {this,"gTowersWriteKey" ,"L1_gFexTriggerTowers", "Write gFexEDM Trigger Tower container"}
StoreGateSvc_t m_evtStore
 Pointer to StoreGate (event store by default)
StoreGateSvc_t m_detStore
 Pointer to StoreGate (detector store by default)
std::vector< SG::VarHandleKeyArray * > m_vhka
bool m_varHandleArraysDeclared

Detailed Description

The gFEXSim class defines the structure of the gFEX Its purpose is:

  • to emulate the steps taken in processing data for gFEX in hardware and firmware
  • It will need to interact with gTowers and produce the gTOBs. It will be created and handed data by gFEXSysSim

Definition at line 38 of file gFEXSim.h.

Member Typedef Documentation

◆ StoreGateSvc_t

typedef ServiceHandle<StoreGateSvc> AthCommonDataStore< AthCommonMsg< AlgTool > >::StoreGateSvc_t
privateinherited

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ gFEXSim()

LVL1::gFEXSim::gFEXSim ( const std::string & type,
const std::string & name,
const IInterface * parent )

Constructors.

Definition at line 21 of file gFEXSim.cxx.

21 :
22 AthAlgTool(type,name,parent)
23 {
24 declareInterface<IgFEXSim>(this);
25 }
AthAlgTool()
Default constructor:

◆ ~gFEXSim()

LVL1::gFEXSim::~gFEXSim ( )
virtual

Destructor.

Definition at line 42 of file gFEXSim.cxx.

42 {
43 }

Member Function Documentation

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareGaudiProperty ( Gaudi::Property< T, V, H > & hndl,
const SG::VarHandleKeyType &  )
inlineprivateinherited

specialization for handling Gaudi::Property<SG::VarHandleKey>

Definition at line 156 of file AthCommonDataStore.h.

158 {
160 hndl.value(),
161 hndl.documentation());
162
163 }
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)

◆ declareProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareProperty ( Gaudi::Property< T, V, H > & t)
inlineinherited

Definition at line 145 of file AthCommonDataStore.h.

145 {
146 typedef typename SG::HandleClassifier<T>::type htype;
148 }
Gaudi::Details::PropertyBase & declareGaudiProperty(Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
specialization for handling Gaudi::Property<SG::VarHandleKey>

◆ detStore()

const ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::detStore ( ) const
inlineinherited

The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 95 of file AthCommonDataStore.h.

◆ evtStore()

ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::evtStore ( )
inlineinherited

The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 85 of file AthCommonDataStore.h.

◆ execute()

void LVL1::gFEXSim::execute ( )
overridevirtual

Implements LVL1::IgFEXSim.

Definition at line 55 of file gFEXSim.cxx.

55 {
56
57 }

◆ executegFEXSim()

StatusCode LVL1::gFEXSim::executegFEXSim ( const gTowersIDs & tmp,
gFEXOutputCollection * gFEXOutputs )
overridevirtual

Implements LVL1::IgFEXSim.

Definition at line 59 of file gFEXSim.cxx.

59 {
60
61 // Container to save gTowers
62 SG::WriteHandle<xAOD::gFexTowerContainer> gTowersContainer(m_gTowersWriteKey);
63 ATH_CHECK(gTowersContainer.record(std::make_unique<xAOD::gFexTowerContainer>(), std::make_unique<xAOD::gFexTowerAuxContainer>()));
64 ATH_MSG_DEBUG("Recorded gFexTriggerTower container with key " << gTowersContainer.key());
65
66 int rows = tmp_gTowersIDs_subset.size();
67 int cols = tmp_gTowersIDs_subset[0].size();
68
69 std::copy(&tmp_gTowersIDs_subset[0][0], &tmp_gTowersIDs_subset[0][0]+(rows*cols),&m_gTowersIDs[0][0]);
70
71 gTowersType Atwr = {{{0}}};
72 gTowersType Btwr = {{{0}}};
73 gTowersType Ctwr = {{{0}}};
74
75 gTowersType Atwr50 = {{{0}}};
76 gTowersType Btwr50 = {{{0}}};
77 gTowersType Ctwr50 = {{{0}}};
78
79 gTowersType Asat = {{{0}}};
80 gTowersType Bsat = {{{0}}};
81 gTowersType Csat = {{{0}}};
82
83
84 //FPGA A----------------------------------------------------------------------------------------------------------------------------------------------
85 gTowersCentral tmp_gTowersIDs_subset_centralFPGA;
86 memset(&tmp_gTowersIDs_subset_centralFPGA, 0, sizeof tmp_gTowersIDs_subset_centralFPGA);
87 for (int myrow = 0; myrow<FEXAlgoSpaceDefs::centralNphi; myrow++){
88 for (int mycol = 0; mycol<12; mycol++){
89 tmp_gTowersIDs_subset_centralFPGA[myrow][mycol] = tmp_gTowersIDs_subset[myrow][mycol+8];
90 }
91 }
92 ATH_CHECK(m_gFEXFPGA_Tool->init(0));
93 m_gFEXFPGA_Tool->FillgTowerEDMCentral(gTowersContainer, tmp_gTowersIDs_subset_centralFPGA, Atwr, Atwr50, Asat);
94 m_gFEXFPGA_Tool->reset();
95
96 //FPGA A----------------------------------------------------------------------------------------------------------------------------------------------
97
98 //FPGA B----------------------------------------------------------------------------------------------------------------------------------------------
99 gTowersCentral tmp_gTowersIDs_subset_centralFPGA_B;
100 memset(&tmp_gTowersIDs_subset_centralFPGA_B, 0, sizeof tmp_gTowersIDs_subset_centralFPGA_B);
101 for (int myrow = 0; myrow<FEXAlgoSpaceDefs::centralNphi; myrow++){
102 for (int mycol = 0; mycol<12; mycol++){
103 tmp_gTowersIDs_subset_centralFPGA_B[myrow][mycol] = tmp_gTowersIDs_subset[myrow][mycol+20];
104 }
105 }
106 ATH_CHECK(m_gFEXFPGA_Tool->init(1));
107 m_gFEXFPGA_Tool->FillgTowerEDMCentral(gTowersContainer, tmp_gTowersIDs_subset_centralFPGA_B, Btwr, Btwr50, Bsat);
108 m_gFEXFPGA_Tool->reset();
109
110 //FPGA B----------------------------------------------------------------------------------------------------------------------------------------------
111
112
113 //FPGA C ----------------------------------------------------------------------------------------------------------------------------------------------
114
115 // C-N
116 //Use a matrix with 32 rows, even if FPGA-N (negative) also deals with regions of 16 bins in phi (those connected to FCAL).
117 //We have 4 columns with 32 rows and 4 columns with 16 rows for each FPGA-C.
118 //So we use a matrix 32x8 but we fill only half of it in the region 3.3<|eta|<4.8.
119 gTowersForward tmp_gTowersIDs_subset_forwardFPGA_N;
120 memset(&tmp_gTowersIDs_subset_forwardFPGA_N, 0, sizeof tmp_gTowersIDs_subset_forwardFPGA_N);
121 for (int myrow = 0; myrow<FEXAlgoSpaceDefs::forwardNphi; myrow++){
122 for (int mycol = 0; mycol<4; mycol++){
123 tmp_gTowersIDs_subset_forwardFPGA_N[myrow][mycol] = tmp_gTowersIDs_subset[myrow][mycol];
124 }
125 }
126 for (int myrow = 0; myrow<FEXAlgoSpaceDefs::centralNphi; myrow++){
127 for (int mycol = 4; mycol<8; mycol++){
128 tmp_gTowersIDs_subset_forwardFPGA_N[myrow][mycol] = tmp_gTowersIDs_subset[myrow][mycol];
129 }
130 }
131
132 // C-P
133 //Use a matrix with 32 rows, even if FPGA-C (positive) also deals with regions of 16 bins in phi (those connected to FCAL).
134 //We have 4 columns with 32 rows and 4 columns with 16 rows for each FPGA-C.
135 //So we use a matrix 32x8 but we fill only half of it in the region 3.3<|eta|<4.8.
136 gTowersForward tmp_gTowersIDs_subset_forwardFPGA_P;
137 memset(&tmp_gTowersIDs_subset_forwardFPGA_P, 0, sizeof tmp_gTowersIDs_subset_forwardFPGA_P);
138 for (int myrow = 0; myrow<FEXAlgoSpaceDefs::centralNphi; myrow++){
139 for (int mycol = 0; mycol<4; mycol++){
140 tmp_gTowersIDs_subset_forwardFPGA_P[myrow][mycol] = tmp_gTowersIDs_subset[myrow][mycol+32];
141 }
142 }
143 for (int myrow = 0; myrow<FEXAlgoSpaceDefs::forwardNphi; myrow++){
144 for (int mycol = 4; mycol<8; mycol++){
145 tmp_gTowersIDs_subset_forwardFPGA_P[myrow][mycol] = tmp_gTowersIDs_subset[myrow][mycol+32];
146 }
147 }
148
149 ATH_CHECK(m_gFEXFPGA_Tool->init(2));
150 m_gFEXFPGA_Tool->FillgTowerEDMForward(gTowersContainer, tmp_gTowersIDs_subset_forwardFPGA_N, tmp_gTowersIDs_subset_forwardFPGA_P, Ctwr, Ctwr50, Csat);
151 m_gFEXFPGA_Tool->reset();
152
153
154 //FPGA C----------------------------------------------------------------------------------------------------------------------------------------------
155
156 // Retrieve the L1 menu configuration
157 SG::ReadHandle<TrigConf::L1Menu> l1Menu (m_l1MenuKey/*, ctx*/);
158 ATH_CHECK(l1Menu.isValid());
159
160 //Parameters related to gLJ (large-R jet objects - gJet)
161 auto & thr_gLJ = l1Menu->thrExtraInfo().gLJ();
162 int gLJ_seedThrA = 0;
163 int gLJ_seedThrB = 0;
164 int gLJ_seedThrC = 0;
165 gLJ_seedThrA = thr_gLJ.seedThrCounts('A'); //defined in GeV by default
166 gLJ_seedThrB = thr_gLJ.seedThrCounts('B'); //defined in GeV by default
167 gLJ_seedThrC = thr_gLJ.seedThrCounts('C'); //defined in GeV by default
168
169 int gLJ_ptMinToTopoCounts1 = 0;
170 int gLJ_ptMinToTopoCounts2 = 0;
171 gLJ_ptMinToTopoCounts1 = thr_gLJ.ptMinToTopoCounts(1);
172 gLJ_ptMinToTopoCounts2 = thr_gLJ.ptMinToTopoCounts(2);
173 float gLJ_rhoMaxA = 0;
174 float gLJ_rhoMaxB = 0;
175 float gLJ_rhoMaxC = 0;
176
177 gLJ_rhoMaxA = (thr_gLJ.rhoTowerMax('A')*1000)/50;//Values are given in GeV, need to be converted with 50MeV scale to be used in PU calculation
178 gLJ_rhoMaxB = (thr_gLJ.rhoTowerMax('B')*1000)/50;//Values are given in GeV, need to be converted with 50MeV scale to be used in PU calculation
179 gLJ_rhoMaxC = (thr_gLJ.rhoTowerMax('C')*1000)/50;//Values are given in GeV, need to be converted with 50MeV scale to be used in PU calculation
180
181
182 //Parameters related to gJ (small-R jet objects - gBlock)
183 auto & thr_gJ = l1Menu->thrExtraInfo().gJ();
184 int gJ_ptMinToTopoCounts1 = 0;
185 int gJ_ptMinToTopoCounts2 = 0;
186 gJ_ptMinToTopoCounts1 = thr_gJ.ptMinToTopoCounts(1);
187 gJ_ptMinToTopoCounts2 = thr_gJ.ptMinToTopoCounts(2);
188
189
190 int pucA = 0;
191 int pucB = 0;
192 int pucC = 0;
193 //note that jetThreshold is not a configurable parameter in firmware, it is used to check that jet values are positive
194 int jetThreshold = FEXAlgoSpaceDefs::jetThr; //this threshold is set by the online software
195
196 if (FEXAlgoSpaceDefs::ENABLE_PUC == true){
197 m_gFEXJetAlgoTool->pileUpCalculation(Atwr50, gLJ_rhoMaxA, 1, pucA);
198 m_gFEXJetAlgoTool->pileUpCalculation(Btwr50, gLJ_rhoMaxB, 1, pucB);
199 m_gFEXJetAlgoTool->pileUpCalculation(Ctwr50, gLJ_rhoMaxC, 1, pucC);
200 }
201
202
203
204 // The output TOBs, to be filled by the gFEXJetAlgoTool
205 std::array<uint32_t, 7> ATOB1_dat = {0};
206 std::array<uint32_t, 7> ATOB2_dat = {0};
207 std::array<uint32_t, 7> BTOB1_dat = {0};
208 std::array<uint32_t, 7> BTOB2_dat = {0};
209 std::array<uint32_t, 7> CTOB1_dat = {0};
210 std::array<uint32_t, 7> CTOB2_dat = {0};
211
212
213 // Pass the energy matrices to the algo tool, and run the algorithms
214 auto tobs_v = m_gFEXJetAlgoTool->largeRfinder(Atwr, Btwr, Ctwr, Asat, Bsat, Csat, pucA, pucB, pucC,
215 gLJ_seedThrA, gLJ_seedThrB, gLJ_seedThrC, gJ_ptMinToTopoCounts1, gJ_ptMinToTopoCounts2,
216 jetThreshold, gLJ_ptMinToTopoCounts1, gLJ_ptMinToTopoCounts2,
217 ATOB1_dat, ATOB2_dat,
218 BTOB1_dat, BTOB2_dat,
219 CTOB1_dat, CTOB2_dat);
220
221 m_gRhoTobWords.resize(3);
222 m_gBlockTobWords.resize(12);
223 m_gJetTobWords.resize(6);
224
225 m_gRhoTobWords[0] = ATOB2_dat[0];//Pile up correction A
226 m_gRhoTobWords[1] = BTOB2_dat[0];//Pile up correction B
227 m_gRhoTobWords[2] = CTOB2_dat[0];//Pile up correction C
228
229 //Placing the gBlock TOBs into a dedicated array
230 m_gBlockTobWords[0] = ATOB1_dat[1];//leading gBlock in FPGA A, eta bins (0--5)
231 m_gBlockTobWords[1] = ATOB2_dat[1];//leading gBlock in FPGA A, eta bins (6--11)
232 m_gBlockTobWords[2] = BTOB1_dat[1];//leading gBlock in FPGA B, eta bins (0--5)
233 m_gBlockTobWords[3] = BTOB2_dat[1];//leading gBlock in FPGA B, eta bins (6--11)
234
235 m_gBlockTobWords[4] = ATOB1_dat[2];//subleading gBlock in FPGA A, eta bins (0--5)
236 m_gBlockTobWords[5] = ATOB2_dat[2];//subleading gBlock in FPGA A, eta bins (6--11)
237 m_gBlockTobWords[6] = BTOB1_dat[2];//subleading gBlock in FPGA B, eta bins (0--5)
238 m_gBlockTobWords[7] = BTOB2_dat[2];//subleading gBlock in FPGA B, eta bins (6--11)
239
240 m_gBlockTobWords[8] = CTOB1_dat[1];//leading gBlock in FPGA C, eta negative
241 m_gBlockTobWords[9] = CTOB2_dat[1];//leading gBlock in FPGA C, eta positive
242 m_gBlockTobWords[10] = CTOB1_dat[2];//sub-leading gBlock in FPGA C, eta negative
243 m_gBlockTobWords[11] = CTOB2_dat[2];//sub-leading gBlock in FPGA C, eta positive
244
245 //Placing the gJet TOBs into a dedicated array
246 m_gJetTobWords[0] = ATOB1_dat[3];//leading gJet in FPGA A, eta bins (0--5)
247 m_gJetTobWords[1] = ATOB2_dat[3];//leading gJet in FPGA A, eta bins (6--11)
248 m_gJetTobWords[2] = BTOB1_dat[3];//leading gJet in FPGA B, eta bins (0--5)
249 m_gJetTobWords[3] = BTOB2_dat[3];//leading gJet in FPGA B, eta bins (6--11)
250 m_gJetTobWords[4] = CTOB1_dat[3];//leading gJet in FPGA C negative
251 m_gJetTobWords[5] = CTOB2_dat[3];//leading gJet in FPGA C positive
252
253
254 // Use the gFEXJetAlgoTool
255 std::array<int32_t, 4> outJwojTOB = {0};
256 std::array<uint32_t, 4> outAltMetTOB = {0};
257
258 //Parameters related to gXE (MET objects, both JwoJ and alternative MET calculation)
259 auto & thr_gXE = l1Menu->thrExtraInfo().gXE();
260 int gXE_seedThrA = 0;
261 int gXE_seedThrB = 0;
262 int gXE_seedThrC = 0;
263 gXE_seedThrA = thr_gXE.seedThr('A'); //defined in counts
264 gXE_seedThrB = thr_gXE.seedThr('B'); //defined in counts
265 gXE_seedThrC = thr_gXE.seedThr('C'); //defined in counts
266
267
268 int aFPGA_A = thr_gXE.JWOJ_param('A','a');// 1003
269 int bFPGA_A = thr_gXE.JWOJ_param('A','b');// 409
270 int aFPGA_B = thr_gXE.JWOJ_param('B','a');// 1003
271 int bFPGA_B = thr_gXE.JWOJ_param('B','b');// 409
272 int aFPGA_C = thr_gXE.JWOJ_param('C','a');// 1003
273 int bFPGA_C = thr_gXE.JWOJ_param('C','b');// 409
274
275 //Set constants for JwoJ and run the algorithm
276 m_gFEXJwoJAlgoTool->setAlgoConstant(aFPGA_A, bFPGA_A,
277 aFPGA_B, bFPGA_B,
278 aFPGA_C, bFPGA_C,
279 gXE_seedThrA, gXE_seedThrB, gXE_seedThrC);
280
281 auto global_tobs = m_gFEXJwoJAlgoTool->jwojAlgo(Atwr, Btwr, Ctwr, outJwojTOB);
282
283 m_gScalarEJwojTobWords.resize(1);
287
288
289 //Placing the global TOBs into a dedicated array
290 m_gScalarEJwojTobWords[0] = outJwojTOB[0];//
291 m_gMETComponentsJwojTobWords[0] = outJwojTOB[1];//
292 m_gMHTComponentsJwojTobWords[0] = outJwojTOB[2];//
293 m_gMSTComponentsJwojTobWords[0] = outJwojTOB[3];//
294
295
296 //Set constants for noise cut and rho+RMS and run the algorithms
297 std::vector<int> thr_A (12, 0);//To be retrieved from COOL database in the future
298 std::vector<int> thr_B (12, 0);//To be retrieved from COOL database in the future
299
300 m_gFEXaltMetAlgoTool->setAlgoConstant(std::move(thr_A) , std::move(thr_B), 10000/200);
301
302 m_gFEXaltMetAlgoTool->altMetAlgo(Atwr, Btwr, outAltMetTOB);
303
307 m_gScalarERmsTobWords.resize(1);
308
309
310 //Placing the global TOBs into a dedicated array
311 m_gMETComponentsNoiseCutTobWords[0] = outAltMetTOB[0];//
312 m_gMETComponentsRmsTobWords[0] = outAltMetTOB[1];//
313 m_gScalarENoiseCutTobWords[0] = outAltMetTOB[2];//
314 m_gScalarERmsTobWords[0] = outAltMetTOB[3];//
315
316 for (int i = 0; i <14; i++){
317 gFEXOutputs->addJetTob(tobs_v[i]->getWord());
318 gFEXOutputs->addValueJet("EtaJet", tobs_v[i]->getEta());
319 gFEXOutputs->addValueJet("PhiJet", tobs_v[i]->getPhi());
320 gFEXOutputs->addValueJet("ETJet", tobs_v[i]->getET());
321 gFEXOutputs->addValueJet("StatusJet", tobs_v[i]->getStatus());
322 gFEXOutputs->addValueJet("TobIDJet", tobs_v[i]->getTobID());
323 gFEXOutputs->fillJet();
324
325 }
326
327 for (int i = 0; i <4; i++){
328 gFEXOutputs->addGlobalTob(global_tobs[i]->getWord());
329 gFEXOutputs->addValueGlobal("GlobalQuantity1", global_tobs[i]->getQuantity1());
330 gFEXOutputs->addValueGlobal("GlobalQuantity2", global_tobs[i]->getQuantity2());
331 gFEXOutputs->addValueGlobal("SaturationGlobal", global_tobs[i]->getSaturation());
332 gFEXOutputs->addValueGlobal("TobIDGlobal", global_tobs[i]->getTobID());
333 gFEXOutputs->addValueGlobal("GlobalStatus1", global_tobs[i]->getStatus1());
334 gFEXOutputs->addValueGlobal("GlobalStatus2", global_tobs[i]->getStatus2());
335 gFEXOutputs->fillGlobal();
336
337 }
338
339 return StatusCode::SUCCESS;
340
341}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_DEBUG(x)
static constexpr int centralNphi
static constexpr int jetThr
static constexpr int forwardNphi
static constexpr bool ENABLE_PUC
ToolHandle< IgFEXFPGA > m_gFEXFPGA_Tool
Definition gFEXSim.h:111
std::vector< uint32_t > m_gMETComponentsJwojTobWords
Definition gFEXSim.h:95
gTowersIDs m_gTowersIDs
Internal data.
Definition gFEXSim.h:83
std::vector< int32_t > m_gScalarEJwojTobWords
Definition gFEXSim.h:93
SG::ReadHandleKey< TrigConf::L1Menu > m_l1MenuKey
Definition gFEXSim.h:119
ToolHandle< IgFEXJwoJAlgo > m_gFEXJwoJAlgoTool
Definition gFEXSim.h:115
SG::WriteHandleKey< xAOD::gFexTowerContainer > m_gTowersWriteKey
Definition gFEXSim.h:121
std::vector< uint32_t > m_gMHTComponentsJwojTobWords
Definition gFEXSim.h:97
std::vector< uint32_t > m_gMSTComponentsJwojTobWords
Definition gFEXSim.h:99
std::vector< uint32_t > m_gScalarERmsTobWords
Definition gFEXSim.h:107
std::vector< uint32_t > m_gScalarENoiseCutTobWords
Definition gFEXSim.h:105
std::vector< uint32_t > m_gMETComponentsRmsTobWords
Definition gFEXSim.h:103
std::vector< uint32_t > m_gRhoTobWords
Definition gFEXSim.h:87
std::vector< uint32_t > m_gJetTobWords
Definition gFEXSim.h:91
ToolHandle< IgFEXJetAlgo > m_gFEXJetAlgoTool
Definition gFEXSim.h:113
ToolHandle< IgFEXaltMetAlgo > m_gFEXaltMetAlgoTool
Definition gFEXSim.h:117
std::vector< uint32_t > m_gMETComponentsNoiseCutTobWords
Definition gFEXSim.h:101
std::vector< uint32_t > m_gBlockTobWords
Definition gFEXSim.h:89
float getPhi(const xAOD::TrackParticle &p)
Accessor utility function for getting the value of phi.
std::array< std::array< int, 8 >, 32 > gTowersForward
std::array< std::array< int, 12 >, 32 > gTowersCentral
std::array< std::array< int, 12 >, 32 > gTowersType
Definition IgFEXFPGA.h:25

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::extraDeps_update_handler ( Gaudi::Details::PropertyBase & ExtraDeps)
protectedinherited

Add StoreName to extra input/output deps as needed.

use the logic of the VarHandleKey to parse the DataObjID keys supplied via the ExtraInputs and ExtraOuputs Properties to add the StoreName if it's not explicitly given

◆ getgBlockTOBs()

std::vector< uint32_t > LVL1::gFEXSim::getgBlockTOBs ( ) const
overridevirtual

Implements LVL1::IgFEXSim.

Definition at line 349 of file gFEXSim.cxx.

350{
351 return m_gBlockTobWords;
352}

◆ getgJetTOBs()

std::vector< uint32_t > LVL1::gFEXSim::getgJetTOBs ( ) const
overridevirtual

Implements LVL1::IgFEXSim.

Definition at line 354 of file gFEXSim.cxx.

355{
356 return m_gJetTobWords;
357}

◆ getgMETComponentsJwojTOBs()

std::vector< uint32_t > LVL1::gFEXSim::getgMETComponentsJwojTOBs ( ) const
overridevirtual

Implements LVL1::IgFEXSim.

Definition at line 364 of file gFEXSim.cxx.

365{
367}

◆ getgMETComponentsNoiseCutTOBs()

std::vector< uint32_t > LVL1::gFEXSim::getgMETComponentsNoiseCutTOBs ( ) const
overridevirtual

Implements LVL1::IgFEXSim.

Definition at line 379 of file gFEXSim.cxx.

380{
382}

◆ getgMETComponentsRmsTOBs()

std::vector< uint32_t > LVL1::gFEXSim::getgMETComponentsRmsTOBs ( ) const
overridevirtual

Implements LVL1::IgFEXSim.

Definition at line 384 of file gFEXSim.cxx.

385{
387}

◆ getgMHTComponentsJwojTOBs()

std::vector< uint32_t > LVL1::gFEXSim::getgMHTComponentsJwojTOBs ( ) const
overridevirtual

Implements LVL1::IgFEXSim.

Definition at line 369 of file gFEXSim.cxx.

370{
372}

◆ getgMSTComponentsJwojTOBs()

std::vector< uint32_t > LVL1::gFEXSim::getgMSTComponentsJwojTOBs ( ) const
overridevirtual

Implements LVL1::IgFEXSim.

Definition at line 374 of file gFEXSim.cxx.

375{
377}

◆ getgRhoTOBs()

std::vector< uint32_t > LVL1::gFEXSim::getgRhoTOBs ( ) const
overridevirtual

Implements LVL1::IgFEXSim.

Definition at line 344 of file gFEXSim.cxx.

345{
346 return m_gRhoTobWords;
347}

◆ getgScalarEJwojTOBs()

std::vector< int32_t > LVL1::gFEXSim::getgScalarEJwojTOBs ( ) const
overridevirtual

Implements LVL1::IgFEXSim.

Definition at line 359 of file gFEXSim.cxx.

360{
362}

◆ getgScalarENoiseCutTOBs()

std::vector< uint32_t > LVL1::gFEXSim::getgScalarENoiseCutTOBs ( ) const
overridevirtual

Implements LVL1::IgFEXSim.

Definition at line 389 of file gFEXSim.cxx.

390{
392}

◆ getgScalarERmsTOBs()

std::vector< uint32_t > LVL1::gFEXSim::getgScalarERmsTOBs ( ) const
overridevirtual

Implements LVL1::IgFEXSim.

Definition at line 394 of file gFEXSim.cxx.

395{
397}

◆ initialize()

StatusCode LVL1::gFEXSim::initialize ( )
overridevirtual

Definition at line 45 of file gFEXSim.cxx.

45 {
46 ATH_CHECK( m_gFEXFPGA_Tool.retrieve() );
47 ATH_CHECK( m_gFEXJetAlgoTool.retrieve() );
48 ATH_CHECK( m_gFEXJwoJAlgoTool.retrieve() );
49 ATH_CHECK( m_gFEXaltMetAlgoTool.retrieve() );
50 ATH_CHECK(m_l1MenuKey.initialize());
51 ATH_CHECK(m_gTowersWriteKey.initialize());
52 return StatusCode::SUCCESS;
53 }

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::inputHandles ( ) const
overridevirtualinherited

Return this algorithm's input handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ interfaceID()

const InterfaceID & LVL1::IgFEXSim::interfaceID ( )
inlinestaticinherited

Definition at line 59 of file IgFEXSim.h.

60 {
61 return IID_IgFEXSim;
62 }
static const InterfaceID IID_IgFEXSim("LVL1::IgFEXSim", 1, 0)

◆ msg()

MsgStream & AthCommonMsg< AlgTool >::msg ( ) const
inlineinherited

Definition at line 24 of file AthCommonMsg.h.

24 {
25 return this->msgStream();
26 }

◆ msgLvl()

bool AthCommonMsg< AlgTool >::msgLvl ( const MSG::Level lvl) const
inlineinherited

Definition at line 30 of file AthCommonMsg.h.

30 {
31 return this->msgLevel(lvl);
32 }

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::outputHandles ( ) const
overridevirtualinherited

Return this algorithm's output handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ renounce()

std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > AthCommonDataStore< AthCommonMsg< AlgTool > >::renounce ( T & h)
inlineprotectedinherited

Definition at line 380 of file AthCommonDataStore.h.

381 {
382 h.renounce();
384 }
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)

◆ renounceArray()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ reset()

void LVL1::gFEXSim::reset ( )
overridevirtual

Implements LVL1::IgFEXSim.

Definition at line 28 of file gFEXSim.cxx.

29 {
30 int rows = m_gTowersIDs.size();
31 int cols = m_gTowersIDs[0].size();
32
33 for (int i=0; i<rows; i++){
34 for (int j=0; j<cols; j++){
35 m_gTowersIDs[i][j] = 0;
36 }
37 }
38
39 }

◆ sysInitialize()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysInitialize ( )
overridevirtualinherited

Perform system initialization for an algorithm.

We override this to declare all the elements of handle key arrays at the end of initialization. See comments on updateVHKA.

Reimplemented in asg::AsgMetadataTool, AthCheckedComponent< AthAlgTool >, AthCheckedComponent<::AthAlgTool >, and DerivationFramework::CfAthAlgTool.

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysStart ( )
overridevirtualinherited

Handle START transition.

We override this in order to make sure that conditions handle keys can cache a pointer to the conditions container.

◆ updateVHKA()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::updateVHKA ( Gaudi::Details::PropertyBase & )
inlineinherited

Definition at line 308 of file AthCommonDataStore.h.

308 {
309 // debug() << "updateVHKA for property " << p.name() << " " << p.toString()
310 // << " size: " << m_vhka.size() << endmsg;
311 for (auto &a : m_vhka) {
313 for (auto k : keys) {
314 k->setOwner(this);
315 }
316 }
317 }
std::vector< SG::VarHandleKeyArray * > m_vhka

Member Data Documentation

◆ m_detStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_evtStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_evtStore
privateinherited

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_gBlockTobWords

std::vector<uint32_t> LVL1::gFEXSim::m_gBlockTobWords
private

Definition at line 89 of file gFEXSim.h.

◆ m_gFEXaltMetAlgoTool

ToolHandle<IgFEXaltMetAlgo> LVL1::gFEXSim::m_gFEXaltMetAlgoTool {this, "gFEXaltMetAlgoTool", "LVL1::gFEXaltMetAlgo", "Tool that runs the gFEX noise cut and rho+RMS algorithms for MET"}
private

Definition at line 117 of file gFEXSim.h.

117{this, "gFEXaltMetAlgoTool", "LVL1::gFEXaltMetAlgo", "Tool that runs the gFEX noise cut and rho+RMS algorithms for MET"};

◆ m_gFEXFPGA_Tool

ToolHandle<IgFEXFPGA> LVL1::gFEXSim::m_gFEXFPGA_Tool {this, "gFEXFPGATool", "LVL1::gFEXFPGA", "Tool that simulates the FPGA hardware"}
private

Definition at line 111 of file gFEXSim.h.

111{this, "gFEXFPGATool", "LVL1::gFEXFPGA", "Tool that simulates the FPGA hardware"};

◆ m_gFEXJetAlgoTool

ToolHandle<IgFEXJetAlgo> LVL1::gFEXSim::m_gFEXJetAlgoTool {this, "gFEXJetAlgoTool", "LVL1::gFEXJetAlgo", "Tool that runs the gFEX jet algorithm"}
private

Definition at line 113 of file gFEXSim.h.

113{this, "gFEXJetAlgoTool", "LVL1::gFEXJetAlgo", "Tool that runs the gFEX jet algorithm"};

◆ m_gFEXJwoJAlgoTool

ToolHandle<IgFEXJwoJAlgo> LVL1::gFEXSim::m_gFEXJwoJAlgoTool {this, "gFEXJwoJAlgoTool", "LVL1::gFEXJwoJAlgo", "Tool that runs the gFEX Jets without Jets algorithm"}
private

Definition at line 115 of file gFEXSim.h.

115{this, "gFEXJwoJAlgoTool", "LVL1::gFEXJwoJAlgo", "Tool that runs the gFEX Jets without Jets algorithm"};

◆ m_gJetTobWords

std::vector<uint32_t> LVL1::gFEXSim::m_gJetTobWords
private

Definition at line 91 of file gFEXSim.h.

◆ m_gMETComponentsJwojTobWords

std::vector<uint32_t> LVL1::gFEXSim::m_gMETComponentsJwojTobWords
private

Definition at line 95 of file gFEXSim.h.

◆ m_gMETComponentsNoiseCutTobWords

std::vector<uint32_t> LVL1::gFEXSim::m_gMETComponentsNoiseCutTobWords
private

Definition at line 101 of file gFEXSim.h.

◆ m_gMETComponentsRmsTobWords

std::vector<uint32_t> LVL1::gFEXSim::m_gMETComponentsRmsTobWords
private

Definition at line 103 of file gFEXSim.h.

◆ m_gMHTComponentsJwojTobWords

std::vector<uint32_t> LVL1::gFEXSim::m_gMHTComponentsJwojTobWords
private

Definition at line 97 of file gFEXSim.h.

◆ m_gMSTComponentsJwojTobWords

std::vector<uint32_t> LVL1::gFEXSim::m_gMSTComponentsJwojTobWords
private

Definition at line 99 of file gFEXSim.h.

◆ m_gRhoTobWords

std::vector<uint32_t> LVL1::gFEXSim::m_gRhoTobWords
private

Definition at line 87 of file gFEXSim.h.

◆ m_gScalarEJwojTobWords

std::vector<int32_t> LVL1::gFEXSim::m_gScalarEJwojTobWords
private

Definition at line 93 of file gFEXSim.h.

◆ m_gScalarENoiseCutTobWords

std::vector<uint32_t> LVL1::gFEXSim::m_gScalarENoiseCutTobWords
private

Definition at line 105 of file gFEXSim.h.

◆ m_gScalarERmsTobWords

std::vector<uint32_t> LVL1::gFEXSim::m_gScalarERmsTobWords
private

Definition at line 107 of file gFEXSim.h.

◆ m_gTowersIDs

gTowersIDs LVL1::gFEXSim::m_gTowersIDs
private

Internal data.

Definition at line 83 of file gFEXSim.h.

◆ m_gTowersWriteKey

SG::WriteHandleKey< xAOD::gFexTowerContainer > LVL1::gFEXSim::m_gTowersWriteKey {this,"gTowersWriteKey" ,"L1_gFexTriggerTowers", "Write gFexEDM Trigger Tower container"}
private

Definition at line 121 of file gFEXSim.h.

121{this,"gTowersWriteKey" ,"L1_gFexTriggerTowers", "Write gFexEDM Trigger Tower container"};

◆ m_l1MenuKey

SG::ReadHandleKey<TrigConf::L1Menu> LVL1::gFEXSim::m_l1MenuKey {this, "L1TriggerMenu", "DetectorStore+L1TriggerMenu","Name of the L1Menu object to read configuration from"}
private

Definition at line 119 of file gFEXSim.h.

119{this, "L1TriggerMenu", "DetectorStore+L1TriggerMenu","Name of the L1Menu object to read configuration from"};

◆ m_sCellsCollection

CaloCellContainer LVL1::gFEXSim::m_sCellsCollection
private

Definition at line 85 of file gFEXSim.h.

◆ m_varHandleArraysDeclared

bool AthCommonDataStore< AthCommonMsg< AlgTool > >::m_varHandleArraysDeclared
privateinherited

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

std::vector<SG::VarHandleKeyArray*> AthCommonDataStore< AthCommonMsg< AlgTool > >::m_vhka
privateinherited

Definition at line 398 of file AthCommonDataStore.h.


The documentation for this class was generated from the following files: