#include <StepNtuple.h>
|
| | StepNtuple (const MSG::Level lvl=MSG::INFO) |
| | Constructor with message level argument for AthMessaging.
|
| virtual void | BeginOfEventAction (const G4Event *) override |
| | the hooks for G4 UA handling
|
| virtual void | EndOfEventAction (const G4Event *) override |
| virtual void | UserSteppingAction (const G4Step *) override |
| virtual void | BeginOfRunAction (const G4Run *) 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.
|
|
| void | initMessaging () const |
| | Initialize our message level and MessageSvc.
|
|
| std::vector< stepdata > | eventSteps |
| | holds data extracted from steps
|
| NTuple::Item< long > | m_nsteps |
| | handles for ntuple writing
|
| NTuple::Array< float > | m_pdgcode |
| NTuple::Array< float > | m_step_x |
| NTuple::Array< float > | m_step_y |
| NTuple::Array< float > | m_step_z |
| NTuple::Array< float > | m_time |
| NTuple::Array< float > | m_dep |
| 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)
|
Definition at line 21 of file StepNtuple.h.
◆ StepNtuple()
| G4UA::StepNtuple::StepNtuple |
( |
const MSG::Level | lvl = MSG::INFO | ) |
|
Constructor with message level argument for AthMessaging.
Definition at line 22 of file StepNtuple.cxx.
23 :
AthMessaging(Gaudi::svcLocator()->service< IMessageSvc >(
"MessageSvc" ),
24 "StepNtuple")
25 {
27 }
void setLevel(MSG::Level lvl)
Change the current logging level.
AthMessaging()
Default constructor:
◆ BeginOfEventAction()
| void G4UA::StepNtuple::BeginOfEventAction |
( |
const G4Event * | | ) |
|
|
overridevirtual |
the hooks for G4 UA handling
Definition at line 29 of file StepNtuple.cxx.
30 {
32 }
std::vector< stepdata > eventSteps
holds data extracted from steps
◆ BeginOfRunAction()
| void G4UA::StepNtuple::BeginOfRunAction |
( |
const G4Run * | | ) |
|
|
overridevirtual |
Definition at line 72 of file StepNtuple.cxx.
73 {
74 NTupleFilePtr file1(
ntupleSvc(),
"/NTUPLES/FILE1");
75
76 SmartDataPtr<NTuple::Directory>
77 ntdir(
ntupleSvc(),
"/NTUPLES/FILE1/StepNtuple");
78 if ( !ntdir ) {
79
80 ntdir =
ntupleSvc()->createDirectory(file1,
"StepNtuple");
81 }
82 if ( ! ntdir ) {
84 }
85
86 NTuplePtr
nt(
ntupleSvc(),
"/NTUPLES/FILE1/StepNtuple/10");
87
88
89 if ( !nt ) {
90 nt =
ntupleSvc()->book( ntdir.ptr(), 10, CLID_ColumnWiseTuple,
"G4 Step Ntuple" );
91
92 if ( nt ) {
94
95
96 if(
nt->addItem (
"NSteps",
m_nsteps,0, 50000).isFailure() ||
104 }
105 }
106 }
107 }
NTuple::Array< float > m_pdgcode
NTuple::Array< float > m_step_y
NTuple::Item< long > m_nsteps
handles for ntuple writing
NTuple::Array< float > m_time
NTuple::Array< float > m_step_x
NTuple::Array< float > m_step_z
NTuple::Array< float > m_dep
◆ EndOfEventAction()
| void G4UA::StepNtuple::EndOfEventAction |
( |
const G4Event * | | ) |
|
|
overridevirtual |
Definition at line 34 of file StepNtuple.cxx.
35 {
38
46 }
47
48 if(!
ntupleSvc()->writeRecord(
"/NTUPLES/FILE1/StepNtuple/10").isSuccess()){
50 }
51 }
◆ 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
43 if (
m_lvl == MSG::NIL) {
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 163 of file AthMessaging.h.
164{
166 if (!ms) {
170 }
171
174}
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 178 of file AthMessaging.h.
179{
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
-
| lvl | The message level to test against |
- Returns
- boolean Indicating if messages at given level will be printed
- Return values
-
| true | Messages at level "lvl" will be printed |
Definition at line 151 of file AthMessaging.h.
152{
155 return true;
156 } else {
157 return false;
158 }
159}
◆ setLevel()
| void AthMessaging::setLevel |
( |
MSG::Level | lvl | ) |
|
|
inherited |
◆ UserSteppingAction()
| void G4UA::StepNtuple::UserSteppingAction |
( |
const G4Step * | aStep | ) |
|
|
overridevirtual |
Definition at line 53 of file StepNtuple.cxx.
54 {
56
58
59 theInfo.
dep=aStep->GetTotalEnergyDeposit();
60 theInfo.time=aStep->GetPreStepPoint()->GetGlobalTime();
61 theInfo.code=aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
62 G4ThreeVector
pos=aStep->GetPreStepPoint()->GetPosition();
66
69 }
70 }
#define ATH_MSG_VERBOSE(x)
simple struct to hold step information
◆ ATLAS_THREAD_SAFE
| std::atomic_flag m_initialized AthMessaging::ATLAS_THREAD_SAFE = ATOMIC_FLAG_INIT |
|
mutableprivateinherited |
◆ eventSteps
| std::vector<stepdata> G4UA::StepNtuple::eventSteps |
|
private |
holds data extracted from steps
Definition at line 44 of file StepNtuple.h.
◆ m_dep
| NTuple::Array<float> G4UA::StepNtuple::m_dep |
|
private |
◆ m_imsg
| std::atomic<IMessageSvc*> AthMessaging::m_imsg { nullptr } |
|
mutableprivateinherited |
◆ m_lvl
| std::atomic<MSG::Level> AthMessaging::m_lvl { MSG::NIL } |
|
mutableprivateinherited |
◆ 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 |
◆ m_nsteps
| NTuple::Item<long> G4UA::StepNtuple::m_nsteps |
|
private |
handles for ntuple writing
Definition at line 47 of file StepNtuple.h.
◆ m_pdgcode
| NTuple::Array<float> G4UA::StepNtuple::m_pdgcode |
|
private |
◆ m_step_x
| NTuple::Array<float> G4UA::StepNtuple::m_step_x |
|
private |
◆ m_step_y
| NTuple::Array<float> G4UA::StepNtuple::m_step_y |
|
private |
◆ m_step_z
| NTuple::Array<float> G4UA::StepNtuple::m_step_z |
|
private |
◆ m_time
| NTuple::Array<float> G4UA::StepNtuple::m_time |
|
private |
The documentation for this class was generated from the following files: