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>
bool AddHit (Args &&... args)
 Templated method to stuff a single hit into the sensitive detector class.

Private Member Functions

TrackRecordCollectiongetTrackRecordCollection () const

Private Attributes

std::string m_outputCollectionName
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 26 of file TrackFastSimSD.cxx.

27 : G4VSensitiveDetector( name )
28 , m_outputCollectionName( outputCollectionName )
29 , m_SD_type( SD_type )
30{
31}
std::string m_outputCollectionName

◆ ~TrackFastSimSD()

TrackFastSimSD::~TrackFastSimSD ( )
inline

Definition at line 27 of file TrackFastSimSD.h.

27{}

Member Function Documentation

◆ AddHit()

template<class... Args>
bool 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.

39 {
41 m_trackRecordCollection->Emplace(std::forward<Args>(args)...);
42 return true;
43 }
44 return false;
45 }
TrackRecordCollection * m_trackRecordCollection

◆ getTrackRecordCollection()

TrackRecordCollection * TrackFastSimSD::getTrackRecordCollection ( ) const
private

Definition at line 132 of file TrackFastSimSD.cxx.

133{
134 auto* eventManager = G4EventManager::GetEventManager();
135 if (!eventManager) {
136 return nullptr;
137 }
138
139 auto* eventInfo =
140 dynamic_cast<AtlasG4EventUserInfo*>(eventManager->GetUserInformation());
141 if (!eventInfo) {
142 return nullptr;
143 }
144
145 std::shared_ptr<HitCollectionMap> hitCollections = eventInfo->GetHitCollectionMap();
146 return hitCollections ? hitCollections->Find<TrackRecordCollection>(m_outputCollectionName)
147 : nullptr;
148}
AtlasHitsVector< TrackRecord > TrackRecordCollection

◆ Initialize()

void TrackFastSimSD::Initialize ( G4HCofThisEvent * )
finaloverride

Definition at line 34 of file TrackFastSimSD.cxx.

35{
37}
TrackRecordCollection * getTrackRecordCollection() const

◆ ProcessHits()

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

Definition at line 39 of file TrackFastSimSD.cxx.

40{
41 if (!m_SD_type) { return false; }
42 const G4StepPoint *preStep=aStep->GetPreStepPoint();
43 const G4VPhysicalVolume *preVol=preStep->GetPhysicalVolume();
44 const G4StepPoint *postStep=aStep->GetPostStepPoint();
45 const G4VPhysicalVolume *postVol=postStep->GetPhysicalVolume();
46
47 if (preVol==postVol) { return false; }
48
49 const G4Track *track(aStep->GetTrack());
50
51 // PDG Code
52 int pdgcode = (track->GetDefinition())?track->GetDefinition()->GetPDGEncoding():0;
53 if (track->GetDefinition() == G4Geantino::Definition() ) pdgcode=999;
54 if (track->GetDefinition() == G4ChargedGeantino::Definition() ) pdgcode=998;
55
56 // Position and Momentum
57 G4ThreeVector pos=postStep->GetPosition();
58 G4ThreeVector mom=postStep->GetMomentum();
59
60 if(m_SD_type==2)
61 {
62 // need to get the local position
63 const G4TouchableHistory* touchHist = static_cast<const G4TouchableHistory*>(preStep->GetTouchable());
64 const G4AffineTransform trans = track->GetTouchable()->GetHistory()->GetTopTransform(); // trans from global to local
65 G4ThreeVector localPos=trans.TransformPoint(pos);
66 G4ThreeVector localMom=mom;
67 trans.ApplyAxisTransform(localMom); // rotate momentum as needed
68 const G4VSolid *shape= touchHist->GetSolid();
69 const G4ThreeVector normal=shape->SurfaceNormal(localPos);
70 // check if particle is entering. if it is extiting don't record it.
71 if(normal.dot(localPos)>=0.) return false;
72 }
73
74 // Energy
75 const double ener=postStep->GetTotalEnergy();
76
77 // Time
78 const double time=postStep->GetGlobalTime();
79
80 // Barcode
81 TrackHelper trHelp(track);
82 const int barcode = trHelp.GetBarcode(); // FIXME barcode based
83 const int id = trHelp.GetUniqueID();
84 const int status = trHelp.GetStatus();
85
86 //create the TimedTrackRecord
87 return AddHit(pdgcode,
88 status,
89 ener,
90 mom,
91 pos,
92 time,
93 barcode, // FIXME barcode based
94 id,
95 preVol->GetName());
96}
bool AddHit(Args &&... args)
Templated method to stuff a single hit into the sensitive detector class.
time(flags, cells_name, *args, **kw)
const std::string barcode
status
Definition merge.py:16

◆ WriteTrack()

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

Definition at line 98 of file TrackFastSimSD.cxx.

99{
100 if (!track) { G4cout << "ERROR: the track pointer was zero" << G4endl; return; }
101
102 G4VPhysicalVolume *preVol=track->GetVolume();
103
104 const int pdgcode = (track->GetDefinition())?track->GetDefinition()->GetPDGEncoding():0;
105
106 const G4ThreeVector pos = originPos?track->GetVertexPosition():track->GetPosition();
107 const double ener=originMom?(track->GetVertexKineticEnergy()+track->GetDynamicParticle()->GetMass()):track->GetTotalEnergy();
108 G4ThreeVector mom = track->GetMomentum();
109 if (originMom){
110 double mommag = std::sqrt(std::pow(ener,2)-std::pow(track->GetDynamicParticle()->GetMass(),2));
111 mom = track->GetVertexMomentumDirection()*mommag;
112 }
113
114 const double time=track->GetGlobalTime();
115 TrackHelper trHelp(track);
116 const int barcode = trHelp.GetBarcode(); // FIXME barcode based
117 const int id = trHelp.GetUniqueID();
118 const int status = trHelp.GetStatus();
119
120 //create the TimedTrackRecord
121 AddHit(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_outputCollectionName

std::string TrackFastSimSD::m_outputCollectionName
private

Definition at line 50 of file TrackFastSimSD.h.

◆ m_SD_type

int TrackFastSimSD::m_SD_type
private

Definition at line 53 of file TrackFastSimSD.h.

◆ m_trackRecordCollection

TrackRecordCollection* TrackFastSimSD::m_trackRecordCollection {}
private

Definition at line 52 of file TrackFastSimSD.h.

52{};

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