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

Extract the raw hists from info in SG. More...

#include <FPGATrackSimSGToRawHitsTool.h>

Inheritance diagram for FPGATrackSimSGToRawHitsTool:

Public Member Functions

 FPGATrackSimSGToRawHitsTool (const std::string &, const std::string &, const IInterface *)
virtual ~FPGATrackSimSGToRawHitsTool ()
virtual StatusCode initialize () override
virtual StatusCode finalize () override
virtual StatusCode readData (FPGATrackSimEventInputHeader *header, const EventContext &eventContext) const override
 This function get from the SG the inner detector raw hits and prepares them for FPGATrackSim simulation.

Private Types

typedef std::map< Identifier, int > HitIndexMap

Private Member Functions

StatusCode readRawSilicon (FPGATrackSimEventInputHeader *header, HitIndexMap &hitIndexMap, const EventContext &eventContext) const
StatusCode readTruthTracks (std::vector< FPGATrackSimTruthTrack > &truth, const EventContext &eventContext) const
StatusCode readOfflineTracks (std::vector< FPGATrackSimOfflineTrack > &Track, const EventContext &eventContext) const
StatusCode readOfflineClusters (std::vector< FPGATrackSimCluster > &Clusters, const EventContext &eventContext) const
StatusCode readPixelSimulation (FPGATrackSimEventInputHeader *header, HitIndexMap &hitIndexMap, unsigned int &hitIndex, const EventContext &eventContext) const
StatusCode readStripSimulation (FPGATrackSimEventInputHeader *header, HitIndexMap &hitIndexMap, unsigned int &hitIndex, const EventContext &eventContext) const
StatusCode dumpPixelClusters (HitIndexMap &pixelClusterIndexMap, const EventContext &eventContext) const
const HepMcParticleLinkgetTruthInformation (InDetSimDataCollection::const_iterator &iter, FPGATrackSimInputUtils::ParentBitmask &parentMask) const

Private Attributes

ToolHandle< Trk::ITruthToTrackm_truthToTrack {this, "TruthToTrackTool", "Trk::TruthToTrack/InDetTruthToTrack" }
 tool to create track parameters from a gen particle
ToolHandle< Trk::IExtrapolatorm_extrapolator {this, "Extrapolator", "Trk::Extrapolator/AtlasExtrapolator"}
 ToolHandle for Extrapolator.
SG::ReadCondHandleKey< InDet::BeamSpotDatam_beamSpotKey { this, "BeamSpotKey", "BeamSpotData", "SG key for beam spot" }
SG::ReadHandleKey< xAOD::EventInfom_eventInfoKey { this, "EventInfo", "EventInfo" }
SG::ReadHandleKey< InDet::SiClusterContainer > m_pixelClusterContainerKey { this, "pixelClustersName", "ITkPixelClusters" }
SG::ReadHandleKey< InDet::SiClusterContainer > m_sctClusterContainerKey { this, "SCT_ClustersName", "ITkStripClusters" }
SG::ReadHandleKey< xAOD::TrackParticleContainerm_offlineTracksKey { this, "OfflineTracks", "InDetTrackParticles"}
SG::ReadHandleKey< McEventCollectionm_mcCollectionKey { this, "McTruth", "TruthEvent" }
SG::ReadHandleKey< InDetSimDataCollectionm_pixelSDOKey { this, "PixelSDO", "ITkPixelSDO_Map" }
SG::ReadHandleKey< InDetSimDataCollectionm_stripSDOKey { this, "StripSDO", "ITkStripSDO_Map" }
SG::ReadHandleKey< PixelRDO_Containerm_pixelRDOKey { this, "PixelRDO", "ITkPixelRDOs" }
SG::ReadHandleKey< SCT_RDO_Containerm_stripRDOKey { this, "StripRDO", "ITkStripRDOs" }
Gaudi::Property< std::string > m_tracksTruthName { this, "OfflineName", "InDetTrackParticles", "name of offline tracks collection" }
Gaudi::Property< bool > m_dumpHitsOnTracks { this, "dumpHitsOnTracks", false }
Gaudi::Property< bool > m_dumpTruthIntersections { this, "dumpTruthIntersections", false }
Gaudi::Property< bool > m_readOfflineClusters { this, "ReadOfflineClusters", true, "flag to enable the offline cluster save" }
Gaudi::Property< bool > m_readTruthTracks { this, "ReadTruthTracks", true, "flag to enable the truth tracking save" }
Gaudi::Property< bool > m_readOfflineTracks { this, "ReadOfflineTracks", true, "flag to enable the offline tracking save" }
Gaudi::Property< bool > m_UseNominalOrigin { this, "UseNominalOrigin", false, "if true truth values are always with respect to (0,0,0)" }
Gaudi::Property< double > m_maxEta { this, "maxEta", 5.0 }
Gaudi::Property< double > m_minPt { this, "minPt", .8*CLHEP::GeV }
Gaudi::Property< bool > m_doMultiTruth { this, "doMultiTruth", true }
const PixelIDm_pixelId = nullptr
const SCT_IDm_sctId = nullptr
const InDetDD::SiDetectorManagerm_PIX_mgr = nullptr
const InDetDD::SiDetectorManagerm_SCT_mgr = nullptr
const HepPDT::ParticleDataTable * m_particleDataTable = nullptr

Detailed Description

Extract the raw hists from info in SG.

Definition at line 45 of file FPGATrackSimSGToRawHitsTool.h.

Member Typedef Documentation

◆ HitIndexMap

typedef std::map<Identifier, int> FPGATrackSimSGToRawHitsTool::HitIndexMap
private

Definition at line 95 of file FPGATrackSimSGToRawHitsTool.h.

Constructor & Destructor Documentation

◆ FPGATrackSimSGToRawHitsTool()

FPGATrackSimSGToRawHitsTool::FPGATrackSimSGToRawHitsTool ( const std::string & algname,
const std::string & name,
const IInterface * ifc )

Definition at line 47 of file FPGATrackSimSGToRawHitsTool.cxx.

47 :
48 base_class(algname, name, ifc)
49{}

◆ ~FPGATrackSimSGToRawHitsTool()

virtual FPGATrackSimSGToRawHitsTool::~FPGATrackSimSGToRawHitsTool ( )
inlinevirtual

Definition at line 49 of file FPGATrackSimSGToRawHitsTool.h.

49{ ; }

Member Function Documentation

◆ dumpPixelClusters()

StatusCode FPGATrackSimSGToRawHitsTool::dumpPixelClusters ( HitIndexMap & pixelClusterIndexMap,
const EventContext & eventContext ) const
private

Definition at line 516 of file FPGATrackSimSGToRawHitsTool.cxx.

