ATLAS Offline Software
Loading...
Searching...
No Matches
RecordingEnvelope.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 "RecordingEnvelope.h"
7
8// Athena headers
11
12//G4 headers
13#include "G4PhysicalVolumeStore.hh"
14#include "G4Step.hh"
15#include "G4TransportationManager.hh"
16#include "G4VPhysicalVolume.hh"
17#include "G4VSolid.hh"
18
19
20RecordingEnvelope::RecordingEnvelope(const std::string& envelopeVolumeName, const std::string& trackRecordCollectionName):
21 m_level(-1),
22 m_logicalVolume(nullptr),
23 m_envelopeVolumeName(envelopeVolumeName),
24 m_trackRecordCollectionName(trackRecordCollectionName)
25{
26
27}
28
33
34bool RecordingEnvelope::checkDaughters(const G4LogicalVolume *possibleParent, const G4VPhysicalVolume *thePhysicalVolume, int& level) const
35{
36 if (possibleParent->IsDaughter(thePhysicalVolume))
37 {
38 ++level;
39 return true;
40 }
41 // Otherwise Loop over the daughters
42 const G4int nDaughters(possibleParent->GetNoDaughters());
43 for(G4int daughter(0); daughter<nDaughters; ++daughter)
44 {
45 const G4VPhysicalVolume* daughterPhysVol = possibleParent->GetDaughter(daughter);
46 const G4LogicalVolume* daughterLogVol = daughterPhysVol->GetLogicalVolume();
47 if (this->checkDaughters(daughterLogVol, thePhysicalVolume,level))
48 {
49 ++level;
50 return true;
51 }
52 }
53 return false;
54}
56{
57 const G4VPhysicalVolume *thePhysicalVolume = G4PhysicalVolumeStore::GetInstance()->GetVolume(m_envelopeVolumeName,false);
58 m_logicalVolume=thePhysicalVolume->GetLogicalVolume();
59 const G4LogicalVolume * logicalWorld = G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume()->GetLogicalVolume();
60 m_level=0;
61 return this->checkDaughters(logicalWorld, thePhysicalVolume,m_level);
62}
63
65{
66 m_trackRecordCollection = trackRecordCollection;
67}
68void RecordingEnvelope::AddTrackRecord(const G4Step* aStep)
69{
70 const std::string pname = aStep->GetTrack()->GetDefinition()->GetParticleName();
71 const int pdgcode = (pname=="geantino") ? 999 : aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
72 TrackHelper trHelp(aStep->GetTrack());
73 const int barcode = trHelp.GetBarcode();
74 const int status = trHelp.GetStatus();
75 const int id = trHelp.GetUniqueID();
76
77 G4StepPoint *postStep=aStep->GetPostStepPoint();
78 G4ThreeVector pos=postStep->GetPosition();
79 G4ThreeVector mom=postStep->GetMomentum();
80 const double ener=postStep->GetTotalEnergy();
81 const double time=postStep->GetGlobalTime();
82
83 G4StepPoint *preStep=aStep->GetPreStepPoint();
84 G4VPhysicalVolume *preVol=preStep->GetPhysicalVolume();
85
86 m_trackRecordCollection->Emplace(pdgcode,status,ener,mom,pos,time,barcode,id,preVol->GetName());
87}
AtlasHitsVector< TrackRecord > TrackRecordCollection
G4LogicalVolume * m_logicalVolume
Pointer to the G4LogicalVolume used by this recording envelope.
int m_level
The number of levels beneath the world that the G4LogicalVolume called m_envelopeVolumeName lies.
const std::string m_envelopeVolumeName
Name of the G4LogicalVolume used by this recording envelope.
TrackRecordCollection * m_trackRecordCollection
TrackRecordCollection used by this envelope.
std::string m_trackRecordCollectionName
void AddTrackRecord(const G4Step *aStep)
Add a TrackRecord to the TrackRecordCollection owned by this recording envelope based on the informat...
bool checkDaughters(const G4LogicalVolume *possibleParent, const G4VPhysicalVolume *thePhysicalVolume, int &level) const
Recursively called method used to hunt for the G4LogicalVolume associated with this recording envelop...
RecordingEnvelope(const std::string &envelopeVolumeName, const std::string &trackRecordCollectionName)
Constructor.
bool Initialize()
Finds the pointer to the G4LogicalVolume called m_envelopeVolumeName and the number of levels beneath...
~RecordingEnvelope()
Destructor.
void BeginOfEvent(TrackRecordCollection *)
Called at the start of each G4 event.
int GetBarcode() const
int GetUniqueID() const
int GetStatus() const