ATLAS Offline Software
Loading...
Searching...
No Matches
TrigFTF_GNN_Layer Class Reference

#include <GNN_Geometry.h>

Collaboration diagram for TrigFTF_GNN_Layer:

Public Member Functions

 TrigFTF_GNN_Layer (const TrigInDetSiLayer &, float, int)
 ~TrigFTF_GNN_Layer ()
int getEtaBin (float, float) const
float getMinBinRadius (int) const
float getMaxBinRadius (int) const
int num_bins () const
bool verifyBin (const TrigFTF_GNN_Layer *, int, int, float, float) const

Public Attributes

const TrigInDetSiLayerm_layer
std::vector< int > m_bins
std::vector< float > m_minRadius
std::vector< float > m_maxRadius
std::vector< float > m_minBinCoord
std::vector< float > m_maxBinCoord
float m_minEta
float m_maxEta
float m_etaBin

Protected Attributes

float m_etaBinWidth
float m_r1
float m_z1
float m_r2
float m_z2
int m_nBins

Detailed Description

Definition at line 16 of file GNN_Geometry.h.

Constructor & Destructor Documentation

◆ TrigFTF_GNN_Layer()

TrigFTF_GNN_Layer::TrigFTF_GNN_Layer ( const TrigInDetSiLayer & ls,
float ew,
int bin0 )

Definition at line 16 of file GNN_Geometry.cxx.