516 {
517 unsigned int pixelClusterIndex = 0;
518 auto pixelSDOHandle = SG::makeHandle(m_pixelSDOKey, eventContext);
519 auto pixelClusterContainerHandle = SG::makeHandle(m_pixelClusterContainerKey, eventContext);
520 // Dump pixel clusters. They're in m_pixelContainer
521 for (const InDet::SiClusterCollection* pixelClusterCollection : *pixelClusterContainerHandle) {
522 if (pixelClusterCollection == nullptr) {
523 ATH_MSG_DEBUG("pixelClusterCollection not available!");
524 continue;
525 }
526
527 for (const InDet::SiCluster* cluster : *pixelClusterCollection) {
528 Identifier theId = cluster->identify();
529 // if there is simulation truth available, try to retrieve the "most likely" barcode for this pixel cluster.
530 FPGATrackSimInputUtils::ParentBitmask parentMask; // FIXME set, but not used
531 if (!m_pixelSDOKey.empty()) {
532 for (const Identifier& rdoId : cluster->rdoList()) {
533 const InDetDD::SiDetectorElement* sielement = m_PIX_mgr->getDetectorElement(rdoId);
534 assert(sielement);
535 InDetDD::SiCellId cellID = sielement->cellIdFromIdentifier(rdoId);
536
537 const int nCells = sielement->numberOfConnectedCells(cellID);
538 InDetSimDataCollection::const_iterator iter(pixelSDOHandle->find(rdoId));
539 // this might be the ganged pixel copy.
540 if (nCells > 1 && iter == pixelSDOHandle->end()) {
541 InDetDD::SiReadoutCellId SiRC(m_pixelId->phi_index(rdoId), m_pixelId->eta_index(rdoId));
542 for (int ii = 0; ii < nCells && iter == pixelSDOHandle->end(); ++ii) {
543 iter = pixelSDOHandle->find(sielement->identifierFromCellId(sielement->design().connectedCell(SiRC, ii)));
544 }
545 } // end search for correct ganged pixel
546 // if SDO found for this pixel, associate the particle. otherwise leave unassociated.
547 if (iter != pixelSDOHandle->end()) { (void) getTruthInformation(iter, parentMask); } // FIXME not used??
548 } // if we have pixel sdo's available
549 }
550 pixelClusterIndexMap[theId] = pixelClusterIndex;
551 pixelClusterIndex++;
552 } // End loop over pixel clusters
553 } // End loop over pixel cluster collection
554
555 return StatusCode::SUCCESS;
556}
#define ATH_MSG_DEBUG(x)
SG::ReadHandleKey< InDet::SiClusterContainer > m_pixelClusterContainerKey
const HepMcParticleLink * getTruthInformation(InDetSimDataCollection::const_iterator &iter, FPGATrackSimInputUtils::ParentBitmask &parentMask) const
const InDetDD::SiDetectorManager * m_PIX_mgr
SG::ReadHandleKey< InDetSimDataCollection > m_pixelSDOKey
virtual SiCellId connectedCell(const SiReadoutCellId &readoutId, int number) const =0
readout id -> id of connected diodes.
virtual SiCellId cellIdFromIdentifier(const Identifier &identifier) const override final
SiCellId from Identifier.
virtual const SiDetectorDesign & design() const override final
access to the local description (inline):
virtual Identifier identifierFromCellId(const SiCellId &cellId) const override final
Identifier <-> SiCellId (ie strip number or pixel eta_index,phi_index) Identifier from SiCellId (ie s...
int numberOfConnectedCells(const SiCellId cellId) const
Test if readout cell has more than one diode associated with it.
std::bitset< NBITS > ParentBitmask
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
setRawEt setRawPhi nCells

◆ finalize()

StatusCode FPGATrackSimSGToRawHitsTool::finalize ( )
overridevirtual

Definition at line 85 of file FPGATrackSimSGToRawHitsTool.cxx.

85 {
86 return StatusCode::SUCCESS;
87}

◆ getTruthInformation()

const HepMcParticleLink * FPGATrackSimSGToRawHitsTool::getTruthInformation ( InDetSimDataCollection::const_iterator & iter,
FPGATrackSimInputUtils::ParentBitmask & parentMask ) const
private

Definition at line 850 of file FPGATrackSimSGToRawHitsTool.cxx.

851 {
852 const HepMcParticleLink* bestTruthLink{};
853 const InDetSimData& sdo(iter->second);
854 const std::vector<InDetSimData::Deposit>& deposits(sdo.getdeposits());
855 float bestPt{-999.f};
856 for (const InDetSimData::Deposit& dep : deposits) {
857
858 const HepMcParticleLink& particleLink = dep.first;
859 // RDO's without SDO's are delta rays or detector noise.
860 if (!particleLink.isValid()) { continue; }
861 const float genEta = particleLink->momentum().pseudoRapidity();
862 const float genPt = particleLink->momentum().perp(); // MeV
863 // reject unstable particles
864 if (!MC::isStable(particleLink.cptr())) { continue; }
865 // reject secondaries and low pT (<400 MeV) pileup
866 if (HepMC::is_simulation_particle(particleLink.cptr()) || particleLink.barcode() == 0 /*HepMC::no_truth_link(particleLink)*/) { continue; } // FIXME
867 // reject far forward particles
868 if (std::fabs(genEta) > m_maxEta) { continue; }
869 // "bestTruthLink" links to the highest pt particle
870 if (bestPt < genPt) {
871 bestPt = genPt;
872 bestTruthLink = &particleLink;
873 }
874 #ifdef HEPMC3
875 parentMask |= FPGATrackSimInputUtils::construct_truth_bitmap(std::shared_ptr<const HepMC3::GenParticle>(particleLink.cptr()));
876 #else
877 parentMask |= FPGATrackSimInputUtils::construct_truth_bitmap(particleLink.cptr());
878 #endif
879 // check SDO
880 } // end for each contributing particle
881 return bestTruthLink;
882}
std::pair< HepMcParticleLink, float > Deposit
const ParentBitmask construct_truth_bitmap(HepMC::ConstGenParticlePtr p)
bool is_simulation_particle(const T &p)
Method to establish if a particle (or barcode) was created during the simulation (TODO update to be s...
bool isStable(const T &p)
Identify if the particle is stable, i.e. has not decayed.

◆ initialize()

StatusCode FPGATrackSimSGToRawHitsTool::initialize ( )
overridevirtual

Definition at line 51 of file FPGATrackSimSGToRawHitsTool.cxx.

51 {
52
53 ATH_MSG_DEBUG("FPGATrackSimSGToRawHitsTool::initialize()");
54
55 if(!m_truthToTrack.empty() ) ATH_CHECK(m_truthToTrack.retrieve());
56 if(!m_extrapolator.empty()) ATH_CHECK(m_extrapolator.retrieve());
57 ATH_CHECK(m_beamSpotKey.initialize());
58
59 SmartIF<IPartPropSvc> partPropSvc{service("PartPropSvc")};
60 ATH_CHECK(partPropSvc.isValid());
61 m_particleDataTable = partPropSvc->PDT();
62
63 ATH_CHECK(detStore()->retrieve(m_PIX_mgr, "ITkPixel"));
64 ATH_CHECK(detStore()->retrieve(m_pixelId, "PixelID"));
65 ATH_CHECK(detStore()->retrieve(m_SCT_mgr, "ITkStrip"));
66 ATH_CHECK(detStore()->retrieve(m_sctId, "SCT_ID"));
67
68 ATH_CHECK(m_eventInfoKey.initialize());
71
73
77 ATH_CHECK(m_pixelRDOKey.initialize());
78 ATH_CHECK(m_stripRDOKey.initialize());
79
80 ATH_MSG_DEBUG("Initialization complete");
81 return StatusCode::SUCCESS;
82}
#define ATH_CHECK
Evaluate an expression and check for errors.
ToolHandle< Trk::ITruthToTrack > m_truthToTrack
tool to create track parameters from a gen particle
const InDetDD::SiDetectorManager * m_SCT_mgr
SG::ReadHandleKey< InDetSimDataCollection > m_stripSDOKey
SG::ReadHandleKey< SCT_RDO_Container > m_stripRDOKey
SG::ReadHandleKey< xAOD::EventInfo > m_eventInfoKey
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_offlineTracksKey
ToolHandle< Trk::IExtrapolator > m_extrapolator
ToolHandle for Extrapolator.
SG::ReadCondHandleKey< InDet::BeamSpotData > m_beamSpotKey
SG::ReadHandleKey< InDet::SiClusterContainer > m_sctClusterContainerKey
SG::ReadHandleKey< PixelRDO_Container > m_pixelRDOKey
SG::ReadHandleKey< McEventCollection > m_mcCollectionKey
const HepPDT::ParticleDataTable * m_particleDataTable

◆ readData()

StatusCode FPGATrackSimSGToRawHitsTool::readData ( FPGATrackSimEventInputHeader * header,
const EventContext & eventContext ) const
overridevirtual

This function get from the SG the inner detector raw hits and prepares them for FPGATrackSim simulation.

Definition at line 92 of file FPGATrackSimSGToRawHitsTool.cxx.

93{
94
95 auto eventInfo = SG::makeHandle(m_eventInfoKey, eventContext);
96 //Filled to variable / start event
97 FPGATrackSimEventInfo event_info;
98 event_info.setRunNumber(eventInfo->runNumber());
99 event_info.setEventNumber(eventInfo->eventNumber());
100 event_info.setLB(eventInfo->lumiBlock());
101 event_info.setBCID(eventInfo->bcid());
102 event_info.setaverageInteractionsPerCrossing(eventInfo->averageInteractionsPerCrossing());
103 event_info.setactualInteractionsPerCrossing(eventInfo->actualInteractionsPerCrossing());
104 event_info.setextendedLevel1ID(eventInfo->extendedLevel1ID());
105 event_info.setlevel1TriggerType(eventInfo->level1TriggerType());
106 eventHeader->newEvent(event_info);
107
108 HitIndexMap hitIndexMap; // keep running index event-unique to each hit
109 HitIndexMap pixelClusterIndexMap;
110 // get pixel and sct cluster containers
111 // dump raw silicon data
112 ATH_MSG_DEBUG("Dump raw silicon data");
113 ATH_CHECK(readRawSilicon(eventHeader, hitIndexMap, eventContext));
114 FPGATrackSimOptionalEventInfo optional;
116 std::vector <FPGATrackSimCluster> clusters;
117 ATH_CHECK(readOfflineClusters(clusters, eventContext));
118 for (const auto& cluster : clusters) optional.addOfflineCluster(cluster);
119 ATH_MSG_DEBUG("Saved " << optional.nOfflineClusters() << " offline clusters");
120 ATH_CHECK(dumpPixelClusters(pixelClusterIndexMap, eventContext));
121 }
122 if (m_readTruthTracks) {
123 std::vector <FPGATrackSimTruthTrack> truth;
124 ATH_CHECK(readTruthTracks(truth, eventContext));
125 for (const FPGATrackSimTruthTrack& trk : truth) optional.addTruthTrack(trk);
126 ATH_MSG_DEBUG("Saved " << optional.nTruthTracks() << " truth tracks");
127 }
128 std::vector <FPGATrackSimOfflineTrack> offline;
130 ATH_CHECK(readOfflineTracks(offline, eventContext));
131 for (const FPGATrackSimOfflineTrack& trk : offline) optional.addOfflineTrack(trk);
132 ATH_MSG_DEBUG("Saved " << optional.nOfflineTracks() << " offline tracks");
133 }
134 eventHeader->setOptional(optional);
135 ATH_MSG_DEBUG(*eventHeader);
136 ATH_MSG_DEBUG("End of execute()");
137 return StatusCode::SUCCESS;
138}
void setaverageInteractionsPerCrossing(const int &val)
void setBCID(const int &val)
void setLB(const int &val)
void setRunNumber(const unsigned long &val)
void setlevel1TriggerType(const unsigned int &val)
void setactualInteractionsPerCrossing(const int &val)
void setEventNumber(const unsigned long &val)
void setextendedLevel1ID(const unsigned int &val)
void addOfflineCluster(const FPGATrackSimCluster &c) const
void addOfflineTrack(const FPGATrackSimOfflineTrack &t) const
void addTruthTrack(const FPGATrackSimTruthTrack &t) const
Gaudi::Property< bool > m_readOfflineTracks
StatusCode readOfflineClusters(std::vector< FPGATrackSimCluster > &Clusters, const EventContext &eventContext) const
StatusCode readRawSilicon(FPGATrackSimEventInputHeader *header, HitIndexMap &hitIndexMap, const EventContext &eventContext) const
StatusCode readOfflineTracks(std::vector< FPGATrackSimOfflineTrack > &Track, const EventContext &eventContext) const
Gaudi::Property< bool > m_readOfflineClusters
std::map< Identifier, int > HitIndexMap
Gaudi::Property< bool > m_readTruthTracks
StatusCode readTruthTracks(std::vector< FPGATrackSimTruthTrack > &truth, const EventContext &eventContext) const
StatusCode dumpPixelClusters(HitIndexMap &pixelClusterIndexMap, const EventContext &eventContext) const

