ATLAS Offline Software
Loading...
Searching...
No Matches
ActsTrk::MaterialTrackRecorder Class Reference

Collects the G4 steps and writes out RecordedMaterialTrackCollection to a store gate. More...

#include <MaterialTrackRecorder.h>

Inheritance diagram for ActsTrk::MaterialTrackRecorder:
Collaboration diagram for ActsTrk::MaterialTrackRecorder:

Classes

struct  Config

Public Member Functions

 MaterialTrackRecorder (const Config &config)
virtual void BeginOfEventAction (const G4Event *) override
virtual void EndOfEventAction (const G4Event *) override
virtual void BeginOfRunAction (const G4Run *) override
virtual void UserSteppingAction (const G4Step *) override
bool msgLvl (const MSG::Level lvl) const
 Test the output level.
MsgStream & msg () const
 The standard message stream.
MsgStream & msg (const MSG::Level lvl) const
 The standard message stream.
void setLevel (MSG::Level lvl)
 Change the current logging level.

Private Member Functions

Acts::Vector3 convertPosition (const G4ThreeVector &g4vec)
 Useful fuction for position coversion.
Acts::Vector3 convertDirection (const G4ThreeVector &g4vec)
 Useful fuction for direction coversion.
void initMessaging () const
 Initialize our message level and MessageSvc.

Private Attributes

RecordedMaterialTrackCollectionm_rmtCollection {nullptr}
 Pointer to the collection, non-owning.
Config m_cfg
 The config class.
std::string m_nm
 Message source name.
boost::thread_specific_ptr< MsgStream > m_msg_tls
 MsgStream instance (a std::cout like with print-out levels)
std::atomic< IMessageSvc * > m_imsg { nullptr }
 MessageSvc pointer.
std::atomic< MSG::Level > m_lvl { MSG::NIL }
 Current logging level.
std::atomic_flag m_initialized ATLAS_THREAD_SAFE = ATOMIC_FLAG_INIT
 Messaging initialized (initMessaging)

Static Private Attributes

static constexpr double s_convertLength = Acts::UnitConstants::mm / CLHEP::mm
 Conversion constants.
static constexpr double s_convertDensity

Detailed Description

Collects the G4 steps and writes out RecordedMaterialTrackCollection to a store gate.

Handles the MaterialTrackRecorder G4UA.

It writes out a MaterialTrack which is usually generated from MaterialTrackRecorder G4UA

Definition at line 34 of file MaterialTrackRecorder.h.

Constructor & Destructor Documentation

◆ MaterialTrackRecorder()

ActsTrk::MaterialTrackRecorder::MaterialTrackRecorder ( const Config & config)

Definition at line 20 of file MaterialTrackRecorder.cxx.

20 :
21 AthMessaging(Gaudi::svcLocator()->service< IMessageSvc >( "MessageSvc" ),"MaterialTrackRecorder"),
22 m_cfg(config)
23 {}
AthMessaging()
Default constructor:

Member Function Documentation

◆ BeginOfEventAction()

void ActsTrk::MaterialTrackRecorder::BeginOfEventAction ( const G4Event * )
overridevirtual

Definition at line 25 of file MaterialTrackRecorder.cxx.

26 {
27 ATH_MSG_DEBUG(" BeginOfEventAction");
28
29 // Get the event context for the write handle below
30 const EventContext& ctx = Gaudi::Hive::currentContext();
31
32 SG::WriteHandle<RecordedMaterialTrackCollection> materialTracks(m_cfg.materialTrackCollectionName, ctx);
33
34 // Only record if it does NOT already exist
35 if (!materialTracks.isPresent()) {
36 auto coll = std::make_unique<RecordedMaterialTrackCollection>();
37 if (materialTracks.record(std::move(coll)).isFailure()) {
38 ATH_MSG_ERROR("Failed to record RecordedMaterialTrackCollection with key "
39 << m_cfg.materialTrackCollectionName);
40 m_rmtCollection = nullptr;
41 return;
42 }
43 } else {
44 // Sanity check (should normally not happen unless someone else records it)
45 ATH_MSG_DEBUG("Collection already present: " << m_cfg.materialTrackCollectionName);
46 }
47
48 m_rmtCollection = materialTracks.ptr();
49
50 if (!m_rmtCollection) {
51 ATH_MSG_ERROR("Recorded collection pointer is null right after record/present for key "
52 << m_cfg.materialTrackCollectionName);
53 }
54
55 ATH_MSG_DEBUG("ctx event=" << ctx.eventID().event_number()
56 << " slot=" << ctx.slot() << " key=" << m_cfg.materialTrackCollectionName);
57 ATH_MSG_DEBUG("recorded? present=" << materialTracks.isPresent() << " valid=" << materialTracks.isValid());
58 }
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
RecordedMaterialTrackCollection * m_rmtCollection
Pointer to the collection, non-owning.

◆ BeginOfRunAction()

void ActsTrk::MaterialTrackRecorder::BeginOfRunAction ( const G4Run * )
overridevirtual

