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

#include <TrackFastSimSD.h>

Inheritance diagram for TrackFastSimSD:
Collaboration diagram for TrackFastSimSD:

Public Member Functions

 TrackFastSimSD (const std::string &name, const std::string &outputCollectionName, const int SD_type=0)
 ~TrackFastSimSD ()
void Initialize (G4HCofThisEvent *) override final
G4bool ProcessHits (G4Step *, G4TouchableHistory *) override final
void WriteTrack (const G4Track *, const bool, const bool)
template<class... Args>
void AddHit (Args &&... args)
 Templated method to stuff a single hit into the sensitive detector class.

Private Attributes

SG::WriteHandle< TrackRecordCollectionm_trackRecordCollection
int m_SD_type

Detailed Description

Definition at line 23 of file TrackFastSimSD.h.

Constructor & Destructor Documentation

◆ TrackFastSimSD()

TrackFastSimSD::TrackFastSimSD ( const std::string & name,
const std::string & outputCollectionName,
const int SD_type = 0 )

Definition at line 22 of file TrackFastSimSD.cxx.

23 : G4VSensitiveDetector( name )
24 , m_trackRecordCollection( outputCollectionName )
25 , m_SD_type( SD_type )
26{
27}
SG::WriteHandle< TrackRecordCollection > m_trackRecordCollection

◆ ~TrackFastSimSD()

TrackFastSimSD::~TrackFastSimSD ( )
inline

Definition at line 27 of file TrackFastSimSD.h.

27{}

Member Function Documentation

◆ AddHit()

template<class... Args>
void TrackFastSimSD::AddHit ( Args &&... args)
inline

Templated method to stuff a single hit into the sensitive detector class.

This could get rather tricky, but the idea is to allow fast simulations to use the very same SD classes as the standard simulation.

Definition at line 38 of file TrackFastSimSD.h.

◆ Initialize()

void TrackFastSimSD::Initialize ( G4HCofThisEvent * )
finaloverride

Definition at line 30 of file TrackFastSimSD.cxx.

31{
32 if (!m_trackRecordCollection.isValid()) m_trackRecordCollection = std::make_unique<TrackRecordCollection>(m_trackRecordCollection.name());
33}

◆ ProcessHits()

G4bool TrackFastSimSD::ProcessHits ( G4Step * aStep,
G4TouchableHistory *  )
finaloverride

Definition at line 35 of file TrackFastSimSD.cxx.

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 = static_cast<const 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
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}
time(flags, cells_name, *args, **kw)
int barcode(const T *p)
Definition Barcode.h:16
status
Definition merge.py:16

◆ WriteTrack()

void TrackFastSimSD::WriteTrack ( const G4Track * track,
const bool originPos,
const bool originMom )

Definition at line 97 of file TrackFastSimSD.cxx.

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
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}

Member Data Documentation

◆ m_SD_type

int TrackFastSimSD::m_SD_type
private

Definition at line 44 of file TrackFastSimSD.h.

◆ m_trackRecordCollection

SG::WriteHandle<TrackRecordCollection> TrackFastSimSD::m_trackRecordCollection
private

Definition at line 42 of file TrackFastSimSD.h.


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