◆ readOfflineClusters()

StatusCode FPGATrackSimSGToRawHitsTool::readOfflineClusters ( std::vector< FPGATrackSimCluster > & Clusters,
const EventContext & eventContext ) const
private

Definition at line 559 of file FPGATrackSimSGToRawHitsTool.cxx.

561{
562
563 //Lets do the Pixel clusters first
564 //Loopover the pixel clusters and convert them into a FPGATrackSimCluster for storage
565 // Dump pixel clusters. They're in m_pixelContainer
566 auto pixelSDOHandle = SG::makeHandle(m_pixelSDOKey, eventContext);
567 auto pixelClusterContainerHandler = SG::makeHandle(m_pixelClusterContainerKey, eventContext);
568 for (const InDet::SiClusterCollection* pixelClusterCollection : *pixelClusterContainerHandler) {
569 if (pixelClusterCollection == nullptr) {
570 ATH_MSG_DEBUG("pixelClusterCollection not available!");
571 continue;
572 }
573 const int size = pixelClusterCollection->size();
574 ATH_MSG_DEBUG("PixelClusterCollection found with " << size << " clusters");
575 for (const InDet::SiCluster* cluster : *pixelClusterCollection) {
576
577 // if there is simulation truth available, try to retrieve the "most likely" barcode for this pixel cluster.
579 const HepMcParticleLink* bestTruthLink{};
580 if (!m_pixelSDOKey.empty()) {
581 for (const Identifier& rdoId : cluster->rdoList()) {
582 const InDetDD::SiDetectorElement* sielement = m_PIX_mgr->getDetectorElement(rdoId);
583 assert(sielement);
584 InDetDD::SiCellId cellID = sielement->cellIdFromIdentifier(rdoId);
585 const int nCells = sielement->numberOfConnectedCells(cellID);
586 InDetSimDataCollection::const_iterator iter(pixelSDOHandle->find(rdoId));
587 // this might be the ganged pixel copy.
588 if (nCells > 1 && iter == pixelSDOHandle->end()) {
589 InDetDD::SiReadoutCellId SiRC(m_pixelId->phi_index(rdoId), m_pixelId->eta_index(rdoId));
590 for (int ii = 0; ii < nCells && iter == pixelSDOHandle->end(); ++ii) {
591 iter = pixelSDOHandle->find(sielement->identifierFromCellId(sielement->design().connectedCell(SiRC, ii)));
592 }
593 } // end search for correct ganged pixel
594 // if SDO found for this pixel, associate the particle. otherwise leave unassociated.
595 if (iter != pixelSDOHandle->end()) { bestTruthLink = getTruthInformation(iter, parentMask); }
596 } // if we have pixel sdo's available
597 }
598 HepMC::ConstGenParticlePtr bestParent = (bestTruthLink) ? bestTruthLink->cptr() : nullptr;
599
600 Identifier theID = cluster->identify();
601 //cluster object to be written out
602 FPGATrackSimCluster clusterOut;
603 //Rawhit object to represent the cluster
604 FPGATrackSimHit clusterEquiv;
605 //Lets get the information of this pixel cluster
606 const InDetDD::SiDetectorElement* sielement = m_PIX_mgr->getDetectorElement(theID);
607 assert(sielement);
608 const InDetDD::SiLocalPosition localPos = sielement->rawLocalPositionOfCell(theID);
609 const Amg::Vector3D globalPos(sielement->globalPosition(localPos));
610 clusterEquiv.setHitType(HitType::clustered);
611 clusterEquiv.setX(globalPos.x());
612 clusterEquiv.setY(globalPos.y());
613 clusterEquiv.setZ(globalPos.z());
614 clusterEquiv.setDetType(SiliconTech::pixel);
615 clusterEquiv.setIdentifierHash(sielement->identifyHash());
616 clusterEquiv.setIdentifier(sielement->identify().get_identifier32().get_compact());
617
618 int barrel_ec = m_pixelId->barrel_ec(theID);
619 if (barrel_ec == 0)
621 else if (barrel_ec == 2)
623 else if (barrel_ec == -2)
625
626 clusterEquiv.setLayerDisk(m_pixelId->layer_disk(theID));
627 clusterEquiv.setPhiModule(m_pixelId->phi_module(theID));
628 clusterEquiv.setEtaModule(m_pixelId->eta_module(theID));
629 clusterEquiv.setPhiIndex(m_pixelId->phi_index(theID));
630 clusterEquiv.setEtaIndex(m_pixelId->eta_index(theID));
631 clusterEquiv.setPhiCoord(localPos.xPhi());
632 clusterEquiv.setEtaCoord(localPos.xEta());
633
634 clusterEquiv.setPhiWidth(cluster->width().colRow()[1]);
635 clusterEquiv.setEtaWidth(cluster->width().colRow()[0]);
636 //Save the truth here as the MultiTruth object is only transient
637 if (bestParent) {
638 clusterEquiv.setEventIndex(bestTruthLink->eventIndex());
639 clusterEquiv.setBarcode(bestTruthLink->barcode()); // FIXME barcode-based
640 clusterEquiv.setUniqueID(bestTruthLink->id());
641 }
642 else {
643 clusterEquiv.setEventIndex(std::numeric_limits<long>::max());
644 clusterEquiv.setBarcode(std::numeric_limits<HepMcParticleLink::barcode_type>::max());
645 clusterEquiv.setUniqueID(std::numeric_limits<HepMcParticleLink::barcode_type>::max());
646 }
647
648 clusterEquiv.setBarcodePt(static_cast<unsigned long>(std::ceil(bestParent ? bestParent->momentum().perp() : 0.)));
649 clusterEquiv.setParentageMask(parentMask.to_ulong());
650 clusterOut.setClusterEquiv(clusterEquiv);
651 clusters.push_back(clusterOut);
652 }
653 }
654
655 //Now lets do the strip clusters
656 //Loopover the pixel clusters and convert them into a FPGATrackSimCluster for storage
657 // Dump pixel clusters. They're in m_pixelContainer
658 auto stripSDOHandle = SG::makeHandle(m_stripSDOKey, eventContext);
659 ATH_MSG_DEBUG("Found SCT SDO Map");
660 auto stripRDOHandle = SG::makeHandle(m_stripRDOKey, eventContext);
661
662 for (const InDetRawDataCollection<SCT_RDORawData>* SCT_Collection : *stripRDOHandle) {
663 if (SCT_Collection == nullptr) { continue; }
664 for (const SCT_RDORawData* sctRawData : *SCT_Collection) {
665 const Identifier rdoId = sctRawData->identify();
666 // get the det element from the det element collection
667 const InDetDD::SiDetectorElement* sielement = m_SCT_mgr->getDetectorElement(rdoId);
668 const InDetDD::SiDetectorDesign& design = dynamic_cast<const InDetDD::SiDetectorDesign&>(sielement->design());
669 const InDetDD::SiLocalPosition localPos = design.localPositionOfCell(m_sctId->strip(rdoId));
670 const Amg::Vector3D gPos = sielement->globalPosition(localPos);
671 // if there is simulation truth available, try to retrieve the
672 // "most likely" barcode for this strip.
674 const HepMcParticleLink* bestTruthLink{};
675 if (!m_stripSDOKey.empty()) {
676 InDetSimDataCollection::const_iterator iter(stripSDOHandle->find(rdoId));
677 // if SDO found for this pixel, associate the particle
678 if (iter != stripSDOHandle->end()) { bestTruthLink = getTruthInformation(iter, parentMask); }
679 } // end if sct truth available
680 HepMC::ConstGenParticlePtr bestParent = (bestTruthLink) ? bestTruthLink->cptr() : nullptr;
681
682 // push back the hit information to DataInput for HitList , copy from RawInput.cxx
683 FPGATrackSimCluster clusterOut;
684 FPGATrackSimHit clusterEquiv;
685 clusterEquiv.setHitType(HitType::clustered);
686 clusterEquiv.setX(gPos.x());
687 clusterEquiv.setY(gPos.y());
688 clusterEquiv.setZ(gPos.z());
689 clusterEquiv.setDetType(SiliconTech::strip);
690 clusterEquiv.setIdentifierHash(sielement->identifyHash());
691 clusterEquiv.setIdentifier(sielement->identify().get_identifier32().get_compact());
692
693 int barrel_ec = m_sctId->barrel_ec(rdoId);
694 if (barrel_ec == 0)
696 else if (barrel_ec == 2)
698 else if (barrel_ec == -2)
700
701 clusterEquiv.setLayerDisk(m_sctId->layer_disk(rdoId));
702 clusterEquiv.setPhiModule(m_sctId->phi_module(rdoId));
703 clusterEquiv.setEtaModule(m_sctId->eta_module(rdoId));
704 clusterEquiv.setPhiIndex(m_sctId->strip(rdoId));
705 clusterEquiv.setEtaIndex(m_sctId->row(rdoId));
706 clusterEquiv.setPhiCoord(localPos.xPhi());
707 clusterEquiv.setEtaCoord(localPos.xEta());
708 clusterEquiv.setSide(m_sctId->side(rdoId));
709 //I think this is the strip "cluster" width
710 clusterEquiv.setPhiWidth(sctRawData->getGroupSize());
711 //Save the truth here as the MultiTruth object is only transient
712 if (bestParent) {
713 clusterEquiv.setEventIndex(bestTruthLink->eventIndex());
714 clusterEquiv.setBarcode(bestTruthLink->barcode()); // FIXME barcode-based
715 clusterEquiv.setUniqueID(bestTruthLink->id());
716 }
717 else {
718 clusterEquiv.setEventIndex(std::numeric_limits<long>::max());
719 clusterEquiv.setBarcode(std::numeric_limits<HepMcParticleLink::barcode_type>::max());
720 clusterEquiv.setUniqueID(std::numeric_limits<HepMcParticleLink::barcode_type>::max());
721 }
722
723 clusterEquiv.setBarcodePt(static_cast<unsigned long>(std::ceil(bestParent ? bestParent->momentum().perp() : 0.)));
724 clusterEquiv.setParentageMask(parentMask.to_ulong());
725 clusterOut.setClusterEquiv(clusterEquiv);
726 clusters.push_back(clusterOut);
727 } // end for each RDO in the strip collection
728 } // end for each strip RDO collection
729 // dump all RDO's and SDO's for a given event, for debugging purposes
730
731 return StatusCode::SUCCESS;
732}
void setClusterEquiv(const FPGATrackSimHit &input)
void setPhiModule(unsigned v)
void setIdentifierHash(unsigned v)
void setEtaIndex(unsigned v)
void setEventIndex(long v)
void setPhiIndex(unsigned v)
void setHitType(HitType type)
void setPhiCoord(float v)
void setIdentifier(unsigned int v)
void setZ(float v)
void setBarcode(const HepMcParticleLink::barcode_type &v)
void setX(float v)
void setParentageMask(unsigned long v)
void setBarcodePt(float v)
void setY(float v)
void setLayerDisk(unsigned v)
void setEtaModule(int v)
void setSide(unsigned v)
void setEtaCoord(float v)
void setEtaWidth(unsigned v)
void setPhiWidth(unsigned v)
void setUniqueID(const HepMcParticleLink::barcode_type &v)
void setDetectorZone(DetectorZone detZone)
void setDetType(SiliconTech detType)
value_type get_compact() const
Get the compact id.
Identifier32 get_identifier32() const
Get the 32-bit version Identifier, will be invalid if >32 bits needed.
virtual SiLocalPosition localPositionOfCell(const SiCellId &cellId) const =0
readout or diode id -> position.
double xPhi() const
position along phi direction:
double xEta() const
position along eta direction:
virtual IdentifierHash identifyHash() const override final
identifier hash (inline)
HepGeom::Point3D< double > globalPosition(const HepGeom::Point3D< double > &localPos) const
transform a reconstruction local position into a global position (inline):
virtual Identifier identify() const override final
identifier of this detector element (inline)
Amg::Vector2D rawLocalPositionOfCell(const SiCellId &cellId) const
Returns position (center) of cell.
Eigen::Matrix< double, 3, 1 > Vector3D
const GenParticle * ConstGenParticlePtr
Definition GenParticle.h:38

