ATLAS Offline Software
Loading...
Searching...
No Matches
MuonTruthSegmentCreationAlg.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6
8#include "Identifier/Identifier.h"
16#include <cstddef>
17#include <memory>
19namespace {
20
21 // Only reject muons from light quark deays
22 const std::set<int> bad_origins{
23 ParticleOrigin ::LightMeson, // 23
24 ParticleOrigin ::StrangeMeson, // 24
25 ParticleOrigin ::LightBaryon, // 30
26 ParticleOrigin ::StrangeBaryon, // 31
27 ParticleOrigin ::PionDecay, // 34
28 ParticleOrigin ::NucReact, // 41
29 ParticleOrigin ::PiZero, // 42
30 };
31
32} // namespace
33
34namespace Muon {
35 using namespace MuonStationIndex;
36 // Initialize method:
38 ATH_CHECK(m_muonTruth.initialize());
39
41 ATH_CHECK(m_SDO_TruthNames.initialize());
43
44 ATH_CHECK(m_detMgrKey.initialize());
45 ATH_CHECK(m_idHelperSvc.retrieve());
46
47 ATH_CHECK(m_truthOriginKey.initialize());
48 if(m_idHelperSvc->hasRPC()) {m_truthHitsKeyArray.emplace_back(m_muonTruth, "truthRpcHits");}
49 if(m_idHelperSvc->hasTGC()) {m_truthHitsKeyArray.emplace_back(m_muonTruth, "truthTgcHits");}
50 if(m_idHelperSvc->hasMDT()) {m_truthHitsKeyArray.emplace_back(m_muonTruth, "truthMdtHits");}
51 if(m_idHelperSvc->hasMM()) {m_truthHitsKeyArray.emplace_back(m_muonTruth, "truthMMHits");}
52 if(m_idHelperSvc->hasSTGC()) {m_truthHitsKeyArray.emplace_back(m_muonTruth, "truthStgcHits");}
53 if(m_idHelperSvc->hasCSC()) {m_truthHitsKeyArray.emplace_back(m_muonTruth, "truthCscHits");}
54 ATH_CHECK(m_truthHitsKeyArray.initialize());
55 return StatusCode::SUCCESS;
56 }
57
58 // Execute method:
59 StatusCode MuonTruthSegmentCreationAlg::execute(const EventContext& ctx) const {
60
61 const xAOD::TruthParticleContainer* muonTruthContainer{nullptr};
62 ATH_CHECK(SG::get(muonTruthContainer, m_muonTruth, ctx));
63
65 ATH_CHECK(truthOrigin.isPresent());
66
67 // create output container
69 ATH_CHECK(segmentContainer.record(std::make_unique<xAOD::MuonSegmentContainer>(),
70 std::make_unique<xAOD::MuonSegmentAuxContainer>()));
71 ATH_MSG_DEBUG("Recorded MuonSegmentContainer with key: " << segmentContainer.name());
72
73 size_t itr = 0;
74 for (const xAOD::TruthParticle* truthParticle : *muonTruthContainer) {
75
76 const int iOrigin = truthOrigin(*truthParticle);
77 bool goodMuon = bad_origins.find(iOrigin) == bad_origins.end();
78
79 // create segments
80 if (goodMuon) {
81 ElementLink<xAOD::TruthParticleContainer> truthLink(*muonTruthContainer, itr);
82 truthLink.toPersistent();
83
84 ChamberIdMap ids{};
85 ATH_CHECK(fillChamberIdMap(ctx, *truthParticle, ids));
86 ATH_CHECK(createSegments(ctx, truthLink, ids, *segmentContainer));
87 }
88 itr++;
89 }
90
91 ATH_MSG_DEBUG("Registered " << segmentContainer->size() << " truth muon segments ");
92
93 return StatusCode::SUCCESS;
94 }
95
96 StatusCode MuonTruthSegmentCreationAlg::fillChamberIdMap(const EventContext& ctx,
97 const xAOD::TruthParticle& truthParticle,
98 ChamberIdMap& ids) const{
99
100 for (SG::ReadDecorHandle<xAOD::TruthParticleContainer, std::vector<unsigned long long>>& hitCollection : m_truthHitsKeyArray.makeHandles(ctx)){
101 for (const unsigned long long& hit_compID : hitCollection(truthParticle)){
102 const Identifier id{hit_compID};
103 if (m_idHelperSvc->isTgc(id)) { // TGCS should be added to both EIL and EIS
104 if (m_idHelperSvc->phiIndex(id) == PhiIndex::T4) {
105 ids[ChIndex::EIS].push_back(id);
106 ids[ChIndex::EIL].push_back(id);
107 } else {
108 ids[ChIndex::EMS].push_back(id);
109 ids[ChIndex::EML].push_back(id);
110 }
111 } else {
112 ids[m_idHelperSvc->chamberIndex(id) ].push_back(id);
113 }
114 }
115 }
116 return StatusCode::SUCCESS;
117 }
118
119 StatusCode MuonTruthSegmentCreationAlg::createSegments(const EventContext& ctx,
121 const ChamberIdMap& ids,
122 xAOD::MuonSegmentContainer& segmentContainer) const {
123
124 const MuonGM::MuonDetectorManager* detMgr{nullptr};
125 ATH_CHECK(SG::get(detMgr, m_detMgrKey, ctx));
126
127 constexpr unsigned techMax = toInt(TechnologyIndex::TechnologyIndexMax);
128 std::array<const MuonSimDataCollection*, techMax> sdoCollections{};
129 bool useSDO = !m_CSC_SDO_TruthNames.empty();
131 const MuonSimDataCollection* coll{nullptr};
132 ATH_CHECK(SG::get(coll, k, ctx));
133 if (coll->empty()) {
134 continue;
135 }
136 Identifier id = coll->begin()->first;
137 sdoCollections[toInt(m_idHelperSvc->technologyIndex(id))] = coll;
138 useSDO = true;
139 }
140 // loop over chamber layers
141 for (const auto& [chIdx, assocIds] : ids) {
142 // skip empty layers
143 Amg::Vector3D firstPos{Amg::Vector3D::Zero()}, secondPos{Amg::Vector3D::Zero()};
144 bool firstPosSet{false}, secondPosSet{false};
145 Identifier chId{};
146 int index = -1;
147 uint8_t nprecLayers{0}, nphiLayers{0}, ntrigEtaLayers{0};
148 std::unordered_set<int> phiLayers{}, etaLayers{}, precLayers{};
149 ATH_MSG_DEBUG(" new chamber layer " << Muon::MuonStationIndex::chName(chIdx) << " hits " << assocIds.size());
150 // loop over hits
151 for (const auto& id : assocIds) {
152 ATH_MSG_VERBOSE(" hit " << m_idHelperSvc->toString(id));
153 bool measPhi = m_idHelperSvc->measuresPhi(id);
154 bool isCsc = m_idHelperSvc->isCsc(id);
155 bool isMM = m_idHelperSvc->isMM(id);
156 bool isTrig = m_idHelperSvc->isTrigger(id);
157 bool isEndcap = m_idHelperSvc->isEndcap(id);
158 if (measPhi) {
159 phiLayers.insert(m_idHelperSvc->gasGap(id));
160 } else {
161 if (!isTrig) {
162 if (!chId.is_valid()) chId = id; // use first precision hit in list
163 if (isCsc || isMM) {
164 precLayers.insert(m_idHelperSvc->gasGap(id));
165 } else {
166 int iid = 10 * m_idHelperSvc->mdtIdHelper().multilayer(id) + m_idHelperSvc->mdtIdHelper().tubeLayer(id);
167 precLayers.insert(iid);
168 // ATH_MSG_VERBOSE("iid " << iid << " precLayers size " << precLayers.size() );
169 }
170 } else {
171 etaLayers.insert(m_idHelperSvc->gasGap(id));
172 }
173 }
174 // use SDO to look-up truth position of the hit
175 if (!useSDO) {
176 continue;
177 }
178 Amg::Vector3D gpos{Amg::Vector3D::Zero()};
179 if (!isCsc) {
180 bool ok = false;
181 TechnologyIndex techIdx = m_idHelperSvc->technologyIndex(id);
182 if (sdoCollections[toInt(techIdx)]) {
183 auto pos = sdoCollections[toInt(techIdx)]->find(id);
184 if (pos != sdoCollections[toInt(techIdx)]->end()) {
185 gpos = pos->second.globalPosition();
186 if (gpos.perp() > 0.1) ok = true; // sanity check
187 }
188 }
189 // look up successful, calculate
190 if (!ok) continue;
191
192 // small comparison function
193 auto isSmaller = [isEndcap](const Amg::Vector3D& p1, const Amg::Vector3D& p2) {
194 if (isEndcap)
195 return std::abs(p1.z()) < std::abs(p2.z());
196 else
197 return p1.perp() < p2.perp();
198 };
199 if (!firstPosSet) {
200 firstPos = gpos;
201 firstPosSet = true;
202 } else if (!secondPosSet) {
203 secondPos = gpos;
204 secondPosSet = true;
205 if (isSmaller(secondPos, firstPos)) std::swap(firstPos, secondPos);
206 } else {
207 // update position if we can increase the distance between the two positions
208 if (isSmaller(gpos, firstPos))
209 firstPos = gpos;
210 else if (isSmaller(secondPos, gpos))
211 secondPos = gpos;
212 }
213 } else {
214 SG::ReadHandle cscCollection(m_CSC_SDO_TruthNames, ctx);
215 ATH_CHECK(cscCollection.isPresent());
216 auto pos = cscCollection->find(id);
217 if (pos == cscCollection->end()) {
218 continue;
219 }
220 const MuonGM::CscReadoutElement* descriptor = detMgr->getCscReadoutElement(id);
221 ATH_MSG_DEBUG("found csc sdo with " << pos->second.getdeposits().size() << " deposits");
222 Amg::Vector3D locpos(0, pos->second.getdeposits()[0].second.ypos(), pos->second.getdeposits()[0].second.zpos());
223 gpos = descriptor->localToGlobalCoords(locpos, m_idHelperSvc->cscIdHelper().elementID(id));
224 ATH_MSG_DEBUG("got CSC global position " << gpos);
225 if (!firstPosSet) {
226 firstPos = gpos;
227 firstPosSet = true;
228 } else if (!secondPosSet) {
229 secondPos = gpos;
230 secondPosSet = true;
231 if (secondPos.perp() < firstPos.perp()) std::swap(firstPos, secondPos);
232 } else {
233 if (gpos.perp() < firstPos.perp())
234 firstPos = gpos;
235 else if (secondPos.perp() < gpos.perp())
236 secondPos = gpos;
237 }
238 }
239 }
240 if (precLayers.size() > 2) {
241 if (!phiLayers.empty()) nphiLayers = phiLayers.size();
242 ntrigEtaLayers = etaLayers.size();
243 nprecLayers = precLayers.size();
244 ATH_MSG_DEBUG(" total counts: precision " << static_cast<int>(nprecLayers) << " phi layers " << static_cast<int>(nphiLayers)
245 << " eta trig layers " << static_cast<int>(ntrigEtaLayers)
246 << " associated reco muon " << index << " unique ID " << HepMC::uniqueID(*truthLink)
247 << " truthLink " << truthLink);
248 xAOD::MuonSegment* segment = segmentContainer.push_back(std::make_unique<xAOD::MuonSegment>());
249
250 segment->setNHits(nprecLayers, nphiLayers, ntrigEtaLayers);
252 truthParticleLinkAcc("truthParticleLink");
253 truthParticleLinkAcc(*segment) = truthLink;
254 if (chId.is_valid()) { //we should always enter here if we have precision measurements
255 int eta = m_idHelperSvc->stationEta(chId);
256 int sector = m_idHelperSvc->sector(chId);
257 MuonStationIndex::TechnologyIndex technology = m_idHelperSvc->technologyIndex(chId);
258 MuonStationIndex::ChIndex chIndex = m_idHelperSvc->chamberIndex(chId);
259 segment->setIdentifier(sector, chIndex, eta, technology);
260 }
261 if (firstPosSet && secondPosSet) {
262 Amg::Vector3D gpos = (firstPos + secondPos) / 2.;
263 Amg::Vector3D gdir = (firstPos - secondPos).unit();
264 ATH_MSG_DEBUG(" got position : r " << gpos.perp() << " z " << gpos.z() << " and direction: theta " << gdir.theta()
265 << " phi " << gdir.phi());
266 segment->setPosition(gpos.x(), gpos.y(), gpos.z());
267 segment->setDirection(gdir.x(), gdir.y(), gdir.z());
268 }
269 }
270 }
271 return StatusCode::SUCCESS;
272 }
273
274
275
276} // namespace Muon
Scalar eta() const
pseudorapidity method
const PlainObject unit() const
This is a plugin that makes Eigen look like CLHEP & defines some convenience methods.
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_DEBUG(x)
Handle class for reading a decoration on an object.
value_type push_back(value_type pElem)
Add an element to the end of the collection.
bool is_valid() const
Check if id is in a valid state.
Amg::Vector3D localToGlobalCoords(const Amg::Vector3D &x, const Identifier &id) const
localToGlobalCoords and Transf connect the Gas Gap Frame (defined as a Sensitive Detector) to the Glo...
The MuonDetectorManager stores the transient representation of the Muon Spectrometer geometry and pro...
const CscReadoutElement * getCscReadoutElement(const Identifier &id) const
access via extended identifier (requires unpacking)
SG::ReadDecorHandleKey< xAOD::TruthParticleContainer > m_truthOriginKey
virtual StatusCode execute(const EventContext &ctx) const override
StatusCode fillChamberIdMap(const EventContext &ctx, const xAOD::TruthParticle &truthParticle, ChamberIdMap &ids) const
This function uses the 6 vectors, contained in.
ServiceHandle< Muon::IMuonIdHelperSvc > m_idHelperSvc
Handle for the muonIdHelper service.
SG::ReadDecorHandleKeyArray< xAOD::TruthParticleContainer, std::vector< unsigned long long > > m_truthHitsKeyArray
Keys of the truth muon decorations that we need to read to (re-)fill the chamberIdMap.
std::map< Muon::MuonStationIndex::ChIndex, std::vector< Identifier > > ChamberIdMap
This map contains all the hits corresponding to truth muons classified by chamber layer that recorded...
SG::ReadCondHandleKey< MuonGM::MuonDetectorManager > m_detMgrKey
MuonDetectorManager from the conditions store.
SG::WriteHandleKey< xAOD::MuonSegmentContainer > m_muonTruthSegmentContainerName
Key for segment container that will be populated with segments.
SG::ReadHandleKey< CscSimDataCollection > m_CSC_SDO_TruthNames
StatusCode createSegments(const EventContext &ctx, const ElementLink< xAOD::TruthParticleContainer > &truthLink, const ChamberIdMap &ids, xAOD::MuonSegmentContainer &segmentContainer) const
This function performs, for each truth muon, the actual segment creation and stores segments into a n...
SG::ReadHandleKey< xAOD::TruthParticleContainer > m_muonTruth
Key for the truth muon container and muon origin decoration.
SG::ReadHandleKeyArray< MuonSimDataCollection > m_SDO_TruthNames
Keys for all çontainers of muon hit simulation data, classified by detector technology.
Helper class to provide type-safe access to aux data.
Handle class for reading a decoration on an object.
Property holding a SG store/key/clid from which a ReadHandle is made.
bool isPresent() const
Is the referenced object present in SG?
const std::string & name() const
Return the StoreGate ID for the referenced object.
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
void setDirection(float px, float py, float pz)
Sets the direction.
void setNHits(int nPrecisionHits, int nPhiLayers, int nTrigEtaLayers)
Set the number of hits/layers.
void setIdentifier(int sector, ::Muon::MuonStationIndex::ChIndex chamberIndex, int etaIndex, ::Muon::MuonStationIndex::TechnologyIndex technology)
Set the identifier.
void setPosition(float x, float y, float z)
Sets the global position.
Eigen::Matrix< double, 3, 1 > Vector3D
int uniqueID(const T &p)
ChIndex chIndex(const std::string &index)
convert ChIndex name string to enum
TechnologyIndex
enum to classify the different layers in the muon spectrometer
constexpr int toInt(const EnumType enumVal)
const std::string & chName(ChIndex index)
convert ChIndex into a string
ChIndex
enum to classify the different chamber layers in the muon spectrometer
NRpcCablingAlg reads raw condition data and writes derived condition data to the condition store.
const T * get(const ReadCondHandleKey< T > &key, const EventContext &ctx)
Convenience function to retrieve an object given a ReadCondHandleKey.
Definition index.py:1
void swap(ElementLinkVector< DOBJ > &lhs, ElementLinkVector< DOBJ > &rhs)
MuonSegmentContainer_v1 MuonSegmentContainer
Definition of the current "MuonSegment container version".
TruthParticle_v1 TruthParticle
Typedef to implementation.
MuonSegment_v1 MuonSegment
Reference the current persistent version:
TruthParticleContainer_v1 TruthParticleContainer
Declare the latest version of the truth particle container.