Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FastTrack.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// G4FastTrack.hh
31//
32// Description:
33// Keeps the current track information and special features
34// for Parameterised Simulation Models.
35//
36// History:
37// Oct 97: Verderi && MoraDeFreitas - First Implementation.
38//
39//---------------------------------------------------------------
40
41#ifndef G4FastTrack_h
42#define G4FastTrack_h
43
44#include "G4AffineTransform.hh"
45#include "G4LogicalVolume.hh"
46#include "G4Navigator.hh"
47#include "G4Region.hh"
48#include "G4Track.hh"
49#include "G4VSolid.hh"
50
51//---------------------------
52// For possible future needs:
53//---------------------------
55
56//-------------------------------------------
57//
58// G4FastTrack class
59//
60//-------------------------------------------
61
62// Class Description:
63// The G4FastTrack provides you access to the current G4Track,
64// gives simple access to envelope related features (G4Region,
65// G4LogicalVolume, G4VSolid, G4AffineTransform references between
66// the global and the envelope local coordinates systems) and
67// simple access to the position, momentum expressed in the
68// envelope coordinate system. Using those quantities and the
69// G4VSolid methods, you can for example easily check how far you
70// are from the envelope boundary.
71//
72
74{
75 public: // without description
76 //------------------------
77 // Constructor/Destructor
78 //------------------------
79 // Only one Constructor. By default the envelope can
80 // be placed n-Times. If the user is sure that it'll be
81 // placed just one time, the IsUnique flag should be set
82 // TRUE to avoid the G4AffineTransform re-calculations each
83 // time we reach the envelope.
84 G4FastTrack(G4Envelope* anEnvelope, G4bool IsUnique);
85 ~G4FastTrack() = default;
86
87 //------------------------------------------------------------
88 // The fast simulation manager uses the SetCurrentTrack
89 // method to setup the current G4FastTrack object
90 //------------------------------------------------------------
91 void SetCurrentTrack(const G4Track&, const G4Navigator* a = nullptr);
92
93 //------------------------------------------------------------
94 // The fast simulation manager uses the OnTheBoundaryButExiting
95 // method to test if the particle is leaving the envelope.
96 //------------------------------------------------------------
98
99 //----------------------------------
100 // Informations useful to the user :
101 // General public get functions.
102 //----------------------------------
103
104 // Returns the current G4Track.
105 const G4Track* GetPrimaryTrack() const;
106
107 // Returns the Envelope G4Region pointer.
108 G4Envelope* GetEnvelope() const;
109
110 // Returns the Envelope G4LogicalVolume pointer.
112
113 // Returns the Envelope G4VPhysicalVolume pointer.
115
116 // Returns the Envelope G4VSolid pointer.
117 G4VSolid* GetEnvelopeSolid() const;
118
119 //-----------------------------------
120 // Primary track informations in the
121 // Envelope coordinate system.
122 //-----------------------------------
123
124 // Returns the particle position in envelope coordinates.
126
127 // Returns the particle momentum in envelope coordinates.
129
130 // Returns the particle direction in envelope coordinates.
132
133 // Returns the particle polarization in envelope coordinates.
135
136 //------------------------------------
137 // 3D transformation of the envelope:
138 //------------------------------------
139
140 // Returns the envelope Global -> Local G4AffineTransform
142
143 // Returns the envelope Local -> Global G4AffineTransform
145
146 private:
147 //-----------------
148 // Private members
149 //-----------------
150 // Current G4Track pointer
151 const G4Track* fTrack{nullptr};
152
153 //------------------------------------------------
154 // Records the Affine/InverseAffine transformation
155 // of the envelope.
156 //------------------------------------------------
157 void FRecordsAffineTransformation(const G4Navigator*);
158 G4bool fAffineTransformationDefined{false};
159 G4Envelope* fEnvelope;
160 G4bool fIsUnique;
161 G4LogicalVolume* fEnvelopeLogicalVolume{nullptr};
162 G4VPhysicalVolume* fEnvelopePhysicalVolume{nullptr};
163 G4VSolid* fEnvelopeSolid{nullptr};
164 G4ThreeVector fLocalTrackPosition, fLocalTrackMomentum, fLocalTrackDirection,
165 fLocalTrackPolarization;
166 G4AffineTransform fAffineTransformation, fInverseAffineTransformation;
167};
168
169// -----------------
170// -- Inline methods
171// -----------------
172
174{
175 return fEnvelope;
176}
177
179{
180 return fEnvelopeLogicalVolume;
181}
182
184{
185 return fEnvelopePhysicalVolume;
186}
187
189{
190 return fEnvelopeSolid;
191}
192
194{
195 return fTrack;
196}
197
199{
200 return fLocalTrackPosition;
201}
202
204{
205 return fLocalTrackMomentum;
206}
207
209{
210 return fLocalTrackDirection;
211}
212
214{
215 return fLocalTrackPolarization;
216}
217
219{
220 return &fAffineTransformation;
221}
222
224{
225 return &fInverseAffineTransformation;
226}
227
229{
230 // tests if particle are on the boundary and leaving.
233 == 0.;
234}
235
236#endif
bool G4bool
Definition G4Types.hh:86
G4ThreeVector GetPrimaryTrackLocalPosition() const
G4Envelope * GetEnvelope() const
const G4Track * GetPrimaryTrack() const
G4ThreeVector GetPrimaryTrackLocalPolarization() const
void SetCurrentTrack(const G4Track &, const G4Navigator *a=nullptr)
const G4AffineTransform * GetInverseAffineTransformation() const
G4VPhysicalVolume * GetEnvelopePhysicalVolume() const
G4ThreeVector GetPrimaryTrackLocalDirection() const
G4ThreeVector GetPrimaryTrackLocalMomentum() const
G4FastTrack(G4Envelope *anEnvelope, G4bool IsUnique)
G4VSolid * GetEnvelopeSolid() const
G4bool OnTheBoundaryButExiting() const
G4LogicalVolume * GetEnvelopeLogicalVolume() const
~G4FastTrack()=default
const G4AffineTransform * GetAffineTransformation() const
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0