Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4SteppingVerboseWithUnits Class Reference

#include <G4SteppingVerboseWithUnits.hh>

+ Inheritance diagram for G4SteppingVerboseWithUnits:

Public Member Functions

 G4SteppingVerboseWithUnits (G4int precision=4)
 
 ~G4SteppingVerboseWithUnits () override
 
G4VSteppingVerboseClone () override
 
void SetManager (G4SteppingManager *const) override
 
void TrackingStarted () override
 
void StepInfo () override
 
void AtRestDoItInvoked () override
 
void AlongStepDoItAllDone () override
 
void PostStepDoItAllDone () override
 
void AlongStepDoItOneByOne () override
 
void PostStepDoItOneByOne () override
 
void DPSLStarted () override
 
void DPSLUserLimit () override
 
void DPSLPostStep () override
 
void DPSLAlongStep () override
 
void VerboseTrack () override
 
void VerboseParticleChange () override
 
void ShowStep () const override
 
- Public Member Functions inherited from G4SteppingVerbose
 G4SteppingVerbose ()=default
 
 ~G4SteppingVerbose () override=default
 
void NewStep () override
 
- Public Member Functions inherited from G4VSteppingVerbose
virtual ~G4VSteppingVerbose ()
 
void CopyState ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4SteppingVerbose
static void UseBestUnit (G4int prec=4)
 
static G4int BestUnitPrecision ()
 
- Static Public Member Functions inherited from G4VSteppingVerbose
static void SetInstance (G4VSteppingVerbose *Instance)
 
static G4VSteppingVerboseGetInstance ()
 
static G4VSteppingVerboseGetMasterInstance ()
 
static G4int GetSilent ()
 
static void SetSilent (G4int fSilent)
 
static G4int GetSilentStepInfo ()
 
static void SetSilentStepInfo (G4int fSilent)
 
- Protected Types inherited from G4VSteppingVerbose
using G4SelectedAtRestDoItVector = std::vector<G4int>
 
using G4SelectedAlongStepDoItVector = std::vector<G4int>
 
using G4SelectedPostStepDoItVector = std::vector<G4int>
 
- Protected Member Functions inherited from G4VSteppingVerbose
 G4VSteppingVerbose ()
 
- Protected Attributes inherited from G4VSteppingVerbose
G4SteppingManagerfManager = nullptr
 
G4UserSteppingActionfUserSteppingAction = nullptr
 
G4double PhysicalStep = 0.0
 
G4double GeometricalStep = 0.0
 
G4double CorrectedStep = 0.0
 
G4bool PreStepPointIsGeom = false
 
G4bool FirstStep = false
 
G4StepStatus fStepStatus = fUndefined
 
G4double TempInitVelocity = 0.0
 
G4double TempVelocity = 0.0
 
G4double Mass = 0.0
 
G4double sumEnergyChange = 0.0
 
G4VParticleChangefParticleChange = nullptr
 
G4TrackfTrack = nullptr
 
G4TrackVectorfSecondary = nullptr
 
G4StepfStep = nullptr
 
G4StepPointfPreStepPoint = nullptr
 
G4StepPointfPostStepPoint = nullptr
 
G4VPhysicalVolumefCurrentVolume = nullptr
 
G4VSensitiveDetectorfSensitive = nullptr
 
G4VProcessfCurrentProcess = nullptr
 
G4ProcessVectorfAtRestDoItVector = nullptr
 
G4ProcessVectorfAlongStepDoItVector = nullptr
 
G4ProcessVectorfPostStepDoItVector = nullptr
 
G4ProcessVectorfAtRestGetPhysIntVector = nullptr
 
G4ProcessVectorfAlongStepGetPhysIntVector = nullptr
 
G4ProcessVectorfPostStepGetPhysIntVector = nullptr
 
std::size_t MAXofAtRestLoops = 0
 
std::size_t MAXofAlongStepLoops = 0
 
std::size_t MAXofPostStepLoops = 0
 
G4double currentMinimumStep = 0.0
 
G4double numberOfInteractionLengthLeft = 0.0
 
std::size_t fAtRestDoItProcTriggered = 0
 
std::size_t fAlongStepDoItProcTriggered = 0
 
std::size_t fPostStepDoItProcTriggered = 0
 
G4int fN2ndariesAtRestDoIt = 0
 
G4int fN2ndariesAlongStepDoIt = 0
 
G4int fN2ndariesPostStepDoIt = 0
 
G4NavigatorfNavigator = nullptr
 
G4int verboseLevel = 0
 
G4SelectedAtRestDoItVectorfSelectedAtRestDoItVector = nullptr
 
G4SelectedAlongStepDoItVectorfSelectedAlongStepDoItVector = nullptr
 
G4SelectedPostStepDoItVectorfSelectedPostStepDoItVector = nullptr
 
G4double fPreviousStepSize = 0.0
 
G4TouchableHandle fTouchableHandle
 
G4SteppingControl StepControlFlag = NormalCondition
 
G4double physIntLength = 0.0
 
G4ForceCondition fCondition = InActivated
 
G4GPILSelection fGPILSelection = NotCandidateForSelection
 
- Static Protected Attributes inherited from G4VSteppingVerbose
static G4ThreadLocal G4VSteppingVerbosefInstance = nullptr
 
