ATLAS Offline Software
Loading...
Searching...
No Matches
MuonCalibExtendedTrack.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7#include <sstream>
8
16
17namespace MuonCalib {
18 const std::vector<std::shared_ptr<MuonCalibExtendedSegment>>& MuonCalibExtendedTrack::associatedSegments() const {
20 }
21
26 m_summary.isTrack = true;
28
29 MuonCalibTrackSummary::ChamberHitSummary* currentChamberSummary = nullptr;
30
31 for (const MuonCalibTrack_E::CalibHitPtr& hit : track.hits()) {
32 int type = hit->type();
33 // count scatters
34 if (type > 10) {
35 ++m_summary.nscatters;
36 continue;
37 }
38
39 const MuonFixedId& id = hit->identify();
40 if (!id.isValid()) {
41 if (type == 3)
42 ++m_summary.npseudo;
43 else {
44 // these are the ID hits
45 if (hit->position().perp() < 200.) {
46 ++m_summary.npixel;
47 } else if (hit->position().perp() < 540.) {
48 ++m_summary.nsct;
49 } else if (hit->position().perp() < 1100.) {
50 if (std::abs(hit->position().z()) < 750.)
51 ++m_summary.ntrtBarrel;
52 else
53 ++m_summary.ntrtEndcap;
54 }
55 }
56 continue;
57 }
58
59 MuonFixedId chId = idManip.chamberIdentifier(id);
60
61 bool measuresPhi = idManip.measuresPhi(id);
62 bool isMdt = id.is_mdt();
63 bool isFirst = isMdt ? id.mdtMultilayer() == 1 : !measuresPhi;
64
65 // check whether first chamber or new chamber
66 if (!currentChamberSummary || currentChamberSummary->chId != chId) {
67 m_summary.chamberHitSummary.push_back(MuonCalibTrackSummary::ChamberHitSummary(chId));
68 currentChamberSummary = &m_summary.chamberHitSummary.back();
69 }
70
72 isFirst ? currentChamberSummary->etaProjection() : currentChamberSummary->phiProjection();
73
75 MuonFixedId chamberId = idManip.chamberIdentifier(id);
76
77 if (type == 6) {
78 ++m_summary.nholes;
79 ++proj.nholes;
80 } else if (type == 4) {
81 ++m_summary.noutliers;
82 ++proj.noutliers;
83 } else if (type == 5) {
84 ++m_summary.noutliers;
85 ++proj.ndeltas;
86 } else {
87 IdHitMap::iterator hitPos = m_hitIdMap.find(id);
88 if (hitPos == m_hitIdMap.end()) m_hitIdMap[id] = hit;
89 m_hitsPerChamber[chamberId].emplace_back(hit);
90
91 ++proj.nhits;
92 // hits counts
93 ++m_summary.nhits;
94 if ((*hit).position().y() < 0.)
95 ++m_summary.nhitsLowerHemisphere;
96 else
97 ++m_summary.nhitsUpperHemisphere;
98
99 if (id.is_mdt()) {
100 ++m_summary.nmdtHits;
101 m_mdtHitsPerStationIndex[stationIndex].emplace_back(hit);
102 if (idManip.isEndcap(id)) {
103 if (id.eta() < 0.)
104 m_summary.hasEndcapA = true;
105 else
106 m_summary.hasEndcapC = true;
107 if (stationIndex == MuonFixedIdManipulator::EIA || stationIndex == MuonFixedIdManipulator::EIC ||
108 stationIndex == MuonFixedIdManipulator::EMA || stationIndex == MuonFixedIdManipulator::EMC)
109 m_summary.hasEndcapLayersWithTGC = true;
110 } else {
111 m_summary.hasBarrel = true;
112 if (stationIndex == MuonFixedIdManipulator::BM || stationIndex == MuonFixedIdManipulator::BO)
113 m_summary.hasBarrelLayersWithRPC = true;
114 }
115
116 m_summary.precisionStationLayers.insert(stationIndex);
117 if (currentChamberSummary->nMdtHitsMl1() > 0 && currentChamberSummary->nMdtHitsMl2() &&
118 currentChamberSummary->nhits() > 3) {
119 m_summary.goodPrecisionStationLayers.insert(stationIndex);
120 }
121 } else if (id.is_rpc()) {
122 if (id.rpcMeasuresPhi()) {
123 ++m_summary.nrpcPhiHits;
124 } else {
125 ++m_summary.nrpcEtaHits;
126 }
127 m_rpcHitsPerStationIndex[stationIndex].emplace_back(hit);
128 m_summary.hasBarrel = true;
129
131 if (measuresPhi) m_summary.phiStationLayers.insert(phiStationIndex);
132 if (currentChamberSummary->netaHits() > 0 && currentChamberSummary->nphiHits() > 0) {
133 m_summary.phiEtaStationLayers.insert(phiStationIndex);
134 }
135 } else if (id.is_tgc()) {
136 if (id.tgcIsStrip())
137 ++m_summary.ntgcPhiHits;
138 else
139 ++m_summary.ntgcEtaHits;
140 m_tgcHitsPerStationIndex[stationIndex].emplace_back(hit);
141
142 if (id.eta() < 0.)
143 m_summary.hasEndcapA = true;
144 else
145 m_summary.hasEndcapC = true;
146 m_summary.hasEndcapLayersWithTGC = true;
147
149 if (measuresPhi) m_summary.phiStationLayers.insert(phiStationIndex);
150 if (currentChamberSummary->netaHits() > 0 && currentChamberSummary->nphiHits() > 0) {
151 m_summary.phiEtaStationLayers.insert(phiStationIndex);
152 }
153 } else if (id.is_csc()) {
154 if (id.cscMeasuresPhi())
155 ++m_summary.ncscPhiHits;
156 else
157 ++m_summary.ncscEtaHits;
158 m_cscHitsPerStationIndex[stationIndex].emplace_back(hit);
159 if (id.eta() < 0.)
160 m_summary.hasEndcapA = true;
161 else
162 m_summary.hasEndcapC = true;
163
165 if (measuresPhi) { m_summary.phiStationLayers.insert(phiStationIndex); }
166 m_summary.precisionStationLayers.insert(stationIndex);
167 if (currentChamberSummary->netaHits() > 2 && currentChamberSummary->nphiHits() > 2) {
168 m_summary.goodPrecisionStationLayers.insert(stationIndex);
169 m_summary.phiEtaStationLayers.insert(phiStationIndex);
170 }
171 }
172 }
173 }
174
175 const StationIndexHitsMap& mdtStHitMap = mdtStationIndexHitsMap();
176 const StationIndexHitsMap& rpcStHitMap = rpcStationIndexHitsMap();
177 const StationIndexHitsMap& tgcStHitMap = tgcStationIndexHitsMap();
178 const StationIndexHitsMap& cscStHitMap = cscStationIndexHitsMap();
179
180 MuonCalibSimpleHoleSearch holeSearch;
181 double tolerance = 100.;
182 MuonCalibSimpleHoleSearch::ResultVec intersectedLayers = holeSearch.intersectsWithGeometry(position(), direction(), tolerance);
183 MuonCalibSimpleHoleSearch::ResultIt rit = intersectedLayers.begin();
184 MuonCalibSimpleHoleSearch::ResultIt rit_end = intersectedLayers.end();
185 for (; rit != rit_end; ++rit) {
186 m_intersectedLayers[rit->stationLayerId] = *rit;
187 const MuonFixedId& stId = rit->stationLayerId;
188 MuonFixedIdManipulator::StationIndex stIndex = idManip.stationLayerIndex(stId);
189 bool isEndcap = idManip.isEndcap(stId);
190
191 StationIndexHitsMap::const_iterator mdtIt = mdtStHitMap.find(stIndex);
192 unsigned int nmdtHits = mdtIt != mdtStHitMap.end() ? mdtIt->second.size() : 0;
193
194 StationIndexHitsMap::const_iterator rpcIt = rpcStHitMap.find(stIndex);
195 unsigned int nrpcHits = rpcIt != rpcStHitMap.end() ? rpcIt->second.size() : 0;
196
197 StationIndexHitsMap::const_iterator tgcIt = tgcStHitMap.find(stIndex);
198 unsigned int ntgcHits = tgcIt != tgcStHitMap.end() ? tgcIt->second.size() : 0;
199
200 StationIndexHitsMap::const_iterator cscIt = cscStHitMap.find(stIndex);
201 unsigned int ncscHits = cscIt != cscStHitMap.end() ? cscIt->second.size() : 0;
202
203 // flag intersected layers without any hits
204 unsigned int ntrig = isEndcap ? ntgcHits : nrpcHits;
205 unsigned int nprec = nmdtHits + ncscHits;
206 if (nprec + ntrig == 0) m_intersectedLayersWithoutHits.insert(stIndex);
207
208 // flag layers where we would expect RPC hits but they are not on the track
209 if (!isEndcap && stIndex != MuonFixedIdManipulator::BI && stIndex != MuonFixedIdManipulator::BE && nrpcHits == 0)
210 m_intersectedRpcLayerWithoutHits.insert(stIndex);
211
212 // flag layers where we would expect TGC hits but they are not on the track
213 if (isEndcap && stIndex != MuonFixedIdManipulator::EIA && stIndex != MuonFixedIdManipulator::EIC && ntgcHits == 0)
214 m_intersectedTgcLayerWithoutHits.insert(stIndex);
215 }
216 }
217
218 std::string MuonCalibExtendedTrack::dump() const {
219 std::ostringstream sout;
220 sout << dumpPars() << std::endl << dumpSummary() << dumpIntersects() << std::endl;
221 return sout.str();
222 }
223
225 std::ostringstream sout;
226 double sign = qOverP() < 0. ? -1 : 1;
227 sout << "Author " << std::setw(4) << author() << std::setprecision(4) << " chi2 " << chi2() << " ndof " << ndof() << " r "
228 << (int)position().perp() << " z " << (int)z0() << std::setprecision(5) << " phi " << phi() << " theta " << theta()
229 << std::setw(6) << " q*mom " << (int)p() * sign << " pt " << std::setw(5) << (int)pt() << " association: segments "
230 << m_associatedSegments.size() << " tracks " << m_associatedTracks.size();
231 return sout.str();
232 }
233
235 std::ostringstream sout;
236 sout << " " << m_summary.dump();
237 return sout.str();
238 }
239
241 MuonFixedIdPrinter printer{};
243
244 const StationIndexHitsMap& mdtStHitMap = mdtStationIndexHitsMap();
245 const StationIndexHitsMap& rpcStHitMap = rpcStationIndexHitsMap();
246 const StationIndexHitsMap& tgcStHitMap = tgcStationIndexHitsMap();
247 const StationIndexHitsMap& cscStHitMap = cscStationIndexHitsMap();
248
249 std::ostringstream sout;
250 sout << " Intersected layers: " << m_intersectedLayers.size() << " without: any hits " << m_intersectedLayersWithoutHits.size()
251 << " rpc hits " << m_intersectedRpcLayerWithoutHits.size() << " tgc hits " << m_intersectedTgcLayerWithoutHits.size()
252 << std::endl;
253 unsigned int nentries = m_intersectedLayers.size();
254 unsigned int currentEntry = 0;
255 StationIntersectedLayerMap::const_iterator it = m_intersectedLayers.begin();
256 StationIntersectedLayerMap::const_iterator it_end = m_intersectedLayers.end();
257 sout.setf(std::ios::left);
258 for (; it != it_end; ++it) {
259 const MuonFixedId& stId = it->first;
261 bool isEndcap = idManip.isEndcap(stId);
262
263 StationIndexHitsMap::const_iterator mdtIt = mdtStHitMap.find(stIndex);
264 unsigned int nmdtHits = mdtIt != mdtStHitMap.end() ? mdtIt->second.size() : 0;
265
266 StationIndexHitsMap::const_iterator rpcIt = rpcStHitMap.find(stIndex);
267 unsigned int nrpcHits = rpcIt != rpcStHitMap.end() ? rpcIt->second.size() : 0;
268
269 StationIndexHitsMap::const_iterator tgcIt = tgcStHitMap.find(stIndex);
270 unsigned int ntgcHits = tgcIt != tgcStHitMap.end() ? tgcIt->second.size() : 0;
271
272 StationIndexHitsMap::const_iterator cscIt = cscStHitMap.find(stIndex);
273 unsigned int ncscHits = cscIt != cscStHitMap.end() ? cscIt->second.size() : 0;
274
275 unsigned int ntrig = isEndcap ? ntgcHits : nrpcHits;
276 unsigned int nprec = nmdtHits + ncscHits;
277
278 sout << " " << std::setw(5) << printer.stationLayerIdentifier(stId);
279 if (nprec + ntrig != 0)
280 sout << " presicion hits " << std::setw(3) << nprec << " trigger hits " << std::setw(3) << ntrig;
281 else
282 sout << std::setw(36) << " no hits in layer";
283 sout << " intersect position " << it->second.intersectPosition;
284
285 // increase count before check to allow equals check
286 ++currentEntry;
287 if (currentEntry != nentries) sout << std::endl;
288 }
289 return sout.str();
290 }
291
293 const IdHitsMap& chHitMap = track.hitsPerChamberMap();
295
297
298 std::set<MuonFixedId> sharedChambers;
299
300 // loop over stations and check whether they are also present in other list
301 for (const auto& rit : m_hitsPerChamber) {
302 IdHitsMap::const_iterator pos = chHitMap.find(rit.first);
303 if (pos != chHitMap.end()) {
304 // shared chamber add to list
305 sharedChambers.insert(rit.first);
306
307 std::set<MuonFixedId> foundIds, sharedEtaLayers, sharedPhiLayers, firstEtaLayers, firstPhiLayers, secondEtaLayers,
308 secondPhiLayers;
309 for (const CalibHitE_Ptr& calib_hit : rit.second) {
310 const MuonFixedId& id = calib_hit->identify();
311 bool measuresPhi = manip.measuresPhi(id);
312 MuonFixedId layerId = manip.moduleIdentifier(id, true);
313
314 MuonCalibTrack_E::HitVector::const_iterator hit =
315 std::find_if(pos->second.begin(), pos->second.end(),
316 [&id](const MuonCalibTrack_E::CalibHitPtr& hit) { return id == hit->identify(); });
317 if (hit != pos->second.end()) {
318 const MuonCalibTrack_E::CalibHitPtr test_hit = (*hit);
319 foundIds.insert(id);
320
321 if (id.is_mdt()) {
323 if (!(std::abs(calib_hit->driftRadius()) > 2. && std::abs(test_hit->driftRadius()) > 2. &&
324 calib_hit->driftRadius() * test_hit->driftRadius() < 0.))
325 sharedEtaLayers.insert(id);
326 } else {
327 if (measuresPhi)
328 sharedPhiLayers.insert(layerId);
329 else
330 sharedEtaLayers.insert(layerId);
331 }
332 } else {
333 if (!id.is_mdt() && foundIds.count(layerId)) continue;
334 if (measuresPhi)
335 firstPhiLayers.insert(layerId);
336 else
337 firstEtaLayers.insert(layerId);
338 }
339 }
340
341 for (const MuonCalibTrack_E::CalibHitPtr& it2 : pos->second) {
342 MuonFixedId id = it2->identify().is_mdt() ? it2->identify() : manip.moduleIdentifier(it2->identify(), true);
343 if (!foundIds.count(id)) {
344 bool measuresPhi = manip.measuresPhi(id);
345 if (measuresPhi)
346 secondPhiLayers.insert(id);
347 else
348 secondEtaLayers.insert(id);
349 }
350 }
351 MuonCalib::MuonCalibExtendedTrackOverlap::TechnologyOverlap *phi_or{nullptr}, *eta_or{nullptr};
352 if (rit.first.is_mdt()) {
353 eta_or = &overlap.mdt;
354 } else if (rit.first.is_rpc()) {
355 phi_or = &overlap.rpcPhi;
356 eta_or = &overlap.rpcEta;
357 } else if (rit.first.is_tgc()) {
358 phi_or = &overlap.tgcPhi;
359 eta_or = &overlap.tgcEta;
360 } else if (rit.first.is_csc()) {
361 phi_or = &overlap.cscPhi;
362 eta_or = &overlap.cscEta;
363 }
364 if (eta_or) {
365 (*eta_or).shared += sharedEtaLayers.size();
366 ;
367 (*eta_or).first += firstEtaLayers.size();
368 (*eta_or).second += secondEtaLayers.size();
369 }
370 if (phi_or) {
371 (*phi_or).shared += sharedPhiLayers.size();
372 (*phi_or).first += firstPhiLayers.size();
373 (*phi_or).second += secondPhiLayers.size();
374 }
375 }
377 else {
378 std::set<MuonFixedId> foundIds;
379 for (const MuonCalibTrack_E::CalibHitPtr& it1 : rit.second) {
380 if (rit.first.is_mdt()) {
381 ++overlap.mdt.first;
382 } else {
383 MuonFixedId id = manip.moduleIdentifier(it1->identify(), true);
384 if (foundIds.count(id)) continue;
385 foundIds.insert(id);
386 bool measuresPhi = manip.measuresPhi(it1->identify());
387
388 if (rit.first.is_rpc()) {
389 if (!measuresPhi)
390 ++overlap.rpcEta.first;
391 else
392 ++overlap.rpcPhi.first;
393 } else if (rit.first.is_tgc()) {
394 if (!measuresPhi)
395 ++overlap.tgcEta.first;
396 else
397 ++overlap.tgcPhi.first;
398 } else if (rit.first.is_csc()) {
399 if (!measuresPhi)
400 ++overlap.cscEta.first;
401 else
402 ++overlap.cscPhi.first;
403 }
404 }
405 }
406 }
407 }
408
409 // loop over stations and check whether they are also present in other list
410 for (const auto& rit : chHitMap) {
411 // skip already handled chambers
412 if (sharedChambers.count(rit.first)) continue;
413 std::set<MuonFixedId> foundIds;
414 for (const MuonCalibTrack_E::CalibHitPtr& it1 : rit.second) {
415 if (rit.first.is_mdt()) {
416 ++overlap.mdt.second;
417 } else {
418 MuonFixedId id = manip.moduleIdentifier(it1->identify(), true);
419 if (foundIds.count(id)) continue;
420 foundIds.insert(id);
421 bool measuresPhi = manip.measuresPhi(it1->identify());
422 if (rit.first.is_rpc()) {
423 if (!measuresPhi)
424 ++overlap.rpcEta.second;
425 else
426 ++overlap.rpcPhi.second;
427 } else if (rit.first.is_tgc()) {
428 if (!measuresPhi)
429 ++overlap.tgcEta.second;
430 else
431 ++overlap.tgcPhi.second;
432 } else if (rit.first.is_csc()) {
433 if (!measuresPhi)
434 ++overlap.cscEta.second;
435 else
436 ++overlap.cscPhi.second;
437 }
438 }
439 }
440 }
441 return overlap;
442 }
443
445
447
450 return false;
451 }
452} // namespace MuonCalib
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
Scalar theta() const
theta method
bool isValid(const T &p)
Av: we implement here an ATLAS-sepcific convention: all particles which are 99xxxxx are fine.
Definition AtlasPID.h:878
if(febId1==febId2)
int sign(int a)
A segment plus everything one can dream of knowing about it.
int pdgCode() const
returns trackparameter d0 at IP
std::vector< std::shared_ptr< MuonCalibExtendedTrack > > m_associatedTracks
list of tracks associated with this track
std::string dumpIntersects() const
dump track intersects to string
StationIndexHitsMap m_mdtHitsPerStationIndex
hit information per station layer index
std::string dumpPars() const
dump track parameters to string
IdHitMap m_hitIdMap
hit information per station
MuonCalibTrackSummary m_summary
track summary
std::string dumpSummary() const
dump track summary to string
bool isIDConfirmed() const
get associated track for the give author, returns zero if not found
IdHitsMap m_hitsPerChamber
hit information per station
const StationIndexHitsMap & tgcStationIndexHitsMap() const
access to hits per station layer index (tgc)
const StationIndexHitsMap & cscStationIndexHitsMap() const
access to hits per station layer index (csc)
void addSegment(MuonCalibExtendedSegment *seg)
access to list of the tracks that are associated to this track *‍/ const std::vector<std::shared_ptr<...
StationIntersectedLayerMap m_intersectedLayers
map with all station layers intersected by track
StationIndexSet m_intersectedTgcLayerWithoutHits
set with all tgc chamber identifiers of layers intersected by the track without hits
const StationIndexHitsMap & mdtStationIndexHitsMap() const
access to hits per station layer index (mdt)
MuonCalibExtendedTrackOverlap calculateHitOverlap(const MuonCalibExtendedTrack &track) const
calculate hit overlap between two tracks
MuonCalibExtendedTrack(const MuonCalibTrack_E &track, int pdgCode=0, int barCode=0) ATLAS_CTORDTOR_NOT_THREAD_SAFE
Constructor taking input track.
const std::vector< std::shared_ptr< MuonCalibExtendedSegment > > & associatedSegments() const
access to list of the segment that are associated to this track
std::vector< std::shared_ptr< MuonCalibExtendedSegment > > m_associatedSegments
list of segments associated with this track
const StationIndexHitsMap & rpcStationIndexHitsMap() const
access to hits per station layer index (rpc)
StationIndexSet m_intersectedRpcLayerWithoutHits
set with all rpc chamber identifiers of layers intersected by the track without hits
std::string dump() const
dump all information to string
StationIndexSet m_intersectedLayersWithoutHits
set with all station layers intersected by the track without hits
ResultVec intersectsWithGeometry(const Amg::Vector3D &parPos, const Amg::Vector3D &parDir, double tolerance=1e9)
int author() const
returns the author
int ndof() const
returns the number of degrees of freedom
const Amg::Vector3D & position() const
position of perigee of track
MuonCalibTrack_E()=default
default constructor
float pt() const
returns pt
float qOverP() const
returns trackparameter q/p
float z0() const
returns trackparameter z0
float p() const
returns momentum
std::shared_ptr< const MuonCalibHit_E > CalibHitPtr
StationIndex
enum defining station layers
PhiStationIndex phiStationLayerIndex(const MuonFixedId &id) const
return phi station layer index for a give identifier
MuonFixedId chamberIdentifier(const MuonFixedId &id) const
returns chamber Identifier for the give id, same as stationIdentifier RPC: includes doubletR
StationIndex stationLayerIndex(const MuonFixedId &id) const
returns station layer index for a give identifier
PhiStationIndex
enum defining trigger phi layers
bool isEndcap(const MuonFixedId &id) const
returns whether this is a phi measurement
MuonFixedId moduleIdentifier(const MuonFixedId &id, bool includeMeasuresPhi=false) const
returns layer Identifier for the give id MDT: station name/eta/phi/ml/lay RPC: station name/eta/phi/d...
bool measuresPhi(const MuonFixedId &id) const
returns whether this is a phi measurement
std::string stationLayerIdentifier(const MuonFixedId &id) const
prints a station layer identifier for the give identifier
Implements fixed identifiers not dependent upon Athena Identifier for internal use in the calibration...
Definition MuonFixedId.h:50
double chi2(TH1 *h0, TH1 *h1)
CscCalcPed - algorithm that finds the Cathode Strip Chamber pedestals from an RDO.
std::shared_ptr< const MuonCalibHit_E > CalibHitE_Ptr
hit information per station layer index
std::map< MuonFixedIdManipulator::StationIndex, std::vector< CalibHitE_Ptr > > StationIndexHitsMap
std::map< MuonFixedId, std::vector< CalibHitE_Ptr > > IdHitsMap
hit information per station