◆ readOfflineTracks()

StatusCode FPGATrackSimSGToRawHitsTool::readOfflineTracks ( std::vector< FPGATrackSimOfflineTrack > & Track,
const EventContext & eventContext ) const
private

Definition at line 141 of file FPGATrackSimSGToRawHitsTool.cxx.

142{
143 auto offlineTracksHandle = SG::makeHandle(m_offlineTracksKey, eventContext);
144 ATH_MSG_DEBUG("read Offline tracks, size= " << offlineTracksHandle->size());
145
146 int iTrk = -1;
147 for (const xAOD::TrackParticle* trackParticle : *offlineTracksHandle) {
148 iTrk++;
149 FPGATrackSimOfflineTrack tmpOfflineTrack;
150 tmpOfflineTrack.setQOverPt(trackParticle->pt() > 0 ? trackParticle->charge() / trackParticle->pt() : 0);
151 tmpOfflineTrack.setEta(trackParticle->eta());
152 tmpOfflineTrack.setPhi(trackParticle->phi());
153 tmpOfflineTrack.setD0(trackParticle->d0());
154 tmpOfflineTrack.setZ0(trackParticle->z0());
155
156 const Trk::TrackStates* trackStates = trackParticle->track()->trackStateOnSurfaces();
157 if (trackStates == nullptr) {
158 ATH_MSG_ERROR("missing trackStatesOnSurface");
159 return StatusCode::FAILURE;
160 }
161 for (const Trk::TrackStateOnSurface* tsos : *trackStates) {
162 if (tsos == nullptr) continue;
164 const Trk::MeasurementBase* measurement = tsos->measurementOnTrack();
165 if (tsos->trackParameters() != nullptr &&
166 tsos->trackParameters()->associatedSurface().associatedDetectorElement() != nullptr &&
167 tsos->trackParameters()->associatedSurface().associatedDetectorElement()->identify() != 0
168 ) {
169 const Trk::RIO_OnTrack* hit = dynamic_cast <const Trk::RIO_OnTrack*>(measurement);
170 const Identifier& hitId = hit->identify();
171 FPGATrackSimOfflineHit tmpOfflineHit;
172 if (m_pixelId->is_pixel(hitId)) {
173 tmpOfflineHit.setIsPixel(true);
174 tmpOfflineHit.setIsBarrel(m_pixelId->is_barrel(hitId));
175
176 const InDetDD::SiDetectorElement* sielement = m_PIX_mgr->getDetectorElement(hitId);
177 tmpOfflineHit.setClusterID(sielement->identifyHash());
178 tmpOfflineHit.setTrackNumber(iTrk);
179 tmpOfflineHit.setLayer(m_pixelId->layer_disk(hitId));
180 tmpOfflineHit.setLocX((float)measurement->localParameters()[Trk::locX]);
181 tmpOfflineHit.setLocY((float)measurement->localParameters()[Trk::locY]);
182 }
183 else if (m_sctId->is_sct(hitId)) {
184 tmpOfflineHit.setIsPixel(false);
185 tmpOfflineHit.setIsBarrel(m_sctId->is_barrel(hitId));
186 const InDetDD::SiDetectorElement* sielement = m_SCT_mgr->getDetectorElement(hitId);
187 tmpOfflineHit.setClusterID(sielement->identifyHash());
188 tmpOfflineHit.setTrackNumber(iTrk);
189 tmpOfflineHit.setLayer(m_sctId->layer_disk(hitId));
190 tmpOfflineHit.setLocX(((float)measurement->localParameters()[Trk::locX]));
191 tmpOfflineHit.setLocY(-99999.9);
192 }
193 tmpOfflineTrack.addHit(tmpOfflineHit);
194 }
195 }
196 }
197 offline.push_back(tmpOfflineTrack);
198 }//end of loop over tracks
199
200
201 return StatusCode::SUCCESS;
202}
#define ATH_MSG_ERROR(x)
void addHit(const FPGATrackSimOfflineHit &s)
const LocalParameters & localParameters() const
Interface method to get the LocalParameters.
Identifier identify() const
return the identifier -extends MeasurementBase
@ Measurement
This is a measurement, and will at least contain a Trk::MeasurementBase.
DataVector< const Trk::TrackStateOnSurface > TrackStates
@ locY
local cartesian
Definition ParamDefs.h:38
@ locX
Definition ParamDefs.h:37
TrackParticle_v1 TrackParticle
Reference the current persistent version:

