ATLAS Offline Software
Loading...
Searching...
No Matches
G4UA::StepHistogram Class Reference

#include <StepHistogram.h>

Inheritance diagram for G4UA::StepHistogram:
Collaboration diagram for G4UA::StepHistogram:

Classes

struct  Config
struct  Report
 this holds all the data from individual threads that needs to be merged at EoR More...

Public Types

typedef std::map< G4String, TH1 * > HistoMap_t
typedef std::map< G4String, HistoMap_tHistoMapMap_t

Public Member Functions

virtual void UserSteppingAction (const G4Step *) override
 the hooks for G4 UA handling
 StepHistogram (const Config &)
 ctor
const ReportgetReport () const
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

void InitializeFillHistogram2D (HistoMapMap_t &hMapMap, const char *suffix, const G4String &pdgId, const G4String &vol, int nbinsx, double xmin, double xmax, int nbinsy, double ymin, double ymax, double valuex, double valuey, double weight)
void InitializeFillHistogram (HistoMapMap_t &hMapMap, const char *suffix, const G4String &pdgId, const G4String &vol, int nbins, double xmin, double xmax, double value, double weight)
void InitializeFillHistogram (HistoMapMap_t &hMapMap, const char *suffix, const G4String &pdgId, const G4String &vol, int nbins, double *edges, double value, double weight)
void initMessaging () const
 Initialize our message level and MessageSvc.

Private Attributes

Report m_report
Config m_config
 configuration data
float m_initialKineticEnergyOfStep
G4String m_initialVolume
G4String m_initialMaterial
G4String m_initialProcess
int m_trackID = 0
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)

Detailed Description

Definition at line 29 of file StepHistogram.h.

Member Typedef Documentation

◆ HistoMap_t

typedef std::map<G4String, TH1*> G4UA::StepHistogram::HistoMap_t

Definition at line 39 of file StepHistogram.h.

◆ HistoMapMap_t

typedef std::map<G4String, HistoMap_t> G4UA::StepHistogram::HistoMapMap_t

Definition at line 40 of file StepHistogram.h.

Constructor & Destructor Documentation

◆ StepHistogram()

G4UA::StepHistogram::StepHistogram ( const Config & config)

ctor

Definition at line 30 of file StepHistogram.cxx.

30 :
31 AthMessaging(Gaudi::svcLocator()->service<IMessageSvc>("MessageSvc"),"StepHistogram"),
32 m_config(config),
37 {}
AthMessaging()
Default constructor:
Config m_config
configuration data

Member Function Documentation

◆ getReport()

const Report & G4UA::StepHistogram::getReport ( ) const
inline

Definition at line 116 of file StepHistogram.h.

116{ return m_report; }

◆ InitializeFillHistogram() [1/2]

void G4UA::StepHistogram::InitializeFillHistogram ( HistoMapMap_t & hMapMap,
const char * suffix,
const G4String & pdgId,
const G4String & vol,
int nbins,
double * edges,
double value,
double weight )
private

Definition at line 256 of file StepHistogram.cxx.

259 {
260 const auto& [pStringMapPair, inserted] = hMapMap.try_emplace(vol,HistoMap_t());
261 HistoMap_t &hMap = pStringMapPair->second;
262 if ( hMap.find(particleName) == hMap.end() ) {
263 // initialize histogram if not yet exist
264 std::ostringstream stringStream;
265 stringStream << vol << "_" << particleName << "_" << suffix;
266 hMap[particleName] = new TH1F(stringStream.str().c_str(), stringStream.str().c_str(), nbins, edges);
267 }
268 hMap[particleName]->Fill(value, weight);
269 }
std::map< G4String, TH1 * > HistoMap_t
std::string particleName(const G4Step *theStep)
TODO.
TH1F(name, title, nxbins, bins_par2, bins_par3=None, path='', **kwargs)

◆ InitializeFillHistogram() [2/2]

void G4UA::StepHistogram::InitializeFillHistogram ( HistoMapMap_t & hMapMap,
const char * suffix,
const G4String & pdgId,
const G4String & vol,
int nbins,
double xmin,
double xmax,
double value,
double weight )
private

Definition at line 241 of file StepHistogram.cxx.

