Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4Step Class Reference

#include <G4Step.hh>

Public Member Functions

 G4Step ()
 
 ~G4Step ()
 
 G4Step (const G4Step &)
 
G4Stepoperator= (const G4Step &)
 
G4TrackGetTrack () const
 
void SetTrack (G4Track *value)
 
G4StepPointGetPreStepPoint () const
 
void SetPreStepPoint (G4StepPoint *value)
 
G4StepPointGetPostStepPoint () const
 
void SetPostStepPoint (G4StepPoint *value)
 
G4double GetStepLength () const
 
void SetStepLength (G4double value)
 
G4double GetTotalEnergyDeposit () const
 
void SetTotalEnergyDeposit (G4double value)
 
G4double GetNonIonizingEnergyDeposit () const
 
void SetNonIonizingEnergyDeposit (G4double value)
 
G4SteppingControl GetControlFlag () const
 
void SetControlFlag (G4SteppingControl StepControlFlag)
 
void AddTotalEnergyDeposit (G4double value)
 
void ResetTotalEnergyDeposit ()
 
void AddNonIonizingEnergyDeposit (G4double value)
 
void ResetNonIonizingEnergyDeposit ()
 
G4bool IsFirstStepInVolume () const
 
G4bool IsLastStepInVolume () const
 
void SetFirstStepFlag ()
 
void ClearFirstStepFlag ()
 
void SetLastStepFlag ()
 
void ClearLastStepFlag ()
 
G4ThreeVector GetDeltaPosition () const
 
G4double GetDeltaTime () const
 
G4ThreeVector GetDeltaMomentum () const
 
G4double GetDeltaEnergy () const
 
void InitializeStep (G4Track *aValue)
 
void UpdateTrack ()
 
void CopyPostToPreStepPoint ()
 
G4PolylineCreatePolyline () const
 
G4int GetNumberOfSecondariesInCurrentStep () const
 
const std::vector< const
G4Track * > * 
GetSecondaryInCurrentStep () const
 
const G4TrackVectorGetSecondary () const
 
G4TrackVectorGetfSecondary ()
 
G4TrackVectorNewSecondaryVector ()
 
void DeleteSecondaryVector ()
 
void SetSecondary (G4TrackVector *value)
 
void SetPointerToVectorOfAuxiliaryPoints (std::vector< G4ThreeVector > *theNewVectorPointer)
 
std::vector< G4ThreeVector > * GetPointerToVectorOfAuxiliaryPoints () const
 

Protected Attributes

G4double fTotalEnergyDeposit
 
G4double fNonIonizingEnergyDeposit
 

Detailed Description

Definition at line 76 of file G4Step.hh.

Constructor & Destructor Documentation

G4Step::G4Step ( )

Definition at line 53 of file G4Step.cc.

56  fStepLength(0.), fpTrack(0),
57  fpSteppingControlFlag(NormalCondition),
58  fFirstStepInVolume(false),
59  fLastStepInVolume(false),
60  fSecondary(0),
61  nSecondaryByLastStep(0), secondaryInCurrentStep(),
62  fpVectorOfAuxiliaryPointsPointer(0)
63 {
64  fpPreStepPoint = new G4StepPoint();
65  fpPostStepPoint = new G4StepPoint();
66 
67  secondaryInCurrentStep = new std::vector<CT>;
68 }
G4double fTotalEnergyDeposit
Definition: G4Step.hh:176
G4double fNonIonizingEnergyDeposit
Definition: G4Step.hh:179
G4Step::~G4Step ( )

Definition at line 71 of file G4Step.cc.

73 {
74  delete fpPreStepPoint;
75  delete fpPostStepPoint;
76 
77  secondaryInCurrentStep->clear();
78  delete secondaryInCurrentStep;
79 
80  if (fSecondary !=0 ) {
81  fSecondary->clear();
82  delete fSecondary;
83  }
84 }
G4Step::G4Step ( const G4Step right)

Definition at line 89 of file G4Step.cc.

