Geant4-11
G4FastStep.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//
28//
29//---------------------------------------------------------------
30//
31// G4FastStep.hh
32//
33// Description:
34// The G4FastStep class insures a friendly interface
35// to manage the primary/secondaries final state for
36// Fast Simulation Models. This includes final states of parent
37// particle (normalized direction of the momentum, energy, etc) and
38// secondary particles generated by the parameterisation.
39//
40// The G4FastStep class acts also as the G4ParticleChange
41// for the Fast Simulation Process. So it inherites from
42// the G4VParticleChange class and redefines the four virtual
43// methods :
44//
45// virtual G4Step* UpdateStepForAtRest(G4Step* Step);
46// virtual G4Step* UpdateStepForAlongStep(G4Step* Step);
47// virtual G4Step* UpdateStepForPostStep(G4Step* Step);
48// virtual void Initialize(const G4Track&);
49//
50// History:
51// Oct 97: Verderi && MoraDeFreitas - First Implementation.
52// Dec 97: Verderi - ForceSteppingHitInvocation(),
53// Set/GetTotalEnergyDeposited() methods.
54// Apr 98: MoraDeFreitas - G4FastStep becomes the G4ParticleChange
55// for the Fast Simulation Process.
56// Nov 04: Verderi - Add ProposeXXX methods. SetXXX ones are kept
57// for backward compatibility.
58//
59//---------------------------------------------------------------
60
61
62#ifndef G4FastStep_h
63#define G4FastStep_h
64
65#include "globals.hh"
66#include "G4ios.hh"
67#include "G4ThreeVector.hh"
68#include "G4ParticleMomentum.hh"
70#include "G4VParticleChange.hh"
71#include "G4FastTrack.hh"
72
73//-------------------------------------------
74//
75// G4FastStep class
76//
77//-------------------------------------------
78
79// Class Description:
80// The final state of the particles after parameterisation has to be returned through a G4FastStep
81// reference. This final state is described as "requests" the tracking will apply after your
82// parameterisation has been invoked.
83//
84// To facilitate the developers work, changes of position/normalized direction of the
85// momentum/polarization can be specified in the local coordinate system of the envelope or in the
86// global one.
87// The default is local system coordinates.
88//
89
91{
92public: // with Description
93 void KillPrimaryTrack();
94 // Set the kinetic energy of the primary to zero, and set the "fStopAndKill" signal
95 // used by the stepping.
96
97 // -- Methods used to change the position, normalized direction of
98 // the momentum, time etc... of the primary.
99 // .. space and time:
101 G4bool localCoordinates = true);
102 // Set the primary track final position.
104 G4bool localCoordinates = true);
105 // Set the primary track final position -- maintained for backward compatibility.
106
107
109 // Set the primary track final time.
111 // Set the primary track final time -- maintained for backward compatibility.
112
113
115 // Set the primary final track Proper Time.
117 // Set the primary final track Proper Time -- maintained for backward compatibility.
118
119
120 // .. dynamics:
122 G4bool localCoordinates = true);
123 // Be careful: the Track Final Momentum means the normalized direction
124 // of the momentum!
126 G4bool localCoordinates = true);
127 // Set the primary track final momentum -- maintained for backward compatibility. Same as ProposePrimaryTrackMomentumDirection(...)
128
129
131 // Set the primary track final kinetic energy.
133 // Set the primary track final kinetic energy-- maintained for backward compatibility.
134
135
137 const G4ThreeVector &,
138 G4bool localCoordinates
139 = true);
140 // Set the primary track final kinetic energy and direction.
142 const G4ThreeVector &,
143 G4bool localCoordinates
144 = true);
145 // Set the primary track final kinetic energy and direction -- maintained for backward compatibility.
146
147
148
150 G4bool localCoordinates = true);
151 // Set the primary track final polarization.
153 G4bool localCoordinates = true);
154 // Set the primary track final polarization.
155
156
158 // Set the true path length of the primary track during the step.
160 // Set the true path length of the primary track during the step -- maintained for backward compatibility.
161
163 // Set the weight applied for event biasing mechanism.
165 // Set the weight applied for event biasing mechanism -- kept for backward compatibility.
166
167 // ------------------------------
168 // -- Management of secondaries:
169 // ------------------------------
170
171 // ----------------------------------------------------
172 // -- The creation of secondaries is Done in two steps:
173 // -- 1) Give the total number of secondaries
174 // -- that the FastStep returns
175 // -- to the tracking using:
176 // -- SetNumberOfSecondaryTracks()
177 // --
178 // -- 2) Invoke the CreateSecondaryTrack() method
179 // -- to create one secondary at each time.
180 // ----------------------------------------------------
181
182 // -- Total Number of secondaries to be created,
183 // -- (to be called first)
185 // Set the total number of secondaries that will be created.
186
187 // -- Number of secondaries effectively stored:
188 // -- (incremented at each CreateSecondaryTrack()
189 // -- call)
191 // Returns the number of secondaries effectively stored.
192
193 // -- Create a secondary: the arguments are:
194 // -- * G4DynamicsParticle: see header file, many constructors exist
195 // -- (allow to set particle type + energy +
196 // -- the normalized direction of momentum...)
197 // -- * G4ThreeVector : Polarization (not in G4ParticleChange constructor)
198 // -- * G4ThreeVector : Position
199 // -- * G4double : Time
200 // -- * G4bool : says if Position/Momentum are given in the
201 // -- local coordinate system (true by default)
202 // -- Returned value: pointer to the track created.
206 G4double,
207 G4bool localCoordinates=true);
208 // Create a secondary. The arguments are:
209 //
210 // G4DynamicsParticle: see the G4DynamicsParticle reference, many constructors exist
211 // (allow to set particle type + energy + the normalized direction of
212 // momentum...);
213 // G4ThreeVector : Polarization;
214 // G4ThreeVector : Position;
215 // G4double : Time;
216 // G4bool : says if Position/Momentum are given in the local envelope coordinate
217 // system (true by default).
218 //
219 // Returned value: pointer to the track created.
220 //
221
222 //-- Create a secondary: the difference with he above declaration
223 //-- is that the Polarization is not given and is assumed already set
224 //-- in the G4DynamicParticle.
225 //-- Returned value: pointer to the track created
228 G4double,
229 G4bool localCoordinates=true);
230 // Create a secondary. The difference with he above declaration is that the Polarization is not
231 // given and is assumed already set in the G4DynamicParticle.
232 //
233 // Returned value: pointer to the track created
234
235
236
238 // Returns a pointer on the i-th secondary track created.
239
240 //------------------------------------------------
241 //
242 // Total energy deposit in the "fast Step"
243 // (a default should be provided in future,
244 // which can be:
245 // delta energy of primary -
246 // energy of the secondaries)
247 // This allow the user to Store a consistent
248 // information in the G4Trajectory.
249 //
250 //------------------------------------------------
252 // Set the total energy deposited.
254 // Set the total energy deposited -- kept for backward compatibility.
255 // It should be the delta energy of primary less the energy of the secondaries.
256
258 // Returns the total energy deposited.
259
261 // Control of the stepping manager Hit invocation.
262 //
263 // In a usual parameterisation, the control of the hits production is under the user
264 // responsability in his G4VFastSimulationModel (he generally produces several hits at once.)
265 //
266 // However, in the particular case the G4FastSimulation user's model acts as the physics
267 // replacement only (ie replaces all the ***DoIt() and leads to the construction of a meaningful
268 // G4Step), the user can delegate to the G4SteppingManager the responsability to invoke
269 // the Hit()method of the current sensitive if any.
270 //
271 // By default, the G4SteppingManager is asked to NOT invoke this Hit() method when parameterisation
272 // is invoked.
273 //
274
275
276public: // Without description
277 //=======================================================
278 // Implementation section and kernel interfaces.
279 //=======================================================
280 //------------------------
281 // Constructor/Destructor
282 //------------------------
283 G4FastStep();
284 virtual ~G4FastStep();
285
286 // equal/unequal operator
287 G4bool operator==(const G4FastStep &right) const;
288 G4bool operator!=(const G4FastStep &right) const;
289
290protected:
291 // hide copy constructor and assignment operator as protected
292 G4FastStep (const G4FastStep &right);
293 G4FastStep & operator= (const G4FastStep &right);
294
295public:
296 // ===============================================
297 // Stepping interface.
298 // ===============================================
299 // --- the following methods are for updating G4Step -----
300 // Return the pointer to the G4Step after updating the Step information
301 // by using final state information of the track given by a Model.
302 //
303 // The Fast Simulation Mechanism doesn't change the track's final
304 // state on the AlongDoIt loop, so the default one all we need.
305 //virtual G4Step* UpdateStepForAlongStep(G4Step* Step);
306
309
310 // A Model gives the final state of the particle
311 // based on information of G4FastTrack. So the
312 // Initialize method is an interface to the
313 // G4FastSimulationManager to Initialize the
314 // G4FastStep.
315
316 void Initialize(const G4FastTrack&);
317
318 // for Debug
319 void DumpInfo() const;
320 G4bool CheckIt(const G4Track&);
321
322private:
323 //===================================================
324 // Private Internal methods (implementation).
325 //===================================================
326
327 // G4FastStep should never be Initialized in this way
328 // but we must define it to avoid compiler warnings.
329 void Initialize(const G4Track&);
330
331 // -- Utility functions --
332 //--- methods to keep information of the final state--
333 // IMPORTANT NOTE: Although the name of the class and methods are
334 // "Change", what it stores (and returns in get) are the "FINAL"
335 // values of the Position, the normalized direction of Momentum,
336 // etc.
337
338 // Set theMomentumChange vector: it is the final unitary momentum
339 // direction.
341 void SetMomentumChange(const G4ThreeVector& Pfinal);
342
343 //=====================================================
344 // Data members.
345 //=====================================================
346 // theMomentumChange is the vector containing the final momentum
347 // direction after the invoked process. The application of the change
348 // of the momentum direction of the particle is not Done here.
349 // The responsibility to apply the change is up the entity
350 // which invoked the process.
352
353 // The changed (final) polarization of a given particle.
355
356 // The final kinetic energy of the current particle.
358
359 // The changed (final) position of a given particle.
361
362 // The changed (final) global time of a given particle.
364
365 // The changed (final) proper time of a given particle.
367
368 // The reference G4FastTrack
369 const G4FastTrack* fFastTrack = nullptr;
370
371 // weight for event biasing mechanism:
373};
374
375//*******************************************************************
376//
377// Inline functions
378//
379//*******************************************************************
380
381#include "G4FastStep.icc"
382
383#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
void SetTotalEnergyDeposited(G4double anEnergyPart)
G4FastStep & operator=(const G4FastStep &right)
Definition: G4FastStep.cc:322
void SetPrimaryTrackFinalKineticEnergyAndDirection(G4double, const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:166
G4Track * GetSecondaryTrack(G4int)
G4Step * UpdateStepForPostStep(G4Step *Step)
Definition: G4FastStep.cc:364
G4ThreeVector thePolarizationChange
Definition: G4FastStep.hh:354
void SetPrimaryTrackFinalProperTime(G4double)
virtual ~G4FastStep()
Definition: G4FastStep.cc:306
void KillPrimaryTrack()
Definition: G4FastStep.cc:87
void SetPrimaryTrackFinalPosition(const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:113
G4bool operator==(const G4FastStep &right) const
Definition: G4FastStep.cc:350
G4bool operator!=(const G4FastStep &right) const
Definition: G4FastStep.cc:355
G4double theEnergyChange
Definition: G4FastStep.hh:357
void SetPrimaryTrackFinalTime(G4double)
void SetPrimaryTrackFinalPolarization(const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:192
void SetPrimaryTrackFinalEventBiasingWeight(G4double)
void SetMomentumChange(const G4ThreeVector &Pfinal)
void DumpInfo() const
Definition: G4FastStep.cc:437
G4double theProperTimeChange
Definition: G4FastStep.hh:366
void ProposePrimaryTrackFinalKineticEnergyAndDirection(G4double, const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:150
void SetPrimaryTrackFinalMomentum(const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:139
void SetNumberOfSecondaryTracks(G4int)
G4bool CheckIt(const G4Track &)
Definition: G4FastStep.cc:460
void SetPrimaryTrackPathLength(G4double)
void ProposeTotalEnergyDeposited(G4double anEnergyPart)
void ProposePrimaryTrackFinalPosition(const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:98
G4double theWeightChange
Definition: G4FastStep.hh:372
void ProposePrimaryTrackFinalProperTime(G4double)
void ProposePrimaryTrackFinalTime(G4double)
G4ParticleMomentum theMomentumChange
Definition: G4FastStep.hh:351
void SetPrimaryTrackFinalKineticEnergy(G4double)
G4int GetNumberOfSecondaryTracks()
G4ThreeVector thePositionChange
Definition: G4FastStep.hh:360
void SetMomentumChange(G4double Px, G4double Py, G4double Pz)
const G4FastTrack * fFastTrack
Definition: G4FastStep.hh:369
void ForceSteppingHitInvocation()
G4double theTimeChange
Definition: G4FastStep.hh:363
void ProposePrimaryTrackFinalPolarization(const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:178
G4double GetTotalEnergyDeposited() const
void ProposePrimaryTrackFinalKineticEnergy(G4double)
void ProposePrimaryTrackFinalMomentumDirection(const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:124
G4Step * UpdateStepForAtRest(G4Step *Step)
Definition: G4FastStep.cc:401
void ProposePrimaryTrackPathLength(G4double)
void ProposePrimaryTrackFinalEventBiasingWeight(G4double)
G4Track * CreateSecondaryTrack(const G4DynamicParticle &, G4ThreeVector, G4ThreeVector, G4double, G4bool localCoordinates=true)
Definition: G4FastStep.cc:202
void Initialize(const G4FastTrack &)
Definition: G4FastStep.cc:53
Definition: G4Step.hh:62