ATLAS Offline Software
Loading...
Searching...
No Matches
gFEXJwoJAlgo.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// gFEXJwoJAlg - Jets without jets algorithm for gFEX
6// -------------------
7// begin : 21 12 2023
8// email : cecilia.tosciri@cern.ch
9//***************************************************************************
10
11
12
13#include "gFEXJwoJAlgo.h"
15
16#include <cmath> //for std::sqrt
17
18namespace LVL1 {
19
20 // default constructor for persistency
21gFEXJwoJAlgo::gFEXJwoJAlgo(const std::string& type, const std::string& name, const IInterface* parent):
22 AthAlgTool(type, name, parent)
23 {
24 declareInterface<IgFEXJwoJAlgo>(this);
25 }
26
27
29
30 return StatusCode::SUCCESS;
31
32}
33
34
35void gFEXJwoJAlgo::setAlgoConstant(int aFPGA_A, int bFPGA_A,
36 int aFPGA_B, int bFPGA_B,
37 int aFPGA_C, int bFPGA_C,
38 int gXE_seedThrA, int gXE_seedThrB, int gXE_seedThrC) {
39 m_aFPGA_A = aFPGA_A;
40 m_bFPGA_A = bFPGA_A;
41 m_aFPGA_B = aFPGA_B;
42 m_bFPGA_B = bFPGA_B;
43 m_aFPGA_C = aFPGA_C;
44 m_bFPGA_C = bFPGA_C;
45 m_gBlockthresholdA = gXE_seedThrA;
46 m_gBlockthresholdB = gXE_seedThrB;
47 m_gBlockthresholdC = gXE_seedThrC;
48}
49
50std::vector<std::unique_ptr<gFEXJwoJTOB>> gFEXJwoJAlgo::jwojAlgo(const gTowersType& Atwr,const gTowersType& Btwr, const gTowersType& Ctwr,
51 std::array<int32_t, 4> & outTOB) const {
52
53
54 // input towers have 200 MeV LSB
55 constexpr bool SumETfast = true;
56
57 // find gBlocks
58 gTowersType AgBlk;
59 gTowersType Ascaled;
60
61 gTowersType BgBlk;
62 gTowersType Bscaled;
63
64 gTowersType CgBlk;
65 gTowersType Cscaled;
66
67 gTowersType hasSeed;
68
69 gBlockAB(Atwr, AgBlk, hasSeed, m_gBlockthresholdA);
70 gBlockAB(Btwr, BgBlk, hasSeed, m_gBlockthresholdB);
71 gBlockAB(Ctwr, CgBlk, hasSeed, m_gBlockthresholdC);
72
73
74 // switch to 10 bit number
75 // DMS -- do we eventaully need to check for overflows here?
76
77 for(int irow = 0; irow<FEXAlgoSpaceDefs::ABCrows; irow++){
78 for(int jcolumn = 0; jcolumn<FEXAlgoSpaceDefs::ABcolumns; jcolumn++){
79 Ascaled[irow][jcolumn] = Atwr[irow][jcolumn] >> 2;
80 AgBlk[irow][jcolumn] = AgBlk[irow][jcolumn] >> 2;
81
82 Bscaled[irow][jcolumn] = Btwr[irow][jcolumn] >> 2;
83 BgBlk[irow][jcolumn] = BgBlk[irow][jcolumn] >> 2;
84
85 Cscaled[irow][jcolumn] = Ctwr[irow][jcolumn] >> 2;
86 CgBlk[irow][jcolumn] = CgBlk[irow][jcolumn] >> 2;
87
88 }
89 }
90
91
92 //FPGA A observables
93 int A_MHT_x = 0x0;
94 int A_MHT_y = 0x0;
95 int A_MST_x = 0x0;
96 int A_MST_y = 0x0;
97 int A_MET_x = 0x0;
98 int A_MET_y = 0x0;
99
100 int A_eth = 0x0;
101 int A_ets = 0x0;
102 int A_etw = 0x0;
103
104 //FPGA B observables
105 int B_MHT_x = 0x0;
106 int B_MHT_y = 0x0;
107 int B_MST_x = 0x0;
108 int B_MST_y = 0x0;
109 int B_MET_x = 0x0;
110 int B_MET_y = 0x0;
111
112 int B_eth = 0x0;
113 int B_ets = 0x0;
114 int B_etw = 0x0;
115
116 //FPGA C observables
117 int C_MHT_x = 0x0;
118 int C_MHT_y = 0x0;
119 int C_MST_x = 0x0;
120 int C_MST_y = 0x0;
121 int C_MET_x = 0x0;
122 int C_MET_y = 0x0;
123
124 int C_eth = 0x0;
125 int C_ets = 0x0;
126 int C_etw = 0x0;
127
128 //Global observables
129 int MHT_x = 0x0;
130 int MHT_y = 0x0;
131 int MST_x = 0x0;
132 int MST_y = 0x0;
133 int MET_x = 0x0;
134 int MET_y = 0x0;
135
136 int ETH = 0x0;
137 int ETS = 0x0;
138 int ETW = 0x0;
139
140 int total_sumEt = 0x0;
141 int MET = 0x0;
142
143
144 // will need to hard code etFPGA ,a's and b's
145 int etBprime =0;
146
147 metFPGA(0, Ascaled, AgBlk, m_gBlockthresholdA, m_aFPGA_A, m_bFPGA_A, A_MHT_x, A_MHT_y, A_MST_x, A_MST_y, A_MET_x, A_MET_y);
148 if (SumETfast) etFastFPGA(0, Ascaled, AgBlk, m_gBlockthresholdA, m_aFPGA_A, etBprime, A_eth, A_ets, A_etw);
149 else etFPGA(0, Ascaled, AgBlk, m_gBlockthresholdA, m_aFPGA_A, etBprime, A_eth, A_ets, A_etw);
150
151 metFPGA(1, Bscaled, BgBlk, m_gBlockthresholdB, m_aFPGA_B, m_bFPGA_B, B_MHT_x, B_MHT_y, B_MST_x, B_MST_y, B_MET_x, B_MET_y);
152 if (SumETfast) etFastFPGA(1, Bscaled, BgBlk, m_gBlockthresholdB, m_aFPGA_B, etBprime, B_eth, B_ets, B_etw);
153 else etFPGA(1, Bscaled, BgBlk, m_gBlockthresholdB, m_aFPGA_B, etBprime, B_eth, B_ets, B_etw);
154
155 metFPGA(2, Cscaled, CgBlk, m_gBlockthresholdC, m_aFPGA_C, m_bFPGA_C, C_MHT_x, C_MHT_y, C_MST_x, C_MST_y, C_MET_x, C_MET_y);
156 if (SumETfast) etFastFPGA(2, Cscaled, CgBlk, m_gBlockthresholdC, m_aFPGA_C, etBprime, C_eth, C_ets, C_etw);
157 else etFPGA(2, Cscaled, CgBlk, m_gBlockthresholdC, m_aFPGA_C, etBprime, C_eth, C_ets, C_etw);
158
159 metTotal(A_MHT_x, A_MHT_y, B_MHT_x, B_MHT_y, C_MHT_x, C_MHT_y, MHT_x, MHT_y);
160 metTotal(A_MST_x, A_MST_y, B_MST_x, B_MST_y, C_MST_x, C_MST_y, MST_x, MST_y);
161 metTotal(A_MET_x, A_MET_y, B_MET_x, B_MET_y, C_MET_x, C_MET_y, MET_x, MET_y);
162
163 etTotal(A_eth, B_eth, C_eth, ETH);
164 etTotal(A_ets, B_ets, C_ets, ETS);
165 etTotal(A_etw, B_etw, C_etw, ETW);
166 total_sumEt = ETW;
167
168 // components should all be less than 12 bits at this point with 200 MeV LSB
169 int MET2 = MET_x * MET_x + MET_y * MET_y;
170
171 if (MET2 > 0x0FFFFFF) {
172 MET = 0x000FFF;
173 } else {
174 // repeat the byte stream converter calculation here -- not what the hardware actually does
175 MET = std::sqrt(MET2);
176
177
178 // best guess at current hardware. Note that this is computed in the bytestream converter
179 // take most signficant 12 bits
180 //int MET12 = MET2 >> 12;
181 // simulate the look up -- only 6 most signficant bits currently set -- to be checked
182 //MET = ( (int)(std::sqrt(MET12)) << 6) & 0x00000FC0 ;
183 }
184
185
186 //Define a vector to be filled with all the TOBs of one event
187 std::vector<std::unique_ptr<gFEXJwoJTOB>> tobs_v;
188 tobs_v.resize(4);
189
190
191 // fill in TOBs
192 // The order of the TOBs is given according to the TOB ID (TODO: check how it's done in fw)
193
194 // First TOB is (MET, SumEt)
195 outTOB[0] = (total_sumEt & 0x00000FFF) << 0; //set the Quantity2 to the corresponding slot (LSB)
196 outTOB[0] = outTOB[0] | (MET & 0x00000FFF) << 12;//Quantity 1 (in bit number 12)
197 if (total_sumEt != 0) outTOB[0] = outTOB[0] | 0x00000001 << 24;//Status bit for Quantity 2 (0 if quantity is null)
198 if (MET != 0) outTOB[0] = outTOB[0] | 0x00000001 << 25;//Status bit for Quantity 1 (0 if quantity is null)
199 outTOB[0] = outTOB[0] | (1 & 0x0000001F) << 26;//TOB ID is 1 for scalar values (5 bits starting at 26)
200
201 // std::cout << "DMS MET " << std::hex << MET << " total_sumEt " << total_sumEt << std::endl << std::dec;
202
203// Second TOB is (MET_x, MET_y)
204 outTOB[1] = (MET_y & 0x00000FFF) << 0; //set the Quantity2 to the corresponding slot (LSB)
205 outTOB[1] = outTOB[1] | (MET_x & 0x00000FFF) << 12;//Quantity 1 (in bit number 12)
206 if (MET_y != 0) outTOB[1] = outTOB[1] | 0x00000001 << 24;//Status bit for Quantity 2 (0 if quantity is null)
207 if (MET_x != 0) outTOB[1] = outTOB[1] | 0x00000001 << 25;//Status bit for Quantity 1 (0 if quantity is null)
208 outTOB[1] = outTOB[1] | (2 & 0x0000001F) << 26;//TOB ID is 2 for MET_x, MET_y (5 bits starting at 26)
209
210// Third TOB is hard components (MHT_x, MHT_y)
211 outTOB[2] = (MHT_y & 0x00000FFF) << 0; //set the Quantity2 to the corresponding slot (LSB)
212 outTOB[2] = outTOB[2] | (MHT_x & 0x00000FFF) << 12;//Quantity 1 (in bit number 12)
213 if (MHT_y != 0) outTOB[2] = outTOB[2] | 0x00000001 << 24;//Status bit for Quantity 2 (0 if quantity is null)
214 if (MHT_x != 0) outTOB[2] = outTOB[2] | 0x00000001 << 25;//Status bit for Quantity 1 (0 if quantity is null)
215 outTOB[2] = outTOB[2] | (3 & 0x0000001F) << 26;//TOB ID is 3 for hard components (5 bits starting at 26)
216
217 // Fourth TOB is hard components (MST_x, MST_y)
218 outTOB[3] = (MST_y & 0x00000FFF) << 0; //set the Quantity2 to the corresponding slot (LSB)
219 outTOB[3] = outTOB[3] | (MST_x & 0x00000FFF) << 12;//Quantity 1 (in bit number 12)
220 if (MST_y != 0) outTOB[3] = outTOB[3] | 0x00000001 << 24;//Status bit for Quantity 2 (0 if quantity is null)
221 if (MST_x != 0) outTOB[3] = outTOB[3] | 0x00000001 << 25;//Status bit for Quantity 1 (0 if quantity is null)
222 outTOB[3] = outTOB[3] | (4 & 0x0000001F) << 26;//TOB ID is 4 for soft components (5 bits starting at 26)
223
224
225 tobs_v[0] = std::make_unique<gFEXJwoJTOB>();
226 tobs_v[0]->setWord(outTOB[0]);
227 tobs_v[0]->setQuantity1(MET);
228 tobs_v[0]->setQuantity2(total_sumEt);
229 tobs_v[0]->setSaturation(0); //Always 0 for now, need a threshold?
230 tobs_v[0]->setTobID(1);
231 if( MET != 0 ) tobs_v[0]->setStatus1(1);
232 else tobs_v[0]->setStatus1(0);
233 if(total_sumEt!= 0) tobs_v[0]->setStatus2(1);
234 else tobs_v[0]->setStatus2(0);
235
236 tobs_v[1] = std::make_unique<gFEXJwoJTOB>();
237 tobs_v[1]->setWord(outTOB[1]);
238 tobs_v[1]->setQuantity1(MET_x);
239 tobs_v[1]->setQuantity2(MET_y);
240 tobs_v[1]->setSaturation(0); //Always 0 for now, need a threshold?
241 tobs_v[1]->setTobID(2);
242 if( MET_x != 0 ) tobs_v[1]->setStatus1(1);
243 else tobs_v[1]->setStatus1(0);
244 if(MET_y!= 0) tobs_v[1]->setStatus2(1);
245 else tobs_v[1]->setStatus2(0);
246
247 tobs_v[2] = std::make_unique<gFEXJwoJTOB>();
248 tobs_v[2]->setWord(outTOB[2]);
249 tobs_v[2]->setQuantity1(MHT_x);
250 tobs_v[2]->setQuantity2(MHT_y);
251 tobs_v[2]->setSaturation(0); //Always 0 for now, need a threshold?
252 tobs_v[2]->setTobID(3);
253 if( MHT_x != 0 ) tobs_v[2]->setStatus1(1);
254 else tobs_v[2]->setStatus1(0);
255 if(MHT_y!= 0) tobs_v[2]->setStatus2(1);
256 else tobs_v[2]->setStatus2(0);
257
258 tobs_v[3] = std::make_unique<gFEXJwoJTOB>();
259 tobs_v[3]->setWord(outTOB[3]);
260 tobs_v[3]->setQuantity1(MST_x);
261 tobs_v[3]->setQuantity2(MST_y);
262 tobs_v[3]->setSaturation(0); //Always 0 for now, need a threshold?
263 tobs_v[3]->setTobID(4);
264 if( MST_x != 0 ) tobs_v[3]->setStatus1(1);
265 else tobs_v[2]->setStatus1(0);
266 if(MST_y!= 0) tobs_v[3]->setStatus2(1);
267 else tobs_v[3]->setStatus2(0);
268
269
270 return tobs_v;
271
272}
273
274
275
276void gFEXJwoJAlgo::gBlockAB(const gTowersType & twrs, gTowersType & gBlkSum, gTowersType & hasSeed, int seedThreshold) const {
277
278 int rows = twrs.size();
279 int cols = twrs[0].size();
280 for( int irow = 0; irow < rows; irow++ ){
281 for(int jcolumn = 0; jcolumn<cols; jcolumn++){
282 // zero jet sum here
283 gBlkSum[irow][jcolumn] = 0;
284 int krowUp = (irow + 1)%32;
285 int krowDn = (irow - 1 +32)%32;
286 if( (jcolumn == 0) || (jcolumn == 6) ) {
287 //left edge case
288 gBlkSum[irow][jcolumn] =
289 twrs[irow][jcolumn] + twrs[krowUp][jcolumn] + twrs[krowDn][jcolumn] +
290 twrs[irow][jcolumn+1] + twrs[krowUp][jcolumn+1] + twrs[krowDn][jcolumn+1];
291 } else if( (jcolumn == 5) || (jcolumn == 11)) {
292 // right edge case
293 gBlkSum[irow][jcolumn] =
294 twrs[irow][jcolumn] + twrs[krowUp][jcolumn] + twrs[krowDn][jcolumn] +
295 twrs[irow][jcolumn-1] + twrs[krowUp][jcolumn-1] + twrs[krowDn][jcolumn-1];
296 } else{
297 // normal case
298 gBlkSum[irow][jcolumn] =
299 twrs[irow][jcolumn] + twrs[krowUp][jcolumn] + twrs[krowDn][jcolumn] +
300 twrs[irow][jcolumn-1] + twrs[krowUp][jcolumn-1] + twrs[krowDn][jcolumn-1] +
301 twrs[irow][jcolumn+1] + twrs[krowUp][jcolumn+1] + twrs[krowDn][jcolumn+1];
302 }
303
304 if( gBlkSum[irow][jcolumn] > seedThreshold) {
305 hasSeed[irow][jcolumn] = 1;
306 } else {
307 hasSeed[irow][jcolumn] = 0;
308 }
309
310 if ( gBlkSum[irow][jcolumn] < 0 )
311 gBlkSum[irow][jcolumn] = 0;
312
313 // was bits 11+3 downto 3, now is 11 downto 0
314 if ( gBlkSum[irow][jcolumn] > FEXAlgoSpaceDefs::gBlockMax ) {
315 gBlkSum[irow][jcolumn] = FEXAlgoSpaceDefs::gBlockMax;
316 }
317 }
318 }
319}
320
321
322void gFEXJwoJAlgo::metFPGA(int FPGAnum, const gTowersType& twrs,
323 const gTowersType & gBlkSum, int gBlockthreshold,
324 int aFPGA, int bFPGA,
325 int & MHT_x, int & MHT_y,
326 int & MST_x, int & MST_y,
327 int & MET_x, int & MET_y) const {
328
329 gBlockthreshold = gBlockthreshold * 200 / 800; //gBlockthreshold is provided in counts with a resolution of 200 MeV, but here needs to be applied with a resolution of 800 GeV
330 // in the RTL code these are 19+ 5 = 24 bits
331 int64_t h_tx_hi = 0;
332 int64_t h_ty_hi = 0;
333 int64_t h_tx_lw = 0;
334 int64_t h_ty_lw = 0;
335
336 int64_t e_tx_hi = 0;
337 int64_t e_ty_hi = 0;
338 int64_t e_tx_lw = 0;
339 int64_t e_ty_lw = 0;
340
341 for( int irow = 0; irow < FEXAlgoSpaceDefs::ABCrows; irow++ ){
342 for(int jcolumn = 6; jcolumn<12; jcolumn++){
343 if( FPGAnum == 2){
344 int frow = 2*(irow/2) + 1;
345
346 if(gBlkSum[irow][jcolumn] > gBlockthreshold){
347 h_tx_hi += (twrs[irow][jcolumn])*(cosLUT(frow, 5));
348 h_ty_hi += (twrs[irow][jcolumn])*(sinLUT(frow, 5));
349 } else {
350 e_tx_hi += (twrs[irow][jcolumn])*(cosLUT(frow, 5));
351 e_ty_hi += (twrs[irow][jcolumn])*(sinLUT(frow, 5));
352 }
353
354 } else {
355
356 if(gBlkSum[irow][jcolumn] > gBlockthreshold){
357 h_tx_hi += (twrs[irow][jcolumn])*(cosLUT(irow, 5));
358 h_ty_hi += (twrs[irow][jcolumn])*(sinLUT(irow, 5));
359 } else {
360 e_tx_hi += (twrs[irow][jcolumn])*(cosLUT(irow, 5));
361 e_ty_hi += (twrs[irow][jcolumn])*(sinLUT(irow, 5));
362 }
363 }
364 }
365
366 for(int jcolumn = 0; jcolumn<6; jcolumn++){
367 if( FPGAnum == 2){
368 int frow = 2*(irow/2) + 1;
369
370 if(gBlkSum[irow][jcolumn] > gBlockthreshold){
371 h_tx_lw += (twrs[irow][jcolumn])*(cosLUT(frow, 5));
372 h_ty_lw += (twrs[irow][jcolumn])*(sinLUT(frow, 5));
373 } else{
374 e_tx_lw += (twrs[irow][jcolumn])*(cosLUT(frow, 5));
375 e_ty_lw += (twrs[irow][jcolumn])*(sinLUT(frow, 5));
376 }
377 } else {
378
379 if(gBlkSum[irow][jcolumn] > gBlockthreshold){
380 h_tx_lw += (twrs[irow][jcolumn])*(cosLUT(irow, 5));
381 h_ty_lw += (twrs[irow][jcolumn])*(sinLUT(irow, 5));
382 } else {
383 e_tx_lw += (twrs[irow][jcolumn])*(cosLUT(irow, 5));
384 e_ty_lw += (twrs[irow][jcolumn])*(sinLUT(irow, 5));
385 }
386 }
387 }
388 }
389
390 // According to https://gitlab.cern.ch/atlas-l1calo/gfex/firmware/-/issues/406#note_5662344
391 // there is no division -- so LSB is indeed 25 MeV
392
393 //but later changed to 200 MeV so divide by 8
394
395 long int fMHT_x = (h_tx_hi + h_tx_lw) ;
396 long int fMHT_y = (h_ty_hi + h_ty_lw) ;
397 long int fMST_x = (e_tx_hi + e_tx_lw) ;
398 long int fMST_y = (e_ty_hi + e_ty_lw) ;
399
400 MHT_x = (h_tx_hi + h_tx_lw) >> 3;
401 MHT_y = (h_ty_hi + h_ty_lw) >> 3;
402 MST_x = (e_tx_hi + e_tx_lw) >> 3;
403 MST_y = (e_ty_hi + e_ty_lw) >> 3;
404
405 // a and b coffecients are 10 bits
406 // multiplication has an addtional 2^10
407 // constant JWJ_OW : integer := 35;--Out width
408 // values are 35 bits long and top 16 bits are taken -- so divide by 2^19
409 // 2^10/2^19 = 1/2^9 = 1/512
410
411 long int fMET_x = ( aFPGA * (fMHT_x) + bFPGA * (fMST_x) ) >> 13 ;
412 long int fMET_y = ( aFPGA * (fMHT_y) + bFPGA * (fMST_y) ) >> 13 ;
413
414 MET_x = fMET_x;
415 MET_y = fMET_y;
416
417}
418
419void gFEXJwoJAlgo::etFPGA(int FPGAnum, const gTowersType& twrs, gTowersType &gBlkSum,
420 int gBlockthreshold, int A, int B, int &eth, int &ets, int &etw) const {
421
422 gBlockthreshold = gBlockthreshold * 200 / 800; //gBlockthreshold is provided in counts with a resolution of 200 MeV, but here needs to be applied with a resolution of 800 GeV
423
424 int64_t ethard_hi = 0;
425 int64_t etsoft_hi = 0;
426 int64_t ethard_lo = 0;
427 int64_t etsoft_lo = 0;
428
429 int64_t ethard = 0.0;
430 int64_t etsoft = 0.0;
431
432 int multiplicitiveFactor = 0;
433
434 if(FPGAnum < 2 ) {
435 multiplicitiveFactor = cosLUT(0, 5);
436 } else{
437 multiplicitiveFactor = cosLUT(1, 5);
438 }
439
440// firmware treats upper and lower columns differnetly
441
442 for( int irow = 0; irow < FEXAlgoSpaceDefs::ABCrows; irow++ ){
443 for(int jcolumn = 0; jcolumn<6; jcolumn++){
444 if(gBlkSum[irow][jcolumn] > gBlockthreshold){
445 ethard_lo = ethard_lo + twrs[irow][jcolumn]*multiplicitiveFactor;
446 } else {
447 etsoft_lo = etsoft_lo + twrs[irow][jcolumn]*multiplicitiveFactor;
448 }
449 }
450 }
451
452 for( int irow = 0; irow < FEXAlgoSpaceDefs::ABCrows; irow++ ){
453 for(int jcolumn = 6; jcolumn<12; jcolumn++){
454 if(gBlkSum[irow][jcolumn] > gBlockthreshold){
455 ethard_hi = ethard_hi + twrs[irow][jcolumn]*multiplicitiveFactor;
456 } else {
457 etsoft_hi = etsoft_hi + twrs[irow][jcolumn]*multiplicitiveFactor;
458 }
459 }
460 }
461
462 ethard = ethard_hi + ethard_lo;
463
464 etsoft = etsoft_hi + etsoft_lo;
465
466
467 int64_t etsum_hi = ethard_hi*A + etsoft_hi*B ;
468 if ( etsum_hi < 0 ) etsum_hi = 0;
469
470 int64_t etsum_lo = ethard_lo*A + etsoft_lo*B ;
471 if ( etsum_lo < 0 ) etsum_lo = 0;
472
473 int64_t etsum = etsum_hi + etsum_lo;
474
475
476 // convert 200 MeV LSB here
477 eth = ethard>>3;
478 ets = etsoft>>3;
479 etw = (etsum >>13 ) ;
480
481 if( etw < 0 ) etw = 0;
482 // max value is 15 bits with 800 MeV LSB -- so 17 bits here
483 if( etw > 0X001FFFF ) etw = 0X001FFFF ;
484
485
486 if(msgLvl(MSG::DEBUG)) {
487 std::cout << "DMS FPGA gTEJWOJ " << std::hex << FPGAnum << "et sum hard " << eth << "etsum soft" << ets << " A " << A << " B " << B << " weighted term " << etw << std::endl << std::dec;
488 }
489}
490
491void gFEXJwoJAlgo::etFastFPGA(int FPGAnum, const gTowersType& twrs, gTowersType &gBlkSum,
492 int gBlockthreshold, int A, int B, int &eth, int &ets, int &etw) const {
493
494 gBlockthreshold = gBlockthreshold * 200 / 800; //gBlockthreshold is provided in counts with a resolution of 200 MeV, but here needs to be applied with a resolution of 800 GeV
495 int64_t ethard_hi = 0;
496 int64_t etsoft_hi = 0;
497 int64_t ethard_lo = 0;
498 int64_t etsoft_lo = 0;
499
500 int64_t ethard = 0.0;
501 int64_t etsoft = 0.0;
502
503 // firmware treats upper and lower columns differently
504 for( int irow = 0; irow < FEXAlgoSpaceDefs::ABCrows; irow++ ){
505 for(int jcolumn = 0; jcolumn<6; jcolumn++){
506 if(gBlkSum[irow][jcolumn] > gBlockthreshold){
507 ethard_lo = ethard_lo + twrs[irow][jcolumn];
508 }
509 else {
510 etsoft_lo = etsoft_lo + twrs[irow][jcolumn];
511 }
512 }
513 }
514
515 for( int irow = 0; irow < FEXAlgoSpaceDefs::ABCrows; irow++ ){
516 for(int jcolumn = 6; jcolumn<12; jcolumn++){
517 if(gBlkSum[irow][jcolumn] > gBlockthreshold){
518 ethard_hi = ethard_hi + twrs[irow][jcolumn];
519 }
520 else {
521 etsoft_hi = etsoft_hi + twrs[irow][jcolumn];
522 }
523 }
524 }
525
526 ethard = ethard_hi + ethard_lo;
527 etsoft = etsoft_hi + etsoft_lo;
528
529 // convert 200 MeV LSB here
530 eth = ethard;
531 ets = etsoft; // Keep for reference -- not used in etFast
532 etw = ethard; // For EtFast the weighted term is just the hard term
533
534 // 16 bits signed, set max and min
535 if( etw < -32768 ) etw = -32768;
536 if( etw > 32767 ) etw = 32767;
537
538 if(msgLvl(MSG::DEBUG)) {
539 std::cout << "DMS FPGA gTEJWOJ " << std::hex << FPGAnum << "et sum hard " << eth << "etsum soft" << ets << " A " << A << " B " << B << " weighted term " << etw << std::endl << std::dec;
540 }
541}
542
543void gFEXJwoJAlgo::metTotal(int A_MET_x, int A_MET_y,
544 int B_MET_x, int B_MET_y,
545 int C_MET_x, int C_MET_y,
546 int & MET_x, int & MET_y) const {
547
548
549 MET_x = A_MET_x + B_MET_x + C_MET_x;
550 MET_y = A_MET_y + B_MET_y+ C_MET_y;
551
552 // Truncation of the result, as the individual quantities are 16 bits, while the TOB field is 12 bits
553 // MET_x = MET_x >> 4;
554 // MET_y = MET_y >> 4;
555
556 if (MET_x < -0x000800) MET_x = -0x000800; //-2048
557 if (MET_y < -0x000800) MET_y = -0x000800; //-2048
558
559 if (MET_x > 0x0007FF) MET_x = 0x0007FF; //2047
560 if (MET_y > 0x0007FF) MET_y = 0x0007FF; //2047
561
562}
563
565 int B_ET,
566 int C_ET,
567 int & ET ) const {
568
569 // leave at 200 MeV scale
570 if (A_ET < 0 ) A_ET = 0;
571 if (B_ET < 0 ) B_ET = 0;
572 if (C_ET < 0 ) C_ET = 0;
573
574 ET = (A_ET + B_ET + C_ET) ;
575
576 // main value of ET is always positive
577 if( ET > 0x0000FFF) ET = 0x0000FFF;
578
579}
580
581//----------------------------------------------------------------------------------
582// bitwise simulation of sine LUT in firmware
583//----------------------------------------------------------------------------------
584float gFEXJwoJAlgo::sinLUT(unsigned int phiIDX, unsigned int aw) const
585{
586 float c = static_cast<float>(phiIDX)/std::pow(2,aw);
587 float rad = (2*M_PI) *c;
588 float rsin = std::sin(rad);
589 int isin = std::round(rsin*(std::pow(2,aw) - 1));
590 return isin;
591
592}
593
594//----------------------------------------------------------------------------------
595// bitwise simulation cosine LUT in firmware
596//----------------------------------------------------------------------------------
597float gFEXJwoJAlgo::cosLUT(unsigned int phiIDX, unsigned int aw) const
598{
599 float c = static_cast<float>(phiIDX)/std::pow(2,aw);
600 float rad = (2*M_PI) *c;
601 float rcos = std::cos(rad);
602 int icos = std::round(rcos*(std::pow(2,aw) - 1));
603 return icos;
604}
605
606} // namespace LVL1
#define M_PI
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
bool msgLvl(const MSG::Level lvl) const
static constexpr int ABCrows
static constexpr int gBlockMax
static constexpr int ABcolumns
void etFastFPGA(int FPGAnum, const gTowersType &twrs, gTowersType &gBlkSum, int gBlockthreshold, int A, int B, int &eth, int &ets, int &etw) const
void metFPGA(int FPGAnum, const gTowersType &twrs, const gTowersType &gBlkSum, int gBlockthreshold, int aFPGA, int bFPGA, int &MHT_x, int &MHT_y, int &MST_x, int &MST_y, int &MET_x, int &MET_y) const
virtual std::vector< std::unique_ptr< gFEXJwoJTOB > > jwojAlgo(const gTowersType &Atwr, const gTowersType &Btwr, const gTowersType &Ctwr, std::array< int32_t, 4 > &outTOB) const override
void etFPGA(int FPGAnum, const gTowersType &twrs, gTowersType &gBlkSum, int gBlockthreshold, int A, int B, int &eth, int &ets, int &etw) const
virtual StatusCode initialize() override
standard Athena-Algorithm method
float cosLUT(unsigned int phiIDX, unsigned int aw) const
void metTotal(int A_MET_x, int A_MET_y, int B_MET_x, int B_MET_y, int C_MET_x, int C_MET_y, int &MET_x, int &MET_y) const
float sinLUT(unsigned int phiIDX, unsigned int aw) const
void gBlockAB(const gTowersType &twrs, gTowersType &gBlkSum, gTowersType &hasSeed, int seedThreshold) const
gFEXJwoJAlgo(const std::string &type, const std::string &name, const IInterface *parent)
Constructors.
void etTotal(int A_ET, int B_ET, int C_ET, int &ET) const
virtual void setAlgoConstant(int aFPGA_A, int bFPGA_A, int aFPGA_B, int bFPGA_B, int aFPGA_C, int bFPGA_C, int gXE_seedThrA, int gXE_seedThrB, int gXE_seedThrC) override
eFexTowerBuilder creates xAOD::eFexTowerContainer from supercells (LATOME) and triggerTowers (TREX) i...
std::array< std::array< int, 12 >, 32 > gTowersType
Definition IgFEXFPGA.h:25
Definition MET.py:1
hold the test vectors and ease the comparison