static G4VSteppingVerbosefMasterInstance = nullptr
 
static G4TRACKING_DLL G4ThreadLocal G4int Silent = 0
 
static G4TRACKING_DLL G4ThreadLocal G4int SilentStepInfo = 0
 

Detailed Description

Definition at line 48 of file G4SteppingVerboseWithUnits.hh.

Constructor & Destructor Documentation

◆ G4SteppingVerboseWithUnits()

G4SteppingVerboseWithUnits::G4SteppingVerboseWithUnits ( G4int precision = 4)

Definition at line 46 of file G4SteppingVerboseWithUnits.cc.

46: fprec(prec) {}

Referenced by Clone().

◆ ~G4SteppingVerboseWithUnits()

G4SteppingVerboseWithUnits::~G4SteppingVerboseWithUnits ( )
override

Definition at line 50 of file G4SteppingVerboseWithUnits.cc.

50{ delete fmessenger; }

Member Function Documentation

◆ AlongStepDoItAllDone()

void G4SteppingVerboseWithUnits::AlongStepDoItAllDone ( )
overridevirtual

Reimplemented from G4SteppingVerbose.

Definition at line 206 of file G4SteppingVerboseWithUnits.cc.

207{
208 G4VProcess* ptProcManager;
209
210 CopyState();
211
212 if (verboseLevel >= 3) {
213 G4cout << G4endl;
214 G4cout << " >>AlongStepDoIt (after all invocations):" << G4endl;
215 G4cout << " ++List of invoked processes " << G4endl;
216
217 for (std::size_t ci = 0; ci < MAXofAlongStepLoops; ++ci) {
218 ptProcManager = (*fAlongStepDoItVector)((G4int)ci);
219 G4cout << " " << ci + 1 << ") ";
220 if (ptProcManager != nullptr) {
221 G4cout << ptProcManager->GetProcessName() << G4endl;
222 }
223 }
224
225 ShowStep();
226 G4cout << G4endl;
227 G4cout << " ++List of secondaries generated "
228 << "(x,y,z,kE,t,PID):"
229 << " No. of secondaries = " << (*fSecondary).size() << G4endl;
230
231 if (! (*fSecondary).empty()) {
232 for (auto& lp1 : *fSecondary) {
233 G4cout << " " << std::setw(9) << G4BestUnit(lp1->GetPosition().x(), "Length") << " "
234 << std::setw(9) << G4BestUnit(lp1->GetPosition().y(), "Length") << " "
235 << std::setw(9) << G4BestUnit(lp1->GetPosition().z(), "Length") << " "
236 << std::setw(9) << G4BestUnit(lp1->GetKineticEnergy(), "Energy") << " "
237 << std::setw(9) << G4BestUnit(lp1->GetGlobalTime(), "Time") << " " << std::setw(18)
238 << lp1->GetDefinition()->GetParticleName() << G4endl;
239 }
240 }
241 }
242}
#define G4BestUnit(a, b)
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
const G4String & GetProcessName() const
G4TrackVector * fSecondary

◆ AlongStepDoItOneByOne()

void G4SteppingVerboseWithUnits::AlongStepDoItOneByOne ( )
overridevirtual

Reimplemented from G4SteppingVerbose.

Definition at line 377 of file G4SteppingVerboseWithUnits.cc.

378{
379 CopyState();
380
381 if (verboseLevel >= 4) {
382 G4cout << G4endl;
383 G4cout << " >>AlongStepDoIt (process by process): "
384 << " Process Name = " << fCurrentProcess->GetProcessName() << G4endl;
385
386 ShowStep();
387 G4cout << " "
388 << "!Note! Safety of PostStep is only valid "
389 << "after all DoIt invocations." << G4endl;
390
392 G4cout << G4endl;
393
394 G4cout << " ++List of secondaries generated "
395 << "(x,y,z,kE,t,PID):"
396 << " No. of secondaries = " << fN2ndariesAlongStepDoIt << G4endl;
397
398 if (fN2ndariesAlongStepDoIt > 0) {
399 for (std::size_t lp1 = (*fSecondary).size() - fN2ndariesAlongStepDoIt;
400 lp1 < (*fSecondary).size(); ++lp1)
401 {
402 G4cout << " " << std::setw(9)
403 << G4BestUnit((*fSecondary)[lp1]->GetPosition().x(), "Length") << " " << std::setw(9)
404 << G4BestUnit((*fSecondary)[lp1]->GetPosition().y(), "Length") << " " << std::setw(9)
405 << G4BestUnit((*fSecondary)[lp1]->GetPosition().z(), "Length") << " " << std::setw(9)
406 << G4BestUnit((*fSecondary)[lp1]->GetKineticEnergy(), "Energy") << " "
407 << std::setw(9) << G4BestUnit((*fSecondary)[lp1]->GetGlobalTime(), "Time") << " "
408 << std::setw(18) << (*fSecondary)[lp1]->GetDefinition()->GetParticleName() << G4endl;
409 }
410 }
411 }
412}

◆ AtRestDoItInvoked()

void G4SteppingVerboseWithUnits::AtRestDoItInvoked ( )
overridevirtual

Reimplemented from G4SteppingVerbose.

Definition at line 156 of file G4SteppingVerboseWithUnits.cc.

157{
158 G4VProcess* ptProcManager;
159 CopyState();
160
161 if (verboseLevel >= 3) {
162 G4int npt = 0;
163 G4cout << " **List of AtRestDoIt invoked:" << G4endl;
164 for (std::size_t np = 0; np < MAXofAtRestLoops; ++np) {
165 std::size_t npGPIL = MAXofAtRestLoops - np - 1;
167 ++npt;
168 ptProcManager = (*fAtRestDoItVector)[(G4int)np];
169 G4cout << " # " << npt << " : " << ptProcManager->GetProcessName() << " (Forced)"
170 << G4endl;
171 }
173 ++npt;
174 ptProcManager = (*fAtRestDoItVector)[(G4int)np];
175 G4cout << " # " << npt << " : " << ptProcManager->GetProcessName() << G4endl;
176 }
177 }
178
179 G4cout << " Generated secondaries = " << fN2ndariesAtRestDoIt << G4endl;
180
181 if (fN2ndariesAtRestDoIt > 0) {
182 G4cout << " -- List of secondaries generated : "
183 << "(x,y,z,kE,t,PID) --" << G4endl;
184 for (std::size_t lp1 = (*fSecondary).size() - fN2ndariesAtRestDoIt;
185 lp1 < (*fSecondary).size(); ++lp1)
186 {
187 G4cout << " " << std::setw(9)
188 << G4BestUnit((*fSecondary)[lp1]->GetPosition().x(), "Length") << " " << std::setw(9)
189 << G4BestUnit((*fSecondary)[lp1]->GetPosition().y(), "Length") << " " << std::setw(9)
190 << G4BestUnit((*fSecondary)[lp1]->GetPosition().z(), "Length") << " " << std::setw(9)
191 << G4BestUnit((*fSecondary)[lp1]->GetKineticEnergy(), "Energy") << " "
192 << std::setw(9) << G4BestUnit((*fSecondary)[lp1]->GetGlobalTime(), "Time") << " "
193 << std::setw(18) << (*fSecondary)[lp1]->GetDefinition()->GetParticleName() << G4endl;
194 }
195 }
196 }
197
198 if (verboseLevel >= 4) {
199 ShowStep();
200 G4cout << G4endl;
201 }
202}
@ NotForced
@ Forced
G4SelectedAtRestDoItVector * fSelectedAtRestDoItVector

◆ Clone()

G4VSteppingVerbose * G4SteppingVerboseWithUnits::Clone ( )
inlineoverridevirtual

Reimplemented from G4SteppingVerbose.

Definition at line 54 of file G4SteppingVerboseWithUnits.hh.

54{ return new G4SteppingVerboseWithUnits; }

◆ DPSLAlongStep()

void G4SteppingVerboseWithUnits::DPSLAlongStep ( )
overridevirtual

Reimplemented from G4SteppingVerbose.

Definition at line 355 of file G4SteppingVerboseWithUnits.cc.

356{
357 CopyState();
358
359 if (verboseLevel > 5) {
360 G4cout << " ++ProposedStep(AlongStep) = " << std::setw(9)
361 << G4BestUnit(physIntLength, "Length")
362 << " : ProcName = " << fCurrentProcess->GetProcessName() << " (";
364 G4cout << "CandidateForSelection)" << G4endl;
365 }
367 G4cout << "NotCandidateForSelection)" << G4endl;
368 }
369 else {
370 G4cout << "?!?)" << G4endl;
371 }
372 }
373}
@ CandidateForSelection
@ NotCandidateForSelection
G4GPILSelection fGPILSelection

◆ DPSLPostStep()

void G4SteppingVerboseWithUnits::DPSLPostStep ( )
overridevirtual

Reimplemented from G4SteppingVerbose.

Definition at line 327 of file G4SteppingVerboseWithUnits.cc.

328{
329 CopyState();
330
331 if (verboseLevel > 5) {
332 G4cout << " ++ProposedStep(PostStep ) = " << std::setw(9)
333 << G4BestUnit(physIntLength, "Length")
334 << " : ProcName = " << fCurrentProcess->GetProcessName() << " (";
336 G4cout << "ExclusivelyForced)" << G4endl;
337 }
338 else if (fCondition == StronglyForced) {
339 G4cout << "StronglyForced)" << G4endl;
340 }
341 else if (fCondition == Conditionally) {
342 G4cout << "Conditionally)" << G4endl;
343 }
344 else if (fCondition == Forced) {
345 G4cout << "Forced)" << G4endl;
346 }
347 else {
348 G4cout << "No ForceCondition)" << G4endl;
349 }
350 }
351}
@ StronglyForced
@ Conditionally
@ ExclusivelyForced
G4ForceCondition fCondition

◆ DPSLStarted()

void G4SteppingVerboseWithUnits::DPSLStarted ( )
overridevirtual

Reimplemented from G4SteppingVerbose.

Definition at line 301 of file G4SteppingVerboseWithUnits.cc.

302{
303 CopyState();
304
305 if (verboseLevel > 5) {
306 G4cout << G4endl << " >>DefinePhysicalStepLength (List of proposed StepLengths): " << G4endl;
307 }
308}

◆ DPSLUserLimit()

void G4SteppingVerboseWithUnits::DPSLUserLimit ( )
overridevirtual

Reimplemented from G4SteppingVerbose.

Definition at line 312 of file G4SteppingVerboseWithUnits.cc.

313{
314 CopyState();
315
316 if (verboseLevel > 5) {
317 G4cout << G4endl << G4endl;
318 G4cout << "=== Defined Physical Step Length (DPSL)" << G4endl;
319 G4cout << " ++ProposedStep(UserLimit) = " << std::setw(9)
320 << G4BestUnit(physIntLength, "Length")
321 << " : ProcName = User defined maximum allowed Step" << G4endl;
322 }
323}

◆ PostStepDoItAllDone()

void G4SteppingVerboseWithUnits::PostStepDoItAllDone ( )
overridevirtual

Reimplemented from G4SteppingVerbose.

Definition at line 246 of file G4SteppingVerboseWithUnits.cc.

247{
248 G4VProcess* ptProcManager;
249
250 CopyState();
251
252 if ((static_cast<int>(fStepStatus == fPostStepDoItProc) | static_cast<int>(fCondition == Forced) |
253 static_cast<int>(fCondition == Conditionally) |
254 static_cast<int>(fCondition == ExclusivelyForced) |
255 static_cast<int>(fCondition == StronglyForced)) != 0)
256 {
257 if (verboseLevel >= 3) {
258 G4int npt = 0;
259 G4cout << G4endl;
260 G4cout << " **PostStepDoIt (after all invocations):" << G4endl;
261 G4cout << " ++List of invoked processes " << G4endl;
262
263 for (std::size_t np = 0; np < MAXofPostStepLoops; ++np) {
264 std::size_t npGPIL = MAXofPostStepLoops - np - 1;
266 ++npt;
267 ptProcManager = (*fPostStepDoItVector)[(G4int)np];
268 G4cout << " " << npt << ") " << ptProcManager->GetProcessName() << " (Forced)"
269 << G4endl;
270 }
272 ++npt;
273 ptProcManager = (*fPostStepDoItVector)[(G4int)np];
274 G4cout << " " << npt << ") " << ptProcManager->GetProcessName() << G4endl;
275 }
276 }
277
278 ShowStep();
279 G4cout << G4endl;
280 G4cout << " ++List of secondaries generated "
281 << "(x,y,z,kE,t,PID):"
282 << " No. of secondaries = " << (*fSecondary).size() << G4endl;
283 G4cout << " [Note]Secondaries from AlongStepDoIt included." << G4endl;
284
285 if (! (*fSecondary).empty()) {
286 for (auto& lp1 : *fSecondary) {
287 G4cout << " " << std::setw(9) << G4BestUnit(lp1->GetPosition().x(), "Length") << " "
288 << std::setw(9) << G4BestUnit(lp1->GetPosition().y(), "Length") << " "
289 << std::setw(9) << G4BestUnit(lp1->GetPosition().z(), "Length") << " "
290 << std::setw(9) << G4BestUnit(lp1->GetKineticEnergy(), "Energy") << " "
291 << std::setw(9) << G4BestUnit(lp1->GetGlobalTime(), "Time") << " " << std::setw(18)
292 << lp1->GetDefinition()->GetParticleName() << G4endl;
293 }
294 }
295 }
296 }
297}
@ fPostStepDoItProc
G4SelectedPostStepDoItVector * fSelectedPostStepDoItVector

◆ PostStepDoItOneByOne()

void G4SteppingVerboseWithUnits::PostStepDoItOneByOne ( )
overridevirtual

Reimplemented from G4SteppingVerbose.

Definition at line 416 of file G4SteppingVerboseWithUnits.cc.

417{
418 CopyState();
419
420 if (verboseLevel >= 4) {
421 G4cout << G4endl;
422 G4cout << " >>PostStepDoIt (process by process): "
423 << " Process Name = " << fCurrentProcess->GetProcessName() << G4endl;
424
425 ShowStep();
426 G4cout << G4endl;
428 G4cout << G4endl;
429
430 G4cout << " ++List of secondaries generated "
431 << "(x,y,z,kE,t,PID):"
432 << " No. of secondaries = " << fN2ndariesPostStepDoIt << G4endl;
433
434 if (fN2ndariesPostStepDoIt > 0) {
435 for (std::size_t lp1 = (*fSecondary).size() - fN2ndariesPostStepDoIt;
436 lp1 < (*fSecondary).size(); ++lp1)
437 {
438 G4cout << " " << std::setw(9)
439 << G4BestUnit((*fSecondary)[lp1]->GetPosition().x(), "Length") << " " << std::setw(9)
440 << G4BestUnit((*fSecondary)[lp1]->GetPosition().y(), "Length") << " " << std::setw(9)
441 << G4BestUnit((*fSecondary)[lp1]->GetPosition().z(), "Length") << " " << std::setw(9)
442 << G4BestUnit((*fSecondary)[lp1]->GetKineticEnergy(), "Energy") << " "
443 << std::setw(9) << G4BestUnit((*fSecondary)[lp1]->GetGlobalTime(), "Time") << " "
444 << std::setw(18) << (*fSecondary)[lp1]->GetDefinition()->GetParticleName() << G4endl;
445 }
446 }
447 }
448}

◆ SetManager()

void G4SteppingVerboseWithUnits::SetManager ( G4SteppingManager * const fMan)
overridevirtual

Reimplemented from G4VSteppingVerbose.

Definition at line 54 of file G4SteppingVerboseWithUnits.cc.

55{
56 fManager = fMan;
57 fmessenger = new G4GenericMessenger(this, "/tracking/", "precision of verbose output");
58 auto& cmd =
59 fmessenger->DeclareProperty("setVerbosePrecision", fprec, "set precision of verbose output");
61}
@ G4State_Idle
@ G4State_PreInit
Command & DeclareProperty(const G4String &name, const G4AnyType &variable, const G4String &doc="")
G4SteppingManager * fManager
Command & SetStates(G4ApplicationState s0)

◆ ShowStep()

void G4SteppingVerboseWithUnits::ShowStep ( ) const
overridevirtual

Reimplemented from G4SteppingVerbose.

Definition at line 551 of file G4SteppingVerboseWithUnits.cc.

552{
553 // Show header
554 G4cout << G4endl;
555 G4cout << " ++G4Step Information " << G4endl;
556 G4long oldprc = G4cout.precision(fprec);
557
558 // Show G4Step specific information
559 G4cout << " Address of G4Track : " << fStep->GetTrack() << G4endl;
560 G4cout << " Step Length : "
561 << G4BestUnit(fStep->GetTrack()->GetStepLength(), "Length") << G4endl;
562 G4cout << " Energy Deposit : " << G4BestUnit(fStep->GetTotalEnergyDeposit(), "Energy")
563 << G4endl;
564
565 // Show G4StepPoint specific information
566 G4cout << " -------------------------------------------------------"
567 << "----------------" << G4endl;
568 G4cout << " StepPoint Information " << std::setw(20) << "PreStep" << std::setw(20)
569 << "PostStep" << G4endl;
570 G4cout << " -------------------------------------------------------"
571 << "----------------" << G4endl;
572 G4cout << " Position - x : " << std::setw(17)
573 << G4BestUnit(fStep->GetPreStepPoint()->GetPosition().x(), "Length") << std::setw(17)
574 << G4BestUnit(fStep->GetPostStepPoint()->GetPosition().x(), "Length") << G4endl;
575 G4cout << " Position - y : " << std::setw(17)
576 << G4BestUnit(fStep->GetPreStepPoint()->GetPosition().y(), "Length") << std::setw(17)
577 << G4BestUnit(fStep->GetPostStepPoint()->GetPosition().y(), "Length") << G4endl;
578 G4cout << " Position - z : " << std::setw(17)
579 << G4BestUnit(fStep->GetPreStepPoint()->GetPosition().z(), "Length") << std::setw(17)
580 << G4BestUnit(fStep->GetPostStepPoint()->GetPosition().z(), "Length") << G4endl;
581 G4cout << " Global Time : " << std::setw(17)
582 << G4BestUnit(fStep->GetPreStepPoint()->GetGlobalTime(), "Time") << std::setw(17)
584 G4cout << " Local Time : " << std::setw(17)
585 << G4BestUnit(fStep->GetPreStepPoint()->GetLocalTime(), "Time") << std::setw(17)
587 G4cout << " Proper Time : " << std::setw(17)
588 << G4BestUnit(fStep->GetPreStepPoint()->GetProperTime(), "Time") << std::setw(17)
590 G4cout << " Momentum Direct - x : " << std::setw(20)
591 << fStep->GetPreStepPoint()->GetMomentumDirection().x() << std::setw(20)
593 G4cout << " Momentum Direct - y : " << std::setw(20)
594 << fStep->GetPreStepPoint()->GetMomentumDirection().y() << std::setw(20)
596 G4cout << " Momentum Direct - z : " << std::setw(20)
597 << fStep->GetPreStepPoint()->GetMomentumDirection().z() << std::setw(20)
599 G4cout << " Momentum - x : " << std::setw(14)
600 << G4BestUnit(fStep->GetPreStepPoint()->GetMomentum().x(), "Momentum") << std::setw(14)
601 << G4BestUnit(fStep->GetPostStepPoint()->GetMomentum().x(), "Momentum") << G4endl;
602 G4cout << " Momentum - y : " << std::setw(14)
603 << G4BestUnit(fStep->GetPreStepPoint()->GetMomentum().y(), "Momentum") << std::setw(14)
604 << G4BestUnit(fStep->GetPostStepPoint()->GetMomentum().y(), "Momentum") << G4endl;
605 G4cout << " Momentum - z : " << std::setw(14)
606 << G4BestUnit(fStep->GetPreStepPoint()->GetMomentum().z(), "Momentum") << std::setw(14)
607 << G4BestUnit(fStep->GetPostStepPoint()->GetMomentum().z(), "Momentum") << G4endl;
608 G4cout << " Total Energy : " << std::setw(16)
609 << G4BestUnit(fStep->GetPreStepPoint()->GetTotalEnergy(), "Energy") << std::setw(16)
611 G4cout << " Kinetic Energy : " << std::setw(16)
612 << G4BestUnit(fStep->GetPreStepPoint()->GetKineticEnergy(), "Energy") << std::setw(16)
614 G4cout << " Velocity : " << std::setw(14)
615 << G4BestUnit(fStep->GetPreStepPoint()->GetVelocity(), "Velocity") << std::setw(14)
616 << G4BestUnit(fStep->GetPostStepPoint()->GetVelocity(), "Velocity") << G4endl;
617 G4cout << " Volume Name : " << std::setw(20)
619 G4String volName = "OutOfWorld";
620 if (fStep->GetPostStepPoint()->GetPhysicalVolume() != nullptr)
622 G4cout << std::setw(20) << volName << G4endl;
623 G4cout << " Safety : " << std::setw(17)
624 << G4BestUnit(fStep->GetPreStepPoint()->GetSafety(), "Length") << std::setw(17)
625 << G4BestUnit(fStep->GetPostStepPoint()->GetSafety(), "Length") << G4endl;
626 G4cout << " Polarization - x : " << std::setw(20)
627 << fStep->GetPreStepPoint()->GetPolarization().x() << std::setw(20)
629 G4cout << " Polarization - y : " << std::setw(20)
630 << fStep->GetPreStepPoint()->GetPolarization().y() << std::setw(20)
632 G4cout << " Polarization - Z : " << std::setw(20)
633 << fStep->GetPreStepPoint()->GetPolarization().z() << std::setw(20)
635 G4cout << " Weight : " << std::setw(20)
636 << fStep->GetPreStepPoint()->GetWeight() << std::setw(20)
638 G4cout << " Step Status : ";
640 if (tStepStatus == fGeomBoundary) {
641 G4cout << std::setw(20) << "Geom Limit";
642 }
643 else if (tStepStatus == fAlongStepDoItProc) {
644 G4cout << std::setw(20) << "AlongStep Proc.";
645 }
646 else if (tStepStatus == fPostStepDoItProc) {
647 G4cout << std::setw(20) << "PostStep Proc";
648 }
649 else if (tStepStatus == fAtRestDoItProc) {
650 G4cout << std::setw(20) << "AtRest Proc";
651 }
652 else if (tStepStatus == fUndefined) {
653 G4cout << std::setw(20) << "Undefined";
654 }
655
656 tStepStatus = fStep->GetPostStepPoint()->GetStepStatus();
657 if (tStepStatus == fGeomBoundary) {
658 G4cout << std::setw(20) << "Geom Limit";
659 }
660 else if (tStepStatus == fAlongStepDoItProc) {
661 G4cout << std::setw(20) << "AlongStep Proc.";
662 }
663 else if (tStepStatus == fPostStepDoItProc) {
664 G4cout << std::setw(20) << "PostStep Proc";
665 }
666 else if (tStepStatus == fAtRestDoItProc) {
667 G4cout << std::setw(20) << "AtRest Proc";
668 }
669 else if (tStepStatus == fUndefined) {
670 G4cout << std::setw(20) << "Undefined";
671 }
672
673 G4cout << G4endl;
674 G4cout << " Process defined Step: ";
675 if (fStep->GetPreStepPoint()->GetProcessDefinedStep() == nullptr) {
676 G4cout << std::setw(20) << "Undefined";
677 }
678 else {
680 }
681 if (fStep->GetPostStepPoint()->GetProcessDefinedStep() == nullptr) {
682 G4cout << std::setw(20) << "Undefined";
683 }
684 else {
686 }
687 G4cout.precision(oldprc);
688
689 G4cout << G4endl;
690 G4cout << " -------------------------------------------------------"
691 << "----------------" << G4endl;
692}
G4StepStatus
@ fGeomBoundary
@ fUndefined
@ fAtRestDoItProc
@ fAlongStepDoItProc
long G4long
Definition G4Types.hh:87
double z() const
double x() const
double y() const
G4double GetTotalEnergy() const
G4StepStatus GetStepStatus() const
G4double GetVelocity() const
G4double GetProperTime() const
G4double GetGlobalTime() const
G4double GetSafety() const
const G4VProcess * GetProcessDefinedStep() const
G4ThreeVector GetMomentum() const
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetLocalTime() const
G4VPhysicalVolume * GetPhysicalVolume() const
const G4ThreeVector & GetPolarization() const
G4double GetKineticEnergy() const
G4double GetWeight() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetTotalEnergyDeposit() const
G4StepPoint * GetPostStepPoint() const
G4double GetStepLength() const
const G4String & GetName() const

