Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4VParticleChange.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$
28//
29//
30// ------------------------------------------------------------
31// GEANT 4 class header file
32//
33//
34// ------------------------------------------------------------
35// Implemented for the new scheme 23 Mar. 1998 H.Kurahige
36//
37// Class Description
38// This class is the abstract class for ParticleChange.
39//-
40// The ParticleChange class ontains the results after invocation
41// of a physics process. This includes final states of parent
42// particle (momentum, energy, etc) and secondary particles generated
43// by the interaction.
44// The tracking assumes that all the values of energy and
45// momentum are in global reference system, therefore all the
46// needed Lorentz transformations must have been already Done
47// when filling the data-members of this class.
48//-
49//-
50// This abstract class has following four virtual methods
51// virtual G4Step* UpdateStepForAtRest(G4Step* Step);
52// virtual G4Step* UpdateStepForAlongStep(G4Step* Step);
53// virtual G4Step* UpdateStepForPostStep(G4Step* Step);
54// virtual void Initialize(const G4Track&);
55// The UpdateStep methods return the pointer to the G4Step
56// after updating the given Step information by using final state
57// information of the track given by a physics process.
58// User must add methods to keep the final state information
59// in his derived class as well as implement UpdateStep methods
60// which he want to use.
61//-
62// The Initialize methods is provided to refresh the final
63// state information and should be called by each process
64// at the beginning of DoIt.
65//
66// ------------------------------------------------------------
67// Implement Event Biasing Scheme 9 Nov.,98 H.Kurashige
68// add CheckIt 13 Apr.,99 H.Kurashige
69// add accuracy leveles 5 May, 99 H.Kurashige
70// add check secondaries 11 June, 03 H.Kurashige
71// add new methods of ProposeXXX 08 May, 04 H.Kurashige
72// remove obsolete methods of SetXXX 19 Sep, 04 H.Kurashige
73// add flag for first/last step in volume 30 Oct. 2006 H.Kurashige
74// add nonIonizingEnergyLoss 26 Mar 2007 H.Kurashige
75// modify/fix bugs related to weight 17 Sep. 2011 H.Kurashige
76// fix bugs related to weight 29 Apr. 2012 H.Kurashige
77//
78
79#ifndef G4VParticleChange_h
80#define G4VParticleChange_h 1
81
82#include "globals.hh"
83#include "G4ios.hh"
84#include <cmath>
85
86class G4Track;
87class G4Step;
88
89#include "G4TrackFastVector.hh"
90#include "G4TrackStatus.hh"
91#include "G4SteppingControl.hh"
92
93
95{
96 public:
97 // default constructor
99
100 // destructor
101 virtual ~G4VParticleChange();
102
103 // equal/unequal operator
104 G4bool operator==(const G4VParticleChange &right) const;
105 G4bool operator!=(const G4VParticleChange &right) const;
106 // "equal" means that teo objects have the same pointer.
107
108 protected:
109 // hide copy constructor and assignment operaor as protected
112
113 public: // with description
114 // --- the following methods are for updating G4Step -----
115 virtual G4Step* UpdateStepForAtRest(G4Step* Step);
116 virtual G4Step* UpdateStepForAlongStep(G4Step* Step);
117 virtual G4Step* UpdateStepForPostStep(G4Step* Step);
118 // Return the pointer to the G4Step after updating the Step information
119 // by using final state information of the track given by a physics
120 // process
121
122 protected: // with description
124 // Update the G4Step specific attributes
125 // (i.e. SteppingControl, LocalEnergyDeposit, and TrueStepLength)
126
127
128 public: // with description
129 virtual void Initialize(const G4Track&);
130 // This methods will be called by each process at the beginning of DoIt
131 // if necessary.
132
133 protected:
139
143 // ------------------------------------------------------
144
145 public: // with description
146 //---- the following methods are for TruePathLength ----
148 void ProposeTrueStepLength(G4double truePathLength);
149 // Get/Propose theTrueStepLength
150
151 //---- the following methods are for LocalEnergyDeposit ----
154 // Get/Propose the locally deposited energy
155
156 //---- the following methods are for nonIonizingEnergyDeposit ----
159 // Get/Propose the non-ionizing deposited energy
160
161 //---- the following methods are for TrackStatus -----
164 // Get/Propose the final TrackStatus of the current particle.
165 // ------------------------------------------------------
166
167 //---- the following methods are for managements of SteppingControl --
170 // Set/Propose a flag to control stepping manager behavier
171 // ------------------------------------------------------
172
173 //---- the following methods are for managements of initial/last step
178
179 //---- the following methods are for managements of secondaries --
180 void Clear();
181 // Clear the contents of this objects
182 // This method should be called after the Tracking(Stepping)
183 // manager removes all secondaries in theListOfSecondaries
184
185 void SetNumberOfSecondaries(G4int totSecondaries);
186 // SetNumberOfSecondaries must be called just before AddSecondary()
187 // in order to secure memory space for theListOfSecondaries
188 // This method resets theNumberOfSecondaries to 0
189 // (that will be incremented at every AddSecondary() call).
190
192 // Returns the number of secondaries current stored in
193 // G4TrackFastVector.
194
195 G4Track* GetSecondary(G4int anIndex) const;
196 // Returns the pointer to the generated secondary particle
197 // which is specified by an Index.
198
199 void AddSecondary(G4Track* aSecondary);
200 // Add a secondary particle to theListOfSecondaries.
201 // ------------------------------------------------------
202
205 // Get weight of the parent (i.e. current) track
206 void ProposeWeight(G4double finalWeight);
207 void ProposeParentWeight(G4double finalWeight);
208 // Propse new weight of the parent (i.e. current) track
209 // As for AlongStepDoIt, the parent weight will be set
210 // in accumulated manner
211 // i.e.) If two processes propose weight of W1 and W2 respectively
212 // for the track with initial weight of W0
213 // the final weight is set to
214 // (W1/W0) * (W2/W0) * W0
215
218 // In default (fSecondaryWeightByProcess flag is false),
219 // the weight of secondary tracks will be set to
220 // the parent weight
221 // If fSecondaryWeightByProcess flag is true,
222 // the weight of secondary tracks will not be changed
223 // by the ParticleChange
224 // (i.e. the process determine the secodary weight)
225 // NOTE:
226 // Make sure that only one processe in AlongStepDoIt
227 // proposes the parent weight,
228 // If several processes in AlongStepDoIt proposes
229 // the parent weight and add secondaties with
230 // fSecondaryWeightByProcess is set to false,
231 // secondary weights may be wrong
232
235 // Obsolete
236
237 virtual void DumpInfo() const;
238 // Print out information
239
240 void SetVerboseLevel(G4int vLevel);
242
243 protected:
244
246 // The vector of secondaries.
247
249 // The total number of secondaries produced by each process.
250
252 // TheSizeOftheListOfSecondaries;
253
255 // The changed (final) track status of a given particle.
256
258 // a flag to control stepping manager behavior
259
261 // It represents the part of the energy lost for discrete
262 // or semi-continuous processes which is due to secondaries
263 // not generated because they would have been below their cut
264 // threshold.
265 // The sum of the locally deposited energy + the delta-energy
266 // coming from the continuous processes gives the
267 // total energy loss localized in the current Step.
268
270 // non-ionizing energu deposit is defined as
271 // a part of local energy deposit, which does not cause
272 // ionization of atoms
273
275 // The value of "True" Step Length
276
277
280 // flag for initial/last step
281
283 // Weight ofparent track
285 // flags for Weight ofparent track
287 // flag for setting weight of secondaries
288
290 // global time of the parent.
291 // This is used only for checking
292
294 // The Verbose level
295
296 public: // with description
297 // CheckIt method is provided for debug
298 virtual G4bool CheckIt(const G4Track&);
299
300 // CheckIt method is activated
301 // if debug flag is set and 'G4VERBOSE' is defined
305
306 protected:
307 // CheckSecondary method is provided for debug
309
312
313 protected:
315
316 // accuracy levels
319
320
321};
322
323#include "G4Step.hh"
324#include "G4Track.hh"
325#include "G4VParticleChange.icc"
326
327#endif
G4SteppingControl
G4TrackStatus
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
Definition: G4Step.hh:78
void InitializeTrueStepLength(const G4Track &)
G4bool GetFirstStepInVolume() const
void InitializeSteppingControl(const G4Track &)
G4double GetNonIonizingEnergyDeposit() const
G4double GetParentWeight() const
virtual G4bool CheckIt(const G4Track &)
static const G4double accuracyForException
void InitializeParentGlobalTime(const G4Track &)
void InitializeStatusChange(const G4Track &)
void ProposeTrackStatus(G4TrackStatus status)
G4TrackStatus theStatusChange
G4TrackFastVector * theListOfSecondaries
G4SteppingControl theSteppingControlFlag
void SetSecondaryWeightByProcess(G4bool)
static const G4double accuracyForWarning
G4bool IsParentWeightSetByProcess() const
void ProposeLastStepInVolume(G4bool flag)
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
G4SteppingControl GetSteppingControl() const
virtual G4Step * UpdateStepForAlongStep(G4Step *Step)
void ProposeWeight(G4double finalWeight)
void InitializeStepInVolumeFlags(const G4Track &)
virtual G4Step * UpdateStepForAtRest(G4Step *Step)
virtual void Initialize(const G4Track &)
G4double GetLocalEnergyDeposit() const
G4bool GetLastStepInVolume() const
void SetVerboseLevel(G4int vLevel)
G4double theNonIonizingEnergyDeposit
virtual G4Step * UpdateStepForPostStep(G4Step *Step)
G4double GetAccuracyForException() const
G4int GetNumberOfSecondaries() const
void AddSecondary(G4Track *aSecondary)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4bool GetDebugFlag() const
void SetParentWeightByProcess(G4bool)
G4int GetVerboseLevel() const
void InitializeLocalEnergyDeposit(const G4Track &)
G4double GetAccuracyForWarning() const
virtual void DumpInfo() const
G4double GetWeight() const
G4VParticleChange & operator=(const G4VParticleChange &right)
G4bool CheckSecondary(G4Track &)
G4Step * UpdateStepInfo(G4Step *Step)
void ProposeSteppingControl(G4SteppingControl StepControlFlag)
void ProposeParentWeight(G4double finalWeight)
void ProposeFirstStepInVolume(G4bool flag)
G4bool IsSecondaryWeightSetByProcess() const
G4Track * GetSecondary(G4int anIndex) const
void InitializeSecondaries(const G4Track &)
void ProposeTrueStepLength(G4double truePathLength)
void InitializeParentWeight(const G4Track &)
G4bool operator==(const G4VParticleChange &right) const
G4bool operator!=(const G4VParticleChange &right) const
G4TrackStatus GetTrackStatus() const
virtual ~G4VParticleChange()
G4double GetTrueStepLength() const