ATLAS Offline Software
Loading...
Searching...
No Matches
JetEnergyModuleKey.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 JetEnergyModuleKey.cpp - description
6 -------------------
7 begin : Tue Sep 26 2000
8 email : e.moyse@qmw.ac.uk
9 ***************************************************************************/
10
11#ifndef TRIGGERSPACE
12// running in Athena
15// #include "TrigT1CaloUtils/ICoordinate.h"
17
19
20#else
22#include "JetEnergyModuleKey.h"
23#include "TrigT1CaloDefs.h"
24#include "ICoordinate.h"
25#include "JetElementKey.h"
26#include "CoordToHardware.h"
27
28//
29#endif
30
31#include <math.h>
32
33namespace LVL1 {
34
37 if (m_debugModuleKey) std::cout << "JetEnergyModuleKey: m_debugModuleKey output turned on...."<<std::endl;
38}
39
42
44unsigned int JetEnergyModuleKey::jemKey(const xAOD::JetElement* jetElement) {
45 return key(jetElement->phi(), jetElement->eta());
46}
47
48
49
51unsigned int JetEnergyModuleKey::jemKey(unsigned int crate, unsigned int module) {
52 double eta = ((module%8)+0.5)*TrigT1CaloDefs::jemEtaSize - 3.2;
53 int quadrant = crate + ((module>7) ? 2 : 0);
54 double phi = (quadrant+0.5)*TrigT1CaloDefs::jemPhiSize;
55 return key(phi, eta);
56}
57
59std::vector<unsigned int> JetEnergyModuleKey::jeKeys(unsigned int crate, unsigned int module) {
60 std::vector<unsigned int> keys;
61
62 // Protection
63 if (crate > 1 || module > 15) return keys;
64
65 // First get coordinates
66 std::vector<Coordinate> coords = jeCoords(crate, module);
67
68 // Then calculate keys
69 JetElementKey get(0.,0.);
70 for (const Coordinate& c : coords) {
71 keys.push_back(get.jeKey(c));
72 }
73
74 return keys;
75}
76
77
79std::vector<Coordinate> JetEnergyModuleKey::jeCoords(unsigned int crate, unsigned int module) {
80 std::vector<Coordinate> coords;
81
82 // Protection
83 if (crate > 1 || module > 15) return coords;
84
90
91 int top = ( module < 8 ? 0 : 1 );
92 double phiStart = crate*TrigT1CaloDefs::jemPhiSize
95 double etaStart;
96 // Central JEMs are easy
97 if ( (module%8) > 0) {
98 etaStart = (module%8)*TrigT1CaloDefs::jemEtaSize - 3.1;
99 }
100 // -ive eta end
101 else {
103 }
104
105 // First coordinate for this JEM
106 Coordinate startCol(phiStart, etaStart);
107 // Nested incrementation of eta and phi
108 CoordToHardware cth;
109 JetElementKey get(0.,0.);
110 while (cth.jepModule( startCol ) == module) {
111 Coordinate next = get.getCentre(startCol);
112 while (cth.jepCrate( next ) == crate) {
113 coords.push_back(next);
114 next = get.upPhi(next);
115 }
116 startCol = get.rightEta(startCol);
117 if (startCol.eta() == TrigT1CaloDefs::RegionERROREtaCentre) break;
118 }
119
120 return coords;
121}
122
125 double phi=coord.phi();
126 double rowHeight=TrigT1CaloDefs::jemPhiSize/8.0;
127
128 unsigned int rowsPerJEM = 8u;
129 int rowNum=static_cast<int>(phi/rowHeight);
130 rowNum=rowNum%rowsPerJEM;
131 return rowNum;
132}
133
137 Coordinate coord(JE->phi(), JE->eta());
138 return row( coord );
139}
140
141
143double LVL1::JetEnergyModuleKey::rowPhiCoord(unsigned int rowNum, const Coordinate& jemCoord){
144 //not sure how to do this.
145 //should really reuse code, but above is tied to KU stuff. Probably way to go though.
146
147 ICoordinate* iCoord = convertCoordsToIntegers(jemCoord.phi(), jemCoord.eta());
148 int etaBin=0; unsigned int phiBin=0;
149 setBins(iCoord,phiBin,etaBin);
150 // now work out row coord.
151
152 if (rowNum>7){
153 std::cerr << "ERROR!!! JetEnergyModuleKey::rowPhiCoord exceeds allowed range. Row, "<<rowNum
154 <<" will be set to zero"<<std::endl;
155 rowNum=0;
156 }
157 double phi =( (static_cast<double>(phiBin))*TrigT1CaloDefs::jemPhiSize ); //bot of JEM
158 phi+= static_cast<double>( rowNum )*0.2+0.1; //mid of row
159 delete iCoord;
160 if (m_debugModuleKey){
161 if ( (phi<(jemCoord.phi()-0.8))||(phi>(jemCoord.phi()+0.8) ) ){
162 std::cerr << "ERROR!!! JetEnergyModuleKey::rowPhiCoord phi of row is "<<phi
163 <<" in a JEM at ("<<jemCoord.phi()<<", "<<jemCoord.eta()<<")"<<std::endl;
164 }//endif phi sensible
165 }//endif debug
166 return phi;
167}
168
170double JetEnergyModuleKey::dPhi(const Coordinate& /*coord*/) const{
171 return M_PI/2;// one quadrant
172}
173
176 if (region(coord) == MidJEM) return 0.8;
178 return 2.5;//end JEMs cover 2.4-4.9
179 return 0.0;
180}
181
182} //end of ns
183
184
188 double absEta=fabs(coord.eta());
189 if (absEta<=2.4) return MidJEM;
190 if ((coord.eta()>2.4)&&(coord.eta()<=4.9))return RightEndJEM;
191 if ((coord.eta()>=-4.9)&&(coord.eta()<-2.4))return LeftEndJEM;
192 std::cerr << "JetEnergyModuleKey::region ... Unknown region!"<<std::endl
193 << "Coordinate is ("<<coord.phi()<<", "<<coord.eta()<<")"
194 <<std::endl;
195 return JEMRegionERROR;
196}
197
198
199
201void LVL1::JetEnergyModuleKey::setBins(ICoordinate* iCoord, unsigned int& phiBin, int& etaBin){
202 unsigned int iPhiSize=16 ; // =PI/2 * 10
203 phiBin=( iCoord->phi() )/(iPhiSize);
204 int iEtaSize=8;// JEMs are 0.8 in eta, i.e ietaSize=8
205 etaBin=sharpRound2(iCoord->eta(),iEtaSize);
206 // the above calculation fails since end JEMs are larger, so ....
207 etaBin=(etaBin < -4 ? -4 : etaBin); // etaBin=-5
208 etaBin=(etaBin > 4 ? 4 : etaBin); // belongs to the JEM in etaBin=-4
209
210 return;
211}
212
215 switch (region(coord)){
216 case MidJEM :
217 return midJEMEtaCol(coord);
218 case LeftEndJEM :
219 return leftEndJEMEtaCol(coord);
220 case RightEndJEM :
221 return rightEndJEMEtaCol(coord);
222 default:
223 std::cerr << "JetEnergyModuleKey::etaColumn UNKNOWN REGION"<<std::endl;
224 }
225 return static_cast<unsigned int>(JEMRegionERROR);
226}
227
229unsigned int LVL1::JetEnergyModuleKey::jem(const Coordinate & coord) const {
230 double min=3.2; double jemEtaWidth=0.8;
231 unsigned int crateModifier=(phiQuadrant(coord.phi())/2);
232 //=0 for quad 0&1, =1 for quad 2&3
233
234 double temp = (coord.eta()+min)/jemEtaWidth;
235 int jem = static_cast<int>(temp);
236 jem=(jem < 0 ? 0 : jem);
237 jem=(jem > 7 ? 7 : jem);
238
239 jem+=(8*crateModifier);
240 return static_cast<unsigned int>(jem);
241}
242
247// double phiBinWidth=((2*M_PI)/64.0);
248 int abs_ieta=abs(iCoord->eta() );
249 int sign=( iCoord->eta() )/abs_ieta;
250
251 int etaBin=0; unsigned int phiBin=0;
252 setBins(iCoord,phiBin,etaBin);
253
254 double centralPhi=( (static_cast<double>(phiBin))*TrigT1CaloDefs::jemPhiSize )
256 double centralEta=0.0;
257 if ((etaBin>-4)&&(etaBin<4)){
258 centralEta=(static_cast<double>(etaBin)*(TrigT1CaloDefs::jemEtaSize) )
260 }else{
261 centralEta=3.65*sign;//centre of end JEMs
262 }
263 if (m_debugModuleKey){
264 std::cout << "JetEnergyModuleKey: start calcTrigBin"<<std::endl;
265 std::cout << "phi, eta : ("<<m_phi<<", "<<m_eta<<")"<<std::endl;
266 std::cout << "iphi, ieta : ("<<( iCoord->phi() )<<", "<<( iCoord->eta() )<<")"<<std::endl;
267 std::cout << "abs_ieta : ("<<abs_ieta<<" and sign : "<<sign<<std::endl;
268 std::cout << "central : ("<<centralPhi<<", "<<centralEta<<")"
269 << "bin : ("<<phiBin<<","<<etaBin<<")"<<std::endl;
270 }
271
272 Coordinate* centralCoords = new Coordinate(centralPhi, centralEta);
273 if (m_debugModuleKey) std::cout <<" JetEnergyModuleKey : created coord "<<(*centralCoords)<<std::endl;
274 return new BinAndCoord(phiBin,etaBin,centralCoords);
275}
276
279 double etaWidth=0.2;//only true for middle JEMS
280 double JEMEtaSize=0.8;
281// unsigned int columnsPerJEM= 4;
282 int etaBin=( jem(coord)%8)-4;//bins range from -4 to +3
283 double etaMin=(etaBin*JEMEtaSize);
284
285 double temp =(coord.eta()-etaMin)/etaWidth;
286 int colNum=sharpRound(temp);
287
288 colNum=(colNum < 0 ? 0 : colNum);//sorts out any rounding errors
289 colNum=(colNum > 3 ? 3 : colNum);
290
291 return colNum;
292}
293
295 // left hand end first
297 switch (get.jeRegion( coord )){
299 return 3;
300 break;
302 return 2;
303 break;
305 return 1;
306 break;
308 return 0;
309 break;
311 // here to stop compiler warnings
312 break;
314 // here to stop compiler warnings
315 break;
316 }
317 return static_cast<unsigned int>(JEMRegionERROR);
318}
319
322 // left hand end first
324 switch (get.jeRegion( coord )){
326 return 0;
327 break;
329 return 1;
330 break;
332 return 2;
333 break;
335 return 3;
336 break;
338 // here to stop compiler warnings
339 break;
341 // here to stop compiler warnings
342 break;
343 }
344 return static_cast<unsigned int>(JEMRegionERROR);
345}
346
347
348
354unsigned int LVL1::JetEnergyModuleKey::phiQuadrant(const double phi) const {
355 double temp = phi/(M_PI/2);
356 unsigned int quad = static_cast<unsigned int>(temp);
357 quad=(quad > 3 ? 0 : quad);
358
359 return static_cast<unsigned int>(quad);
360}
#define M_PI
double coord
Type of coordination system.
@ top
#define min(a, b)
Definition cfImp.cxx:40
Used to pass data between the methods of the Key Classes: Returns the Eta and Phi bins,...
Definition BinAndCoord.h:35
returns the trigger hardware components associated with a given Coordinate
unsigned int jepModule(const Coordinate &coord)
returns ID of JEP module (i.e.
unsigned int jepCrate(const Coordinate &Coord)
returns ID of JEP Crate that covers this coordinate
double phi() const
return phi
double eta() const
return eta
Used by Key Classes, returns and integer coorginate for the bin Eta-Phi.
Definition ICoordinate.h:26
int phi() const
return phi
int eta() const
return eta
The JetElementKey object provides the key for each JetElement depending on its eta,...
unsigned int col(const Coordinate &coord)
return row of passed coordinate
BinAndCoord * calculateTriggerBin(ICoordinate *iCoord)
converts integer phi, eta coordinates to phi, eta trigger bins, and central coords
double dPhi(const Coordinate &coord) const
height
void setBins(ICoordinate *iCoord, unsigned int &phiBin, int &etaBin)
sets the eta and phi bins
unsigned int phiQuadrant(const double phi) const
returns the quadrant number associated with the phi coordinate, 0 - 90 = 0 90 - 180 = 1 180-270 = 2 2...
std::vector< unsigned int > jeKeys(unsigned int crate, unsigned int module)
calculates keys of all JetElements in given crate and module
unsigned int jemKey(const xAOD::JetElement *jetElement)
returns the key of the passed Coordinate
unsigned int midJEMEtaCol(const Coordinate &coord) const
No descriptions.
unsigned int jem(const Coordinate &coord) const
returns ID of JEP module (i.e.
unsigned int leftEndJEMEtaCol(const Coordinate &coord) const
returns eta row of JEMs 0, or 8
JetEnergyModuleKey()
constructs a JetEnergyModuleKey object
std::vector< Coordinate > jeCoords(unsigned int crate, unsigned int module)
returns coordinates of all JetElements in given crate and module
bool m_debugModuleKey
set to true to turn debugging info on
unsigned int row(const Coordinate &coord) const
returns the phi row of a coord within the JEM that contains it.
double rowPhiCoord(unsigned int row, const Coordinate &jemCoord)
returns the phi coord of the row of the JEM at the coordinates passed
unsigned int rightEndJEMEtaCol(const Coordinate &coord) const
returns eta row of JEMs 7 or 15
double dEta(const Coordinate &coord) const
width
JEMRegion region(const Coordinate &coord) const
region
virtual double eta() const
returns the centre of the TT at eta_coord:
virtual int sign(int temp) const
returns -1 if temp is -ve and +1 if it is +ve.
int sharpRound2(int a, int b) const
divides a/b and returns a number as follows (where +b means +ve b): if 0
virtual double phi() const
returns phi coordinate of centre of relevant trigger tower.
int sharpRound(double a) const
rounds number as follows (-1.0 to 0.0) -> -1, (0.0 to 1.0) -> 0, (1.0 to 2.0)->1 etc.
ICoordinate * convertCoordsToIntegers(double phi, double eta)
converts the coordinates and corrects for overflows etc.
unsigned int key(double phi, double eta)
calculates a map key from passed phi, eta coordinates
Coordinate coord() const
return central coords of current key value.
KeyUtilities()
the constructor is protected so a user can never make a KeyUtilities object
double m_phi
phi coordinate of key
double m_eta
eta coordinate of key
static const double Region5EtaCentre
static const double RegionERROREtaCentre
static const double jemPhiSize
static const double jemEtaSize
static const double Region0Height
float phi() const
get phi (note that for L1Calo phi runs from 0 to 2pi)
float eta() const
get eta
T * get(TKey *tobj)
get a TObject* from a TKey* (why can't a TObject be a TKey?)
Definition hcg.cxx:130
eFexTowerBuilder creates xAOD::eFexTowerContainer from supercells (LATOME) and triggerTowers (TREX) i...
JetElement_v2 JetElement
Define the latest version of the JetElement class.