Referenced by AlongStepDoItAllDone(), AlongStepDoItOneByOne(), AtRestDoItInvoked(), PostStepDoItAllDone(), and PostStepDoItOneByOne().

◆ StepInfo()

void G4SteppingVerboseWithUnits::StepInfo ( )
overridevirtual

Reimplemented from G4SteppingVerbose.

Definition at line 98 of file G4SteppingVerboseWithUnits.cc.

99{
100 CopyState();
101 G4long oldprec = G4cout.precision(fprec);
102
103 if (verboseLevel >= 1) {
104 if (verboseLevel >= 4) VerboseTrack();
105 if (verboseLevel >= 3) {
106 G4cout << G4endl;
107 G4cout << std::setw(5) << "#Step#"
108 << " " << std::setw(fprec + 3) << "X"
109 << " " << std::setw(fprec + 3) << "Y"
110 << " " << std::setw(fprec + 3) << "Z"
111 << " " << std::setw(fprec + 6) << "KineE"
112 << " " << std::setw(fprec + 10) << "dEStep"
113 << " " << std::setw(fprec + 7) << "StepLeng" << std::setw(fprec + 7) << "TrakLeng"
114 << std::setw(10) << "Volume"
115 << " " << std::setw(10) << "Process" << G4endl;
116 }
117
118 G4cout << std::setw(5) << fTrack->GetCurrentStepNumber() << " " << std::setw(fprec + 3)
119 << G4BestUnit(fTrack->GetPosition().x(), "Length") << std::setw(fprec + 3)
120 << G4BestUnit(fTrack->GetPosition().y(), "Length") << std::setw(fprec + 3)
121 << G4BestUnit(fTrack->GetPosition().z(), "Length") << std::setw(fprec + 3)
122 << G4BestUnit(fTrack->GetKineticEnergy(), "Energy") << std::setw(fprec + 7)
123 << G4BestUnit(fStep->GetTotalEnergyDeposit(), "Energy") << std::setw(fprec + 3)
124 << G4BestUnit(fStep->GetStepLength(), "Length") << std::setw(fprec + 3)
125 << G4BestUnit(fTrack->GetTrackLength(), "Length") << std::setw(10)
126 << fTrack->GetVolume()->GetName();
127
129 G4String procName = " UserLimit";
130 if (process != nullptr) procName = process->GetProcessName();
131 if (fStepStatus == fWorldBoundary) procName = "OutOfWorld";
132 G4cout << " " << std::setw(9) << procName;
133 G4cout << G4endl;
134
135 if (verboseLevel == 2) {
136 const std::vector<const G4Track*>* secondary = fStep->GetSecondaryInCurrentStep();
137 if (! secondary->empty()) {
138 G4cout << "\n :----- List of secondaries ----------------" << G4endl;
139 G4cout.precision(4);
140 for (auto lp : *secondary) {
141 G4cout << " " << std::setw(13) << lp->GetDefinition()->GetParticleName()
142 << ": energy =" << std::setw(6) << G4BestUnit(lp->GetKineticEnergy(), "Energy")
143 << " time =" << std::setw(6) << G4BestUnit(lp->GetGlobalTime(), "Time");
144 G4cout << G4endl;
145 }
146
147 G4cout << " :------------------------------------------\n" << G4endl;
148 }
149 }
150 }
151 G4cout.precision(oldprec);
152}
@ fWorldBoundary
const std::vector< const G4Track * > * GetSecondaryInCurrentStep() const
Definition G4Step.cc:202
G4double GetStepLength() const
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
G4double GetTrackLength() const
G4int GetCurrentStepNumber() const
G4double GetKineticEnergy() const