Definition at line 65 of file MaterialTrackRecorder.cxx.

66 {
67 ATH_MSG_DEBUG(" BeginOfRunAction");
68 m_rmtCollection = nullptr;
69 }

◆ convertDirection()

Acts::Vector3 ActsTrk::MaterialTrackRecorder::convertDirection ( const G4ThreeVector & g4vec)
inlineprivate

Useful fuction for direction coversion.

Definition at line 68 of file MaterialTrackRecorder.h.

68 {
69 return Acts::Vector3{g4vec[0], g4vec[1], g4vec[2]};
70 }

◆ convertPosition()

Acts::Vector3 ActsTrk::MaterialTrackRecorder::convertPosition ( const G4ThreeVector & g4vec)
inlineprivate

Useful fuction for position coversion.

Definition at line 63 of file MaterialTrackRecorder.h.

63 {
64 return Acts::Vector3(g4vec[0] * s_convertLength, g4vec[1] * s_convertLength,
65 g4vec[2] * s_convertLength);
66 }
static constexpr double s_convertLength
Conversion constants.

◆ EndOfEventAction()

void ActsTrk::MaterialTrackRecorder::EndOfEventAction ( const G4Event * )
overridevirtual

Definition at line 60 of file MaterialTrackRecorder.cxx.

61 {
62 ATH_MSG_DEBUG(" EndOfEventAction");
63 }

◆ initMessaging()

void AthMessaging::initMessaging ( ) const
privateinherited

Initialize our message level and MessageSvc.

This method should only be called once.

Definition at line 39 of file AthMessaging.cxx.

40{
42 // If user did not set an explicit level, set a default
43 if (m_lvl == MSG::NIL) {
44 m_lvl = m_imsg ?
45 static_cast<MSG::Level>( m_imsg.load()->outputLevel(m_nm) ) :
46 MSG::INFO;
47 }
48}
std::string m_nm
Message source name.
std::atomic< IMessageSvc * > m_imsg
MessageSvc pointer.
std::atomic< MSG::Level > m_lvl
Current logging level.
IMessageSvc * getMessageSvc(bool quiet=false)

◆ msg() [1/2]

MsgStream & AthMessaging::msg ( ) const
inlineinherited

The standard message stream.

Returns a reference to the default message stream May not be invoked before sysInitialize() has been invoked.

Definition at line 167 of file AthMessaging.h.

168{
169 MsgStream* ms = m_msg_tls.get();
170 if (!ms) {
171 if (!m_initialized.test_and_set()) initMessaging();
172 ms = new MsgStream(m_imsg,m_nm);
173 m_msg_tls.reset( ms );
174 }
175
176 ms->setLevel (m_lvl);
177 return *ms;
178}
boost::thread_specific_ptr< MsgStream > m_msg_tls
MsgStream instance (a std::cout like with print-out levels)
void initMessaging() const
Initialize our message level and MessageSvc.

◆ msg() [2/2]

MsgStream & AthMessaging::msg ( const MSG::Level lvl) const
inlineinherited

The standard message stream.

Returns a reference to the default message stream May not be invoked before sysInitialize() has been invoked.

Definition at line 182 of file AthMessaging.h.

183{ return msg() << lvl; }
MsgStream & msg() const
The standard message stream.

◆ msgLvl()

bool AthMessaging::msgLvl ( const MSG::Level lvl) const
inlineinherited

Test the output level.

Parameters
lvlThe message level to test against
Returns
boolean Indicating if messages at given level will be printed
Return values
trueMessages at level "lvl" will be printed

Definition at line 151 of file AthMessaging.h.

152{
153 // If user did not set explicit message level we have to initialize
154 // the messaging and retrieve the default via the MessageSvc.
155 if (m_lvl==MSG::NIL && !m_initialized.test_and_set()) initMessaging();
156
157 if (m_lvl <= lvl) {
158 msg() << lvl;
159 return true;
160 } else {
161 return false;
162 }
163}

◆ setLevel()

void AthMessaging::setLevel ( MSG::Level lvl)
inherited

Change the current logging level.

Use this rather than msg().setLevel() for proper operation with MT.

Definition at line 28 of file AthMessaging.cxx.

29{
30 m_lvl = lvl;
31}

◆ UserSteppingAction()

void ActsTrk::MaterialTrackRecorder::UserSteppingAction ( const G4Step * step)
overridevirtual

Definition at line 71 of file MaterialTrackRecorder.cxx.