16 : m_layer(ls), m_etaBinWidth(ew) {
17
18 if(m_layer.m_type == 0) {//barrel
19 m_r1 = m_layer.m_refCoord;
20 m_r2 = m_layer.m_refCoord;
21 m_z1 = m_layer.m_minBound;
22 m_z2 = m_layer.m_maxBound;
23 }
24 else {//endcap
25 m_r1 = m_layer.m_minBound;
26 m_r2 = m_layer.m_maxBound;
27 m_z1 = m_layer.m_refCoord;
28 m_z2 = m_layer.m_refCoord;
29 }
30
31 float t1 = m_z1/m_r1;
32 float eta1 = -std::log(std::sqrt(1+t1*t1)-t1);
33
34
35 float t2 = m_z2/m_r2;
36 float eta2 = -std::log(std::sqrt(1+t2*t2)-t2);
37
38 m_minEta = eta1;
39 m_maxEta = eta2;
40
41 if(m_maxEta < m_minEta) {
42 m_minEta = eta2;
43 m_maxEta = eta1;
44 }
45
46 m_maxEta += 1e-6;//increasing them slightly to avoid range_check exceptions
47 m_minEta -= 1e-6;
48
49 float deltaEta = m_maxEta - m_minEta;
50
51 int binCounter = bin0;
52
53 if(deltaEta < m_etaBinWidth) {
54 m_nBins = 1;
55 m_bins.push_back(binCounter++);
57 if(m_layer.m_type == 0) {//barrel
58 m_minRadius.push_back(m_layer.m_refCoord - 2.0);
59 m_maxRadius.push_back(m_layer.m_refCoord + 2.0);
60 m_minBinCoord.push_back(m_layer.m_minBound);
61 m_maxBinCoord.push_back(m_layer.m_maxBound);
62 }
63 else {//endcap
64 m_minRadius.push_back(m_layer.m_minBound - 2.0);
65 m_maxRadius.push_back(m_layer.m_maxBound + 2.0);
66 m_minBinCoord.push_back(m_layer.m_minBound);
67 m_maxBinCoord.push_back(m_layer.m_maxBound);
68 }
69 }
70 else {
71 float nB = static_cast<int>(deltaEta/m_etaBinWidth);
72 m_nBins = nB;
73 if(deltaEta - m_etaBinWidth*nB > 0.5*m_etaBinWidth) m_nBins++;
74
76
77 if(m_nBins == 1) {
78 m_bins.push_back(binCounter++);
79 if(m_layer.m_type == 0) {//barrel
80 m_minRadius.push_back(m_layer.m_refCoord - 2.0);
81 m_maxRadius.push_back(m_layer.m_refCoord + 2.0);
82 m_minBinCoord.push_back(m_layer.m_minBound);
83 m_maxBinCoord.push_back(m_layer.m_maxBound);
84 }
85 else {//endcap
86 m_minRadius.push_back(m_layer.m_minBound - 2.0);
87 m_maxRadius.push_back(m_layer.m_maxBound + 2.0);
88 m_minBinCoord.push_back(m_layer.m_minBound);
89 m_maxBinCoord.push_back(m_layer.m_maxBound);
90 }
91 }
92 else {
93
94 float eta = m_minEta+0.5*m_etaBin;
95
96 for(int i=1;i<=m_nBins;i++) {
97
98 m_bins.push_back(binCounter++);
99
100 float e1 = eta - 0.5*m_etaBin;
101 float e2 = eta + 0.5*m_etaBin;
102
103 if(m_layer.m_type == 0) {//barrel
104 m_minRadius.push_back(m_layer.m_refCoord - 2.0);
105 m_maxRadius.push_back(m_layer.m_refCoord + 2.0);
106 float z1 = m_layer.m_refCoord*std::sinh(e1);
107 m_minBinCoord.push_back(z1);
108 float z2 = m_layer.m_refCoord*std::sinh(e2);
109 m_maxBinCoord.push_back(z2);
110 }
111 else {//endcap
112 if (m_layer.m_refCoord > 0) {//for the positive endcap larger eta corresponds to smaller radius
113 std::swap(e1, e2);
114 }
115 float r = m_layer.m_refCoord/std::sinh(e1);
116 m_minBinCoord.push_back(r);
117 m_minRadius.push_back(r - 2.0);
118 r = m_layer.m_refCoord/std::sinh(e2);
119 m_maxBinCoord.push_back(r);
120 m_maxRadius.push_back(r + 2.0);
121
122 }
123
124 eta += m_etaBin;
125 }
126 }
127 }
128}
Scalar eta() const
pseudorapidity method
std::vector< float > m_minBinCoord
const TrigInDetSiLayer & m_layer
std::vector< float > m_maxRadius
std::vector< int > m_bins
std::vector< float > m_minRadius
std::vector< float > m_maxBinCoord
int r
Definition globals.cxx:22
std::vector< ALFA_RawDataContainer_p1 > t2
std::vector< ALFA_RawDataCollection_p1 > t1
double deltaEta(const I4Momentum &p1, const I4Momentum &p2)
Computes efficiently .
Definition P4Helpers.h:66
double e2(const xAOD::CaloCluster &cluster)
return the uncorrected cluster energy in 2nd sampling
double e1(const xAOD::CaloCluster &cluster)
return the uncorrected cluster energy in 1st sampling
void swap(ElementLinkVector< DOBJ > &lhs, ElementLinkVector< DOBJ > &rhs)
setEt setPhi setE277 setWeta2 eta1

◆ ~TrigFTF_GNN_Layer()

TrigFTF_GNN_Layer::~TrigFTF_GNN_Layer ( )

Definition at line 296 of file GNN_Geometry.cxx.

296 {
297 m_bins.clear();
298}

Member Function Documentation

◆ getEtaBin()

int TrigFTF_GNN_Layer::getEtaBin ( float zh,
float rh ) const

Definition at line 265 of file GNN_Geometry.cxx.

265 {
266
267 if(m_bins.size() == 1) return m_bins.at(0);
268
269 float t1 = zh/rh;
270 float eta = -std::log(std::sqrt(1+t1*t1)-t1);
271
272
273 int idx = static_cast<int>((eta-m_minEta)/m_etaBin);
274
275
276 if(idx < 0) idx = 0;
277 if(idx >= static_cast<int>(m_bins.size())) idx = static_cast<int>(m_bins.size())-1;
278
279 return m_bins.at(idx);//index in the global storage
280}

◆ getMaxBinRadius()

float TrigFTF_GNN_Layer::getMaxBinRadius ( int idx) const

Definition at line 289 of file GNN_Geometry.cxx.

289 {
290 if(idx >= static_cast<int>(m_maxRadius.size())) idx = idx-1;
291 if(idx < 0) idx = 0;
292
293 return m_maxRadius.at(idx);
294}

◆ getMinBinRadius()

float TrigFTF_GNN_Layer::getMinBinRadius ( int idx) const

Definition at line 282 of file GNN_Geometry.cxx.

282 {
283 if(idx >= static_cast<int>(m_minRadius.size())) idx = idx-1;
284 if(idx < 0) idx = 0;
285
286 return m_minRadius.at(idx);
287}

◆ num_bins()

int TrigFTF_GNN_Layer::num_bins ( ) const
inline

Definition at line 26 of file GNN_Geometry.h.

26{return m_bins.size();}

◆ verifyBin()

bool TrigFTF_GNN_Layer::verifyBin ( const TrigFTF_GNN_Layer * pL,
int b1,
int b2,
float min_z0,
float max_z0 ) const

Definition at line 130 of file GNN_Geometry.cxx.

130 {
131
132 float z1min = m_minBinCoord.at(b1);
133 float z1max = m_maxBinCoord.at(b1);
134 float r1 = m_layer.m_refCoord;
135
136 const float tol = 5.0;
137
138 if(m_layer.m_type == 0 && pL->m_layer.m_type == 0) {//barrel <- barrel
139
140 float min_b2 = pL->m_minBinCoord.at(b2);
141 float max_b2 = pL->m_maxBinCoord.at(b2);
142
143 float r2 = pL->m_layer.m_refCoord;
144
145 float A = r2/(r2-r1);
146 float B = r1/(r2-r1);
147
148 float z0_min = z1min*A - max_b2*B;
149 float z0_max = z1max*A - min_b2*B;
150
151 if(z0_max < min_z0-tol || z0_min > max_z0+tol) return false;
152
153 return true;
154 }
155
156 if(m_layer.m_type == 0 && pL->m_layer.m_type != 0) {//barrel <- endcap
157
158 float z2 = pL->m_layer.m_refCoord;
159 float r2max = pL->m_maxBinCoord.at(b2);
160 float r2min = pL->m_minBinCoord.at(b2);
161
162 if(r2max <= r1) return false;
163
164 if(r2min <= r1) {
165 r2min = r1 + 1e-3;
166 }
167
168 float z0_max = 0.0;
169 float z0_min = 0.0;
170
171 if(z2 > 0) {
172 z0_max = (z1max*r2max - z2*r1)/(r2max-r1);
173 z0_min = (z1min*r2min - z2*r1)/(r2min-r1);
174 }
175 else {
176 z0_max = (z1max*r2min - z2*r1)/(r2min-r1);
177 z0_min = (z1min*r2max - z2*r1)/(r2max-r1);
178 }
179
180 if(z0_max < min_z0-tol || z0_min > max_z0+tol) return false;
181 return true;
182 }
183
184 if(m_layer.m_type != 0 && pL->m_layer.m_type != 0) {//endcap <- endcap
185
186 float z2 = pL->m_layer.m_refCoord;
187 float z1 = m_layer.m_refCoord;
188 float r2max = pL->m_maxBinCoord.at(b2);
189 float r2min = pL->m_minBinCoord.at(b2);
190 float r1max = m_maxBinCoord.at(b1);
191 float r1min = m_minBinCoord.at(b1);
192
193 if (r1min >= r2max) return false;
194
195 if (z2 > 0) { // positive endcap
196
197 float z0_max = z1 - r1min*(z2-z1)/(r2max-r1min);
198
199 if(z0_max < min_z0-tol) return false;
200
201 if (r2min > r1max) {
202
203 float z0_min = z1 - r1max*(z2-z1)/(r2min-r1max);
204
205 if (z0_min > max_z0+tol) return false;
206 }
207 }
208 else { // negative endcap
209 float z0_min = z1 - r1min*(z2-z1)/(r2max-r1min);
210
211 if(z0_min > max_z0+tol) return false;
212
213 if (r2min > r1max) {
214
215 float z0_max = z1 - r1max*(z2-z1)/(r2min-r1max);
216
217 if (z0_max < min_z0-tol) return false;
218 }
219 }
220 return true;
221 }
222
223 if(m_layer.m_type != 0 && pL->m_layer.m_type == 0) {//endcap <- barrel
224
225 float z1 = m_layer.m_refCoord;
226 float r1max = m_maxBinCoord.at(b1);
227 float r1min = m_minBinCoord.at(b1);
228
229 float z2min = pL->m_minBinCoord.at(b2);
230 float z2max = pL->m_maxBinCoord.at(b2);
231 float r2 = pL->m_layer.m_refCoord;
232
233 if (r2 < r1min) return false;
234
235 // interval 1
236
237 float z0_min = z1 - (z2max - z1)/(r2/r1max - 1);
238 float z0_max = z1 - (z2max - z1)/(r2/r1min - 1);
239
240 if (z0_min > z0_max) std::swap(z0_min, z0_max);
241
242 bool beyond_range = (z0_max < min_z0-tol || z0_min > max_z0+tol);
243
244 if (!beyond_range) return true;
245
246 // interval 2
247
248 z0_min = z1 - (z2min - z1)/(r2/r1max - 1);
249 z0_max = z1 - (z2min - z1)/(r2/r1min - 1);
250
251 if (z0_min > z0_max) std::swap(z0_min, z0_max);
252
253 beyond_range = (z0_max < min_z0-tol || z0_min > max_z0+tol);
254
255 if (!beyond_range) return true;
256
257 return false;
258
259 }
260
261 return true;
262}