◆ TrackingStarted()

void G4SteppingVerboseWithUnits::TrackingStarted ( )
overridevirtual

Reimplemented from G4SteppingVerbose.

Definition at line 65 of file G4SteppingVerboseWithUnits.cc.

66{
67 CopyState();
68 G4long oldprec = G4cout.precision(fprec);
69
70 // Step zero
71 //
72 if (verboseLevel > 0) {
73 G4cout << std::setw(5) << "Step#"
74 << " " << std::setw(fprec + 3) << "X"
75 << " " << std::setw(fprec + 3) << "Y"
76 << " " << std::setw(fprec + 3) << "Z"
77 << " " << std::setw(fprec + 6) << "KineE"
78 << " " << std::setw(fprec + 10) << "dEStep"
79 << " " << std::setw(fprec + 7) << "StepLeng" << std::setw(fprec + 7) << "TrakLeng"
80 << std::setw(10) << "Volume"
81 << " " << std::setw(10) << "Process" << G4endl;
82
83 G4cout << std::setw(5) << fTrack->GetCurrentStepNumber() << " " << std::setw(fprec + 3)
84 << G4BestUnit(fTrack->GetPosition().x(), "Length") << std::setw(fprec + 3)
85 << G4BestUnit(fTrack->GetPosition().y(), "Length") << std::setw(fprec + 3)
86 << G4BestUnit(fTrack->GetPosition().z(), "Length") << std::setw(fprec + 3)
87 << G4BestUnit(fTrack->GetKineticEnergy(), "Energy") << std::setw(fprec + 7)
88 << G4BestUnit(fStep->GetTotalEnergyDeposit(), "Energy") << std::setw(fprec + 3)
89 << G4BestUnit(fStep->GetStepLength(), "Length") << std::setw(fprec + 3)
90 << G4BestUnit(fTrack->GetTrackLength(), "Length") << std::setw(10)
91 << fTrack->GetVolume()->GetName() << std::setw(9) << " initStep" << G4endl;
92 }
93 G4cout.precision(oldprec);
94}

◆ VerboseParticleChange()

void G4SteppingVerboseWithUnits::VerboseParticleChange ( )
overridevirtual

Reimplemented from G4SteppingVerbose.

Definition at line 542 of file G4SteppingVerboseWithUnits.cc.

543{
544 G4cout << G4endl;
545 G4cout << " ++G4ParticleChange Information " << G4endl;
547}
virtual void DumpInfo() const
G4VParticleChange * fParticleChange

Referenced by AlongStepDoItOneByOne(), and PostStepDoItOneByOne().

◆ VerboseTrack()

void G4SteppingVerboseWithUnits::VerboseTrack ( )
overridevirtual

Reimplemented from G4SteppingVerbose.

Definition at line 452 of file G4SteppingVerboseWithUnits.cc.