◆ readPixelSimulation()

StatusCode FPGATrackSimSGToRawHitsTool::readPixelSimulation ( FPGATrackSimEventInputHeader * header,
HitIndexMap & hitIndexMap,
unsigned int & hitIndex,
const EventContext & eventContext ) const
private

Definition at line 224 of file FPGATrackSimSGToRawHitsTool.cxx.

228 {
229
230 auto pixelSDOHandle = SG::makeHandle(m_pixelSDOKey, eventContext);
231 auto pixelRDOHandle = SG::makeHandle(m_pixelRDOKey, eventContext);
232
233 ATH_MSG_DEBUG("Found Pixel SDO Map");
234
235 for (const InDetRawDataCollection<PixelRDORawData>* pixel_rdoCollection : *pixelRDOHandle) {
236 if (pixel_rdoCollection == nullptr) { continue; }
237 // loop on all RDOs
238 for (const PixelRDORawData* pixelRawData : *pixel_rdoCollection) {
239 Identifier rdoId = pixelRawData->identify();
240 // get the det element from the det element collection
241 const InDetDD::SiDetectorElement* sielement = m_PIX_mgr->getDetectorElement(rdoId); assert(sielement);
242
243 Amg::Vector2D localPos = sielement->rawLocalPositionOfCell(rdoId);
244 Amg::Vector3D globalPos = sielement->globalPosition(localPos);
245 InDetDD::SiCellId cellID = sielement->cellIdFromIdentifier(rdoId);
246
247 // update map between pixel identifier and event-unique hit index.
248 // ganged pixels (nCells==2) get two entries.
249 hitIndexMap[rdoId] = hitIndex;
250 const int nCells = sielement->numberOfConnectedCells(cellID);
251 if (nCells == 2) {
252 const InDetDD::SiCellId tmpCell = sielement->connectedCell(cellID, 1);
253 const Identifier tmpId = sielement->identifierFromCellId(tmpCell);
254 hitIndexMap[tmpId] = hitIndex; // add second entry for ganged pixel ID
255 }
256 // if there is simulation truth available, try to retrieve the "most likely" barcode for this pixel.
258 const HepMcParticleLink* bestTruthLink{};
259 if (!m_pixelSDOKey.empty()) {
260 InDetSimDataCollection::const_iterator iter(pixelSDOHandle->find(rdoId));
261 if (nCells > 1 && iter == pixelSDOHandle->end()) {
262 InDetDD::SiReadoutCellId SiRC(m_pixelId->phi_index(rdoId), m_pixelId->eta_index(rdoId));
263 for (int ii = 0; ii < nCells && iter == pixelSDOHandle->end(); ++ii) {
264 iter = pixelSDOHandle->find(sielement->identifierFromCellId(sielement->design().connectedCell(SiRC, ii)));
265 }
266 } // end search for correct ganged pixel
267 // if SDO found for this pixel, associate the particle. otherwise leave unassociated.
268 if (iter != pixelSDOHandle->end()) { bestTruthLink = getTruthInformation(iter, parentMask); }
269 } // end if pixel truth available
270 HepMC::ConstGenParticlePtr bestParent = (bestTruthLink) ? bestTruthLink->cptr() : nullptr;
271 ++hitIndex;
272
273 // push back the hit information to DataInput for HitList
274 FPGATrackSimHit tmpSGhit;
277 tmpSGhit.setIdentifierHash(sielement->identifyHash());
278 tmpSGhit.setIdentifier(sielement->identify().get_identifier32().get_compact());
279 tmpSGhit.setRdoIdentifier(rdoId.get_compact()); // full 64 bit hit identifier
280
281 int barrel_ec = m_pixelId->barrel_ec(rdoId);
282 if (barrel_ec == 0)
284 else if (barrel_ec == 2)
286 else if (barrel_ec == -2)
288
289 tmpSGhit.setLayerDisk(m_pixelId->layer_disk(rdoId));
290 tmpSGhit.setPhiModule(m_pixelId->phi_module(rdoId));
291 tmpSGhit.setEtaModule(m_pixelId->eta_module(rdoId));
292 tmpSGhit.setPhiIndex(m_pixelId->phi_index(rdoId));
293 tmpSGhit.setEtaIndex(m_pixelId->eta_index(rdoId));
294 tmpSGhit.setPhiCoord(localPos[0]);
295 tmpSGhit.setEtaCoord(localPos[1]);
296 tmpSGhit.setEtaWidth(0);
297 tmpSGhit.setPhiWidth(0);
298 tmpSGhit.setX(globalPos[Amg::x]);
299 tmpSGhit.setY(globalPos[Amg::y]);
300 tmpSGhit.setZ(globalPos[Amg::z]);
301 tmpSGhit.setToT(pixelRawData->getToT());
302 tmpSGhit.setisValidForITkHit(true); // Pixel clusters are close enough right now that they all can be considered valid for ITK
303 if (bestParent) {
304 tmpSGhit.setEventIndex(bestTruthLink->eventIndex());
305 tmpSGhit.setBarcode(bestTruthLink->barcode()); // FIXME barcode-based
306 tmpSGhit.setUniqueID(bestTruthLink->id()); // May need fixing when uid will be used.
307 }
308 else {
309 tmpSGhit.setEventIndex(std::numeric_limits<long>::max());
310 tmpSGhit.setBarcode(std::numeric_limits<HepMcParticleLink::barcode_type>::max());
311 tmpSGhit.setUniqueID(std::numeric_limits<HepMcParticleLink::barcode_type>::max());
312 }
313
314 tmpSGhit.setBarcodePt(static_cast<unsigned long>(std::ceil(bestParent ? bestParent->momentum().perp() : 0.)));
315 tmpSGhit.setParentageMask(parentMask.to_ulong());
316
317 if (m_doMultiTruth) {
318 // Add truth
319 FPGATrackSimMultiTruth mt;
320 FPGATrackSimMultiTruth::Barcode uniqueID(tmpSGhit.getEventIndex(), tmpSGhit.getBarcode()); // FIXME barcode-based
321 mt.maximize(uniqueID, tmpSGhit.getBarcodePt()); // FIXME barcode-based
322 tmpSGhit.setTruth(mt);
323 }
324
325 eventHeader->addHit(tmpSGhit);
326 } // end for each RDO in the collection
327 } // for each pixel RDO collection
328
329 return StatusCode::SUCCESS;
330}
long getEventIndex() const
float getBarcodePt() const
void setisValidForITkHit(bool v)
void setRdoIdentifier(Identifier::value_type v)
HepMcParticleLink::barcode_type getBarcode() const
void setToT(unsigned v)
void setTruth(const FPGATrackSimMultiTruth &v)
void maximize(const FPGATrackSimMultiTruth::Barcode &code, const FPGATrackSimMultiTruth::Weight &weight)
std::pair< unsigned long, unsigned long > Barcode
value_type get_compact() const
Get the compact id.
SiCellId connectedCell(const SiCellId cellId, int number) const
Get the cell ids sharing the readout for this cell.
Eigen::Matrix< double, 2, 1 > Vector2D
int uniqueID(const T &p)

◆ readRawSilicon()

StatusCode FPGATrackSimSGToRawHitsTool::readRawSilicon ( FPGATrackSimEventInputHeader * header,
HitIndexMap & hitIndexMap,
const EventContext & eventContext ) const
private

Definition at line 208 of file FPGATrackSimSGToRawHitsTool.cxx.

212{
213 ATH_MSG_DEBUG("read silicon hits");
214 unsigned int hitIndex = 0u;
215
216 ATH_CHECK(readPixelSimulation(eventHeader, hitIndexMap, hitIndex, eventContext));
217 ATH_CHECK(readStripSimulation(eventHeader, hitIndexMap, hitIndex, eventContext));
218
219 return StatusCode::SUCCESS;
220}
StatusCode readPixelSimulation(FPGATrackSimEventInputHeader *header, HitIndexMap &hitIndexMap, unsigned int &hitIndex, const EventContext &eventContext) const
StatusCode readStripSimulation(FPGATrackSimEventInputHeader *header, HitIndexMap &hitIndexMap, unsigned int &hitIndex, const EventContext &eventContext) const
@ u
Enums for curvilinear frames.
Definition ParamDefs.h:77

◆ readStripSimulation()

StatusCode FPGATrackSimSGToRawHitsTool::readStripSimulation ( FPGATrackSimEventInputHeader * header,
HitIndexMap & hitIndexMap,
unsigned int & hitIndex,
const EventContext & eventContext ) const
private

Definition at line 333 of file FPGATrackSimSGToRawHitsTool.cxx.

337 {
338
339 constexpr int MaxChannelinStripRow = 128;
340
341 auto stripSDOHandle = SG::makeHandle(m_stripSDOKey, eventContext);
342 ATH_MSG_DEBUG("Found SCT SDO Map");
343 auto stripRDOHandle = SG::makeHandle(m_stripRDOKey, eventContext);
344 for (const InDetRawDataCollection<SCT_RDORawData>* SCT_Collection : *stripRDOHandle) {
345 if (SCT_Collection == nullptr) { continue; }
346
347 std::map<int, bool> firedStrips;
348 std::map<int, const SCT_RDORawData*> firedStripsToRDO;
349 // Preprocess the SCT collection hits to get information for encoding strip in ITK format
350 // All strips fired read into a map to an overview of full module that should be used to encode
351 // the data into the ITk formatl
352 for (const SCT_RDORawData* sctRawData : *SCT_Collection)
353 {
354 const Identifier rdoId = sctRawData->identify();
355 const int baseLineStrip{m_sctId->strip(rdoId)};
356 for(int i = 0; i < sctRawData->getGroupSize(); i++) {
357 firedStrips[baseLineStrip+ i] = true;
358 firedStripsToRDO[baseLineStrip + i] = sctRawData;
359 }
360 }
361
362 // Loop over the fired hits and encode them in the ITk strips hit map
363 // It find unique hits in the list that can be encoded and don't overlap
364 std::map<int, int> stripEncodingForITK;
365 std::map<int, const SCT_RDORawData* > stripEncodingForITKToRDO;
366 for(const auto& [stripID, fired]: firedStrips)
367 {
368 // Don't use the strip that has been set false.
369 // This will be the case where neighbouring strip will "used up in the cluster"
370 // And then we don't want to re use them
371 if(!fired) continue;
372
373 // Check the next 3 hits if they are there and have a hit in them
374 std::bitset<3> hitMap;
375
376
377 // Get the current chip id of the strip
378 int currChipID = stripID / MaxChannelinStripRow;
379 // Compute the maximum stripID this chip can have
380 int maxStripIDForCurrChip = (currChipID + 1) * MaxChannelinStripRow;
381
382 for(int i = 0; i < 3; i++)
383 {
384 // We don't want to "cluster" strips that are outside the range of this chip
385 if((stripID + 1 + i) >= maxStripIDForCurrChip) continue;
386
387 if(firedStrips.find(stripID + 1 + i) != firedStrips.end())
388 {
389 if(firedStrips.at(stripID + 1 + i))
390 {
391 hitMap[2 - i] = 1;
392 firedStrips[stripID + 1 + i] = false;
393 }
394 else
395 {
396 hitMap[2 - i] = 0;
397 }
398 }
399 }
400
401 // Encode the hit map into a int
402 stripEncodingForITK[stripID] = (int)(hitMap.to_ulong());
403 stripEncodingForITKToRDO[stripID] = firedStripsToRDO[stripID];
404 }
405
406 // Actual creation of the FPGAHit objects
407 for(const auto& [stripID, fired]: firedStrips)
408 {
409 const SCT_RDORawData* sctRawData = firedStripsToRDO[stripID];
410 const Identifier rdoId = sctRawData->identify();
411 // get the det element from the det element collection
412 const InDetDD::SiDetectorElement* sielement = m_SCT_mgr->getDetectorElement(rdoId);
413 const InDetDD::SiDetectorDesign& design = dynamic_cast<const InDetDD::SiDetectorDesign&>(sielement->design());
414
415 InDetDD::SiCellId frontId(stripID);
416 Amg::Vector2D localPos = design.localPositionOfCell(frontId);
417 std::pair<Amg::Vector3D, Amg::Vector3D> endsOfStrip = sielement->endsOfStrip(localPos);
418
419 hitIndexMap[rdoId] = hitIndex;
420 ++hitIndex;
421 // if there is simulation truth available, try to retrieve the
422 // "most likely" barcode for this strip.
424 const HepMcParticleLink* bestTruthLink{};
425 if (!m_stripSDOKey.empty()) {
426 InDetSimDataCollection::const_iterator iter(stripSDOHandle->find(rdoId));
427 // if SDO found for this strip, associate the particle
428 if (iter != stripSDOHandle->end()) { bestTruthLink = getTruthInformation(iter, parentMask); }
429 } // end if sct truth available
430 HepMC::ConstGenParticlePtr bestParent = (bestTruthLink) ? bestTruthLink->cptr() : nullptr;
431 // push back the hit information to DataInput for HitList , copy from RawInput.cxx
432
433 FPGATrackSimHit tmpSGhit;
436 tmpSGhit.setIdentifierHash(sielement->identifyHash());
437 tmpSGhit.setIdentifier(sielement->identify().get_identifier32().get_compact());
438 tmpSGhit.setRdoIdentifier(rdoId.get_compact()); // full 64 bit hit identifier
439
440 int barrel_ec = m_sctId->barrel_ec(rdoId);
441 if (barrel_ec == 0)
443 else if (barrel_ec == 2)
445 else if (barrel_ec == -2)
447
448 tmpSGhit.setLayerDisk(m_sctId->layer_disk(rdoId));
449 tmpSGhit.setPhiModule(m_sctId->phi_module(rdoId));
450 tmpSGhit.setEtaModule(m_sctId->eta_module(rdoId));
451 tmpSGhit.setPhiIndex(stripID);
452 tmpSGhit.setEtaIndex(m_sctId->row(rdoId));
453 tmpSGhit.setPhiCoord(localPos[0]);
454 tmpSGhit.setEtaCoord(localPos[1]);
455 tmpSGhit.setSide(m_sctId->side(rdoId));
456 tmpSGhit.setEtaWidth(0);
457 tmpSGhit.setPhiWidth(1);
458 if (bestParent) {
459 tmpSGhit.setEventIndex(bestTruthLink->eventIndex());
460 tmpSGhit.setBarcode(bestTruthLink->barcode()); // FIXME barcode-based
461 tmpSGhit.setUniqueID(bestTruthLink->id());
462 }
463 else {
464 tmpSGhit.setEventIndex(std::numeric_limits<long>::max());
465 tmpSGhit.setBarcode(std::numeric_limits<HepMcParticleLink::barcode_type>::max());
466 tmpSGhit.setUniqueID(std::numeric_limits<HepMcParticleLink::barcode_type>::max());
467 }
468
469 // If the strip has been identified by the previous for loop as a valid hit that can be encoded into ITk Strip format
470 if(stripEncodingForITK.find(stripID) != stripEncodingForITK.end())
471 {
472 // Each ITK ABC chip reads 128 channels in one row, so we just need to divide the current strip with 128 to get the chip index
473 // for the Strip ID, it is the remainder left after dividing by 128
474 int chipID = stripID / MaxChannelinStripRow;
475 int ITkStripID = stripID % MaxChannelinStripRow;
476
477 // for each ABC chip readout, each reads 256 channels actually. 0-127 corresponds to lower row and then 128-255 corresponds to the
478 // upper. This can be simulated in the code by using the eta module index. Even index are not offset, while odd index, the
479 // strip id is offset by 128
480 // One point to not is that for barrel, the eta module index start at 1, and not zero. Hence a shift of 1 is needed
481 int offset = m_sctId->eta_module(rdoId) % 2;
482 if(m_sctId->barrel_ec(rdoId) == 0) offset = (std::abs(m_sctId->eta_module(rdoId)) - 1) % 2;
483
484 ITkStripID += offset * MaxChannelinStripRow;
485
486 tmpSGhit.setisValidForITkHit(true);
487 tmpSGhit.setStripRowIDForITk(ITkStripID);
488 tmpSGhit.setStripChipIDForITk(chipID);
489 tmpSGhit.setStripHitMapForITk(stripEncodingForITK.at(stripID));
490 }
491
492 tmpSGhit.setBarcodePt(static_cast<unsigned long>(std::ceil(bestParent ? bestParent->momentum().perp() : 0.)));
493 tmpSGhit.setParentageMask(parentMask.to_ulong());
494 tmpSGhit.setX(0.5 * (endsOfStrip.first.x() + endsOfStrip.second.x()));
495 tmpSGhit.setY(0.5 * (endsOfStrip.first.y() + endsOfStrip.second.y()));
496 tmpSGhit.setZ(0.5 * (endsOfStrip.first.z() + endsOfStrip.second.z()));
497
498 if (m_doMultiTruth) {
499 // Add truth
500 FPGATrackSimMultiTruth mt;
501 FPGATrackSimMultiTruth::Barcode uniqueID(tmpSGhit.getEventIndex(), tmpSGhit.getBarcode()); // FIXME barcode-based
502 mt.maximize(uniqueID, tmpSGhit.getBarcodePt()); // FIMXE barcode-based
503 tmpSGhit.setTruth(mt);
504 }
505
506 eventHeader->addHit(tmpSGhit);
507 } // end for each RDO in the strip collection
508 } // end for each strip RDO collection
509 // dump all RDO's and SDO's for a given event, for debugging purposes
510
511 return StatusCode::SUCCESS;
512}
void setStripChipIDForITk(int v)
void setStripHitMapForITk(int v)
void setStripRowIDForITk(int v)
std::pair< Amg::Vector3D, Amg::Vector3D > endsOfStrip(const Amg::Vector2D &position) const
Special method for SCT to retrieve the two ends of a "strip" Returned coordinates are in global frame...
virtual Identifier identify() const override final

