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

This is an internal class, used in the jet trigger. More...

#include <JEMJetAlgorithm.h>

Collaboration diagram for LVL1::JEMJetAlgorithm:

Public Member Functions

 JEMJetAlgorithm (double eta, double phi, const std::map< int, JetInput * > *jiContainer, const TrigConf::L1Menu *l1menu)
 ~JEMJetAlgorithm ()
double eta ()
 Accessors.
double phi ()
 Returns phi coordinate of RoI, using standard ATLAS convention.
int Core ()
 Returns RoI Core ET.
int ET4x4 ()
 Returns 4x4 TT cluster ET.
int ET6x6 ()
 Returns 6x6 TT cluster ET.
int ET8x8 ()
 Returns 8x8 TT cluster ET.
int ETLarge ()
 Returns Large cluster ET.
int ETSmall ()
 Returns Small cluster ET.
bool isEtMax ()
 Does this window pass the local ET maximum condition.
bool isRoI ()
 Does this window pass the local ET maximum condition.
xAOD::JEMTobRoIjemTobRoI ()
 Create JEMTobRoI and return pointers to it.

Private Member Functions

void setRoICoord (double eta, double phi)
 Compute RoI coordinate.
void testEtMax (const std::vector< int > &cores)
 Form all 2x2 clusters within window and test centre is a local ET maximum.
void passesTrigger ()
 Check trigger condition and set ET values if TOB created.

Private Attributes

double m_refEta
double m_refPhi
const TrigConf::L1Menum_l1menu {nullptr}
double m_eta
 Algorithm results.
double m_phi
int m_ET4x4
int m_ET6x6
int m_ET8x8
int m_ETLarge
int m_ETSmall
bool m_EtMax

Static Private Attributes

static const int m_satLarge = 0x3FF
 Algorithm parameters.
static const int m_satSmall = 0x1FF

Detailed Description

This is an internal class, used in the jet trigger.

The JEMJetAlgorithm:

  • identifies JetInputs used in a particular trigger window (4x4 towers)
  • forms sums required for trigger algorithm
  • identifies whether passes the RoI condition (local ET maximum)
  • returns sums (with appropriate precision)
  • returns hit bits and RoI word (if all requirements met)

Definition at line 40 of file JEMJetAlgorithm.h.

Constructor & Destructor Documentation

◆ JEMJetAlgorithm()

LVL1::JEMJetAlgorithm::JEMJetAlgorithm ( double eta,
double phi,
const std::map< int, JetInput * > * jiContainer,
const TrigConf::L1Menu * l1menu )

This could potentially be called with the coordinate being the centre
of an RoI or the centre of a JetElement. We want to ensure that we
do the right thing in either case.

Offsetting by JE size/4 should put the coordinates inside the
reference JetElement whether input was reference JetElement
coordinate or RoI centre coordinate. Can't subtract a fixed amount
as the eta size of JetElements varies, so need to work out how much
to subtract first!
To be extra safe, the code below will then compute the centre
of the reference JetElement.

This should all ensure consistent & safe key generation however the
initial coordinate was specified.

Now loop over JetInputs within window and form clusters. ET8x8 = all elements ET4x4 = central 2x2 ET6x6 = most energetic of 4 clusters Interesting one is the 9 RoI cores. The logic below fills these as an a 9-element std::vector. The ordering of the elements is as follows: 2 5 8 1 4 7 0 3 6 So the RoI "ET maximum" condition is that 4 >= 0-3, 4 > 5-8 If you picture the cores in a 3x3 array, cores[iphi][ieta], then the index in the vector above = iphi + 3*ieta.

We use the fact that core[4] = ET4x4 to safe a little time here

Definition at line 30 of file JEMJetAlgorithm.cxx.