453{
454 CopyState();
455
456 G4cout << G4endl;
457 G4cout << " ++G4Track Information " << G4endl;
458 G4long oldprec = G4cout.precision(fprec);
459
460 G4cout << " -----------------------------------------------" << G4endl;
461 G4cout << " G4Track Information " << std::setw(20) << G4endl;
462 G4cout << " -----------------------------------------------" << G4endl;
463
464 G4cout << " Step number : " << std::setw(20) << fTrack->GetCurrentStepNumber()
465 << G4endl;
466 G4cout << " Position - x : " << std::setw(20)
467 << G4BestUnit(fTrack->GetPosition().x(), "Length") << G4endl;
468 G4cout << " Position - y : " << std::setw(20)
469 << G4BestUnit(fTrack->GetPosition().y(), "Length") << G4endl;
470 G4cout << " Position - z : " << std::setw(20)
471 << G4BestUnit(fTrack->GetPosition().z(), "Length") << G4endl;
472 G4cout << " Global Time : " << std::setw(20)
473 << G4BestUnit(fTrack->GetGlobalTime(), "Time") << G4endl;
474 G4cout << " Local Time : " << std::setw(20)
475 << G4BestUnit(fTrack->GetLocalTime(), "Time") << G4endl;
476 G4cout << " Momentum Direct - x : " << std::setw(20) << fTrack->GetMomentumDirection().x()
477 << G4endl;
478 G4cout << " Momentum Direct - y : " << std::setw(20) << fTrack->GetMomentumDirection().y()
479 << G4endl;
480 G4cout << " Momentum Direct - z : " << std::setw(20) << fTrack->GetMomentumDirection().z()
481 << G4endl;
482 G4cout << " Kinetic Energy : " << std::setw(20)
483 << G4BestUnit(fTrack->GetKineticEnergy(), "Energy") << G4endl;
484 G4cout << " Polarization - x : " << std::setw(20) << fTrack->GetPolarization().x()
485 << G4endl;
486 G4cout << " Polarization - y : " << std::setw(20) << fTrack->GetPolarization().y()
487 << G4endl;
488 G4cout << " Polarization - z : " << std::setw(20) << fTrack->GetPolarization().z()
489 << G4endl;
490 G4cout << " Track Length : " << std::setw(20)
491 << G4BestUnit(fTrack->GetTrackLength(), "Length") << G4endl;
492 G4cout << " Track ID # : " << std::setw(20) << fTrack->GetTrackID() << G4endl;
493 G4cout << " Parent Track ID # : " << std::setw(20) << fTrack->GetParentID() << G4endl;
494 G4cout << " Next Volume : " << std::setw(20);
495 if (fTrack->GetNextVolume() != nullptr)
496 G4cout << fTrack->GetNextVolume()->GetName() << " ";
497 else
498 G4cout << "OutOfWorld"
499 << " ";
500 G4cout << G4endl;
501 G4cout << " Track Status : " << std::setw(20);
502 if (fTrack->GetTrackStatus() == fAlive)
503 G4cout << " Alive";
504 else if (fTrack->GetTrackStatus() == fStopButAlive)
505 G4cout << " StopButAlive";
506 else if (fTrack->GetTrackStatus() == fStopAndKill)
507 G4cout << " StopAndKill";
509 G4cout << " KillTrackAndSecondaries";
510 else if (fTrack->GetTrackStatus() == fSuspend)
511 G4cout << " Suspend";
513 G4cout << " PostponeToNextEvent";
514 G4cout << G4endl;
515 G4cout << " Vertex - x : " << std::setw(20)
516 << G4BestUnit(fTrack->GetVertexPosition().x(), "Length") << G4endl;
517 G4cout << " Vertex - y : " << std::setw(20)
518 << G4BestUnit(fTrack->GetVertexPosition().y(), "Length") << G4endl;
519 G4cout << " Vertex - z : " << std::setw(20)
520 << G4BestUnit(fTrack->GetVertexPosition().z(), "Length") << G4endl;
521 G4cout << " Vertex - Px (MomDir): " << std::setw(20)
523 G4cout << " Vertex - Py (MomDir): " << std::setw(20)
525 G4cout << " Vertex - Pz (MomDir): " << std::setw(20)
527 G4cout << " Vertex - KineE : " << std::setw(20)
529
530 G4cout << " Creator Process : " << std::setw(20);
531 if (fTrack->GetCreatorProcess() == nullptr)
532 G4cout << " Event Generator" << G4endl;
533 else
535
536 G4cout << " -----------------------------------------------" << G4endl;
537 G4cout.precision(oldprec);
538}
@ fKillTrackAndSecondaries
@ fSuspend
@ fAlive
@ fStopAndKill
@ fStopButAlive
@ fPostponeToNextEvent
G4TrackStatus GetTrackStatus() const
G4double GetVertexKineticEnergy() const
G4int GetTrackID() const
const G4VProcess * GetCreatorProcess() const
G4VPhysicalVolume * GetNextVolume() const
const G4ThreeVector & GetVertexMomentumDirection() const
G4double GetGlobalTime() const
const G4ThreeVector & GetVertexPosition() const
G4double GetLocalTime() const
const G4ThreeVector & GetMomentumDirection() const
const G4ThreeVector & GetPolarization() const
G4int GetParentID() const

Referenced by StepInfo().


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