ATLAS Offline Software
Loading...
Searching...
No Matches
TgcRoadDefiner.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#include "TgcRoadDefiner.h"
6#include "MdtRegion.h"
9
11
12#include <cmath>
13
14// --------------------------------------------------------------------------------
15// --------------------------------------------------------------------------------
16
18 const std::string& name,
19 const IInterface* parent):
20 AthAlgTool(type, name, parent)
21{
22}
23
24// --------------------------------------------------------------------------------
25// --------------------------------------------------------------------------------
26
28{
29
30 ATH_CHECK(AthAlgTool::initialize());
31
32 ATH_CHECK(m_regionSelector.retrieve());
33 ATH_MSG_DEBUG("Retrieved the RegionSelector tool ");
34
35 ATH_CHECK(m_tgcFit.retrieve());
36 ATH_MSG_DEBUG("Retrieved service " << m_tgcFit);
37
38 ATH_CHECK(m_idHelperSvc.retrieve());
39
40 return StatusCode::SUCCESS;
41}
42
43// --------------------------------------------------------------------------------
44// --------------------------------------------------------------------------------
45
46StatusCode TrigL2MuonSA::TgcRoadDefiner::defineRoad(const EventContext& ctx,
47 const TrigRoiDescriptor* p_roids,
48 const bool insideOut,
49 const TrigL2MuonSA::TgcHits& tgcHits,
50 TrigL2MuonSA::MuonRoad& muonRoad,
51 TrigL2MuonSA::TgcFitResult& tgcFitResult) const
52{
53 const int N_STATION = 10;
54 const int N_LAYER = 8;
55
56 bool isMiddleFailure = true;
57
58 // Define road by using TGC fit result
59 const double R_WIDTH_DEFAULT = 100;
60 const double R_WIDTH_INNER_NO_HIT = 200;
61
62 const double ZERO_LIMIT = 1e-5;
63
64 muonRoad.isEndcap = true;
65
66 int side;
67 double roiEta;
68 double theta;
69 double aw;
77
78 TrigL2MuonSA::TgcFit::PointArray tgcStripMidPoints; // List of TGC strip middle station points.
79 TrigL2MuonSA::TgcFit::PointArray tgcWireMidPoints; // List of TGC wire middle station points.
80 TrigL2MuonSA::TgcFit::PointArray tgcStripInnPoints; // List of TGC strip inner station points.
81 TrigL2MuonSA::TgcFit::PointArray tgcWireInnPoints; // List of TGC wire inner station points.
82
83 if (tgcHits.size()>0) {
84 // TGC data is properly read
85
86 // Split digits to Strip/Wire points.
87 if( ! prepareTgcPoints(tgcHits, tgcStripInnPoints, tgcWireInnPoints, tgcStripMidPoints, tgcWireMidPoints) ) {
88 ATH_MSG_ERROR("Preparation of Tgc points failed");
89 return StatusCode::FAILURE;
90 }
91
92
93 // Fit lines to TGC middle station
94 isMiddleFailure = false;
95 TgcFit::Status status = m_tgcFit->runTgcMiddle(tgcStripMidPoints, tgcWireMidPoints, tgcFitResult);
96 if (status == TgcFit::FIT_NONE) {
97 ATH_MSG_WARNING("Fit to TGC middle station points failed");
98 isMiddleFailure = true;
99 }
100 else if (status == TgcFit::FIT_POINT) {
101 ATH_MSG_DEBUG("Fit to TGC middle station returns only a point");
102 }
103
104 // Fit lines to TGC inner station
105 status = m_tgcFit->runTgcInner(tgcStripInnPoints, tgcWireInnPoints, tgcFitResult);
106 if (status == TgcFit::FIT_NONE) {
107 ATH_MSG_DEBUG("Fit to TGC inner station points failed");
108 }
109 else if (status == TgcFit::FIT_POINT) {
110 ATH_MSG_DEBUG("Fit to TGC inner station returns only a point");
111 }
112
113 ATH_MSG_DEBUG("tgcFitResult.tgcInn[0/1/2/3]=" << tgcFitResult.tgcInn[0] << "/" << tgcFitResult.tgcInn[1]
114 << "/" << tgcFitResult.tgcInn[2] << "/" << tgcFitResult.tgcInn[3]);
115 ATH_MSG_DEBUG("tgcFitResult.tgcMid1[0/1/2/3]=" << tgcFitResult.tgcMid1[0] << "/" << tgcFitResult.tgcMid1[1]
116 << "/" << tgcFitResult.tgcMid1[2] << "/" << tgcFitResult.tgcMid1[3]);
117 ATH_MSG_DEBUG("tgcFitResult.tgcMid2[0/1/2/3]=" << tgcFitResult.tgcMid2[0] << "/" << tgcFitResult.tgcMid2[1]
118 << "/" << tgcFitResult.tgcMid2[2] << "/" << tgcFitResult.tgcMid2[3]);
119 } else {
120 ATH_MSG_DEBUG("Skip TGC Fit due to zero-tgcHits");
121 }
122
123 if (tgcHits.size()>0 && !isMiddleFailure){
124
125 tgcFitResult.isSuccess = true;
126
127 // PT calculation by using TGC fit result
128 const double PHI_RANGE = 12./(M_PI/8.);
129 side = (tgcFitResult.tgcMid2[3]<=0) ? 0 : 1;
130 double alpha = (*m_ptEndcapLUT)->alpha(tgcFitResult.tgcMid1[3], tgcFitResult.tgcMid1[2],
131 tgcFitResult.tgcMid2[3], tgcFitResult.tgcMid2[2]);
132
133 int Octant = (int)(tgcFitResult.tgcMid1[1] / (M_PI/4.));
134 double PhiInOctant = std::abs(tgcFitResult.tgcMid1[1] - Octant * (M_PI/4.));
135 if (PhiInOctant > (M_PI/8.)) PhiInOctant = (M_PI/4.) - PhiInOctant;
136
137 int phiBin = static_cast<int>(PhiInOctant * PHI_RANGE);
138 int etaBin = static_cast<int>((std::abs(tgcFitResult.tgcMid1[0]) - 1.)/0.05);
139
140 int charge = (tgcFitResult.intercept * tgcFitResult.tgcMid2[3]) < 0.0 ? 0: 1;
141
142 tgcFitResult.tgcPT = (*m_ptEndcapLUT)->lookup(side, charge, PtEndcapLUT::TGCALPHAPOL2, etaBin, phiBin, alpha) / 1000.;
143 if (charge==0) tgcFitResult.tgcPT = -1.*tgcFitResult.tgcPT;
144
145 // Determine phi direction
146 if (std::abs(tgcFitResult.tgcMid1[3])<=ZERO_LIMIT || std::abs(tgcFitResult.tgcMid2[3])<=ZERO_LIMIT) {
147
148 tgcFitResult.isPhiDir = false;
149 if (std::abs(tgcFitResult.tgcMid1[3])>=0.) tgcFitResult.phi = tgcFitResult.tgcMid1[1];
150 if (std::abs(tgcFitResult.tgcMid2[3])>=0.) tgcFitResult.phi = tgcFitResult.tgcMid2[1];
151
152 } else {
153
154 tgcFitResult.isPhiDir = true;
155
156 if( tgcFitResult.tgcMid1[1]*tgcFitResult.tgcMid2[1] < 0
157 && std::abs(tgcFitResult.tgcMid1[1])>M_PI/2. ) {
158
159 double tmp1 = (tgcFitResult.tgcMid1[1]>0)?
160 tgcFitResult.tgcMid1[1] - M_PI : tgcFitResult.tgcMid1[1] + M_PI;
161
162 double tmp2 = (tgcFitResult.tgcMid2[1]>0)?
163 tgcFitResult.tgcMid2[1] - M_PI : tgcFitResult.tgcMid2[1] + M_PI;
164
165 double tmp = (tmp1+tmp2)/2.;
166
167 tgcFitResult.dPhidZ = (std::abs(tgcFitResult.tgcMid2[3]-tgcFitResult.tgcMid1[3]) > ZERO_LIMIT)?
168 (tmp2-tmp1)/std::abs(tgcFitResult.tgcMid2[3]-tgcFitResult.tgcMid1[3]): 0;
169
170 tgcFitResult.phi = (tmp>0.)? tmp - M_PI : tmp + M_PI;
171
172 } else {
173
174 tgcFitResult.dPhidZ = (std::abs(tgcFitResult.tgcMid2[3]-tgcFitResult.tgcMid1[3]) > ZERO_LIMIT)?
175 (tgcFitResult.tgcMid2[1]-tgcFitResult.tgcMid1[1])/std::abs(tgcFitResult.tgcMid2[3]-tgcFitResult.tgcMid1[3]): 0;
176 tgcFitResult.phi = (tgcFitResult.tgcMid2[1]+tgcFitResult.tgcMid1[1])/2.;
177
178 }
179 }
180 float X1 = tgcFitResult.tgcMid1[3] * std::cos(tgcFitResult.tgcMid1[1]);
181 float Y1 = tgcFitResult.tgcMid1[3] * std::sin(tgcFitResult.tgcMid1[1]);
182 float X2 = tgcFitResult.tgcMid2[3] * std::cos(tgcFitResult.tgcMid2[1]);
183 float Y2 = tgcFitResult.tgcMid2[3] * std::sin(tgcFitResult.tgcMid2[1]);
184 if (X1>ZERO_LIMIT && X2>ZERO_LIMIT) tgcFitResult.phiDir = (Y1/X1 + Y2/X2)/2.;
185
186 muonRoad.aw[endcap_middle][0] = tgcFitResult.slope;
187 muonRoad.bw[endcap_middle][0] = tgcFitResult.intercept;
188 muonRoad.aw[endcap_outer][0] = tgcFitResult.slope;
189 muonRoad.bw[endcap_outer][0] = tgcFitResult.intercept;
190 muonRoad.aw[endcap_extra][0] = tgcFitResult.slope;
191 muonRoad.bw[endcap_extra][0] = tgcFitResult.intercept;
192 muonRoad.aw[bee][0] = tgcFitResult.slope;
193 muonRoad.bw[bee][0] = tgcFitResult.intercept;
194 for (int i_layer=0; i_layer<N_LAYER; i_layer++) {
195 muonRoad.rWidth[endcap_middle][i_layer] = R_WIDTH_DEFAULT;
196 muonRoad.rWidth[endcap_outer][i_layer] = R_WIDTH_DEFAULT;
197 muonRoad.rWidth[endcap_extra][i_layer] = R_WIDTH_DEFAULT;
198 muonRoad.rWidth[bee][i_layer] = R_WIDTH_DEFAULT;
199 }
200
201 if( std::abs(tgcFitResult.tgcInn[3]) > ZERO_LIMIT ) {
202 muonRoad.aw[endcap_inner][0] = tgcFitResult.tgcInn[2]/tgcFitResult.tgcInn[3];
203 muonRoad.aw[barrel_inner][0] = tgcFitResult.tgcInn[2]/tgcFitResult.tgcInn[3];
204 muonRoad.aw[csc][0] = tgcFitResult.tgcInn[2]/tgcFitResult.tgcInn[3];
205 for (int i_layer=0; i_layer<N_LAYER; i_layer++) muonRoad.rWidth[endcap_inner][i_layer] = R_WIDTH_DEFAULT;
206 for (int i_layer=0; i_layer<N_LAYER; i_layer++) muonRoad.rWidth[barrel_inner][i_layer] = R_WIDTH_DEFAULT;
207 for (int i_layer=0; i_layer<N_LAYER; i_layer++) muonRoad.rWidth[csc][i_layer] = R_WIDTH_DEFAULT;
208 } else {
209 // use the back extrapolator to retrieve the Etain the Innermost
210
211 double etaMiddle = (tgcFitResult.tgcMid1[3])? tgcFitResult.tgcMid1[0] : tgcFitResult.tgcMid2[0];
212 double phiMiddle = (tgcFitResult.tgcMid1[3])? tgcFitResult.tgcMid1[1] : tgcFitResult.tgcMid2[1];
213 double eta;
214 double sigma_eta;
215 double extrInnerEta = 0;
216
218 muonSA->makePrivateStore();
219 muonSA->setSAddress(-1);
220 muonSA->setPt(tgcFitResult.tgcPT);
221 muonSA->setEtaMS(etaMiddle);
222 muonSA->setPhiMS(phiMiddle);
223 muonSA->setRMS(0.);
224 muonSA->setZMS(0.);
225
226 double phi;
227 double sigma_phi;
228
230 StatusCode sc
231 = (*m_backExtrapolatorTool)->give_eta_phi_at_vertex(muonSA, eta,sigma_eta,phi,sigma_phi,0.);
232 if (sc.isSuccess() ){
233 extrInnerEta = eta;
234 } else {
235 extrInnerEta = etaMiddle;
236 }
237 } else {
238 ATH_MSG_ERROR("Null pointer to ITrigMuonBackExtrapolator");
239 return StatusCode::FAILURE;
240 }
241
242 if (muonSA) delete muonSA;
243
244 double theta = 0.;
245 if (extrInnerEta != 0.) {
246 theta = std::atan(std::exp(-std::abs(extrInnerEta)))*2.;
247 muonRoad.aw[endcap_inner][0] = std::tan(theta)*(std::abs(extrInnerEta)/extrInnerEta);
248 muonRoad.aw[barrel_inner][0] = std::tan(theta)*(std::abs(extrInnerEta)/extrInnerEta);
249 muonRoad.aw[csc][0] = std::tan(theta)*(std::abs(extrInnerEta)/extrInnerEta);
250 } else {
251 muonRoad.aw[endcap_inner][0] = 0;
252 muonRoad.aw[barrel_inner][0] = 0;
253 muonRoad.aw[csc][0] = 0;
254 }
255
256 for (int i_layer=0; i_layer<N_LAYER; i_layer++) muonRoad.rWidth[endcap_inner][i_layer] = R_WIDTH_INNER_NO_HIT;
257 for (int i_layer=0; i_layer<N_LAYER; i_layer++) muonRoad.rWidth[barrel_inner][i_layer] = R_WIDTH_INNER_NO_HIT;
258 for (int i_layer=0; i_layer<N_LAYER; i_layer++) muonRoad.rWidth[csc][i_layer] = R_WIDTH_INNER_NO_HIT;
259
260 }
261
262 muonRoad.bw[endcap_inner][0] = 0.;
263 muonRoad.bw[barrel_inner][0] = 0.;
264 muonRoad.bw[csc][0] = 0.;
265
266 } else {
267 // If no TGC hit are available, estimate the road from RoI
268 // or if inside-out mode, width is tuned based on FTF track extrapolation resolution
269 ATH_MSG_DEBUG("Because no TGC hits are available, estimate the road from RoI");
270
271 roiEta = p_roids->eta();
272 theta = std::atan(std::exp(-std::abs(roiEta)))*2.;
273 aw = (std::abs(roiEta) > ZERO_LIMIT)? std::tan(theta)*(std::abs(roiEta)/roiEta): 0.;
274
275 muonRoad.aw[endcap_inner][0] = aw;
276 muonRoad.bw[endcap_inner][0] = 0;
277 muonRoad.aw[endcap_middle][0] = aw;
278 muonRoad.bw[endcap_middle][0] = 0;
279 muonRoad.aw[endcap_outer][0] = aw;
280 muonRoad.bw[endcap_outer][0] = 0;
281 muonRoad.aw[endcap_extra][0] = aw;
282 muonRoad.bw[endcap_extra][0] = 0;
283 muonRoad.aw[barrel_inner][0] = aw;
284 muonRoad.bw[barrel_inner][0] = 0;
285 muonRoad.aw[bee][0] = aw;
286 muonRoad.bw[bee][0] = 0;
287 muonRoad.aw[csc][0] = aw;
288 muonRoad.bw[csc][0] = 0;
289 for (int i_layer=0; i_layer<N_LAYER; i_layer++) {
290 if(insideOut) {
291 muonRoad.rWidth[endcap_inner][i_layer] = 300;
292 muonRoad.rWidth[endcap_middle][i_layer] = 400;
293 muonRoad.rWidth[endcap_outer][i_layer] = 600;
294 muonRoad.rWidth[endcap_extra][i_layer] = 400;
295 muonRoad.rWidth[barrel_inner][i_layer] = 250;
296 muonRoad.rWidth[bee][i_layer] = 500;
297 muonRoad.rWidth[csc][i_layer] = 200;
298 } else {
299 muonRoad.rWidth[endcap_inner][i_layer] = m_rWidth_TGC_Failed;
300 muonRoad.rWidth[endcap_middle][i_layer] = m_rWidth_TGC_Failed;
301 muonRoad.rWidth[endcap_outer][i_layer] = m_rWidth_TGC_Failed;
302 muonRoad.rWidth[endcap_extra][i_layer] = m_rWidth_TGC_Failed;
303 muonRoad.rWidth[barrel_inner][i_layer] = m_rWidth_TGC_Failed;
304 muonRoad.rWidth[bee][i_layer] = m_rWidth_TGC_Failed;
305 muonRoad.rWidth[csc][i_layer] = m_rWidth_TGC_Failed;
306 }
307 }
308 }
309
310 for(int i=0; i<N_STATION; i++) {
311 muonRoad.aw[i][1] = muonRoad.aw[i][0];
312 muonRoad.bw[i][1] = muonRoad.bw[i][0];
313 }
314
315 //
316 side = 0;
317 float phiMiddle = 0;
318 if( ! isMiddleFailure ) {
319 side = (tgcFitResult.tgcMid1[3]<0.)? 0 : 1;
320 phiMiddle = tgcFitResult.tgcMid1[1];
321 } else {
322 side = (p_roids->eta()<0.)? 0 : 1;
323 phiMiddle = p_roids->phi();
324 }
325 muonRoad.side = side;
326 muonRoad.phiMiddle = phiMiddle;
327 muonRoad.phiRoI = p_roids->phi();
328
329 int sector_trigger = 99;
330 int sector_overlap = 99;
331 int temp_sector=99;
332 float deltaPhi=99;
333 float tempDeltaPhi=99;
334 std::vector<Identifier> stationList;
335 std::vector<IdentifierHash> mdtHashList;
336
337 // get sector_trigger and sector_overlap by using the region selector
338 IdContext context = m_idHelperSvc->mdtIdHelper().module_context();
339
340 double etaMin = p_roids->eta()-.02;
341 double etaMax = p_roids->eta()+.02;
342 double phiMin = muonRoad.phiMiddle-.01;
343 double phiMax = muonRoad.phiMiddle+.01;
344 if(phiMax > M_PI) phiMax -= M_PI*2.;
345 if(phiMin < M_PI*-1) phiMin += M_PI*2.;
346 TrigRoiDescriptor* roi = new TrigRoiDescriptor( p_roids->eta(), etaMin, etaMax, p_roids->phi(), phiMin, phiMax );
347 const IRoiDescriptor* iroi = static_cast<IRoiDescriptor*> (roi);
348 if (iroi) m_regionSelector->lookup(ctx)->HashIDList(*iroi, mdtHashList);
349 else {
350 TrigRoiDescriptor fullscan_roi( true );
351 m_regionSelector->lookup(ctx)->HashIDList(fullscan_roi, mdtHashList);
352 }
353 if(roi) delete roi;
354
355 for(const IdentifierHash& i_hash : mdtHashList ){
356 Identifier id;
357 int convert = m_idHelperSvc->mdtIdHelper().get_id(i_hash, id, &context);
358
359 if(convert!=0) ATH_MSG_ERROR("problem converting hash list to id");
360
361 muonRoad.stationList.push_back(id);
362 std::string name = m_idHelperSvc->mdtIdHelper().stationNameString(m_idHelperSvc->mdtIdHelper().stationName(id));
363 if ( name.substr(0, 1) == "B" ) continue;
364 if ( name.substr(1, 1) != "M" ) continue;
365 int stationPhi = m_idHelperSvc->mdtIdHelper().stationPhi(id);
366 float floatPhi = (stationPhi-1)*M_PI/4;
367 if (name[2]=='S' || name[2]=='E') floatPhi = floatPhi + M_PI/8;
368 tempDeltaPhi = std::abs(floatPhi-muonRoad.phiMiddle);
369 if (phiMiddle<0) tempDeltaPhi = std::abs(floatPhi-muonRoad.phiMiddle-2*M_PI);
370 if(tempDeltaPhi > M_PI) tempDeltaPhi = std::abs(tempDeltaPhi - 2*M_PI);
371
372 int LargeSmall = 0;
373 if(name[2]=='S' || name[2]=='E') LargeSmall = 1;
374 int sector = (stationPhi-1)*2 + LargeSmall;
375 if(sector_trigger == 99)
376 sector_trigger = sector;
377 else if(sector_trigger != sector)
378 sector_overlap = sector;
379
380 if(tempDeltaPhi < deltaPhi){
381 deltaPhi = tempDeltaPhi;
382 temp_sector = sector;
383 muonRoad.LargeSmall=LargeSmall;
384 }
385
386 }
387 if(temp_sector != sector_trigger){
388 sector_overlap = sector_trigger;
389 sector_trigger = temp_sector;
390 }
391
392 muonRoad.MDT_sector_trigger = sector_trigger;
393 muonRoad.MDT_sector_overlap = sector_overlap;
394
395 if (insideOut) {
396 muonRoad.side = (muonRoad.extFtfMiddleEta<0.)? 0 : 1;
397 muonRoad.phiMiddle = muonRoad.extFtfMiddlePhi;
398 muonRoad.phiRoI = p_roids->phi();
399 for (int i_sector=0; i_sector<N_SECTOR; i_sector++) {
400 ATH_MSG_DEBUG("Use aw_ftf and bw_ftf as aw and bw");
401 muonRoad.aw[endcap_inner][i_sector] = muonRoad.aw_ftf[3][0];
402 muonRoad.bw[endcap_inner][i_sector] = muonRoad.bw_ftf[3][0];
403 muonRoad.aw[endcap_middle][i_sector] = muonRoad.aw_ftf[4][0];
404 muonRoad.bw[endcap_middle][i_sector] = muonRoad.bw_ftf[4][0];
405 muonRoad.aw[endcap_outer][i_sector] = muonRoad.aw_ftf[5][0];
406 muonRoad.bw[endcap_outer][i_sector] = muonRoad.bw_ftf[5][0];
407 muonRoad.aw[endcap_extra][i_sector] = muonRoad.aw_ftf[6][0];
408 muonRoad.bw[endcap_extra][i_sector] = muonRoad.bw_ftf[6][0];
409 muonRoad.aw[csc][i_sector] = muonRoad.aw_ftf[7][0];
410 muonRoad.bw[csc][i_sector] = muonRoad.bw_ftf[7][0];
411 muonRoad.aw[bee][i_sector] = muonRoad.aw_ftf[8][0];
412 muonRoad.bw[bee][i_sector] = muonRoad.bw_ftf[8][0];
413 muonRoad.aw[barrel_inner][i_sector] = muonRoad.aw_ftf[0][0];
414 muonRoad.bw[barrel_inner][i_sector] = muonRoad.bw_ftf[0][0];
415 }
416 }
417
418 //
419 return StatusCode::SUCCESS;
420}
421
422// --------------------------------------------------------------------------------
423// --------------------------------------------------------------------------------
424
426 TrigL2MuonSA::TgcFit::PointArray& tgcStripInnPoints,
427 TrigL2MuonSA::TgcFit::PointArray& tgcWireInnPoints,
428 TrigL2MuonSA::TgcFit::PointArray& tgcStripMidPoints,
429 TrigL2MuonSA::TgcFit::PointArray& tgcWireMidPoints) const
430{
431 const double PHI_BOUNDARY = 0.2;
432
433 tgcStripMidPoints.clear();
434 tgcStripInnPoints.clear();
435 tgcWireMidPoints.clear();
436 tgcWireInnPoints.clear();
437
438 // loop over TGC digits.
439 unsigned int iHit;
440 for (iHit = 0; iHit < tgcHits.size(); iHit++)
441 {
442 // Get the digit point.
443 const TrigL2MuonSA::TgcHitData& hit = tgcHits[iHit];
444
445 // reject width=0 hits
446 const double ZERO_LIMIT = 1e-5;
447 if( std::abs(hit.width) < ZERO_LIMIT ) continue;
448
449 double w = 12.0 / hit.width / hit.width;
450 if (hit.isStrip)
451 {
452 w *= hit.r * hit.r;
453 double phi = hit.phi;
454 if( phi < 0 && ( (M_PI+phi)<PHI_BOUNDARY) ) phi += M_PI*2;
455 if ( hit.sta < 3 ) { tgcStripMidPoints.push_back(TgcFit::Point(iHit + 1, hit.sta, hit.z, phi, w)); }
456 else if ( hit.sta ==3 ) { tgcStripInnPoints.push_back(TgcFit::Point(iHit + 1, hit.sta, hit.z, phi, w)); }
457 }
458 else
459 {
460 if ( hit.sta < 3 ) { tgcWireMidPoints.push_back(TgcFit::Point(iHit + 1, hit.sta, hit.z, hit.r, w)); }
461 else if ( hit.sta ==3 ) { tgcWireInnPoints.push_back(TgcFit::Point(iHit + 1, hit.sta, hit.z, hit.r, w)); }
462 }
463 }
464
465 ATH_MSG_DEBUG(", tgcStripMidPoints.size()=" << tgcStripMidPoints.size() <<
466 ", tgcStripInnPoints.size()=" << tgcStripInnPoints.size() <<
467 ", tgcWireMidPoints.size()=" << tgcWireMidPoints.size() <<
468 ", tgcWireInnPoints.size()=" << tgcWireInnPoints.size());
469
470 return true;
471}
#define M_PI
Scalar eta() const
pseudorapidity method
Scalar deltaPhi(const MatrixBase< Derived > &vec) const
Scalar phi() const
phi method
Scalar theta() const
theta method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
double charge(const T &p)
Definition AtlasPID.h:997
static Double_t sc
Athena::TPCnvVers::Current TrigRoiDescriptor
const float ZERO_LIMIT
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
Describes the API of the Region of Ineterest geometry.
This class saves the "context" of an expanded identifier (ExpandedIdentifier) for compact or hash ver...
Definition IdContext.h:26
This is a "hash" representation of an Identifier.
virtual double phi() const override final
Methods to retrieve data members.
virtual double eta() const override final
void makePrivateStore()
Create a new (empty) private store for this object.
std::vector< Identifier > stationList
Definition MuonRoad.h:102
double aw[N_STATION][N_SECTOR]
Definition MuonRoad.h:83
double bw_ftf[N_STATION][N_SECTOR]
Definition MuonRoad.h:95
double aw_ftf[N_STATION][N_SECTOR]
Definition MuonRoad.h:94
double bw[N_STATION][N_SECTOR]
Definition MuonRoad.h:84
double rWidth[N_STATION][N_LAYER]
Definition MuonRoad.h:86
std::vector< Point > PointArray
Definition TgcFit.h:86
const ToolHandle< ITrigMuonBackExtrapolator > * m_backExtrapolatorTool
bool prepareTgcPoints(const TrigL2MuonSA::TgcHits &tgcHits, TrigL2MuonSA::TgcFit::PointArray &tgcStripInnPoints, TrigL2MuonSA::TgcFit::PointArray &tgcWireInnPoints, TrigL2MuonSA::TgcFit::PointArray &tgcStripMidPoints, TrigL2MuonSA::TgcFit::PointArray &tgcWireMidPoints) const
ToolHandle< TgcFit > m_tgcFit
ServiceHandle< Muon::IMuonIdHelperSvc > m_idHelperSvc
TgcRoadDefiner(const std::string &type, const std::string &name, const IInterface *parent)
StatusCode defineRoad(const EventContext &ctx, const TrigRoiDescriptor *p_roids, const bool insideOut, const TrigL2MuonSA::TgcHits &tgcHits, TrigL2MuonSA::MuonRoad &muonRoad, TrigL2MuonSA::TgcFitResult &tgcFitResult) const
ToolHandle< IRegSelTool > m_regionSelector
virtual StatusCode initialize() override
nope - should be used for standalone also, perhaps need to protect the class def bits ifndef XAOD_ANA...
void setEtaMS(float value)
Set the eta at muon spectrometer.
void setRMS(float value)
Set the R at muon spectrometer.
void setSAddress(int value)
Set the station address of the muon.
void setPhiMS(float value)
Set the phi at muon spectrometer.
void setPt(float pt)
Set the transverse momentum ( ) of the muon.
void setZMS(float value)
Set the Z at muon spectrometer.
constexpr int N_STATION
Definition MuonRoad.h:15
std::vector< TgcHitData > TgcHits
Definition TgcData.h:43
constexpr int N_LAYER
Definition MuonRoad.h:17
constexpr int N_SECTOR
Definition MuonRoad.h:16
@ BarrelInner
Inner station in the barrel spectrometer.
@ EndcapOuter
Outer station in the endcap spectrometer.
@ CSC
CSC measurement point.
@ EndcapMiddle
Middle station in the endcap spectrometer.
@ BEE
BEE measurement point.
@ EndcapExtra
Extra station in the endcap spectrometer.
@ EndcapInner
Inner station in the endcap spectrometer.
L2StandAloneMuon_v2 L2StandAloneMuon
Define the latest version of the muon SA class.