93  fStepLength(right.fStepLength),
94  fpTrack(right.fpTrack),
95  fpSteppingControlFlag(right.fpSteppingControlFlag),
96  fFirstStepInVolume(right.fFirstStepInVolume),
97  fLastStepInVolume(right.fLastStepInVolume),
98  nSecondaryByLastStep(right.nSecondaryByLastStep),
99  secondaryInCurrentStep(right.secondaryInCurrentStep),
100  fpVectorOfAuxiliaryPointsPointer(right.fpVectorOfAuxiliaryPointsPointer)
101 {
102  if (right.fpPreStepPoint !=0) {
103  fpPreStepPoint = new G4StepPoint(*(right.fpPreStepPoint));
104  } else {
105  fpPreStepPoint = new G4StepPoint();
106  }
107  if (right.fpPostStepPoint !=0) {
108  fpPostStepPoint = new G4StepPoint(*(right.fpPostStepPoint));
109  } else {
110  fpPostStepPoint = new G4StepPoint();
111  }
112 
113  if (right.fSecondary !=0) {
114  fSecondary = new G4TrackVector(*(right.fSecondary));
115  } else {
116  fSecondary = new G4TrackVector();
117  }
118  secondaryInCurrentStep = new std::vector<CT>;
119 }
G4double fTotalEnergyDeposit
Definition: G4Step.hh:176
G4double fNonIonizingEnergyDeposit
Definition: G4Step.hh:179
std::vector< G4Track * > G4TrackVector

Member Function Documentation

void G4Step::AddNonIonizingEnergyDeposit ( G4double  value)

Here is the caller graph for this function:

void G4Step::AddTotalEnergyDeposit ( G4double  value)

Here is the caller graph for this function:

void G4Step::ClearFirstStepFlag ( )

Here is the caller graph for this function:

void G4Step::ClearLastStepFlag ( )

Here is the caller graph for this function:

void G4Step::CopyPostToPreStepPoint ( )

Here is the caller graph for this function:

G4Polyline* G4Step::CreatePolyline ( ) const
void G4Step::DeleteSecondaryVector ( )

Here is the caller graph for this function:

G4SteppingControl G4Step::GetControlFlag ( ) const

Here is the caller graph for this function:

G4double G4Step::GetDeltaEnergy ( ) const

Definition at line 176 of file G4Step.cc.

178 {
179  static G4ThreadLocal G4bool isFirstTime = true;
180  if (isFirstTime) {
181  isFirstTime = false;
182 #ifdef G4VERBOSE
183  G4Exception( "G4Step::GetDeltaEnergy()","Warning", JustWarning,
184  "This method is obsolete and will be removed soon");
185 #endif
186  }
187 
188  return fpPostStepPoint->GetKineticEnergy()
189  - fpPreStepPoint->GetKineticEnergy();
190 }
#define G4ThreadLocal
Definition: tls.hh:89
bool G4bool
Definition: G4Types.hh:79
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double GetKineticEnergy() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4ThreeVector G4Step::GetDeltaMomentum ( ) const

Definition at line 159 of file G4Step.cc.

161 {
162  static G4ThreadLocal G4bool isFirstTime = true;
163  if (isFirstTime) {
164  isFirstTime = false;
165 #ifdef G4VERBOSE
166  G4Exception( "G4Step::GetDeltaMomentum()","Warning", JustWarning,
167  "This method is obsolete and will be removed soon");
168 #endif
169  }
170 
171  return fpPostStepPoint->GetMomentum()
172  - fpPreStepPoint->GetMomentum();
173 }
G4ThreeVector GetMomentum() const
#define G4ThreadLocal
Definition: tls.hh:89
bool G4bool
Definition: G4Types.hh:79
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

Here is the caller graph for this function:

G4ThreeVector G4Step::GetDeltaPosition ( ) const

Here is the caller graph for this function:

G4double G4Step::GetDeltaTime ( ) const

Here is the caller graph for this function:

G4TrackVector* G4Step::GetfSecondary ( )

Here is the caller graph for this function:

G4double G4Step::GetNonIonizingEnergyDeposit ( ) const

Here is the caller graph for this function:

G4int G4Step::GetNumberOfSecondariesInCurrentStep ( ) const
std::vector<G4ThreeVector>* G4Step::GetPointerToVectorOfAuxiliaryPoints ( ) const
inline

Definition at line 242 of file G4Step.hh.

242  {
243  return fpVectorOfAuxiliaryPointsPointer;
244  }

Here is the caller graph for this function:

G4StepPoint* G4Step::GetPostStepPoint ( ) const
G4StepPoint* G4Step::GetPreStepPoint ( ) const
const G4TrackVector* G4Step::GetSecondary ( ) const

Here is the caller graph for this function:

const std::vector< const G4Track * > * G4Step::GetSecondaryInCurrentStep ( ) const