◆ readTruthTracks()

StatusCode FPGATrackSimSGToRawHitsTool::readTruthTracks ( std::vector< FPGATrackSimTruthTrack > & truth,
const EventContext & eventContext ) const
private

Definition at line 735 of file FPGATrackSimSGToRawHitsTool.cxx.

736{
737 auto simTracksHandle = SG::makeHandle(m_mcCollectionKey, eventContext);
738 ATH_MSG_DEBUG("Dump truth tracks, size " << simTracksHandle->size());
739
740 // dump each truth track
741 for (unsigned int ievt = 0; ievt < simTracksHandle->size(); ++ievt) {
742 const HepMC::GenEvent* genEvent = simTracksHandle->at(ievt);
743 // retrieve the primary interaction vertex here. for now, use the dummy origin.
744 HepGeom::Point3D<double> primaryVtx(0., 0., 0.);
745 // the event should have signal process vertex unless it was generated as single particles.
746 // if it exists, use it for the primary vertex.
748 if (spv) {
749 primaryVtx.set(spv->position().x(),
750 spv->position().y(),
751 spv->position().z());
752 ATH_MSG_DEBUG("using signal process vertex for eventIndex " << ievt << ":"
753 << primaryVtx.x() << "\t" << primaryVtx.y() << "\t" << primaryVtx.z());
754 }
755 for (const auto& particle: *genEvent) {
756 const int pdgcode = particle->pdg_id();
757 // reject generated particles without a production vertex.
758 if (particle->production_vertex() == nullptr) {
759 continue;
760 }
761 // reject neutral or unstable particles
762 const HepPDT::ParticleData* pd = m_particleDataTable->particle(abs(pdgcode));
763 if (pd == nullptr) {
764 continue;
765 }
766 float charge = pd->charge();
767 if (pdgcode < 0) charge *= -1.; // since we took absolute value above
768 if (std::abs(charge) < 0.5) {
769 continue;
770 }
771 if (!MC::isStable(particle)) {
772 continue;
773 }
774 // truth-to-track tool
775 const Amg::Vector3D momentum(particle->momentum().px(), particle->momentum().py(), particle->momentum().pz());
776 const Amg::Vector3D position(particle->production_vertex()->position().x(), particle->production_vertex()->position().y(), particle->production_vertex()->position().z());
777 const Trk::CurvilinearParameters cParameters(position, momentum, charge);
778 Trk::PerigeeSurface persf;
779 if (m_UseNominalOrigin) {
780 Amg::Vector3D origin(0, 0, 0);
781 persf = Trk::PerigeeSurface(origin);
782 }
783 else {
784 SG::ReadCondHandle<InDet::BeamSpotData> beamSpotHandle{ m_beamSpotKey, eventContext };
785 Trk::PerigeeSurface persf(beamSpotHandle->beamPos());
786 }
787 const std::unique_ptr<Trk::TrackParameters> tP = m_extrapolator->extrapolate(eventContext, cParameters, persf, Trk::anyDirection, false);
788 const double track_truth_d0 = tP ? tP->parameters()[Trk::d0] : 999.;
789 const double track_truth_phi = tP ? tP->parameters()[Trk::phi] : 999.;
790 const double track_truth_p = (tP && fabs(tP->parameters()[Trk::qOverP]) > 1.e-8) ?
791 tP->charge() / tP->parameters()[Trk::qOverP] : 10E7;
792 const double track_truth_x0 = tP ? tP->position().x() : 999.;
793 const double track_truth_y0 = tP ? tP->position().y() : 999.;
794 const double track_truth_z0 = tP ? tP->parameters()[Trk::z0] : 999.;
795 const double track_truth_q = tP ? tP->charge() : 0.;
796 const double track_truth_sinphi = tP ? std::sin(tP->parameters()[Trk::phi]) : -1.;
797 const double track_truth_cosphi = tP ? std::cos(tP->parameters()[Trk::phi]) : -1.;
798 const double track_truth_sintheta = tP ? std::sin(tP->parameters()[Trk::theta]) : -1.;
799 const double track_truth_costheta = tP ? std::cos(tP->parameters()[Trk::theta]) : -1.;
800 double truth_d0corr = track_truth_d0 - (primaryVtx.y() * cos(track_truth_phi) - primaryVtx.x() * sin(track_truth_phi));
801 double truth_zvertex = 0.;
802 const HepGeom::Point3D<double> startVertex(particle->production_vertex()->position().x(), particle->production_vertex()->position().y(), particle->production_vertex()->position().z());
803 // categorize particle (prompt, secondary, etc.) based on InDetPerformanceRTT/detector paper criteria.
804 bool isPrimary = true;
805 if (std::abs(truth_d0corr) > 2.) { isPrimary = false; }
806 const int bc = HepMC::barcode(particle); // FIXME update barcode-based syntax
807 const int uid = HepMC::uniqueID(particle);
808 if (HepMC::is_simulation_particle(particle) || bc == 0) { isPrimary = false; } // FIXME update barcode-based syntax
809 if (isPrimary && particle->production_vertex()) {
810 const HepGeom::Point3D<double> startVertex(particle->production_vertex()->position().x(), particle->production_vertex()->position().y(), particle->production_vertex()->position().z());
811 if (std::abs(startVertex.z() - truth_zvertex) > 100.) { isPrimary = false; }
812 if (particle->end_vertex()) {
813 HepGeom::Point3D<double> endVertex(particle->end_vertex()->position().x(), particle->end_vertex()->position().y(), particle->end_vertex()->position().z());
814 if (endVertex.perp() < FPGATrackSim_PT_TRUTHMIN && std::abs(endVertex.z()) < FPGATrackSim_Z_TRUTHMIN) { isPrimary = false; }
815 }
816 }
817 else {
818 isPrimary = false;
819 }
820
821 HepMcParticleLink truthLink2(uid, ievt, HepMcParticleLink::IS_POSITION, HepMcParticleLink::IS_ID);
822
823 FPGATrackSimTruthTrack tmpSGTrack;
824 tmpSGTrack.setVtxX(track_truth_x0);
825 tmpSGTrack.setVtxY(track_truth_y0);
826 tmpSGTrack.setVtxZ(track_truth_z0);
827 tmpSGTrack.setD0(track_truth_d0);
828 tmpSGTrack.setZ0(track_truth_z0);
829 tmpSGTrack.setVtxZ(primaryVtx.z());
830 tmpSGTrack.setQ(track_truth_q);
831 tmpSGTrack.setPX(track_truth_p * (track_truth_cosphi * track_truth_sintheta));
832 tmpSGTrack.setPY(track_truth_p * (track_truth_sinphi * track_truth_sintheta));
833 tmpSGTrack.setPZ(track_truth_p * track_truth_costheta);
834 tmpSGTrack.setPDGCode(pdgcode);
835 tmpSGTrack.setStatus(particle->status());
836
837 tmpSGTrack.setBarcode(truthLink2.barcode());
838 tmpSGTrack.setUniqueID(truthLink2.id());
839 tmpSGTrack.setEventIndex(truthLink2.eventIndex());
840
841 truth.push_back(tmpSGTrack);
842 } // end for each GenParticle in this GenEvent
843 } // end for each GenEvent
844
845
846 return StatusCode::SUCCESS;
847}
double charge(const T &p)
Definition AtlasPID.h:997
Gaudi::Property< bool > m_UseNominalOrigin
void setUniqueID(const HepMcParticleLink::barcode_type &v)
void setBarcode(const HepMcParticleLink::barcode_type &v)
int barcode(const T *p)
Definition Barcode.h:16
GenVertex * signal_process_vertex(const GenEvent *e)
Definition GenEvent.h:625
const HepMC::GenVertex * ConstGenVertexPtr
Definition GenVertex.h:60
constexpr ParticleHypothesis particle[PARTICLEHYPOTHESES]
the array of masses
@ anyDirection
CurvilinearParametersT< TrackParametersDim, Charged, PlaneSurface > CurvilinearParameters
@ theta
Definition ParamDefs.h:66
@ qOverP
perigee
Definition ParamDefs.h:67
@ phi
Definition ParamDefs.h:75
@ d0
Definition ParamDefs.h:63
@ z0
Definition ParamDefs.h:64
@ isPrimary
true if matched track has a hit in first or second pixel layer

Member Data Documentation

◆ m_beamSpotKey

SG::ReadCondHandleKey<InDet::BeamSpotData> FPGATrackSimSGToRawHitsTool::m_beamSpotKey { this, "BeamSpotKey", "BeamSpotData", "SG key for beam spot" }
private

Definition at line 60 of file FPGATrackSimSGToRawHitsTool.h.