72 {
73 ATH_MSG_DEBUG(" UserSteppingAction");
74
75 if (!m_rmtCollection) {
76 ATH_MSG_ERROR("No per-event collection pointer cached (m_coll is null).");
77 return;
78 }
79
80 // Get the material & check if it is present
81 G4Material* material = step->GetPreStepPoint()->GetMaterial();
82 if (material == nullptr) {
83 return;
84 }
85
86 // First check for exclusion
87 std::string materialName = material->GetName();
88 for (const auto& emat : m_cfg.excludeMaterials) {
89 if (emat == materialName) {
90 ATH_MSG_DEBUG("Exclude step in material '" << materialName << "'.");
91 return;
92 }
93 }
94
95 ATH_MSG_DEBUG("Performing a step with step size = "
96 << s_convertLength * step->GetStepLength());
97
98 // Quantities valid for elemental materials and mixtures
99 double X0 = s_convertLength * material->GetRadlen();
100 double L0 = s_convertLength * material->GetNuclearInterLength();
101 double rho = s_convertDensity * material->GetDensity();
102
103 // Get{A,Z} is only meaningful for single-element materials (according to
104 // the Geant4 docs). Need to compute average manually.
105 const G4ElementVector* elements = material->GetElementVector();
106 const G4double* fraction = material->GetFractionVector();
107 std::size_t nElements = material->GetNumberOfElements();
108
109 double Ar = 0.;
110 double Z = 0.;
111 if (nElements == 1) {
112 Ar = material->GetA() / (CLHEP::gram / CLHEP::mole);
113 Z = material->GetZ();
114 } else {
115 for (std::size_t i = 0; i < nElements; i++) {
116 Ar += elements->at(i)->GetA() * fraction[i] / (CLHEP::gram / CLHEP::mole);
117 Z += elements->at(i)->GetZ() * fraction[i];
118 }
119 }
120
121 // Construct passed material slab for the step
122 const auto slab = Acts::MaterialSlab(Acts::Material::fromMassDensity(X0, L0, Ar, Z, rho), s_convertLength * step->GetStepLength());
123
124 // Create the RecordedMaterialSlab
125 Acts::MaterialInteraction mInteraction;
126 mInteraction.position = convertPosition(step->GetPreStepPoint()->GetPosition());
127 mInteraction.direction = convertDirection(step->GetPreStepPoint()->GetMomentum()).normalized();
128 mInteraction.materialSlab = slab;
129 mInteraction.pathCorrection = (step->GetStepLength() / CLHEP::mm);
130
131 G4Track* g4Track = step->GetTrack();
132
134 Acts::Vector3 vertex = convertPosition(g4Track->GetVertexPosition());
135 Acts::Vector3 direction = convertDirection(g4Track->GetMomentumDirection());
136 rmTrack.first = {vertex, direction};
137 rmTrack.second.materialInteractions.push_back(mInteraction);
138 m_rmtCollection->push_back(rmTrack);
139 }
static constexpr double s_convertDensity
Acts::Vector3 convertPosition(const G4ThreeVector &g4vec)
Useful fuction for position coversion.
Acts::Vector3 convertDirection(const G4ThreeVector &g4vec)
Useful fuction for direction coversion.
std::pair< std::pair< Acts::Vector3, Acts::Vector3 >, RecordedMaterial > RecordedMaterialTrack
Recorded material track.

Member Data Documentation

◆ ATLAS_THREAD_SAFE

std::atomic_flag m_initialized AthMessaging::ATLAS_THREAD_SAFE = ATOMIC_FLAG_INIT
mutableprivateinherited

Messaging initialized (initMessaging)

Definition at line 141 of file AthMessaging.h.

◆ m_cfg

Config ActsTrk::MaterialTrackRecorder::m_cfg
private

The config class.

Definition at line 60 of file MaterialTrackRecorder.h.

◆ m_imsg

std::atomic<IMessageSvc*> AthMessaging::m_imsg { nullptr }
mutableprivateinherited

MessageSvc pointer.

Definition at line 135 of file AthMessaging.h.

135{ nullptr };

◆ m_lvl

std::atomic<MSG::Level> AthMessaging::m_lvl { MSG::NIL }
mutableprivateinherited

Current logging level.

Definition at line 138 of file AthMessaging.h.

138{ MSG::NIL };

◆ m_msg_tls

boost::thread_specific_ptr<MsgStream> AthMessaging::m_msg_tls
mutableprivateinherited

MsgStream instance (a std::cout like with print-out levels)

Definition at line 132 of file AthMessaging.h.

◆ m_nm

std::string AthMessaging::m_nm
privateinherited

Message source name.

Definition at line 129 of file AthMessaging.h.

◆ m_rmtCollection

RecordedMaterialTrackCollection* ActsTrk::MaterialTrackRecorder::m_rmtCollection {nullptr}
private

Pointer to the collection, non-owning.

Definition at line 53 of file MaterialTrackRecorder.h.

53{nullptr};

◆ s_convertDensity

double ActsTrk::MaterialTrackRecorder::s_convertDensity
staticconstexprprivate
Initial value:
=
(Acts::UnitConstants::g / Acts::UnitConstants::mm3) / (CLHEP::gram / CLHEP::mm3)

Definition at line 56 of file MaterialTrackRecorder.h.

◆ s_convertLength

double ActsTrk::MaterialTrackRecorder::s_convertLength = Acts::UnitConstants::mm / CLHEP::mm
staticconstexprprivate

Conversion constants.

Definition at line 55 of file MaterialTrackRecorder.h.


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