ATLAS Offline Software
Loading...
Searching...
No Matches
TruthSegmentMaker.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4#include "TruthSegmentMaker.h"
5
6
8
15
18
19#include "GaudiKernel/PhysicalConstants.h"
20
21#include <unordered_map>
22
23namespace{
24 constexpr double c_inv = 1./Gaudi::Units::c_light;
25}
26
27namespace MuonR4{
28 using namespace SegmentFit;
29
31 ATH_CHECK(m_idHelperSvc.retrieve());
32 ATH_CHECK(m_readKeys.initialize());
33 if (m_readKeys.empty()){
34 ATH_MSG_ERROR("No simulated hit containers have been parsed to build the segments from ");
35 return StatusCode::FAILURE;
36 }
37 for (const auto& truthLink : m_readKeys) {
38 m_segLinkKeys.emplace_back(truthLink, m_segLinkKey);
39 }
40 ATH_CHECK(m_segLinkKeys.initialize());
41 ATH_CHECK(m_mdtCalibKey.initialize(m_idHelperSvc->hasMDT()));
42 ATH_CHECK(m_nswUncertKey.initialize(m_idHelperSvc->hasMM() || m_idHelperSvc->hasSTGC()));
43
44 ATH_CHECK(m_segmentKey.initialize());
45 ATH_CHECK(m_eleLinkKey.initialize());
46 ATH_CHECK(m_ptKey.initialize());
47
48 ATH_CHECK(m_locParKey.initialize());
49 ATH_CHECK(m_qKey.initialize());
50
51 ATH_CHECK(m_geoCtxKey.initialize());
52 ATH_CHECK(detStore()->retrieve(m_detMgr));
53 return StatusCode::SUCCESS;
54 }
55 float TruthSegmentMaker::hitUncertainty(const EventContext& ctx, const xAOD::MuonSimHit& hit) const {
56 const Identifier hitId{hit.identify()};
57 switch (const auto techIdx = m_idHelperSvc->technologyIndex(hitId)) {
59 case MDT: {
60 const MuonCalib::MdtCalibDataContainer* calibCont{nullptr};
61 if (!SG::get(calibCont, m_mdtCalibKey, ctx).isSuccess()) {
62 THROW_EXCEPTION("Failed to retrieve Mdt calib constants");
63 }
64 const auto& rtCalib{calibCont->getCalibData(hitId, msgStream())->rtRelation};
65 const double driftTime = rtCalib->tr()->driftTime(hit.localPosition().perp()).value_or(rtCalib->tr()->maxRadius());
66 return rtCalib->rtRes()->resolution(driftTime);
67 } case RPC: {
68 const auto* re = m_detMgr->getRpcReadoutElement(hitId);
69 return re->stripEtaPitch() / std::sqrt(12.);
70 } case TGC: {
71 const auto* re = m_detMgr->getTgcReadoutElement(hitId);
72 const IdentifierHash measHash = re->measurementHash(hitId);
73 const auto& design = re->wireGangLayout(measHash);
74 return design.stripPitch() / std::sqrt(12.);
75 } case STGC:
76 case MM: {
77 const NswErrorCalibData* errorCalibDB{nullptr};
78 if (!SG::get(errorCalibDB, m_nswUncertKey, ctx).isSuccess()) {
79 THROW_EXCEPTION("Failed to retrieve the STGC calibration constants");
80 }
81 NswErrorCalibData::Input errorCalibInput{};
82 errorCalibInput.stripId= hitId;
83 errorCalibInput.locTheta = M_PI - hit.localDirection().theta();
84 if (techIdx == STGC) {
85 errorCalibInput.clusterAuthor = 3; // centroid
86 } else {
87 errorCalibInput.clusterAuthor=66; // cluster time projection method
88 }
89 return errorCalibDB->clusterUncertainty(errorCalibInput);
90 }
91 default:
92 break;
93 }
94 return 0.;
95 }
96 float TruthSegmentMaker::muonPt(const xAOD::MuonSimHit& hit, const Amg::Vector3D& globDir) const{
97 const auto& link = hit.genParticleLink();
98 if (link.isValid()) {
99 return link->momentum().perp();
100 }
101 const float e{hit.kineticEnergy()},m{hit.mass()};
102 return std::sqrt(std::max(e*e - m*m, 0.f)) * std::sin(globDir.theta());
103 }
105 const Identifier& chanId) const {
106 const MuonGMR4::MuonReadoutElement* reEle = m_detMgr->getReadoutElement(chanId);
107 const IdentifierHash trfHash{reEle->detectorType() == ActsTrk::DetectorType::Mdt ?
108 reEle->measurementHash(chanId) :
109 reEle->layerHash(chanId)};
110 return reEle->msSector()->globalToLocalTransform(gctx) * reEle->localToGlobalTransform(gctx, trfHash);
111
112 }
113
114 void TruthSegmentMaker::buildSegmentsFromBkg(const EventContext& ctx,
115 const Amg::Transform3D& locToGlob,
116 const SimHitVec_t& simHits,
117 WriteDecorHolder& out)const {
118
119 ATH_MSG_VERBOSE("Assemble segments from "<<simHits.size()<<" background hits.");
120 std::vector<char> alreadyUsed(simHits.size(), 0);
121 for (std::size_t h = 0 ;h < simHits.size(); ++h) {
122 if (alreadyUsed[h]) {
123 continue;
124 }
125 const auto& [refHit, refPos, refDir] = simHits[h];
127 ATH_MSG_VERBOSE("Try to find other hits on trajectory given by "<<m_idHelperSvc->toString(refHit->identify())
128 <<", "<<Amg::toString(refPos)<<", "<<Amg::toString(refDir)<<", E: "<<refHit->kineticEnergy() / Gaudi::Units::GeV<<" [GeV] "
129 <<", pt: "<<muonPt(*refHit, locToGlob.linear()* refDir)/ Gaudi::Units::GeV <<" [GeV].");
130 std::vector<std::size_t> indicesOnSeg{h};
131 for (std::size_t h1 = h +1; h1 < simHits.size(); ++h1) {
132 if (alreadyUsed[h1]) {
133 continue;
134 }
135 const auto&[testHit, testPos, testDir] = simHits[h1];
136
138 if (std::abs(refHit->kineticEnergy() - testHit->kineticEnergy()) > m_pileUpHitELoss) {
139 continue;
140 }
141 if (MC::charge(refHit) != MC::charge(testHit)){
142 continue;
143 }
144 const double angleDir = Amg::angle(testDir, refDir);
145 const double hitSep = std::abs(Amg::signedDistance(refPos, refDir, testPos,testDir));
146 ATH_MSG_VERBOSE("Test hit "<<m_idHelperSvc->toString(testHit->identify())<<", "
147 <<Amg::toString(testPos)<<", testDir: "<<Amg::toString(testDir)
148 <<", E: "<<testHit->kineticEnergy() / (Gaudi::Units::GeV)<<", angle: "
149 <<(angleDir / Gaudi::Units::deg) <<", distance: "<<hitSep<<".");
150
151 if (std::abs(angleDir) > m_pileUpHitAngleCone || hitSep > m_pileUpHitDistance) {
152 continue;
153 }
154 indicesOnSeg.push_back(h1);
155 }
157 const Amg::Vector3D globPos = locToGlob * refPos;
158 const Amg::Vector3D globDir = locToGlob.linear() * refDir;
159 const Amg::Vector3D perigee = globPos + Amg::intersect<3>(Amg::Vector3D::Zero(), Amg::Vector3D::UnitZ(),
160 globPos, globDir).value_or(0.) * globDir;
161 ATH_MSG_VERBOSE("Found "<<indicesOnSeg.size()<<" matching hits. Closest perigee "
162 <<Amg::toString(perigee)<<", "<<m_idCylinderR<<", "<<m_idCylinderHalfZ);
163 SimHitVec_t hitsOnSeg{};
164 std::ranges::transform(indicesOnSeg,
165 std::back_inserter(hitsOnSeg),
166 [&simHits](const auto idx) {
167 return simHits[idx];
168 });
169 if (perigee.perp() > m_idCylinderR || std::abs(perigee.z()) > m_idCylinderHalfZ ||
170 constructSegmentFromHits(ctx, locToGlob, hitsOnSeg, out)){
171 std::ranges::for_each(indicesOnSeg,
172 [&alreadyUsed](const auto idx) {
173 alreadyUsed[idx] = 1;
174 });
175 }
176 }
177 }
180 const Amg::Transform3D& locToGlob,
181 const SimHitVec_t& simHits,
182 WriteDecorHolder& out) const {
185 auto refHit_itr = std::ranges::min_element(simHits,
186 [](const HitPosTuple_t& a, const HitPosTuple_t& b){
187 return std::abs(std::get<1>(a).z()) <std::abs(std::get<1>(b).z());
188 });
189
190 const auto [simHit, localPos, chamberDir] = (*refHit_itr);
191 ATH_MSG_VERBOSE("Create segement from hit: "<<m_idHelperSvc->toString(simHit->identify())<<
192 " pdgId: "<<simHit->pdgId()<<", energy: "<<simHit->kineticEnergy()
193 <<", genParticle: "<<simHit->genParticleLink().cptr());
194 const Identifier segId{simHit->identify()};
195
197 const double distance = Amg::intersect<3>(localPos, chamberDir, Amg::Vector3D::UnitZ(), 0.).value_or(0.);
198 const Amg::Vector3D chamberPos = localPos + distance*chamberDir;
199
200 const Amg::Vector3D globPos = locToGlob * chamberPos;
201 const Amg::Vector3D globDir = locToGlob.linear() * chamberDir;
202 HitLinkVec_t associatedHits{};
203 unsigned nMdt{0}, nRpcEta{0}, nRpcPhi{0}, nTgcEta{0}, nTgcPhi{0};
204 unsigned nMm{0}, nStgcEta{0}, nStgcPhi{0};
206 float chi2{0.f};
207 for (const auto& [assocMe, pos, dir] : simHits) {
208 chi2 += std::pow(Amg::signedDistance(pos,dir, chamberPos, chamberDir) / hitUncertainty(ctx,*assocMe), 2);
209 const MuonGMR4::MuonReadoutElement* assocRE = m_detMgr->getReadoutElement(assocMe->identify());
210 switch (assocRE->detectorType()) {
212 ++nMdt;
213 break;
215 auto castRE{static_cast<const MuonGMR4::RpcReadoutElement*>(assocRE)};
216 if (castRE->nEtaStrips()) ++nRpcEta;
217 if (castRE->nPhiStrips()) ++nRpcPhi;
218 break;
220 auto castRE{static_cast<const MuonGMR4::TgcReadoutElement*>(assocRE)};
221 const IdentifierHash gapHash = assocRE->measurementHash(assocMe->identify());
222 if (castRE->numStrips(gapHash)){
223 ++nTgcPhi;
224 }
225 if (castRE->numWireGangs(gapHash)) {
226 ++nTgcEta;
227 }
228 break;
230 ++nStgcEta;
231 ++nStgcPhi;
232 break;
234 ++nMm;
235 break;
236 } default:
237 ATH_MSG_WARNING("Csc are not defined "<<m_idHelperSvc->toString(simHit->identify()));
238 }
239 ATH_MSG_VERBOSE("Associate hit "<<m_idHelperSvc->toString(assocMe->identify())
240 <<" pdgId: "<<assocMe->pdgId()<<", energy: "<<assocMe->kineticEnergy()
241 <<", genParticle: "<<assocMe->genParticleLink().cptr()<<", beta: "<<simHit->beta()
242 <<" global time: "<<simHit->globalTime()<<", pos: "<<Amg::toString(pos)
243 <<", dir: "<<Amg::toString(dir));
244 EleLink_t link{*static_cast<const xAOD::MuonSimHitContainer*>(assocMe->container()), assocMe->index()};
245 associatedHits.push_back(std::move(link));
246 }
247 int nPrecisionHits = nMdt + nMm + nStgcEta;
248 int nPhiLayers = nTgcPhi + nRpcPhi + nStgcPhi;
249 // if nMdt + nMm + nStgcEta < 3, do not create a segment
250 if (nPrecisionHits < 3) {
251 return nullptr;
252 }
253
254 xAOD::MuonSegment* truthSegment = out.segments.push_back(std::make_unique<xAOD::MuonSegment>());
255 out.ptDecor(*truthSegment) = muonPt(*simHit, globDir);
256 out.chargeDecor(*truthSegment) = MC::charge(simHit);
257 SegPars_t& locPars{out.paramDecor(*truthSegment)};
258 locPars[Acts::toUnderlying(ParamDefs::x0)] = chamberPos.x();
259 locPars[Acts::toUnderlying(ParamDefs::y0)] = chamberPos.y();
260 constexpr float betaLowLimit = 1.e-6;
261 locPars[Acts::toUnderlying(ParamDefs::t0)] = simHit->globalTime() + distance *c_inv / std::max(simHit->beta(), betaLowLimit);
262 locPars[Acts::toUnderlying(ParamDefs::theta)] = chamberDir.theta();
263 locPars[Acts::toUnderlying(ParamDefs::phi)] = chamberDir.phi();
264
265 truthSegment->setPosition(globPos.x(), globPos.y(), globPos.z());
266 truthSegment->setDirection(globDir.x(), globDir.y(), globDir.z());
267 truthSegment->setT0Error(locPars[Acts::toUnderlying(ParamDefs::t0)], 0.);
268
269 truthSegment->setNHits(nPrecisionHits, nPhiLayers, nTgcEta + nRpcEta);
270 truthSegment->setIdentifier(m_idHelperSvc->sector(segId),
271 m_idHelperSvc->chamberIndex(segId),
272 m_idHelperSvc->stationEta(segId),
273 m_idHelperSvc->technologyIndex(segId));
274 // adding chi2 and ndof (nHits - 5 for 2 position, 2 direction and 1 time)
275 if (nPhiLayers == 0){
276 truthSegment->setFitQuality(chi2, (nPrecisionHits + nTgcEta + nRpcEta - 3));
277 } else {
278 truthSegment->setFitQuality(chi2, (nPrecisionHits + nPhiLayers + nTgcEta + nRpcEta - 5));
279 }
280 out.hitLinkDecor(*truthSegment) = std::move(associatedHits);
281 return truthSegment;
282 }
283
284 StatusCode TruthSegmentMaker::execute(const EventContext& ctx) const {
285 const ActsTrk::GeometryContext* gctx{nullptr};
286 ATH_CHECK(SG::get(gctx, m_geoCtxKey, ctx));
287
288 using HitsPerParticle = std::unordered_map<HepMC::ConstGenParticlePtr, SimHitVec_t>;
289 using HitCollector = std::unordered_map<const MuonGMR4::SpectrometerSector*, HitsPerParticle>;
290 HitCollector hitCollector{};
291
293 const xAOD::MuonSimHitContainer* simHits{nullptr};
294 ATH_CHECK(SG::get(simHits, key, ctx));
295 for (const xAOD::MuonSimHit* simHit : *simHits) {
296 const MuonGMR4::MuonReadoutElement* reElement = m_detMgr->getReadoutElement(simHit->identify());
297 const MuonGMR4::SpectrometerSector* id{reElement->msSector()};
298 auto genLink = simHit->genParticleLink();
299 HepMC::ConstGenParticlePtr genParticle = nullptr;
300 if (genLink.isValid()){
301 genParticle = genLink.cptr();
302 }
303 if ( (!genParticle && (!m_includePileUpHits || simHit->kineticEnergy() < m_pileUpHitMinE
304 || !MC::isMuon(simHit))) || (m_useOnlyMuonHits && !MC::isMuon(simHit))) {
305 ATH_MSG_VERBOSE("Skip hit "<<m_idHelperSvc->toString(simHit->identify())<<
306 " pdgId: "<<simHit->pdgId()<<", energy: "<<simHit->kineticEnergy()
307 <<", genParticle: "<<genParticle);
308 continue;
309 }
310 const Amg::Transform3D toChTrf{toChamber(*gctx, simHit->identify())};
311 hitCollector[id][genParticle].emplace_back(simHit,
312 toChTrf *xAOD::toEigen(simHit->localPosition()),
313 toChTrf.linear()* xAOD::toEigen(simHit->localDirection()));
314 }
315 }
316
317 SG::WriteHandle writeHandle{m_segmentKey, ctx};
318 ATH_CHECK(writeHandle.record(std::make_unique<xAOD::MuonSegmentContainer>(),
319 std::make_unique<xAOD::MuonSegmentAuxContainer>()));
320
321 WriteDecorHolder writerHolder{*writeHandle,*this, ctx};
322
323 for (auto& [chamber, collectedParts] : hitCollector) {
324 const Amg::Transform3D& locToGlob{chamber->localToGlobalTransform(*gctx)};
325 for (auto& [particle, simHits]: collectedParts) {
326 /* Sort hits by local z */
327 std::ranges::stable_sort(simHits,[](const HitPosTuple_t& a, const HitPosTuple_t& b){
328 return std::get<1>(a).z() < std::get<1>(b).z();
329 });
330 if (!particle) {
331 buildSegmentsFromBkg(ctx, locToGlob, simHits, writerHolder);
332 continue;
333 }
334 constructSegmentFromHits(ctx, locToGlob, simHits, writerHolder);
335 }
336 }
337 ATH_CHECK(linkSegmentsToHits(ctx, *writeHandle));
338 ATH_MSG_DEBUG("Constructed "<<writeHandle->size()<<" truth segments in total ");
339 return StatusCode::SUCCESS;
340 }
341 StatusCode TruthSegmentMaker::linkSegmentsToHits(const EventContext& ctx,
342 const xAOD::MuonSegmentContainer& segments) const {
345 std::unordered_map<const SG::AuxVectorData*, DecorHandle_t> handleMap{};
346 for (const auto& decorKey : m_segLinkKeys) {
347 DecorHandle_t decorHandle{decorKey, ctx};
348 if (decorHandle->empty()) {
349 ATH_MSG_DEBUG("Don't setup a decoration handle for "<<decorKey.fullKey());
350 continue;
351 }
352 decorHandle(*decorHandle->front()) = SegLink_t{};
353 handleMap.insert(std::make_pair(decorHandle.cptr(), std::move(decorHandle)));
354 }
355 for (const xAOD::MuonSegment* segment: segments) {
356 SegLink_t segLink{&segments, segment->index()};
357 for (const xAOD::MuonSimHit* simHit : getMatchingSimHits(*segment)) {
358 handleMap.at(simHit->container())(*simHit) = segLink;
359 }
360 }
361 return StatusCode::SUCCESS;
362 }
363}
const boost::regex re(r_e)
#define M_PI
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
ATLAS-specific HepMC functions.
static Double_t a
@ STGC
Definition RegSelEnums.h:39
@ MM
Definition RegSelEnums.h:38
@ RPC
Definition RegSelEnums.h:32
@ MDT
Definition RegSelEnums.h:31
#define z
virtual DetectorType detectorType() const =0
Returns the detector element type.
const ServiceHandle< StoreGateSvc > & detStore() const
Header file for AthHistogramAlgorithm.
This is a "hash" representation of an Identifier.
const MdtFullCalibData * getCalibData(const Identifier &measId, MsgStream &msg) const
Returns the calibration data associated with this station.
MuonReadoutElement is an abstract class representing the geometry of a muon detector.
const Amg::Transform3D & localToGlobalTransform(const ActsTrk::GeometryContext &ctx) const
Returns the transformation from the local coordinate system of the readout element into the global AT...
const SpectrometerSector * msSector() const
Returns the pointer to the envelope volume enclosing all chambers in the sector.
virtual IdentifierHash measurementHash(const Identifier &measId) const =0
The measurement hash is a continous numbering schema of all readout channels described by the specifi...
virtual IdentifierHash layerHash(const Identifier &measId) const =0
The layer hash removes the bits from the IdentifierHash corresponding to the measurement's channel nu...
A spectrometer sector forms the envelope of all chambers that are placed in the same MS sector & laye...
Amg::Transform3D globalToLocalTransform(const ActsTrk::GeometryContext &gctx) const
Returns the global -> local transformation from the ATLAS global.
SG::WriteDecorHandleKey< xAOD::MuonSegmentContainer > m_locParKey
Decoration key of the local parameters.
std::vector< HitPosTuple_t > SimHitVec_t
Gaudi::Property< float > m_pileUpHitMinE
Minimum energy threshold for pile up hits to be converted.
SG::ReadHandleKey< ActsTrk::GeometryContext > m_geoCtxKey
Key to the geometry context.
StatusCode execute(const EventContext &ctx) const override
ServiceHandle< Muon::IMuonIdHelperSvc > m_idHelperSvc
IdHelperSvc to decode the Identifiers.
Gaudi::Property< float > m_pileUpHitAngleCone
Maximum scattering angle between two pile-up hits.
Gaudi::Property< float > m_idCylinderR
ID / ITk cylinder radius.
Amg::Transform3D toChamber(const ActsTrk::GeometryContext &gctx, const Identifier &chanId) const
Returns the transform from the local simHit frame -> chamber frame.
const MuonGMR4::MuonDetectorManager * m_detMgr
Pointer to the muon readout geometry.
StatusCode initialize() override final
SG::WriteDecorHandleKey< xAOD::MuonSegmentContainer > m_eleLinkKey
Decoration key of the associated sim hit links.
SG::ReadCondHandleKey< MuonCalib::MdtCalibDataContainer > m_mdtCalibKey
Data dependency on the Mdt calibration container to calculate the uncertainty.
ElementLink< xAOD::MuonSimHitContainer > EleLink_t
Gaudi::Property< float > m_pileUpHitELoss
Maximum energy loss between two pile-up hits.
std::vector< EleLink_t > HitLinkVec_t
Gaudi::Property< float > m_idCylinderHalfZ
ID / Itk cylinder half length.
SG::WriteDecorHandleKey< xAOD::MuonSegmentContainer > m_qKey
Decoration key of the muon charge.
xAOD::MuonSegment * constructSegmentFromHits(const EventContext &ctx, const Amg::Transform3D &locToGlob, const SimHitVec_t &hits, WriteDecorHolder &writeShip) const
Takes a list of related sim hits and transforms them into a truth segment.
SG::WriteDecorHandleKey< xAOD::MuonSegmentContainer > m_ptKey
Decoration key of the associated particle pt.
Gaudi::Property< bool > m_useOnlyMuonHits
Build segments from muon hits only.
SG::WriteDecorHandleKeyArray< xAOD::MuonSimHitContainer > m_segLinkKeys
Decorate the truth segment link to the simHit.
Gaudi::Property< std::string > m_segLinkKey
Name of the link to the truth segment.
void buildSegmentsFromBkg(const EventContext &ctx, const Amg::Transform3D &locToGlob, const SimHitVec_t &simHits, WriteDecorHolder &writeShip) const
Attempts to assemble truth segments from a list of loose sim hits, i.e., the hits are stemming from a...
std::tuple< const xAOD::MuonSimHit *, Amg::Vector3D, Amg::Vector3D > HitPosTuple_t
Tuple consisting out of pointer to the sim hit and the position & direction expressed in the chamber'...
xAOD::MeasVector< Acts::toUnderlying(SegmentFit::ParamDefs::nPars)> SegPars_t
SG::WriteHandleKey< xAOD::MuonSegmentContainer > m_segmentKey
Key under which the segment Container will be recorded in StoreGate.
SG::ReadCondHandleKey< NswErrorCalibData > m_nswUncertKey
Data dependency on the Nsw calibration container to estimate the uncertaintys.
Gaudi::Property< bool > m_includePileUpHits
Construct segments from pile-up hits without GenParticleLink.
Gaudi::Property< float > m_pileUpHitDistance
Maximum separation between two pile-up hits.
float hitUncertainty(const EventContext &ctx, const xAOD::MuonSimHit &hit) const
Evaluates the hit uncertainty.
float muonPt(const xAOD::MuonSimHit &hit, const Amg::Vector3D &globDir) const
Returns the muon pt from the sim hit.
StatusCode linkSegmentsToHits(const EventContext &ctx, const xAOD::MuonSegmentContainer &segments) const
Establish the link from the simulated hit -> truth segment.
SG::ReadHandleKeyArray< xAOD::MuonSimHitContainer > m_readKeys
List of sim hit containers from which the truth segments shall be retrieved.
double clusterUncertainty(const Input &clustInfo) const
Property holding a SG store/key/clid from which a ReadHandle is made.
Handle class for adding a decoration to an 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 setFitQuality(float chiSquared, float numberDoF)
Set the 'Fit Quality' information.
void setNHits(int nPrecisionHits, int nPhiLayers, int nTrigEtaLayers)
Set the number of hits/layers.
void setT0Error(float t0, float t0Error)
Sets the time error.
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.
ConstVectorMap< 3 > localDirection() const
Returns the local direction of the traversing particle.
float mass() const
Returns the rest-mass of the traversing particle.
Identifier identify() const
Returns the global ATLAS identifier of the SimHit.
ConstVectorMap< 3 > localPosition() const
Returns the local postion of the traversing particle.
const HepMcParticleLink & genParticleLink() const
Returns the link to the HepMC particle producing this hit.
float kineticEnergy() const
Returns the kinetic energy of the traversing particle.
double chi2(TH1 *h0, TH1 *h1)
@ Mm
Maybe not needed in the migration.
@ Tgc
Resitive Plate Chambers.
@ sTgc
Micromegas (NSW).
@ Rpc
Monitored Drift Tubes.
@ Mdt
MuonSpectrometer.
std::optional< double > intersect(const AmgVector(N)&posA, const AmgVector(N)&dirA, const AmgVector(N)&posB, const AmgVector(N)&dirB)
Calculates the point B' along the line B that's closest to a second line A.
std::string toString(const Translation3D &translation, int precision=4)
GeoPrimitvesToStringConverter.
double signedDistance(const Amg::Vector3D &posA, const Amg::Vector3D &dirA, const Amg::Vector3D &posB, const Amg::Vector3D &dirB)
Calculates the signed distance between two lines in 3D space.
double angle(const Amg::Vector3D &v1, const Amg::Vector3D &v2)
calculates the opening angle between two vectors
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 3, 1 > Vector3D
const GenParticle * ConstGenParticlePtr
Definition GenParticle.h:38
bool isMuon(const T &p)
double charge(const T &p)
This header ties the generic definitions in this package.
std::unordered_set< const xAOD::MuonSimHit * > getMatchingSimHits(const xAOD::MuonSegment &segment)
: Returns all sim hits matched to a xAOD::MuonSegment
ElementLink< MuonR4::SegmentContainer > SegLink_t
Abrivation of the link to the reco segment container.
TechnologyIndex
enum to classify the different layers in the muon spectrometer
const T * get(const ReadCondHandleKey< T > &key, const EventContext &ctx)
Convenience function to retrieve an object given a ReadCondHandleKey.
Definition TgcBase.h:6
MuonSimHitContainer_v1 MuonSimHitContainer
Define the version of the pixel cluster container.
MuonSegmentContainer_v1 MuonSegmentContainer
Definition of the current "MuonSegment container version".
MuonSimHit_v1 MuonSimHit
Defined the version of the MuonSimHit.
Definition MuonSimHit.h:12
MuonSegment_v1 MuonSegment
Reference the current persistent version:
Helper struct to ship the write DecorHandles and the reference to the output segment container throug...
Helper struct to be parsed to the object to derive the specific error of the cluster.
double locTheta
Direction of the muon in the local coordinate frame.
Identifier stripId
Identifier of the strip.
uint8_t clusterAuthor
Author of the cluster.
#define THROW_EXCEPTION(MESSAGE)
Definition throwExcept.h:10