244 {
245 const auto& [pStringMapPair, inserted] = hMapMap.try_emplace(vol,HistoMap_t());
246 HistoMap_t &hMap = pStringMapPair->second;
247 if ( hMap.find(particleName) == hMap.end() ) {
248 // initialize histogram if not yet exist
249 std::ostringstream stringStream;
250 stringStream << vol << "_" << particleName << "_" << suffix;
251 hMap[particleName] = new TH1F(stringStream.str().c_str(), stringStream.str().c_str(), nbins, xmin, xmax);
252 }
253 hMap[particleName]->Fill(value, weight);
254 }
double xmax
Definition listroot.cxx:61
double xmin
Definition listroot.cxx:60

◆ InitializeFillHistogram2D()

void G4UA::StepHistogram::InitializeFillHistogram2D ( HistoMapMap_t & hMapMap,
const char * suffix,
const G4String & pdgId,
const G4String & vol,
int nbinsx,
double xmin,
double xmax,
int nbinsy,
double ymin,
double ymax,
double valuex,
double valuey,
double weight )
private

Definition at line 224 of file StepHistogram.cxx.

229 {
230 const auto& [pStringMapPair, inserted] = hMapMap.try_emplace(vol,HistoMap_t());
231 HistoMap_t &hMap = pStringMapPair->second;
232 if ( hMap.find(particleName) == hMap.end() ) {
233 // initialize histogram if not yet exist
234 std::ostringstream stringStream;
235 stringStream << vol << "_" << particleName << "_" << suffix;
236 hMap[particleName] = new TH2F(stringStream.str().c_str(), stringStream.str().c_str(), nbinsx, xmin, xmax, nbinsy, ymin, ymax);
237 }
238 static_cast<TH2*>(hMap[particleName])->Fill(valuex, valuey, weight);
239 }
double ymin
Definition listroot.cxx:63
double ymax
Definition listroot.cxx:64
TH2F(name, title, nxbins, bins_par2, bins_par3, bins_par4, bins_par5=None, bins_par6=None, path='', **kwargs)

◆ 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 163 of file AthMessaging.h.

164{
165 MsgStream* ms = m_msg_tls.get();
166 if (!ms) {
167 if (!m_initialized.test_and_set()) initMessaging();
168 ms = new MsgStream(m_imsg,m_nm);
169 m_msg_tls.reset( ms );
170 }
171
172 ms->setLevel (m_lvl);
173 return *ms;
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
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 (m_lvl <= lvl) {
154 msg() << lvl;
155 return true;
156 } else {
157 return false;
158 }
159}

◆ 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 G4UA::StepHistogram::UserSteppingAction ( const G4Step * aStep)
overridevirtual

the hooks for G4 UA handling

Definition at line 39 of file StepHistogram.cxx.

