ATLAS Offline Software
Loading...
Searching...
No Matches
ISFTrajectory.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5// class header
6#include "ISFTrajectory.h"
7
8// ISF includes
12
13// ISF Geant4 includes
16
17// Athena includes
18//#include "FadsActions/TrackingAction.h"
19
20// MCTruth includes
23#include "MCTruth/TrackHelper.h"
24
25//HepMC includes
27
28#undef _ISFTRAJECTORY_DEBUG_
29
31 : G4Trajectory()
32 , m_truthRecordSvcQuick(nullptr)
33{
34}
35
37 ISF::ITruthSvc* truthSvc)
38 : G4Trajectory(aTrack)
39 , m_truthRecordSvcQuick(truthSvc)
40{
41}
42
47
48void iGeant4::ISFTrajectory::AppendStep(const G4Step* aStep)
49{
50
51 // only use truth service if there are new any secondaries
52 const int numSecondaries = aStep->GetSecondaryInCurrentStep()->size();
53 //const int numSecondaries = aStep->GetNumberOfSecondariesInCurrentStep(); //Once we switch to G4 10.2 or later
54
55 if (numSecondaries) {
56
57 G4Track* track=aStep->GetTrack();
58#ifdef _ISFTRAJECTORY_DEBUG_
59 std::cout << "A new track "
60 << " (trackID " << track->GetTrackID()
61 << "), track pos: "<<track->GetPosition()
62 << ", mom: "<<track->GetMomentum()
63 << ", parentID " << track->GetParentID()
64 << ", numSec="<<numSecondaries<<" is in AppendStep." << std::endl;
65#endif
66 // OK, there was an interaction. look at the track, if it
67 // is not a secondary (i.e. we have a connected tree) we
68 // apply the MC truth machinery...
69 TrackHelper tHelper(track);
70 if (tHelper.IsSecondary()) {
71#ifdef _ISFTRAJECTORY_DEBUG_
72 std::cout<<"is secondary, returning"<<std::endl;
73#endif
74 return;
75 }
76
77 // get base ISFParticle
79 if (!trackInfo) {
80 G4ExceptionDescription description;
81 description << G4String("AppendStep: ") + "No VTrackInformation associated with G4Track (trackID: "
82 << track->GetTrackID() << ", track pos: "<<track->GetPosition() << ", mom: "<<track->GetMomentum()
83 << ", parentID " << track->GetParentID() << ")";
84 G4Exception("iGeant4::ISFTrajectory", "NoVTrackInformation", FatalException, description);
85 return; //The G4Exception call above should abort the job, but Coverity does not seem to pick this up.
86 }
87
88 ISF::ISFParticle* baseIsp = trackInfo->GetBaseISFParticle();
89 if (!baseIsp) {
90 G4ExceptionDescription description;
91 description << G4String("AppendStep: ") + "NULL ISFParticle pointer for current G4Step (trackID "
92 << track->GetTrackID() << ", track pos: "<<track->GetPosition() << ", mom: "<<track->GetMomentum()
93 << ", parentID " << track->GetParentID() << ", numSec="<<numSecondaries << ")";
94 G4Exception("iGeant4::ISFTrajectory", "NoISFParticle", FatalException, description);
95 return; //The G4Exception call above should abort the job, but Coverity does not seem to pick this up.
96 }
97
98 AtlasDetDescr::AtlasRegion geoID = baseIsp->nextGeoID();
99
100 auto* atlasG4EvtUserInfo = ISFG4Helper::getAtlasG4EventUserInfo();
101 if (atlasG4EvtUserInfo->GetCurrentGenParticle() &&
102 atlasG4EvtUserInfo->GetCurrentGenParticle()->end_vertex()) {
103 HepMC::GenParticlePtr currentGenParticle = atlasG4EvtUserInfo->GetCurrentGenParticle();
104 G4ExceptionDescription description;
105 description << G4String("AppendStep: ") + "Currently Traced Particle has an end vertex!\n";
106 description << "G4Track Properties: trackID = " << track->GetTrackID()<< ", Step Number = "<<track->GetCurrentStepNumber() << ", parentID = " << track->GetParentID() << ", TrackStatus = " << track->GetTrackStatus() << "\n";
107 description << "G4Step Properties: number of secondaries in the current step: " << aStep->GetSecondaryInCurrentStep()->size() << "\n";
108 description << "currentGenParticle : " << currentGenParticle << ", barcode: " << HepMC::barcode(currentGenParticle) << "\n";
109 description << "currentGenParticle->end_vertex(): " << currentGenParticle->end_vertex() << ", barcode: " << HepMC::barcode(currentGenParticle->end_vertex()) << "\n";
110 description << "ISFParticle (from TrackInformation): " << *baseIsp;
111 HepMC::GenParticlePtr currentTrackInfoGenParticle = trackInfo->GetCurrentGenParticle();
112 if (currentTrackInfoGenParticle) {
113 description << "currentTrackInfoGenParticle : " << currentTrackInfoGenParticle << ", barcode: " << HepMC::barcode(currentTrackInfoGenParticle) << "\n";
114 if (currentTrackInfoGenParticle->end_vertex()) {
115 description << "currentTrackInfoGenParticle->end_vertex(): " << currentTrackInfoGenParticle->end_vertex() << ", barcode: " << HepMC::barcode(currentTrackInfoGenParticle->end_vertex()) << "\n";
116 }
117 else {
118 description << "currentTrackInfoGenParticle has no end_vertex!\n";
119 }
120 }
121 else {
122 description << " trackInfo->GetHepMCParticle() == nullptr \n";
123 }
124 G4Exception("iGeant4::ISFTrajectory", "EndVertexExists", FatalException, description);
125 }
126 iGeant4::Geant4TruthIncident truth(aStep, *baseIsp, geoID, atlasG4EvtUserInfo);
127
129 m_truthRecordSvcQuick->registerTruthIncident(truth);
130
131 // read the TrackInformation to determine whether the G4Track was
132 // returned to the ISF in this step
133 if ( trackInfo->GetReturnedToISF()==true ) {
134 // make sure that the TruthBinding of the ISFParticle points to the newest
135 // HepMC::GenParticle instance in case it got updated by the
136 // ITruthSvc::registerTruthIncident call above
137 auto currentGenPart = atlasG4EvtUserInfo->GetCurrentGenParticle();
138 baseIsp->getTruthBinding()->setCurrentGenParticle( currentGenPart );
139 int newBarcode = HepMC::barcode(currentGenPart); // FIXME barcode-based
140 baseIsp->setBarcode( newBarcode ); // FIXME barcode-based
141 baseIsp->setStatus( currentGenPart->status() );
142 int id = HepMC::uniqueID(currentGenPart);
143 baseIsp->setId( id );
144 }
145 }
146 else {
147 G4ExceptionDescription description;
148 description << G4String("AppendStep: ") + "m_truthRecordSvcQuick is NULL!";
149 G4Exception("iGeant4::ISFTrajectory", "NoTruthRecordSvc", FatalException, description);
150 }
151 }
152
153 G4Trajectory::AppendStep(aStep);
154}
The generic ISF particle definition,.
Definition ISFParticle.h:42
const TruthBinding * getTruthBinding() const
pointer to the simulation truth - optional, can be 0
void setBarcode(int bc)
set a new barcode
void setId(int id)
set a new unique ID
void setStatus(int a)
AtlasDetDescr::AtlasRegion nextGeoID() const
next geoID the particle will be simulated in
@ class ITruthSvc
Definition ITruthSvc.h:29
void setCurrentGenParticle(HepMC::GenParticlePtr p)
bool IsSecondary() const
ISF_Geant4 specific implementation of the ISF::ITruthIncident.
static VTrackInformation * getISFTrackInfo(const G4Track &aTrack)
return a valid UserInformation object of the G4Track for use within the ISF
static AtlasG4EventUserInfo * getAtlasG4EventUserInfo()
return pointer to current AtlasG4EventUserInfo
void AppendStep(const G4Step *aStep)
ISF::ITruthSvc * m_truthRecordSvcQuick
std::string description
glabal timer - how long have I taken so far?
Definition hcg.cxx:91
AtlasRegion
A simple enum of ATLAS regions and sub-detectors.
Definition AtlasRegion.h:21
int barcode(const T *p)
Definition Barcode.h:16
int uniqueID(const T &p)
GenParticle * GenParticlePtr
Definition GenParticle.h:37