60{ this, "BeamSpotKey", "BeamSpotData", "SG key for beam spot" };

◆ m_doMultiTruth

Gaudi::Property<bool> FPGATrackSimSGToRawHitsTool::m_doMultiTruth { this, "doMultiTruth", true }
private

Definition at line 85 of file FPGATrackSimSGToRawHitsTool.h.

85{ this, "doMultiTruth", true };

◆ m_dumpHitsOnTracks

Gaudi::Property<bool> FPGATrackSimSGToRawHitsTool::m_dumpHitsOnTracks { this, "dumpHitsOnTracks", false }
private

Definition at line 77 of file FPGATrackSimSGToRawHitsTool.h.

77{ this, "dumpHitsOnTracks", false };

◆ m_dumpTruthIntersections

Gaudi::Property<bool> FPGATrackSimSGToRawHitsTool::m_dumpTruthIntersections { this, "dumpTruthIntersections", false }
private

Definition at line 78 of file FPGATrackSimSGToRawHitsTool.h.

78{ this, "dumpTruthIntersections", false };

◆ m_eventInfoKey

SG::ReadHandleKey<xAOD::EventInfo> FPGATrackSimSGToRawHitsTool::m_eventInfoKey { this, "EventInfo", "EventInfo" }
private

Definition at line 61 of file FPGATrackSimSGToRawHitsTool.h.

61{ this, "EventInfo", "EventInfo" };

◆ m_extrapolator

ToolHandle<Trk::IExtrapolator> FPGATrackSimSGToRawHitsTool::m_extrapolator {this, "Extrapolator", "Trk::Extrapolator/AtlasExtrapolator"}
private

ToolHandle for Extrapolator.

Definition at line 58 of file FPGATrackSimSGToRawHitsTool.h.

58{this, "Extrapolator", "Trk::Extrapolator/AtlasExtrapolator"};

◆ m_maxEta

Gaudi::Property<double> FPGATrackSimSGToRawHitsTool::m_maxEta { this, "maxEta", 5.0 }
private

Definition at line 83 of file FPGATrackSimSGToRawHitsTool.h.

83{ this, "maxEta", 5.0 };

◆ m_mcCollectionKey

SG::ReadHandleKey<McEventCollection> FPGATrackSimSGToRawHitsTool::m_mcCollectionKey { this, "McTruth", "TruthEvent" }
private

Definition at line 66 of file FPGATrackSimSGToRawHitsTool.h.

66{ this, "McTruth", "TruthEvent" };

◆ m_minPt

Gaudi::Property<double> FPGATrackSimSGToRawHitsTool::m_minPt { this, "minPt", .8*CLHEP::GeV }
private

Definition at line 84 of file FPGATrackSimSGToRawHitsTool.h.

84{ this, "minPt", .8*CLHEP::GeV };

◆ m_offlineTracksKey

SG::ReadHandleKey<xAOD::TrackParticleContainer> FPGATrackSimSGToRawHitsTool::m_offlineTracksKey { this, "OfflineTracks", "InDetTrackParticles"}
private

Definition at line 65 of file FPGATrackSimSGToRawHitsTool.h.

65{ this, "OfflineTracks", "InDetTrackParticles"};

◆ m_particleDataTable

const HepPDT::ParticleDataTable* FPGATrackSimSGToRawHitsTool::m_particleDataTable = nullptr
private

Definition at line 93 of file FPGATrackSimSGToRawHitsTool.h.

◆ m_PIX_mgr

const InDetDD::SiDetectorManager* FPGATrackSimSGToRawHitsTool::m_PIX_mgr = nullptr
private

Definition at line 91 of file FPGATrackSimSGToRawHitsTool.h.

◆ m_pixelClusterContainerKey

SG::ReadHandleKey<InDet::SiClusterContainer> FPGATrackSimSGToRawHitsTool::m_pixelClusterContainerKey { this, "pixelClustersName", "ITkPixelClusters" }
private

Definition at line 62 of file FPGATrackSimSGToRawHitsTool.h.

62{ this, "pixelClustersName", "ITkPixelClusters" };

◆ m_pixelId

const PixelID* FPGATrackSimSGToRawHitsTool::m_pixelId = nullptr
private

Definition at line 88 of file FPGATrackSimSGToRawHitsTool.h.

◆ m_pixelRDOKey

SG::ReadHandleKey<PixelRDO_Container> FPGATrackSimSGToRawHitsTool::m_pixelRDOKey { this, "PixelRDO", "ITkPixelRDOs" }
private

Definition at line 69 of file FPGATrackSimSGToRawHitsTool.h.

69{ this, "PixelRDO", "ITkPixelRDOs" };

◆ m_pixelSDOKey

SG::ReadHandleKey<InDetSimDataCollection> FPGATrackSimSGToRawHitsTool::m_pixelSDOKey { this, "PixelSDO", "ITkPixelSDO_Map" }
private

Definition at line 67 of file FPGATrackSimSGToRawHitsTool.h.

67{ this, "PixelSDO", "ITkPixelSDO_Map" };

◆ m_readOfflineClusters

Gaudi::Property<bool> FPGATrackSimSGToRawHitsTool::m_readOfflineClusters { this, "ReadOfflineClusters", true, "flag to enable the offline cluster save" }
private

Definition at line 79 of file FPGATrackSimSGToRawHitsTool.h.

79{ this, "ReadOfflineClusters", true, "flag to enable the offline cluster save" };

◆ m_readOfflineTracks

Gaudi::Property<bool> FPGATrackSimSGToRawHitsTool::m_readOfflineTracks { this, "ReadOfflineTracks", true, "flag to enable the offline tracking save" }
private

Definition at line 81 of file FPGATrackSimSGToRawHitsTool.h.

81{ this, "ReadOfflineTracks", true, "flag to enable the offline tracking save" };

◆ m_readTruthTracks

Gaudi::Property<bool> FPGATrackSimSGToRawHitsTool::m_readTruthTracks { this, "ReadTruthTracks", true, "flag to enable the truth tracking save" }
private

Definition at line 80 of file FPGATrackSimSGToRawHitsTool.h.

80{ this, "ReadTruthTracks", true, "flag to enable the truth tracking save" };

◆ m_SCT_mgr

const InDetDD::SiDetectorManager* FPGATrackSimSGToRawHitsTool::m_SCT_mgr = nullptr
private

Definition at line 92 of file FPGATrackSimSGToRawHitsTool.h.

◆ m_sctClusterContainerKey

SG::ReadHandleKey<InDet::SiClusterContainer> FPGATrackSimSGToRawHitsTool::m_sctClusterContainerKey { this, "SCT_ClustersName", "ITkStripClusters" }
private

Definition at line 63 of file FPGATrackSimSGToRawHitsTool.h.

63{ this, "SCT_ClustersName", "ITkStripClusters" };

◆ m_sctId

const SCT_ID* FPGATrackSimSGToRawHitsTool::m_sctId = nullptr
private

Definition at line 89 of file FPGATrackSimSGToRawHitsTool.h.

◆ m_stripRDOKey

SG::ReadHandleKey<SCT_RDO_Container> FPGATrackSimSGToRawHitsTool::m_stripRDOKey { this, "StripRDO", "ITkStripRDOs" }
private

Definition at line 70 of file FPGATrackSimSGToRawHitsTool.h.

70{ this, "StripRDO", "ITkStripRDOs" };

◆ m_stripSDOKey

SG::ReadHandleKey<InDetSimDataCollection> FPGATrackSimSGToRawHitsTool::m_stripSDOKey { this, "StripSDO", "ITkStripSDO_Map" }
private

Definition at line 68 of file FPGATrackSimSGToRawHitsTool.h.

68{ this, "StripSDO", "ITkStripSDO_Map" };

◆ m_tracksTruthName

Gaudi::Property<std::string> FPGATrackSimSGToRawHitsTool::m_tracksTruthName { this, "OfflineName", "InDetTrackParticles", "name of offline tracks collection" }
private

Definition at line 76 of file FPGATrackSimSGToRawHitsTool.h.

76{ this, "OfflineName", "InDetTrackParticles", "name of offline tracks collection" };

◆ m_truthToTrack

ToolHandle<Trk::ITruthToTrack> FPGATrackSimSGToRawHitsTool::m_truthToTrack {this, "TruthToTrackTool", "Trk::TruthToTrack/InDetTruthToTrack" }
private

tool to create track parameters from a gen particle

Definition at line 57 of file FPGATrackSimSGToRawHitsTool.h.

57{this, "TruthToTrackTool", "Trk::TruthToTrack/InDetTruthToTrack" };

◆ m_UseNominalOrigin

Gaudi::Property<bool> FPGATrackSimSGToRawHitsTool::m_UseNominalOrigin { this, "UseNominalOrigin", false, "if true truth values are always with respect to (0,0,0)" }
private

Definition at line 82 of file FPGATrackSimSGToRawHitsTool.h.

82{ this, "UseNominalOrigin", false, "if true truth values are always with respect to (0,0,0)" };

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