39 {
40
41 // track
42 G4Track *tr = aStep->GetTrack();
43
44 // pre-step point
45 G4StepPoint *PreStepPoint = aStep->GetPreStepPoint();
46
47 // post-step point
48 G4StepPoint *PostStepPoint = aStep->GetPostStepPoint();
49
50 // pre-step kinetic energy
51 double stepKinetic = PreStepPoint->GetKineticEnergy();
52
53 // post-step kinetic energy
54 double postStepKinetic = PostStepPoint->GetKineticEnergy();
55
56 // pre-step position
57 const G4ThreeVector& myPos = PreStepPoint->GetPosition();
58
59 // particle name
60 G4String particleName = "nucleus";
61 if (!(tr->GetDefinition()->GetParticleType() == "nucleus"))
62 particleName = G4DebuggingHelpers::ClassifyParticle(tr->GetParticleDefinition());
63
64 // pre-step volume
65 G4String volumeName = PreStepPoint->GetPhysicalVolume()->GetName();
66 volumeName = G4DebuggingHelpers::ClassifyVolume(volumeName);
67
68 // pre-step material
69 G4String materialName = PreStepPoint->GetMaterial()->GetName();
70 materialName = G4DebuggingHelpers::ClassifyMaterial(materialName);
71
72 // process name (uses post-step point)
73 G4String processName = PostStepPoint->GetProcessDefinedStep() ?
74 PostStepPoint->GetProcessDefinedStep()->GetProcessName() : "Unknown";
75
76 // secondaries
77 const std::vector<const G4Track*>* secondaries = aStep->GetSecondaryInCurrentStep();
78
79 // 2D map
80 if (m_config.do2DHistograms) {
81 float DepositedE = aStep->GetTotalEnergyDeposit();
82 InitializeFillHistogram2D(m_report.histoMapMap2D_vol_RZ, "vol_RZ", particleName, volumeName,
83 2000, -20000, 20000, 1000, 0, 5000, myPos.getZ(), myPos.perp(), 1.);
84 InitializeFillHistogram2D(m_report.histoMapMap2D_mat_RZ, "mat_RZ", particleName, materialName,
85 2000, -20000, 20000, 1000, 0, 5000, myPos.getZ(), myPos.perp(), 1.);
86 InitializeFillHistogram2D(m_report.histoMapMap2D_prc_RZ, "prc_RZ", particleName, processName,
87 2000, -20000, 20000, 1000, 0, 5000, myPos.getZ(), myPos.perp(), 1.);
88 InitializeFillHistogram2D(m_report.histoMapMap2D_vol_RZ_E, "vol_RZ_E", particleName, volumeName,
89 2000, -20000, 20000, 1000, 0, 5000, myPos.getZ(), myPos.perp(), DepositedE);
90 InitializeFillHistogram2D(m_report.histoMapMap2D_mat_RZ_E, "mat_RZ_E", particleName, materialName,
91 2000, -20000, 20000, 1000, 0, 5000, myPos.getZ(), myPos.perp(), DepositedE);
92 }
93
94 // step length
95 InitializeFillHistogram(m_report.histoMapMap_vol_stepSize, "vol_stepLength", particleName, volumeName,
96 1000, -12, 4, std::log10(aStep->GetStepLength()), 1.);
97 InitializeFillHistogram(m_report.histoMapMap_mat_stepSize, "mat_stepLength", particleName, materialName,
98 1000, -12, 4, std::log10(aStep->GetStepLength()), 1.);
99 InitializeFillHistogram(m_report.histoMapMap_prc_stepSize, "prc_stepLength", particleName, processName,
100 1000, -12, 4, std::log10(aStep->GetStepLength()), 1.);
101
102 // step pseudorapidity
103 InitializeFillHistogram(m_report.histoMapMap_vol_stepPseudorapidity, "vol_stepPseudorapidity", particleName, volumeName,
104 200, -10, 10, myPos.eta(), 1.);
105 InitializeFillHistogram(m_report.histoMapMap_mat_stepPseudorapidity, "mat_stepPseudorapidity", particleName, materialName,
106 200, -10, 10, myPos.eta(), 1.);
107 InitializeFillHistogram(m_report.histoMapMap_prc_stepPseudorapidity, "prc_stepPseudorapidity", particleName, processName,
108 200, -10, 10, myPos.eta(), 1.);
109
110 // step kinetic energy
111 InitializeFillHistogram(m_report.histoMapMap_vol_stepKineticEnergy, "vol_stepKineticEnergy", particleName, volumeName,
112 1000, -9, 7, std::log10(stepKinetic), 1.);
113 InitializeFillHistogram(m_report.histoMapMap_mat_stepKineticEnergy, "mat_stepKineticEnergy", particleName, materialName,
114 1000, -9, 7, std::log10(stepKinetic), 1.);
115 InitializeFillHistogram(m_report.histoMapMap_prc_stepKineticEnergy, "prc_stepKineticEnergy", particleName, processName,
116 1000, -9, 7, std::log10(stepKinetic), 1.);
117 InitializeFillHistogram(m_report.histoMapMap_stepKinetic, "stepKineticEnergy", particleName, "AllATLAS",
118 1000, -9, 7, std::log10(stepKinetic), 1.);
119
120 // post step kinetic energy
121 InitializeFillHistogram(m_report.histoMapMap_vol_postStepKineticEnergy, "vol_postStepKineticEnergy", particleName, volumeName,
122 1000, -9, 7, std::log10(postStepKinetic), 1.);
123 InitializeFillHistogram(m_report.histoMapMap_mat_postStepKineticEnergy, "mat_postStepKineticEnergy", particleName, materialName,
124 1000, -9, 7, std::log10(postStepKinetic), 1.);
125 InitializeFillHistogram(m_report.histoMapMap_prc_postStepKineticEnergy, "prc_postStepKineticEnergy", particleName, processName,
126 1000, -9, 7, std::log10(postStepKinetic), 1.);
127 InitializeFillHistogram(m_report.histoMapMap_postStepKinetic, "postStepKineticEnergy", particleName, "AllATLAS",
128 1000, -9, 7, std::log10(postStepKinetic), 1.);
129
130 // step energy deposit
131 InitializeFillHistogram(m_report.histoMapMap_vol_stepEnergyDeposit, "vol_stepEnergyDeposit", particleName, volumeName,
132 1000, -11, 3, std::log10(aStep->GetTotalEnergyDeposit()), 1.);
133 InitializeFillHistogram(m_report.histoMapMap_mat_stepEnergyDeposit, "mat_stepEnergyDeposit", particleName, materialName,
134 1000, -11, 3, std::log10(aStep->GetTotalEnergyDeposit()), 1.);
135 InitializeFillHistogram(m_report.histoMapMap_prc_stepEnergyDeposit, "prc_stepEnergyDeposit", particleName, processName,
136 1000, -11, 3, std::log10(aStep->GetTotalEnergyDeposit()), 1.);
137
138 // step non-ionizing energy deposit
139 InitializeFillHistogram(m_report.histoMapMap_vol_stepEnergyNonIonDeposit, "vol_stepEnergyNonIonDeposit", particleName, volumeName,
140 1000, -11, 1, std::log10(aStep->GetNonIonizingEnergyDeposit()), 1.);
141 InitializeFillHistogram(m_report.histoMapMap_mat_stepEnergyNonIonDeposit, "mat_stepEnergyNonIonDeposit", particleName, materialName,
142 1000, -11, 1, std::log10(aStep->GetNonIonizingEnergyDeposit()), 1.);
143 InitializeFillHistogram(m_report.histoMapMap_prc_stepEnergyNonIonDeposit, "prc_stepEnergyNonIonDeposit", particleName, processName,
144 1000, -11, 1, std::log10(aStep->GetNonIonizingEnergyDeposit()), 1.);
145
146 // secondary kinetic energy
147 for (const auto &track : *secondaries) {
148 G4String secondary_particleName = G4DebuggingHelpers::ClassifyParticle(track->GetParticleDefinition());
149 InitializeFillHistogram(m_report.histoMapMap_vol_stepSecondaryKinetic, "vol_stepSecondaryKinetic", secondary_particleName, volumeName,
150 1000, -7, 5, std::log10(track->GetKineticEnergy()), 1.);
151 InitializeFillHistogram(m_report.histoMapMap_mat_stepSecondaryKinetic, "mat_stepSecondaryKinetic", secondary_particleName, materialName,
152 1000, -7, 5, std::log10(track->GetKineticEnergy()), 1.);
153 InitializeFillHistogram(m_report.histoMapMap_prc_stepSecondaryKinetic, "prc_stepSecondaryKinetic", secondary_particleName, processName,
154 1000, -7, 5, std::log10(track->GetKineticEnergy()), 1.);
155 }
156
157 // stop here if 'general' histograms not activated
158 // _______________________________________________
159 if (!m_config.doGeneralHistograms)
160 return;
161
162 // first step (after initial step)
163 if (tr->GetCurrentStepNumber()==1) {
164 // initial kinetic energy
165 m_initialKineticEnergyOfStep = stepKinetic;
166
167 // initial volume/material/processes
168 m_initialVolume = std::move(volumeName);
169 m_initialMaterial = std::move(materialName);
170 m_initialProcess = std::move(processName);
171
172 // save track ID for checking if we later have the same track
173 m_trackID = tr->GetTrackID();
174
175 // initial energy
176 InitializeFillHistogram(m_report.histoMapMap_InitialE, "InitialE", particleName, "AllATLAS",
177 1000, -9, 7, std::log10(m_initialKineticEnergyOfStep), 1.0);
178 InitializeFillHistogram(m_report.histoMapMap_vol_InitialE, "vol_InitialE", particleName, m_initialVolume,
179 1000, -9, 7, std::log10(m_initialKineticEnergyOfStep), 1.0);
180 InitializeFillHistogram(m_report.histoMapMap_mat_InitialE, "mat_InitialE", particleName, m_initialMaterial,
181 1000, -9, 7, std::log10(m_initialKineticEnergyOfStep), 1.0);
182 InitializeFillHistogram(m_report.histoMapMap_prc_InitialE, "prc_InitialE", particleName, m_initialProcess,
183 1000, -9, 7, std::log10(m_initialKineticEnergyOfStep), 1.0);
184 }
185
186 // last step
187 if ( tr->GetTrackStatus() == 2 ) {
188 // assert to check if we have the correct track
189 if (not (tr->GetTrackID() == m_trackID)) {
190 ATH_MSG_ERROR("Track ID changed between the assumed first step and the last.");
191 throw std::exception();
192 }
193 // number of steps
194 int nSteps = tr->GetCurrentStepNumber() + 1;
195 InitializeFillHistogram(m_report.histoMapMap_numberOfSteps, "numberOfSteps", particleName, "AllATLAS",
196 10000, 0.5, 10000.5, nSteps, 1.);
197 InitializeFillHistogram(m_report.histoMapMap_vol_numberOfSteps, "vol_numberOfSteps", particleName, m_initialVolume,
198 10000, 0.5, 10000.5, nSteps, 1.);
199 InitializeFillHistogram(m_report.histoMapMap_mat_numberOfSteps, "mat_numberOfSteps", particleName, m_initialMaterial,
200 10000, 0.5, 10000.5, nSteps, 1.);
201 InitializeFillHistogram(m_report.histoMapMap_prc_numberOfSteps, "prc_numberOfSteps", particleName, m_initialProcess,
202 10000, 0.5, 10000.5, nSteps, 1.);
203 // number of steps vs initial energy
204 InitializeFillHistogram(m_report.histoMapMap_numberOfStepsPerInitialE, "numberOfStepsPerInitialE", particleName, "AllATLAS",
205 1000, -9, 7, std::log10(m_initialKineticEnergyOfStep), nSteps);
206 InitializeFillHistogram(m_report.histoMapMap_vol_numberOfStepsPerInitialE, "vol_numberOfStepsPerInitialE", particleName, m_initialVolume,
207 1000, -9, 7, std::log10(m_initialKineticEnergyOfStep), nSteps);
208 InitializeFillHistogram(m_report.histoMapMap_mat_numberOfStepsPerInitialE, "mat_numberOfStepsPerInitialE", particleName, m_initialMaterial,
209 1000, -9, 7, std::log10(m_initialKineticEnergyOfStep), nSteps);
210 InitializeFillHistogram(m_report.histoMapMap_prc_numberOfStepsPerInitialE, "prc_numberOfStepsPerInitialE", particleName, m_initialProcess,
211 1000, -9, 7, std::log10(m_initialKineticEnergyOfStep), nSteps);
212 // track length vs initial energy
213 InitializeFillHistogram(m_report.histoMapMap_trackLengthPerInitialE, "trackLengthPerInitialE", particleName, "AllATLAS",
214 1000, -9, 7, std::log10(tr->GetTrackLength()), 1.);
215 InitializeFillHistogram(m_report.histoMapMap_vol_trackLengthPerInitialE, "vol_trackLengthPerInitialE", particleName, m_initialVolume,
216 1000, -9, 7, std::log10(tr->GetTrackLength()), 1.);
217 InitializeFillHistogram(m_report.histoMapMap_mat_trackLengthPerInitialE, "mat_trackLengthPerInitialE", particleName, m_initialMaterial,
218 1000, -9, 7, std::log10(tr->GetTrackLength()), 1.);
219 InitializeFillHistogram(m_report.histoMapMap_prc_trackLengthPerInitialE, "prc_trackLengthPerInitialE", particleName, m_initialProcess,
220 1000, -9, 7, std::log10(tr->GetTrackLength()), 1.);
221 }
222 }
#define ATH_MSG_ERROR(x)
void InitializeFillHistogram2D(HistoMapMap_t &hMapMap, const char *suffix, const G4String &pdgId, const G4String &vol, int nbinsx, double xmin, double xmax, int nbinsy, double ymin, double ymax, double valuex, double valuey, double weight)
void InitializeFillHistogram(HistoMapMap_t &hMapMap, const char *suffix, const G4String &pdgId, const G4String &vol, int nbins, double xmin, double xmax, double value, double weight)
const G4String ClassifyMaterial(const G4String &nom)
const G4String ClassifyVolume(const G4String &nom)
const G4String ClassifyParticle(const G4ParticleDefinition *def)

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_config

Config G4UA::StepHistogram::m_config
private

configuration data

Definition at line 123 of file StepHistogram.h.

◆ m_imsg

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

MessageSvc pointer.

Definition at line 135 of file AthMessaging.h.

135{ nullptr };

◆ m_initialKineticEnergyOfStep

float G4UA::StepHistogram::m_initialKineticEnergyOfStep
private

Definition at line 140 of file StepHistogram.h.

◆ m_initialMaterial

G4String G4UA::StepHistogram::m_initialMaterial
private

Definition at line 142 of file StepHistogram.h.

◆ m_initialProcess

G4String G4UA::StepHistogram::m_initialProcess
private

Definition at line 143 of file StepHistogram.h.

◆ m_initialVolume

G4String G4UA::StepHistogram::m_initialVolume
private

Definition at line 141 of file StepHistogram.h.

◆ 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_report

Report G4UA::StepHistogram::m_report
private

Definition at line 120 of file StepHistogram.h.

◆ m_trackID

int G4UA::StepHistogram::m_trackID = 0
private

Definition at line 144 of file StepHistogram.h.


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