Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4Step.hh
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4Step.hh 93028 2015-09-30 16:09:00Z gcosmo $
28 //
29 //
30 //---------------------------------------------------------------
31 //
32 // G4Step.hh
33 //
34 // Class Description:
35 // This class represents the Step of a particle tracked.
36 // It includes information of
37 // 1) List of Step points which compose the Step,
38 // 2) static information of particle which generated the
39 // Step,
40 // 3) trackID and parent particle ID of the Step,
41 // 4) termination condition of the Step,
42 //
43 // Contact:
44 // Questions and comments to this code should be sent to
45 // Katsuya Amako (e-mail: Katsuya.Amako@kek.jp)
46 // Takashi Sasaki (e-mail: Takashi.Sasaki@kek.jp)
47 //
48 // ---------------------------------------------------------------
49 // Modified for the new G4ParticleChange 12 Mar. 1998 H.Kurahige
50 // Correct treatment of touchable in G4Step::UpdateTrack
51 // 12 May. 1998 H.Kurashige
52 // ---------------------------------------------------------------
53 // Separate implementation of inline functions inti G4Step.icc
54 // Add updating mass/charge 6 Oct. 1999 H.Kurashige
55 // add nonIonizingEnergyLoss 26 Mar. 2007 H.Kurashige
56 //
57 // Repository test - Dennis Wright
58 //
59 #ifndef G4Step_h
60 #define G4Step_h 1
61 
62 #include <stdlib.h> // Include from 'system'
63 #include <cmath> // Include from 'system'
64 #include "G4ios.hh" // Include from 'system'
65 #include <iomanip> // Include from 'system'
66 #include "globals.hh" // Include from 'global'
67 #include "G4ThreeVector.hh" // Include from 'global'
68 #include "G4VPhysicalVolume.hh" // Include from 'geometry'
69 #include "G4StepPoint.hh" // Include from 'track'
70 #include "G4StepStatus.hh" // Include from 'track'
71 class G4Polyline; // Forward declaration.
72 class G4Track; // Forward declaration.
73 #include "G4TrackVector.hh" // Include from 'tracking'
74 
76 class G4Step
78 {
79 
80 //--------
81  public:
82 
83 // Constructor/Destrcutor
84  G4Step();
85  ~G4Step();
86 
87 // Copy Counstructor and assignment operator
88  G4Step(const G4Step& );
89  G4Step & operator=(const G4Step &);
90 
91 //--------
92  public: // WIth description
93 
94 // Get/Set functions
95  // currnet track
96  G4Track* GetTrack() const;
97  void SetTrack(G4Track* value);
98 
99  // step points
100  G4StepPoint* GetPreStepPoint() const;
101  void SetPreStepPoint(G4StepPoint* value);
102 
103  G4StepPoint* GetPostStepPoint() const;
104  void SetPostStepPoint(G4StepPoint* value);
105 
106  // step length
107  G4double GetStepLength() const;
108  void SetStepLength(G4double value);
109  // Before the end of the AlongStepDoIt loop,StepLength keeps
110  // the initial value which is determined by the shortest geometrical Step
111  // proposed by a physics process. After finishing the AlongStepDoIt,
112  // it will be set equal to 'StepLength' in G4Step.
113 
114  // total energy deposit
116  void SetTotalEnergyDeposit(G4double value);
117 
118  // total non-ionizing energy deposit
121 
122  // cotrole flag for stepping
124  void SetControlFlag(G4SteppingControl StepControlFlag);
125 
126  // manipulation of total energy deposit
127  void AddTotalEnergyDeposit(G4double value);
129 
130  // manipulation of non-ionizng energy deposit
133 
134 
135  // Get/Set/Clear flag for initial/last step
136  // NOTE: following flags are not used
137  // will be ready in later release
138  G4bool IsFirstStepInVolume() const;
139  G4bool IsLastStepInVolume() const;
140 
141  void SetFirstStepFlag();
142  void ClearFirstStepFlag();
143  void SetLastStepFlag();
144  void ClearLastStepFlag();
145 
146  // difference of position, time, momentum and energy
148  G4double GetDeltaTime() const;
149 
150  // These methods will be deleted
151  // NOTE: use GetTotalEnergyDeposit() to obtain
152  // energy loss in the material
153  //
155  G4double GetDeltaEnergy() const;
156 
157 
158 // Other member functions
159  void InitializeStep( G4Track* aValue );
160  // initiaize contents of G4Step
161 
162  void UpdateTrack( );
163  // update track by using G4Step information
164 
165  void CopyPostToPreStepPoint( );
166  // copy PostStepPoint to PreStepPoint
167 
168  G4Polyline* CreatePolyline () const;
169  // for visualization
170 
171 //-----------
172  protected:
173 //-----------
174 
175 // Member data
177  // Accummulated total energy desposit in the current Step
178 
180  // Accummulated non-ionizing energy desposit in the current Step
181 
182 //---------
183  private:
184 //---------
185 
186 // Member data
187  G4StepPoint* fpPreStepPoint;
188  G4StepPoint* fpPostStepPoint;
189  G4double fStepLength;
190  // Step length which may be updated at each invocation of
191  // AlongStepDoIt and PostStepDoIt
192  G4Track* fpTrack;
193  //
194  G4SteppingControl fpSteppingControlFlag;
195  // A flag to control SteppingManager behavier from process
196 
197  // flag for initial/last step
198  G4bool fFirstStepInVolume;
199  G4bool fLastStepInVolume;
200 
201 // Secondary buckets
202 public:
203  // secodaries in the current step
205 
206  const std::vector<const G4Track*>* GetSecondaryInCurrentStep() const;
207 
208  // NOTE: Secondary bucket of the Step contains
209  // all secondaries during tracking the current track
210  // (i.e. NOT secondaries produced in the current step)
211  // all following methods give same object (i.e. G4TrackVector )
212  // but 2nd one will create bucket in addition
213  const G4TrackVector* GetSecondary() const ;
216 
217  // just delete secondary bucket
218  // NOTE: G4Track objects inside the bucket are not deleted
219  void DeleteSecondaryVector();
220 
221  // Add secondary tracks to the bucket
222  void SetSecondary( G4TrackVector* value);
223 
224 private:
225  // Secondaty bucket implemented by using std::vector of G4Track*
226  G4TrackVector* fSecondary;
227 
228  // number of secondaries which have been created by the last step
229  G4int nSecondaryByLastStep;
230 
231  typedef const G4Track* CT;
232  std::vector<CT>* secondaryInCurrentStep;
233 
234  // Prototyping implementation of smooth representation of curved
235  // trajectories. (jacek 30/10/2002)
236 public:
237  // Auxiliary points are ThreeVectors for now; change to
238  // G4VAuxiliaryPoints or some such (jacek 30/10/2002)
239  void SetPointerToVectorOfAuxiliaryPoints( std::vector<G4ThreeVector>* theNewVectorPointer ) {
240  fpVectorOfAuxiliaryPointsPointer = theNewVectorPointer;
241  }
242  std::vector<G4ThreeVector>* GetPointerToVectorOfAuxiliaryPoints() const {
243  return fpVectorOfAuxiliaryPointsPointer;
244  }
245 private:
246  // Explicity including the word "Pointer" in the name as I keep
247  // forgetting the * (jacek 30/10/2002)
248  std::vector<G4ThreeVector>* fpVectorOfAuxiliaryPointsPointer;
249 
250 };
251 
252 #include "G4Step.icc"
253 
254 
255 #endif
G4double GetDeltaEnergy() const
Definition: G4Step.cc:176
void SetFirstStepFlag()
void SetStepLength(G4double value)
G4double fTotalEnergyDeposit
Definition: G4Step.hh:176
void AddNonIonizingEnergyDeposit(G4double value)
void SetTrack(G4Track *value)
G4double GetStepLength() const
void DeleteSecondaryVector()
G4double fNonIonizingEnergyDeposit
Definition: G4Step.hh:179
void SetLastStepFlag()
G4SteppingControl
G4double GetNonIonizingEnergyDeposit() const
G4SteppingControl GetControlFlag() const
void SetSecondary(G4TrackVector *value)
G4int GetNumberOfSecondariesInCurrentStep() const
void SetPointerToVectorOfAuxiliaryPoints(std::vector< G4ThreeVector > *theNewVectorPointer)
Definition: G4Step.hh:239
void ClearLastStepFlag()
void ClearFirstStepFlag()
int G4int
Definition: G4Types.hh:78
G4StepPoint * GetPreStepPoint() const
void UpdateTrack()
G4TrackVector * NewSecondaryVector()
const XML_Char int const XML_Char * value
Definition: expat.h:331
void SetControlFlag(G4SteppingControl StepControlFlag)
bool G4bool
Definition: G4Types.hh:79
G4Step()
Definition: G4Step.cc:53
void ResetTotalEnergyDeposit()
void SetPreStepPoint(G4StepPoint *value)
G4ThreeVector GetDeltaMomentum() const
Definition: G4Step.cc:159
Definition: G4Step.hh:76
G4double GetDeltaTime() const
const std::vector< const G4Track * > * GetSecondaryInCurrentStep() const
Definition: G4Step.cc:193
G4Step & operator=(const G4Step &)
Definition: G4Step.cc:122
void SetPostStepPoint(G4StepPoint *value)
G4double GetTotalEnergyDeposit() const
void ResetNonIonizingEnergyDeposit()
G4TrackVector * GetfSecondary()
G4Polyline * CreatePolyline() const
std::vector< G4Track * > G4TrackVector
G4StepPoint * GetPostStepPoint() const
std::vector< G4ThreeVector > * GetPointerToVectorOfAuxiliaryPoints() const
Definition: G4Step.hh:242
void AddTotalEnergyDeposit(G4double value)
void SetNonIonizingEnergyDeposit(G4double value)
void SetTotalEnergyDeposit(G4double value)
void InitializeStep(G4Track *aValue)
~G4Step()
Definition: G4Step.cc:71
double G4double
Definition: G4Types.hh:76
G4Track * GetTrack() const
G4ThreeVector GetDeltaPosition() const
G4bool IsFirstStepInVolume() const
G4bool IsLastStepInVolume() const
const G4TrackVector * GetSecondary() const
void CopyPostToPreStepPoint()