Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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#ifndef G4FastStep_h
62#define G4FastStep_h
63
64#include "G4ParticleMomentum.hh"
65#include "G4ThreeVector.hh"
66#include "G4ios.hh"
67#include "globals.hh"
69#include "G4FastTrack.hh"
70#include "G4VParticleChange.hh"
71
72//-------------------------------------------
73//
74// G4FastStep class
75//
76//-------------------------------------------
77
78// Class Description:
79// The final state of the particles after parameterisation has to be returned through a G4FastStep
80// reference. This final state is described as "requests" the tracking will apply after your
81// parameterisation has been invoked.
82//
83// To facilitate the developers work, changes of position/normalized direction of the
84// momentum/polarization can be specified in the local coordinate system of the envelope or in the
85// global one.
86// The default is local system coordinates.
87//
88
90{
91 public: // Without description
92 //------------------------
93 // Constructor/Destructor
94 //------------------------
95 G4FastStep() = default;
96 ~G4FastStep() override = default;
97
98 G4FastStep(const G4FastStep& right) = delete;
99 G4FastStep& operator=(const G4FastStep& right) = delete;
100
101 // Set the kinetic energy of the primary to zero, and set the "fStopAndKill" signal
102 // used by the stepping.
103 void KillPrimaryTrack();
104
105 // -- Methods used to change the position, normalized direction of
106 // the momentum, time etc... of the primary.
107 // .. space and time:
108
109 // Set the primary track final position.
110 void ProposePrimaryTrackFinalPosition(const G4ThreeVector&, G4bool localCoordinates = true);
111
112 // Set the primary track final position -- maintained for backward compatibility.
113 [[deprecated("use ProposePrimaryTrackFinalPosition instead")]]
114 void SetPrimaryTrackFinalPosition(const G4ThreeVector&, G4bool localCoordinates = true);
115
116 // Set the primary track final time.
118
119 // Set the primary track final time -- maintained for backward compatibility.
120 [[deprecated("use ProposePrimaryTrackFinalTime instead")]]
122
123 // Set the primary final track Proper Time.
125
126 // Set the primary final track Proper Time -- maintained for backward compatibility.
127 [[deprecated("use ProposePrimaryTrackProperTime instead")]]
129
130 // .. dynamics:
131
132 // Be careful: the Track Final Momentum means the normalized direction
133 // of the momentum!
135 G4bool localCoordinates = true);
136
137 // Set the primary track final momentum -- maintained for backward compatibility. Same as
138 // ProposePrimaryTrackMomentumDirection(...)
139 [[deprecated("use ProposePrimaryTrackMomentumDirection instead")]]
140 void SetPrimaryTrackFinalMomentum(const G4ThreeVector&, G4bool localCoordinates = true);
141
142 // Set the primary track final kinetic energy.
144
145 // Set the primary track final kinetic energy-- maintained for backward compatibility.
146 [[deprecated("use ProposePrimaryTrackFinalKineticEnergy instead")]]
148
149 // Set the primary track final kinetic energy and direction.
151 G4bool localCoordinates = true);
152
153 // Set the primary track final kinetic energy and direction -- maintained for backward
154 // compatibility.
155 [[deprecated("use ProposePrimaryTrackFinalKineticEnergyAndDirection instead")]]
157 G4bool localCoordinates = true);
158
159 // Set the primary track final polarization.
160 void ProposePrimaryTrackFinalPolarization(const G4ThreeVector&, G4bool localCoordinates = true);
161
162 // Set the primary track final polarization.
163 [[deprecated("use ProposePrimaryTrackFinalPolarization instead")]]
164 void SetPrimaryTrackFinalPolarization(const G4ThreeVector&, G4bool localCoordinates = true);
165
166 // Set the true path length of the primary track during the step.
168
169 // Set the true path length of the primary track during the step -- maintained for backward
170 // compatibility.
171 [[deprecated("use ProposePrimaryTrackPathLength instead")]]
173
174 // Set the weight applied for event biasing mechanism.
176
177 // Set the weight applied for event biasing mechanism -- kept for backward compatibility.
178 [[deprecated("use ProposePrimaryTrackFinalEventBiasingWeight instead")]]
180
181 // ------------------------------
182 // -- Management of secondaries:
183 // ------------------------------
184
185 // ----------------------------------------------------
186 // -- The creation of secondaries is Done in two steps:
187 // -- 1) Give the total number of secondaries
188 // -- that the FastStep returns
189 // -- to the tracking using:
190 // -- SetNumberOfSecondaryTracks()
191 // --
192 // -- 2) Invoke the CreateSecondaryTrack() method
193 // -- to create one secondary at each time.
194 // ----------------------------------------------------
195
196 // Set the total number of secondaries that will be created.
197 // -- Total Number of secondaries to be created,
198 // -- (to be called first)
200
201 // Returns the number of secondaries effectively stored.
202 // -- Number of secondaries effectively stored:
203 // -- (incremented at each CreateSecondaryTrack()
204 // -- call)
206
207 // -- Create a secondary: the arguments are:
208 // -- * G4DynamicsParticle: see header file, many constructors exist
209 // -- (allow to set particle type + energy +
210 // -- the normalized direction of momentum...)
211 // -- * G4ThreeVector : Polarization (not in G4ParticleChange constructor)
212 // -- * G4ThreeVector : Position
213 // -- * G4double : Time
214 // -- * G4bool : says if Position/Momentum are given in the
215 // -- local coordinate system (true by default)
216 // -- Returned value: pointer to the track created.
218 G4bool localCoordinates = true);
219
220 //-- Create a secondary: the difference with he above declaration
221 //-- is that the Polarization is not given and is assumed already set
222 //-- in the G4DynamicParticle.
223 //-- Returned value: pointer to the track created
225 G4bool localCoordinates = true);
226
227 // Returns a pointer on the i-th secondary track created.
229
230 //------------------------------------------------
231 //
232 // Total energy deposit in the "fast Step"
233 // (a default should be provided in future,
234 // which can be:
235 // delta energy of primary -
236 // energy of the secondaries)
237 // This allow the user to Store a consistent
238 // information in the G4Trajectory.
239 //
240 //------------------------------------------------
241 // Set the total energy deposited.
243
244 // Set the total energy deposited -- kept for backward compatibility.
245 // It should be the delta energy of primary less the energy of the secondaries.
246 [[deprecated("use ProposeTotalEnergyDeposited instead")]]
248
249 // Returns the total energy deposited.
251
252 // Control of the stepping manager Hit invocation.
253 //
254 // In a usual parameterisation, the control of the hits production is under the user
255 // responsability in his G4VFastSimulationModel (he generally produces several hits at once.)
256 //
257 // However, in the particular case the G4FastSimulation user's model acts as the physics
258 // replacement only (ie replaces all the ***DoIt() and leads to the construction of a meaningful
259 // G4Step), the user can delegate to the G4SteppingManager the responsability to invoke
260 // the Hit()method of the current sensitive if any.
261 //
262 // By default, the G4SteppingManager is asked to NOT invoke this Hit() method when
263 // parameterisation is invoked.
265
266 // ===============================================
267 // Stepping interface.
268 // ===============================================
269 // --- the following methods are for updating G4Step -----
270 // Return the pointer to the G4Step after updating the Step information
271 // by using final state information of the track given by a Model.
272 //
273 // The Fast Simulation Mechanism doesn't change the track's final
274 // state on the AlongDoIt loop, so the default one all we need.
275 // virtual G4Step* UpdateStepForAlongStep(G4Step* Step);
276
277 G4Step* UpdateStepForAtRest(G4Step* Step) override;
278 G4Step* UpdateStepForPostStep(G4Step* Step) override;
279
280 // A Model gives the final state of the particle
281 // based on information of G4FastTrack. So the
282 // Initialize method is an interface to the
283 // G4FastSimulationManager to Initialize the
284 // G4FastStep.
285
286 void Initialize(const G4FastTrack&);
287
288 // for Debug
289 void DumpInfo() const override;
290 G4bool CheckIt(const G4Track&) override;
291
292 private:
293 //===================================================
294 // Private Internal methods (implementation).
295 //===================================================
296
297 // G4FastStep should never be Initialized in this way
298 // but we must define it to avoid compiler warnings.
299 void Initialize(const G4Track&) override;
300
301 // -- Utility functions --
302 //--- methods to keep information of the final state--
303 // IMPORTANT NOTE: Although the name of the class and methods are
304 // "Change", what it stores (and returns in get) are the "FINAL"
305 // values of the Position, the normalized direction of Momentum,
306 // etc.
307
308 // Set theMomentumChange vector: it is the final unitary momentum
309 // direction.
310 void SetMomentumChange(G4double Px, G4double Py, G4double Pz);
311 void SetMomentumChange(const G4ThreeVector& Pfinal);
312
313 //=====================================================
314 // Data members.
315 //=====================================================
316 // theMomentumChange is the vector containing the final momentum
317 // direction after the invoked process. The application of the change
318 // of the momentum direction of the particle is not Done here.
319 // The responsibility to apply the change is up the entity
320 // which invoked the process.
321 G4ParticleMomentum theMomentumChange;
322
323 // The changed (final) polarization of a given particle.
324 G4ThreeVector thePolarizationChange;
325
326 // The final kinetic energy of the current particle.
327 G4double theEnergyChange = 0.0;
328
329 // The changed (final) position of a given particle.
330 G4ThreeVector thePositionChange;
331
332 // The changed (final) global time of a given particle.
333 G4double theTimeChange = 0.0;
334
335 // The changed (final) proper time of a given particle.
336 G4double theProperTimeChange = 0.0;
337
338 // The reference G4FastTrack
339 const G4FastTrack* fFastTrack = nullptr;
340
341 // weight for event biasing mechanism:
342 G4double theWeightChange = 0.0;
343};
344
345//*******************************************************************
346//
347// Inline functions
348//
349//*******************************************************************
350
351#include "G4FastStep.icc"
352
353#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void SetTotalEnergyDeposited(G4double anEnergyPart)
void SetPrimaryTrackFinalKineticEnergyAndDirection(G4double, const G4ThreeVector &, G4bool localCoordinates=true)
G4Track * GetSecondaryTrack(G4int)
void SetPrimaryTrackFinalProperTime(G4double)
void KillPrimaryTrack()
Definition G4FastStep.cc:87
void SetPrimaryTrackFinalPosition(const G4ThreeVector &, G4bool localCoordinates=true)
G4bool CheckIt(const G4Track &) override
void SetPrimaryTrackFinalTime(G4double)
void SetPrimaryTrackFinalPolarization(const G4ThreeVector &, G4bool localCoordinates=true)
G4Step * UpdateStepForPostStep(G4Step *Step) override
void SetPrimaryTrackFinalEventBiasingWeight(G4double)
void ProposePrimaryTrackFinalKineticEnergyAndDirection(G4double, const G4ThreeVector &, G4bool localCoordinates=true)
void DumpInfo() const override
void SetPrimaryTrackFinalMomentum(const G4ThreeVector &, G4bool localCoordinates=true)
~G4FastStep() override=default
void SetNumberOfSecondaryTracks(G4int)
void SetPrimaryTrackPathLength(G4double)
void ProposeTotalEnergyDeposited(G4double anEnergyPart)
void ProposePrimaryTrackFinalPosition(const G4ThreeVector &, G4bool localCoordinates=true)
Definition G4FastStep.cc:96
G4FastStep(const G4FastStep &right)=delete
G4FastStep & operator=(const G4FastStep &right)=delete
void ProposePrimaryTrackFinalProperTime(G4double)
void ProposePrimaryTrackFinalTime(G4double)
void SetPrimaryTrackFinalKineticEnergy(G4double)
G4int GetNumberOfSecondaryTracks()
G4Step * UpdateStepForAtRest(G4Step *Step) override
void ForceSteppingHitInvocation()
void ProposePrimaryTrackFinalPolarization(const G4ThreeVector &, G4bool localCoordinates=true)
G4double GetTotalEnergyDeposited() const
void ProposePrimaryTrackFinalKineticEnergy(G4double)
void ProposePrimaryTrackFinalMomentumDirection(const G4ThreeVector &, G4bool localCoordinates=true)
void ProposePrimaryTrackPathLength(G4double)
void ProposePrimaryTrackFinalEventBiasingWeight(G4double)
G4Track * CreateSecondaryTrack(const G4DynamicParticle &, G4ThreeVector, G4ThreeVector, G4double, G4bool localCoordinates=true)
void Initialize(const G4FastTrack &)
Definition G4FastStep.cc:53
G4FastStep()=default