Definition at line 193 of file G4Step.cc.

195 {
196  secondaryInCurrentStep->clear();
197  G4int nSecondary = fSecondary->size();
198  for (G4int i=nSecondaryByLastStep; i<nSecondary; i++) {
199  secondaryInCurrentStep->push_back((*fSecondary)[i]);
200  }
201  return secondaryInCurrentStep;
202 }
int G4int
Definition: G4Types.hh:78

Here is the caller graph for this function:

G4double G4Step::GetStepLength ( ) const
G4double G4Step::GetTotalEnergyDeposit ( ) const
G4Track* G4Step::GetTrack ( ) const
void G4Step::InitializeStep ( G4Track aValue)

Here is the caller graph for this function:

G4bool G4Step::IsFirstStepInVolume ( ) const
G4bool G4Step::IsLastStepInVolume ( ) const
G4TrackVector* G4Step::NewSecondaryVector ( )

Here is the caller graph for this function:

G4Step & G4Step::operator= ( const G4Step right)

Definition at line 122 of file G4Step.cc.

124 {
125  if (this != &right){
128  fStepLength = right.fStepLength;
129  fpTrack = right.fpTrack;
130  fpSteppingControlFlag = right.fpSteppingControlFlag;
131  fFirstStepInVolume = right.fFirstStepInVolume;
132  fLastStepInVolume = right.fLastStepInVolume;
133  nSecondaryByLastStep = right.nSecondaryByLastStep;
134  secondaryInCurrentStep = right.secondaryInCurrentStep;
135  fpVectorOfAuxiliaryPointsPointer = right.fpVectorOfAuxiliaryPointsPointer;
136 
137  if (fpPreStepPoint !=0 ) delete fpPreStepPoint;
138  if (right.fpPreStepPoint !=0) {
139  fpPreStepPoint = new G4StepPoint(*(right.fpPreStepPoint));
140  } else {
141  fpPreStepPoint = new G4StepPoint();
142  }
143  if (fpPostStepPoint !=0 ) delete fpPostStepPoint;
144  if (right.fpPostStepPoint !=0) {
145  fpPostStepPoint = new G4StepPoint(*(right.fpPostStepPoint));
146  } else {
147  fpPostStepPoint = new G4StepPoint();
148  }
149  if (right.fSecondary !=0) {
150  fSecondary = new G4TrackVector(*(right.fSecondary));
151  } else {
152  fSecondary = new G4TrackVector();
153  }
154  }
155  return *this;
156 }
G4double fTotalEnergyDeposit
Definition: G4Step.hh:176
G4double fNonIonizingEnergyDeposit
Definition: G4Step.hh:179
std::vector< G4Track * > G4TrackVector
void G4Step::ResetNonIonizingEnergyDeposit ( )
void G4Step::ResetTotalEnergyDeposit ( )

Here is the caller graph for this function:

void G4Step::SetControlFlag ( G4SteppingControl  StepControlFlag)

Here is the caller graph for this function:

void G4Step::SetFirstStepFlag ( )

Here is the caller graph for this function:

void G4Step::SetLastStepFlag ( )

Here is the caller graph for this function:

void G4Step::SetNonIonizingEnergyDeposit ( G4double  value)

Here is the caller graph for this function:

void G4Step::SetPointerToVectorOfAuxiliaryPoints ( std::vector< G4ThreeVector > *  theNewVectorPointer)
inline

Definition at line 239 of file G4Step.hh.

239  {
240  fpVectorOfAuxiliaryPointsPointer = theNewVectorPointer;
241  }

Here is the caller graph for this function:

void G4Step::SetPostStepPoint ( G4StepPoint value)
void G4Step::SetPreStepPoint ( G4StepPoint value)
void G4Step::SetSecondary ( G4TrackVector value)
void G4Step::SetStepLength ( G4double  value)

Here is the caller graph for this function:

void G4Step::SetTotalEnergyDeposit ( G4double  value)

Here is the caller graph for this function:

void G4Step::SetTrack ( G4Track value)

Here is the caller graph for this function:

void G4Step::UpdateTrack ( )

Here is the caller graph for this function:

Member Data Documentation

G4double G4Step::fNonIonizingEnergyDeposit
protected

Definition at line 179 of file G4Step.hh.

G4double G4Step::fTotalEnergyDeposit
protected

Definition at line 176 of file G4Step.hh.


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