Member Data Documentation

◆ m_bins

std::vector<int> TrigFTF_GNN_Layer::m_bins

Definition at line 31 of file GNN_Geometry.h.

◆ m_etaBin

float TrigFTF_GNN_Layer::m_etaBin

Definition at line 37 of file GNN_Geometry.h.

◆ m_etaBinWidth

float TrigFTF_GNN_Layer::m_etaBinWidth
protected

Definition at line 41 of file GNN_Geometry.h.

◆ m_layer

const TrigInDetSiLayer& TrigFTF_GNN_Layer::m_layer

Definition at line 30 of file GNN_Geometry.h.

◆ m_maxBinCoord

std::vector<float> TrigFTF_GNN_Layer::m_maxBinCoord

Definition at line 35 of file GNN_Geometry.h.

◆ m_maxEta

float TrigFTF_GNN_Layer::m_maxEta

Definition at line 37 of file GNN_Geometry.h.

◆ m_maxRadius

std::vector<float> TrigFTF_GNN_Layer::m_maxRadius

Definition at line 33 of file GNN_Geometry.h.

◆ m_minBinCoord

std::vector<float> TrigFTF_GNN_Layer::m_minBinCoord

Definition at line 34 of file GNN_Geometry.h.

◆ m_minEta

float TrigFTF_GNN_Layer::m_minEta

Definition at line 37 of file GNN_Geometry.h.

◆ m_minRadius

std::vector<float> TrigFTF_GNN_Layer::m_minRadius

Definition at line 32 of file GNN_Geometry.h.

◆ m_nBins

int TrigFTF_GNN_Layer::m_nBins
protected

Definition at line 44 of file GNN_Geometry.h.

◆ m_r1

float TrigFTF_GNN_Layer::m_r1
protected

Definition at line 43 of file GNN_Geometry.h.

◆ m_r2

float TrigFTF_GNN_Layer::m_r2
protected

Definition at line 43 of file GNN_Geometry.h.

◆ m_z1

float TrigFTF_GNN_Layer::m_z1
protected

Definition at line 43 of file GNN_Geometry.h.

◆ m_z2

float TrigFTF_GNN_Layer::m_z2
protected

Definition at line 43 of file GNN_Geometry.h.


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