ATLAS Offline Software
TrackFastSimSD.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
7 
8 // Athena headers
9 #include "MCTruth/TrackHelper.h"
10 
11 // Geant4 headers
12 #include "G4ChargedGeantino.hh"
13 #include "G4DynamicParticle.hh"
14 #include "G4Geantino.hh"
15 #include "G4Track.hh"
16 #include "G4VPhysicalVolume.hh"
17 #include "G4VSolid.hh"
18 
19 // STL headers
20 #include <cmath>
21 
22 TrackFastSimSD::TrackFastSimSD(const std::string& name, const std::string& outputCollectionName, const int SD_type)
23  : G4VSensitiveDetector( name )
24  , m_trackRecordCollection( outputCollectionName )
25  , m_SD_type( SD_type )
26 {
27 }
28 
29 // Initialize from G4 - necessary to new the write handle for now
30 void TrackFastSimSD::Initialize(G4HCofThisEvent *)
31 {
32  if (!m_trackRecordCollection.isValid()) m_trackRecordCollection = std::make_unique<TrackRecordCollection>(m_trackRecordCollection.name());
33 }
34 
35 G4bool TrackFastSimSD::ProcessHits(G4Step* aStep,G4TouchableHistory* )
36 {
37  if (!m_SD_type) { return false; }
38  const G4StepPoint *preStep=aStep->GetPreStepPoint();
39  const G4VPhysicalVolume *preVol=preStep->GetPhysicalVolume();
40  const G4StepPoint *postStep=aStep->GetPostStepPoint();
41  const G4VPhysicalVolume *postVol=postStep->GetPhysicalVolume();
42 
43  if (preVol==postVol) { return false; }
44 
45  const G4Track *track(aStep->GetTrack());
46 
47  // PDG Code
48  int pdgcode = (track->GetDefinition())?track->GetDefinition()->GetPDGEncoding():0;
49  if (track->GetDefinition() == G4Geantino::Definition() ) pdgcode=999;
50  if (track->GetDefinition() == G4ChargedGeantino::Definition() ) pdgcode=998;
51 
52  // Position and Momentum
53  G4ThreeVector pos=postStep->GetPosition();
54  G4ThreeVector mom=postStep->GetMomentum();
55 
56  if(m_SD_type==2)
57  {
58  // need to get the local position
59  const G4TouchableHistory* touchHist = (G4TouchableHistory*)preStep->GetTouchable();
60  const G4AffineTransform trans = track->GetTouchable()->GetHistory()->GetTopTransform(); // trans from global to local
61  G4ThreeVector localPos=trans.TransformPoint(pos);
62  G4ThreeVector localMom=mom;
63  trans.ApplyAxisTransform(localMom); // rotate momentum as needed
64  const G4VSolid *shape= touchHist->GetSolid();
65  const G4ThreeVector normal=shape->SurfaceNormal(localPos);
66  // check if particle is entering. if it is extiting don't record it.
67  if(normal.dot(localPos)>=0.) return false;
68  }
69 
70  // Energy
71  const double ener=postStep->GetTotalEnergy();
72 
73  // Time
74  const double time=postStep->GetGlobalTime();
75 
76  // Barcode
77  TrackHelper trHelp(track);
78  const int barcode = trHelp.GetBarcode(); // FIXME barcode based
79  const int id = trHelp.GetUniqueID();
80  const int status = trHelp.GetStatus();
81 
82  //create the TimedTrackRecord
83  m_trackRecordCollection->Emplace(
84  pdgcode,
85  status,
86  ener,
87  mom,
88  pos,
89  time,
90  barcode, // FIXME barcode based
91  id,
92  preVol->GetName());
93 
94  return true;
95 }
96 
97 void TrackFastSimSD::WriteTrack(const G4Track* track, const bool originPos, const bool originMom)
98 {
99  if (!track) { G4cout << "ERROR: the track pointer was zero" << G4endl; return; }
100 
101  G4VPhysicalVolume *preVol=track->GetVolume();
102 
103  const int pdgcode = (track->GetDefinition())?track->GetDefinition()->GetPDGEncoding():0;
104 
105  const G4ThreeVector pos = originPos?track->GetVertexPosition():track->GetPosition();
106  const double ener=originMom?(track->GetVertexKineticEnergy()+track->GetDynamicParticle()->GetMass()):track->GetTotalEnergy();
107  G4ThreeVector mom = track->GetMomentum();
108  if (originMom){
109  double mommag = std::sqrt(std::pow(ener,2)-std::pow(track->GetDynamicParticle()->GetMass(),2));
110  mom = track->GetVertexMomentumDirection()*mommag;
111  }
112 
113  const double time=track->GetGlobalTime();
114  TrackHelper trHelp(track);
115  const int barcode = trHelp.GetBarcode(); // FIXME barcode based
116  const int id = trHelp.GetUniqueID();
117  const int status = trHelp.GetStatus();
118 
119  //create the TimedTrackRecord
120  m_trackRecordCollection->Emplace(
121  pdgcode,
122  status,
123  ener,
124  mom,
125  pos,
126  time,
127  barcode, // FIXME barcode based
128  id,
129  preVol?preVol->GetName():"Unknown");
130 }
131 
TrackFastSimSD::ProcessHits
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override final
Definition: TrackFastSimSD.cxx:35
TrackFastSimSD::m_SD_type
int m_SD_type
Definition: TrackFastSimSD.h:44
TrackFastSimSD.h
TrackHelper::GetBarcode
int GetBarcode() const
Definition: TrackHelper.cxx:35
TrackHelper.h
TrackFastSimSD::TrackFastSimSD
TrackFastSimSD(const std::string &name, const std::string &outputCollectionName, const int SD_type=0)
Definition: TrackFastSimSD.cxx:22
TrackHelper
Definition: TrackHelper.h:14
ParticleGun_EoverP_Config.mom
mom
Definition: ParticleGun_EoverP_Config.py:63
TrackFastSimSD::m_trackRecordCollection
SG::WriteHandle< TrackRecordCollection > m_trackRecordCollection
Definition: TrackFastSimSD.h:42
HepMC::barcode
int barcode(const T *p)
Definition: Barcode.h:16
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:228
TrackFastSimSD::Initialize
void Initialize(G4HCofThisEvent *) override final
Definition: TrackFastSimSD.cxx:30
TrackHelper::GetStatus
int GetStatus() const
Definition: TrackHelper.cxx:47
python.LumiBlobConversion.pos
pos
Definition: LumiBlobConversion.py:18
CaloSwCorrections.time
def time(flags, cells_name, *args, **kw)
Definition: CaloSwCorrections.py:242
TrackFastSimSD::WriteTrack
void WriteTrack(const G4Track *, const bool, const bool)
Definition: TrackFastSimSD.cxx:97
merge.status
status
Definition: merge.py:17
xAOD::track
@ track
Definition: TrackingPrimitives.h:512
TrackHelper::GetUniqueID
int GetUniqueID() const
Definition: TrackHelper.cxx:41
pow
constexpr int pow(int base, int exp) noexcept
Definition: ap_fixedTest.cxx:15