ATLAS Offline Software
Loading...
Searching...
No Matches
gFEXaltMetAlgo.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4//***************************************************************************
5// gFEXaltMetAlgo - Noise cut and Rho+RMS algorithm for gFEX MET
6// -------------------
7// begin : 31 03 2022
8// email : cecilia.tosciri@cern.ch
9//***************************************************************************
10
11#include <cmath>
12#include <vector>
13
14#include "gFEXaltMetAlgo.h"
16#include "L1CaloFEXSim/gTower.h"
17
18namespace LVL1 {
19
20gFEXaltMetAlgo::gFEXaltMetAlgo(const std::string& type, const std::string& name, const IInterface* parent):
21base_class(type, name, parent)
22{}
23
25
26 return StatusCode::SUCCESS;
27
28}
29
30
31void gFEXaltMetAlgo::setAlgoConstant(std::vector<int>&& A_thr,
32 std::vector<int>&& B_thr,
33 const int rhoPlusThr) {
34 m_etaThr[0] = std::move(A_thr);
35 m_etaThr[1] = std::move(B_thr);
36 m_rhoPlusThr = rhoPlusThr;
37}
38
39
40
42 std::array<uint32_t, 4> & outTOB) const {
43
44 //FPGA A observables
45 int A_MET_x_nc = 0x0;
46 int A_MET_y_nc = 0x0;
47 int A_MET_x_rms = 0x0;
48 int A_MET_y_rms = 0x0;
49
50 int A_sumEt_nc = 0x0;
51 int A_sumEt_rms = 0x0;
52
53 //FPGA B observables
54 int B_MET_x_nc = 0x0;
55 int B_MET_y_nc = 0x0;
56 int B_MET_x_rms = 0x0;
57 int B_MET_y_rms = 0x0;
58
59 int B_sumEt_nc = 0x0;
60 int B_sumEt_rms = 0x0;
61
62 //Global observables
63 int MET_x_nc = 0x0;
64 int MET_y_nc = 0x0;
65 int MET_nc = 0x0;
66 int MET_x_rms = 0x0;
67 int MET_y_rms = 0x0;
68 int MET_rms = 0x0;
69
70
71 int total_sumEt_nc = 0x0;
72 int total_sumEt_rms = 0x0;
73
74 metFPGA(Atwr, A_MET_x_nc, A_MET_y_nc, 0);
75 metFPGA(Btwr, B_MET_x_nc, B_MET_y_nc, 1);
76
77 metTotal(A_MET_x_nc, A_MET_y_nc, B_MET_x_nc, B_MET_y_nc, MET_x_nc, MET_y_nc, MET_nc);
78
79 int A_rho{get_rho(Atwr)};
80 int B_rho{get_rho(Btwr)};
81 int A_sigma{3*get_sigma(Atwr)};
82 int B_sigma{3*get_sigma(Btwr)};
83
84 rho_MET(Atwr, A_MET_x_rms, A_MET_y_rms, A_rho, A_sigma);
85 rho_MET(Btwr, B_MET_x_rms, B_MET_y_rms, B_rho, B_sigma);
86
87 metTotal(A_MET_x_rms, A_MET_y_rms, B_MET_x_rms, B_MET_y_rms, MET_x_rms, MET_y_rms, MET_rms);
88
89 A_sumEt_nc = sumEtFPGAnc(Atwr, 0);
90 B_sumEt_nc = sumEtFPGAnc(Btwr, 1);
91 total_sumEt_nc = sumEt(A_sumEt_nc, B_sumEt_nc);
92 total_sumEt_nc = total_sumEt_nc/4;
93
94 A_sumEt_rms = sumEtFPGArms(Atwr, A_sigma);
95 B_sumEt_rms = sumEtFPGArms(Btwr, B_sigma);
96 total_sumEt_rms = sumEt(A_sumEt_rms, B_sumEt_rms);
97 total_sumEt_rms = total_sumEt_rms/4;
98 //Define a vector to be filled with all the TOBs of one event
99
100 //TOB order
101 // 1) MET_x | MET_y <- ncMET
102 // 2) MET_x | MET_y <- rms
103 // 3) MET | sumET <- ncMET
104 // 4) MET | sumET <- rms
105
106 // fill in TOBs
107 // The order of the TOBs is given according to the TOB ID (TODO: check how it's done in fw)
108
109 // First TOB is (MET, SumEt)
110 outTOB[0] = (MET_y_nc& 0x00000FFF) << 0; //set the Quantity2 to the corresponding slot (LSB)
111 outTOB[0] = outTOB[0] | (MET_x_nc & 0x00000FFF) << 12;//Quantity 1 (in bit number 12)
112 if (MET_y_nc != 0) outTOB[0] = outTOB[0] | 0x00000001 << 24;//Status bit for Quantity 2 (0 if quantity is null)
113 if (MET_x_nc != 0) outTOB[0] = outTOB[0] | 0x00000001 << 25;//Status bit for Quantity 1 (0 if quantity is null)
114 outTOB[0] = outTOB[0] | (2 & 0x0000001F) << 26;//TOB ID temporary set to 2 according to JwoJ convention (need updates in EDM)
115
116// Second TOB is (MET_x, MET_y)
117 outTOB[1] = (MET_y_rms& 0x00000FFF) << 0; //set the Quantity2 to the corresponding slot (LSB)
118 outTOB[1] = outTOB[1] | (MET_x_rms & 0x00000FFF) << 12;//Quantity 1 (in bit number 12)
119 if (MET_y_rms != 0) outTOB[1] = outTOB[1] | 0x00000001 << 24;//Status bit for Quantity 2 (0 if quantity is null)
120 if (MET_x_rms != 0) outTOB[1] = outTOB[1] | 0x00000001 << 25;//Status bit for Quantity 1 (0 if quantity is null)
121 outTOB[1] = outTOB[1] | (2 & 0x0000001F) << 26;//TOB ID temporary set to 2 according to JwoJ convention (need updates in EDM)
122
123// Third TOB is hard components (MHT_x, MHT_y)
124 outTOB[2] = (total_sumEt_nc& 0x00000FFF) << 0; //set the Quantity2 to the corresponding slot (LSB)
125 outTOB[2] = outTOB[2] | (MET_nc & 0x00000FFF) << 12;//Quantity 1 (in bit number 12)
126 if (total_sumEt_nc != 0) outTOB[2] = outTOB[2] | 0x00000001 << 24;//Status bit for Quantity 2 (0 if quantity is null)
127 if (MET_nc != 0) outTOB[2] = outTOB[2] | 0x00000001 << 25;//Status bit for Quantity 1 (0 if quantity is null)
128 outTOB[2] = outTOB[2] | (1 & 0x0000001F) << 26;//TOB ID temporary set to 1 according to JwoJ convention (need updates in EDM)
129
130 // Fourth TOB is hard components (MST_x, MST_y)
131 outTOB[3] = (total_sumEt_rms& 0x00000FFF) << 0; //set the Quantity2 to the corresponding slot (LSB)
132 outTOB[3] = outTOB[3] | (MET_rms & 0x00000FFF) << 12;//Quantity 1 (in bit number 12)
133 if (total_sumEt_rms != 0) outTOB[3] = outTOB[3] | 0x00000001 << 24;//Status bit for Quantity 2 (0 if quantity is null)
134 if (MET_rms != 0) outTOB[3] = outTOB[3] | 0x00000001 << 25;//Status bit for Quantity 1 (0 if quantity is null)
135 outTOB[3] = outTOB[3] | (1 & 0x0000001F) << 26;//TOB ID temporary set to 1 according to JwoJ convention (need updates in EDM)
136
137
138}
139
140
141
142void gFEXaltMetAlgo::metFPGA(const gTowersCentral &twrs, int & MET_x, int & MET_y, const unsigned short FPGA_NO) const {
143
144 int rows = twrs.size();
145 int cols = twrs[0].size();
146 for( int irow = 0; irow < rows; irow++ ){
147 for(int jcolumn = 0; jcolumn<cols; jcolumn++){
148 bool filter{twrs[irow][jcolumn] > m_etaThr[FPGA_NO][jcolumn]};
149 MET_x += filter ? (twrs[irow][jcolumn])*cosLUT(irow, 5) : 0;
150 MET_y += filter ? (twrs[irow][jcolumn])*sinLUT(irow, 5) : 0;
151
152 }
153 }
154}
155
156
157inline void gFEXaltMetAlgo::metTotal(const int A_MET_x, const int A_MET_y,
158 const int B_MET_x, const int B_MET_y,
159 int & MET_x, int & MET_y, int & MET) const {
160
161 MET_x = A_MET_x + B_MET_x;
162 MET_y = A_MET_y + B_MET_y;
163
164 if (MET_x < -0x0007FF) MET_x = -0x0007FF;
165 if (MET_y < -0x0007FF) MET_y = -0x0007FF;
166
167 if (MET_x > 0x0007FF) MET_x = 0x0007FF;
168 if (MET_y > 0x0007FF) MET_y = 0x0007FF;
169
170 int MET2 = MET_x * MET_x + MET_y * MET_y;
171
172 if (MET2 > 0x000FFF) MET = 0x000FFF;
173 else if (MET2 < 0) MET = 0x000FFF;
174 else MET = std::sqrt(MET2);
175
176}
177
178//Function to calculate rho for the given set of gtowers
180 const int rows = twrs.size();
181 const int cols = twrs[0].size();
182 const int n{rows*cols};
183 float rho = 0;
184 for(int i = 0; i < rows; i++) {
185 for(int j = 0; j < cols; j++) {
186 rho += twrs[i][j] < m_rhoPlusThr ? twrs[i][j] : 0;
187 }
188 }
189 return rho/n;
190}
191
192//Function calculates standard deviation of the gtowers
194
195 int rows = twrs.size();
196 int cols = twrs[0].size();
197 const int n{rows*cols};
198 int sigma = 0;
199 for(int i = 0; i < rows; ++i) {
200 for(int j = 0; j < cols; ++j) {
201 const int towers{twrs[i][j]};
202 sigma += twrs[i][j] < m_rhoPlusThr ? towers*towers: 0;
203 }
204 }
205
206 return static_cast<int>(std::sqrt(sigma * 1. / n));
207
208}
209
210
211void gFEXaltMetAlgo::rho_MET(const gTowersCentral &twrs, int & MET_x, int & MET_y, const int rho, const int sigma) const {
212
213 int rows = twrs.size();
214 int cols = twrs[0].size();
215 for( int irow = 0; irow < rows; irow++ ){
216 for(int jcolumn = 0; jcolumn < cols; jcolumn++){
217 const int ET_gTower_sub{(twrs[irow][jcolumn] - rho) & 0xFFF};
218 const bool filter{ET_gTower_sub > sigma && !(ET_gTower_sub & 0x800)};
219 MET_x += filter ? (ET_gTower_sub)*cosLUT(irow, 5) : 0;
220 MET_y += filter ? (ET_gTower_sub)*sinLUT(irow, 5) : 0;
221
222 }
223 }
224}
225
226int gFEXaltMetAlgo::sumEtFPGAnc(const gTowersCentral &twrs, const unsigned short FPGA_NO) const {
227
228 int partial_sumEt = 0;
229 const int rows = twrs.size();
230 const int cols = twrs[0].size();
231 for( int irow = 0; irow < rows; irow++ ){
232 for(int jcolumn = 0; jcolumn<cols; jcolumn++){
233 partial_sumEt += twrs[irow][jcolumn] > m_etaThr[FPGA_NO][jcolumn] ? twrs[irow][jcolumn] : 0;
234 }
235 }
236 return partial_sumEt;
237}
238
239int gFEXaltMetAlgo::sumEtFPGArms(const gTowersCentral &twrs, const int sigma) const {
240
241 int partial_sumEt = 0;
242 const int rows = twrs.size();
243 const int cols = twrs[0].size();
244 for(int i{0}; i < rows; ++i){
245 for(int j{0}; j < cols; ++j) {
246 partial_sumEt += twrs[i][j] > sigma ? twrs[i][j] : 0;
247 }
248 }
249 return partial_sumEt;
250}
251
252
253int gFEXaltMetAlgo::sumEt(const int A_sumEt, const int B_sumEt) const {
254
255 return A_sumEt + B_sumEt;
256}
257
258
259//----------------------------------------------------------------------------------
260// bitwise simulation of sine LUT in firmware
261//----------------------------------------------------------------------------------
262float gFEXaltMetAlgo::sinLUT(const unsigned int phiIDX, const unsigned int aw) const
263{
264 float c = ((float)phiIDX)/std::pow(2,aw);
265 float rad = (2*M_PI) *c;
266 float rsin = std::sin(rad);
267 return rsin;
268}
269
270
271//----------------------------------------------------------------------------------
272// bitwise simulation cosine LUT in firmware
273//----------------------------------------------------------------------------------
274float gFEXaltMetAlgo::cosLUT(const unsigned int phiIDX, const unsigned int aw) const
275{
276 float c = ((float)phiIDX)/std::pow(2,aw);
277 float rad = (2*M_PI) *c;
278 float rcos = std::cos(rad);
279 return rcos;
280}
281
282
283} // namespace LVL1
#define M_PI
void metFPGA(const gTowersCentral &twrs, int &MET_x, int &MET_y, const unsigned short FPGA_NO) const
virtual void altMetAlgo(const gTowersCentral &Atwr, const gTowersCentral &Btwr, std::array< uint32_t, 4 > &outTOB) const override
gFEXaltMetAlgo(const std::string &type, const std::string &name, const IInterface *parent)
Constructor.
float cosLUT(const unsigned int phiIDX, const unsigned int aw) const
void metTotal(const int A_MET_x, const int A_MET_y, const int B_MET_x, const int B_MET_y, int &MET_x, int &MET_y, int &MET) const
float sinLUT(const unsigned int phiIDX, const unsigned int aw) const
virtual StatusCode initialize() override
standard Athena-Algorithm method
std::array< std::vector< int >, 2 > m_etaThr
virtual void setAlgoConstant(std::vector< int > &&A_thr, std::vector< int > &&B_thr, const int rhoPlusThr) override
int sumEtFPGArms(const gTowersCentral &twrs, const int sigma) const
int get_sigma(const gTowersCentral &twrs) const
int sumEtFPGAnc(const gTowersCentral &twrs, const unsigned short FPGA_NO) const
void rho_MET(const gTowersCentral &twrs, int &MET_x, int &MET_y, const int rho, const int sigma) const
int sumEt(const int A_sumEt, const int B_sumEt) const
int get_rho(const gTowersCentral &twrs) const
eFexTowerBuilder creates xAOD::eFexTowerContainer from supercells (LATOME) and triggerTowers (TREX) i...
std::array< std::array< int, 12 >, 32 > gTowersCentral
Definition MET.py:1