Geant4 VMC Version 6.6
Loading...
Searching...
No Matches
TG4StackPopper.h
Go to the documentation of this file.
1#ifndef TG4_STACK_POPPER_H
2#define TG4_STACK_POPPER_H
3
4//------------------------------------------------
5// The Geant4 Virtual Monte Carlo package
6// Copyright (C) 2007 - 2014 Ivana Hrivnacova
7// All rights reserved.
8//
9// For the licensing terms see geant4_vmc/LICENSE.
10// Contact: root-vmc@cern.ch
11//-------------------------------------------------
12
17
18#include <G4VProcess.hh>
19
20class TVirtualMCStack;
21
22class G4Track;
23
28
33
39
41{
42 public:
43 TG4StackPopper(const G4String& processName = "stackPopper");
44 virtual ~TG4StackPopper();
45
46 // static access method
47 static TG4StackPopper* Instance();
48
49 // methods
50 virtual G4bool IsApplicable(
51 const G4ParticleDefinition& /*particleDefinition*/);
52
53 virtual G4double PostStepGetPhysicalInteractionLength(const G4Track& track,
54 G4double previousStepSize, G4ForceCondition* condition);
55
56 virtual G4VParticleChange* PostStepDoIt(
57 const G4Track& track, const G4Step& step);
58
59 // No operation in AlongStepDoIt and AtRestDoIt
60
62 const G4Track& /*track*/, G4double /*previousStepSize*/,
63 G4double /*currentMinimumStep*/, G4double& /*proposedSafety*/,
64 G4GPILSelection* /*selection*/)
65 {
66 return -1.0;
67 }
68
70 const G4Track& /*track*/, G4ForceCondition* /*condition*/)
71 {
72 return -1.0;
73 }
74
75 virtual G4VParticleChange* AlongStepDoIt(
76 const G4Track& /*track*/, const G4Step& /*step*/)
77 {
78 return 0;
79 }
80
81 virtual G4VParticleChange* AtRestDoIt(
82 const G4Track& /*track*/, const G4Step& /*step*/)
83 {
84 return 0;
85 }
86
87 void Notify();
88 void Reset();
89 void SetMCStack(TVirtualMCStack* mcStack);
90 void SetDoExclusiveStep(G4TrackStatus trackStatus);
91
92 G4bool HasPoppedTracks() const;
93
94 private:
99
101 static G4ThreadLocal TG4StackPopper* fgInstance;
102
104 TVirtualMCStack* fMCStack;
105
108
114
116 G4TrackStatus fTrackStatus;
117};
118
119// inline methods
120
122{
124 return fgInstance;
125}
126
128 const G4ParticleDefinition& /*particleDefinition*/)
129{
131 return true;
132}
133
134inline void TG4StackPopper::SetMCStack(TVirtualMCStack* mcStack)
135{
137 fMCStack = mcStack;
138}
139
140#endif // TG4_STACK_POPPER_H
The process which pops particles defined by user from the VMC stack and passes them to tracking.
G4bool HasPoppedTracks() const
void SetDoExclusiveStep(G4TrackStatus trackStatus)
static G4ThreadLocal TG4StackPopper * fgInstance
this instance
G4TrackStatus fTrackStatus
The track status to be restored after performing exclusive step.
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
static TG4StackPopper * Instance()
TG4StackPopper(const G4String &processName="stackPopper")
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step)
void SetMCStack(TVirtualMCStack *mcStack)
TG4StackPopper(const TG4StackPopper &right)
Not implemented.
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
TVirtualMCStack * fMCStack
Cached pointer to thread-local VMC stack.
TG4StackPopper & operator=(const TG4StackPopper &right)
Not implemented.
virtual G4bool IsApplicable(const G4ParticleDefinition &)
G4int fNofDoneTracks
the counter for popped tracks
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual ~TG4StackPopper()
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
TG4StackPopperProcessType
Definition of Stack popper process type.
@ fStackPopper