31 :
32 m_l1menu(l1menu),
33 m_ET4x4(0),
34 m_ET6x6(0),
35 m_ET8x8(0),
36 m_ETLarge(0),
37 m_ETSmall(0),
38 m_EtMax(false)
39 //m_debug(false)
40{
56
57 // Offset depends on eta. Need to protect against rounding errors even here
58 JetInputKey get(phi-0.01,eta-0.01);
59 double de = get.dEta()/4.;
60
61 // Get coordinate of centre of this "reference element"
62 Coordinate refCoord = get.getCentre(phi-M_PI/64., eta-de);
63 m_refEta = refCoord.eta();
64 m_refPhi = refCoord.phi();
65
66 // Set coordinate of centre of RoI, starting from reference tower coordinate
68
69 // Get coordinate of bottom-left JetInput in this window
70 int minEtaOffset = -1;
71 Coordinate startCoord = get.lowerLeft(m_refPhi,m_refEta);
72 if (startCoord.eta() == TrigT1CaloDefs::RegionERROREtaCentre) { // already at left edge
73 startCoord = get.downPhi(m_refPhi,m_refEta);
74 minEtaOffset = 0;
75 }
76 double tempEta = startCoord.eta();
77 double startPhi = startCoord.phi();
78
95
96 std::vector<int> et6x6(4);
97 std::vector<int> cores(9);
98 for (int etaOffset=minEtaOffset; etaOffset<=2 ; etaOffset++){
99 Coordinate tempCoord(startPhi,tempEta);
100 for (int phiOffset=-1; phiOffset<=2 ; phiOffset++){
101 int tempKey = get.jeKey(tempCoord);
102 std::map<int, JetInput*>::const_iterator ji = jiContainer->find(tempKey);
103 if (ji != jiContainer->end() ){
104 // get ET once here, rather than repeat function calls
105 int ET = (ji->second)->energy();
106 // 8x8 jet is easy
107 m_ET8x8 += ET;
108 // but there are 4 possible 6x6 clusters
109 if (phiOffset > -1) {
110 if (etaOffset < 2) et6x6[1] += ET;
111 if (etaOffset > -1) et6x6[2] += ET;
112 }
113 if (phiOffset < 2) {
114 if (etaOffset < 2) et6x6[0] += ET;
115 if (etaOffset > -1) et6x6[3] += ET;
116 }
117 // Each JetInput is part of up to 4 RoI core clusters
118
119 if (etaOffset >= 0) {
120 if (phiOffset >= 0) cores[phiOffset+3*etaOffset] += ET;
121 if (phiOffset < 2) cores[phiOffset+3*etaOffset+1] += ET;
122 }
123 if (etaOffset < 2) {
124 if (phiOffset >= 0) cores[phiOffset+3*etaOffset+3] += ET;
125 if (phiOffset < 2) cores[phiOffset+3*etaOffset+4] += ET;
126 }
127 } // end of check that jetinput exists in container
128 tempCoord = get.upPhi(tempCoord); // Increment phi coordinate
129 } // end phi offset loop
130 tempCoord = get.rightEta(tempCoord); // Increment eta coordinate
131 tempEta = tempCoord.eta();
132 if (tempEta == TrigT1CaloDefs::RegionERROREtaCentre) break; // gone outside coverage
133 } // end eta offset loop
134
135 // 4x4 cluster = central RoI core
136 m_ET4x4 = cores[4];
137
138 // find most energetic 6x6 cluster
139 for (int i = 0; i < 4; i++) if (et6x6[i] > m_ET6x6) m_ET6x6 = et6x6[i];
140
141 // Check whether RoI condition is met.
142 testEtMax(cores);
143
144 // test trigger conditions
146
147}
#define M_PI
const TrigConf::L1Menu * m_l1menu
double eta()
Accessors.
double phi()
Returns phi coordinate of RoI, using standard ATLAS convention.
void passesTrigger()
Check trigger condition and set ET values if TOB created.
void setRoICoord(double eta, double phi)
Compute RoI coordinate.
void testEtMax(const std::vector< int > &cores)
Form all 2x2 clusters within window and test centre is a local ET maximum.
static const double RegionERROREtaCentre
T * get(TKey *tobj)
get a TObject* from a TKey* (why can't a TObject be a TKey?)
Definition hcg.cxx:130

◆ ~JEMJetAlgorithm()

LVL1::JEMJetAlgorithm::~JEMJetAlgorithm ( )

Definition at line 149 of file JEMJetAlgorithm.cxx.

149 {
150}

Member Function Documentation

◆ Core()

int LVL1::JEMJetAlgorithm::Core ( )

Returns RoI Core ET.

Definition at line 246 of file JEMJetAlgorithm.cxx.

246 {
247 return m_ET4x4;
248}

◆ ET4x4()

int LVL1::JEMJetAlgorithm::ET4x4 ( )

Returns 4x4 TT cluster ET.

Definition at line 251 of file JEMJetAlgorithm.cxx.

251 {
252 return ( (m_ET4x4 < m_satLarge) ? m_ET4x4 : m_satLarge );
253}
static const int m_satLarge
Algorithm parameters.

◆ ET6x6()

int LVL1::JEMJetAlgorithm::ET6x6 ( )

Returns 6x6 TT cluster ET.

Definition at line 256 of file JEMJetAlgorithm.cxx.

256 {
257 return ( (m_ET6x6 < m_satLarge) ? m_ET6x6 : m_satLarge );
258}

◆ ET8x8()

int LVL1::JEMJetAlgorithm::ET8x8 ( )

Returns 8x8 TT cluster ET.

Definition at line 261 of file JEMJetAlgorithm.cxx.

261 {
262 return ( (m_ET8x8 < m_satLarge) ? m_ET8x8 : m_satLarge );
263}

◆ eta()

double LVL1::JEMJetAlgorithm::eta ( )

Accessors.

Returns eta coordinate of RoI.

Definition at line 286 of file JEMJetAlgorithm.cxx.

286 {
287 return m_eta;
288}
double m_eta
Algorithm results.

◆ ETLarge()

int LVL1::JEMJetAlgorithm::ETLarge ( )

Returns Large cluster ET.

Definition at line 266 of file JEMJetAlgorithm.cxx.

266 {
267 return ( (m_ETLarge < m_satLarge) ? m_ETLarge : m_satLarge );
268}

◆ ETSmall()

int LVL1::JEMJetAlgorithm::ETSmall ( )

Returns Small cluster ET.

Definition at line 271 of file JEMJetAlgorithm.cxx.

271 {
272 return ( (m_ETSmall < m_satSmall) ? m_ETSmall : m_satSmall );
273}
static const int m_satSmall

◆ isEtMax()

bool LVL1::JEMJetAlgorithm::isEtMax ( )

Does this window pass the local ET maximum condition.

Definition at line 276 of file JEMJetAlgorithm.cxx.

276 {
277 return m_EtMax;
278}

◆ isRoI()

bool LVL1::JEMJetAlgorithm::isRoI ( )

Does this window pass the local ET maximum condition.

Definition at line 281 of file JEMJetAlgorithm.cxx.

281 {
282 return ( m_EtMax && (m_ETLarge > 0 || m_ETSmall > 0) );
283}

◆ jemTobRoI()

xAOD::JEMTobRoI * LVL1::JEMJetAlgorithm::jemTobRoI ( )

Create JEMTobRoI and return pointers to it.

Returns a JEMTobRoI object, provided the TOB conditions were met.

Will return a null pointer if object does not pass hypothesis. It is the user's responsibility to check

If not will return a null pointer - user's responsibility to check

Definition at line 298 of file JEMJetAlgorithm.cxx.

298 {
299
300 if (isRoI()) {
301
302 // Need to calculate hardware coordinate
303 Coordinate coord(m_refPhi, m_refEta);
304 CoordToHardware convertor;
305
306 int crate = convertor.jepCrate(coord);
307 int jem = convertor.jepModule(coord);
308 unsigned int jemCoord = convertor.jepLocalCoordinate(coord);
309 int frame = (jemCoord>>2);
310 int lc = jemCoord&3;
311
312 /*
313 JetEnergyModuleKey get();
314 Coordinate tempCoord(m_refPhi, m_refEta);
315 unsigned int quadrant = m_refPhi*2/M_PI;
316 unsigned int row = get.row(tempCoord);
317 unsigned int col = get.col(tempCoord);
318
319 unsigned int crate = quadrant&1;
320 unsigned int jem = get.jem(tempCoord);
321 unsigned int frame = ((col&2)<<1) + (row>>1);
322 unsigned int lc = (col&1) + (row&1)*2;
323 */
324 xAOD::JEMTobRoI* roi = new xAOD::JEMTobRoI();
325 roi->makePrivateStore();
326 roi->initialize(crate, jem, frame, lc, ETLarge(), ETSmall());
327 return roi;
328 }
329
330 return 0;
331
332}
double coord
Type of coordination system.
int ETSmall()
Returns Small cluster ET.
bool isRoI()
Does this window pass the local ET maximum condition.
int ETLarge()
Returns Large cluster ET.
void makePrivateStore()
Create a new (empty) private store for this object.
virtual void initialize(const int crate, const int jem, const int frame, const int location, const int energyLarge, const int energySmall)
JEMTobRoI_v1 JEMTobRoI
Define the latest version of the JEMTobRoI class.

◆ passesTrigger()

void LVL1::JEMJetAlgorithm::passesTrigger ( )
private

Check trigger condition and set ET values if TOB created.

Definition at line 207 of file JEMJetAlgorithm.cxx.

207 {
208
209 // Belt and braces
210 m_ETLarge = 0;
211 m_ETSmall = 0;
212
213 // Don't waste time if it isn't an RoI candidate
214 if (!m_EtMax) return;
215
216 // Does this pass min TOB pT cut?
217 unsigned int sizeSmall{4}; // the size of the small jets (by default 4)
218 unsigned int sizeLarge{8}; // the size of the large jets (by default 8)
219 int threshSmall{0}; // the minimum pT of small jet objects sent to TOPO (in counts, not in GeV)
220 int threshLarge{0}; // the minimum pT of large jet objects sent to TOPO (in counts, not in GeV)
221
222 sizeSmall = 4; // not part of the new menu
223 sizeLarge = 8; // not part of the new menu
224 float scale = m_l1menu->thrExtraInfo().JET().jetScale();
225 threshSmall = m_l1menu->thrExtraInfo().JET().ptMinToTopoSmallWindowCounts()*scale;
226 threshLarge = m_l1menu->thrExtraInfo().JET().ptMinToTopoLargeWindowCounts()*scale;
227
228 int etLarge = m_ET8x8;
229 if (sizeLarge == 6) etLarge = m_ET6x6;
230 else if (sizeLarge == 4) etLarge = m_ET4x4;
231
232 int etSmall = m_ET4x4;
233 if (sizeSmall == 6) etLarge = m_ET6x6;
234 else if (sizeSmall == 8) etLarge = m_ET8x8;
235
236 if (etLarge <= threshLarge && etSmall <= threshSmall) return;
237
238 m_ETLarge = etLarge;
239 m_ETSmall = etSmall;
240
241}

◆ phi()

double LVL1::JEMJetAlgorithm::phi ( )

Returns phi coordinate of RoI, using standard ATLAS convention.

Definition at line 291 of file JEMJetAlgorithm.cxx.

291 {
292 return ( (m_phi <= M_PI) ? m_phi : m_phi - 2.*M_PI);
293}

◆ setRoICoord()

void LVL1::JEMJetAlgorithm::setRoICoord ( double eta,
double phi )
private

Compute RoI coordinate.

Input coordinate should be inside "reference tower" of window. Computes RoI coordinate as centre of 2x2 element RoI core

Definition at line 154 of file JEMJetAlgorithm.cxx.

154 {
155
156 JetInputKey keyLL(phi,eta);
157 // Get coordinate of centre of this "reference element"
158 Coordinate lowerLeft = keyLL.getCentre(phi, eta);
159
160 // Hence find lower-left corner of RoI core
161 double lowerEta = lowerLeft.eta() - ( keyLL.dEta()/2. );
162 double lowerPhi = lowerLeft.phi() - ( keyLL.dPhi()/2. );
163
164 // Get coordinate of opposite corner of RoI core:
165 Coordinate upperRight;
166 if (keyLL.isFCAL(lowerLeft.eta()) && eta > 0) {
167 // Special case: no JE to right, so get centre of element at eta, phi+1
168 upperRight = keyLL.upPhi(phi, eta);
169 }
170 else {
171 // Get centre of element at eta+1, phi+1, i.e. opposite corner of RoI core
172 upperRight = keyLL.upperRight(phi, eta);
173 }
174
175 // Now get upper-right corner of RoI core:
176 JetInputKey keyUR(upperRight);
177 double upperEta = upperRight.eta() + ( keyUR.dEta()/2. );
178 double upperPhi = upperRight.phi() + ( keyUR.dPhi()/2. );
179
180
181 // CoordinateRange object will compute centre, correcting for wrap-around
182 CoordinateRange roi(lowerPhi,upperPhi,lowerEta,upperEta);
183 m_eta = roi.eta();
184 m_phi = roi.phi();
185
186}

◆ testEtMax()

void LVL1::JEMJetAlgorithm::testEtMax ( const std::vector< int > & cores)
private

Form all 2x2 clusters within window and test centre is a local ET maximum.

Input clusters form a list, arranged as
2 5 8
1 4 7
0 3 6
Then 4 = this window's RoI core, and ET max condition is
4 >= 0, 1, 2, 3 && 4 < 5, 6, 7, 8

Definition at line 189 of file JEMJetAlgorithm.cxx.

189 {
190
198
199 // RoI condition test
200 m_EtMax = true;
201 for (int i = 0; i < 4; i++) if (cores[4] < cores[i]) m_EtMax = false;
202 for (int i = 5; i < 9; i++) if (cores[4] <= cores[i]) m_EtMax = false;
203
204}

Member Data Documentation

◆ m_ET4x4

int LVL1::JEMJetAlgorithm::m_ET4x4
private

Definition at line 73 of file JEMJetAlgorithm.h.

◆ m_ET6x6

int LVL1::JEMJetAlgorithm::m_ET6x6
private

Definition at line 74 of file JEMJetAlgorithm.h.

◆ m_ET8x8

int LVL1::JEMJetAlgorithm::m_ET8x8
private

Definition at line 75 of file JEMJetAlgorithm.h.

◆ m_eta

double LVL1::JEMJetAlgorithm::m_eta
private

Algorithm results.

Definition at line 71 of file JEMJetAlgorithm.h.

◆ m_ETLarge

int LVL1::JEMJetAlgorithm::m_ETLarge
private

Definition at line 76 of file JEMJetAlgorithm.h.

◆ m_EtMax

bool LVL1::JEMJetAlgorithm::m_EtMax
private

Definition at line 78 of file JEMJetAlgorithm.h.

◆ m_ETSmall

int LVL1::JEMJetAlgorithm::m_ETSmall
private

Definition at line 77 of file JEMJetAlgorithm.h.

◆ m_l1menu

const TrigConf::L1Menu* LVL1::JEMJetAlgorithm::m_l1menu {nullptr}
private

Definition at line 68 of file JEMJetAlgorithm.h.

68{nullptr};

◆ m_phi

double LVL1::JEMJetAlgorithm::m_phi
private

Definition at line 72 of file JEMJetAlgorithm.h.

◆ m_refEta

double LVL1::JEMJetAlgorithm::m_refEta
private

Definition at line 66 of file JEMJetAlgorithm.h.

◆ m_refPhi

double LVL1::JEMJetAlgorithm::m_refPhi
private

Definition at line 67 of file JEMJetAlgorithm.h.

◆ m_satLarge

const int LVL1::JEMJetAlgorithm::m_satLarge = 0x3FF
staticprivate

Algorithm parameters.

Definition at line 83 of file JEMJetAlgorithm.h.

◆ m_satSmall

const int LVL1::JEMJetAlgorithm::m_satSmall = 0x1FF
staticprivate

Definition at line 84 of file JEMJetAlgorithm.h.


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