Geant4 11.1.1
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
42#ifndef G4FastTrack_h
43#define G4FastTrack_h
44
45#include "G4VSolid.hh"
46#include "G4LogicalVolume.hh"
47#include "G4Region.hh"
48#include "G4AffineTransform.hh"
49#include "G4Track.hh"
50#include "G4Navigator.hh"
51
52//---------------------------
53// For possible future needs:
54//---------------------------
56
57
58//-------------------------------------------
59//
60// G4FastTrack class
61//
62//-------------------------------------------
63
64// Class Description:
65// The G4FastTrack provides you access to the current G4Track,
66// gives simple access to envelope related features (G4Region,
67// G4LogicalVolume, G4VSolid, G4AffineTransform references between
68// the global and the envelope local coordinates systems) and
69// simple access to the position, momentum expressed in the
70// envelope coordinate system. Using those quantities and the
71// G4VSolid methods, you can for example easily check how far you
72// are from the envelope boundary.
73//
74
75
77{
78public: // without description
79 //------------------------
80 // Constructor/Destructor
81 //------------------------
82 // Only one Constructor. By default the envelope can
83 // be placed n-Times. If the user is sure that it'll be
84 // placed just one time, the IsUnique flag should be set
85 // TRUE to avoid the G4AffineTransform re-calculations each
86 // time we reach the envelope.
87 G4FastTrack(G4Envelope *anEnvelope,
88 G4bool IsUnique);
90
91 //------------------------------------------------------------
92 // The fast simulation manager uses the SetCurrentTrack
93 // method to setup the current G4FastTrack object
94 //------------------------------------------------------------
95 void SetCurrentTrack(const G4Track&, const G4Navigator* a = 0);
96
97 //------------------------------------------------------------
98 // The fast simulation manager uses the OnTheBoundaryButExiting
99 // method to test if the particle is leaving the envelope.
100 //------------------------------------------------------------
102
103 //----------------------------------
104 // Informations useful to the user :
105 // General public get functions.
106 //----------------------------------
107
108public: // with Description
109
110 const G4Track* GetPrimaryTrack() const;
111 // Returns the current G4Track.
112
113 G4Envelope* GetEnvelope() const;
114 // Returns the Envelope G4Region pointer.
115
117 // Returns the Envelope G4LogicalVolume pointer.
118
120 // Returns the Envelope G4VPhysicalVolume pointer.
121
122 G4VSolid* GetEnvelopeSolid() const;
123 // Returns the Envelope G4VSolid pointer.
124
125 //-----------------------------------
126 // Primary track informations in the
127 // Envelope coordinate system.
128 //-----------------------------------
129
131 // Returns the particle position in envelope coordinates.
132
134 // Returns the particle momentum in envelope coordinates.
135
137 // Returns the particle direction in envelope coordinates.
138
140 // Returns the particle polarization in envelope coordinates.
141
142 //------------------------------------
143 // 3D transformation of the envelope:
144 //------------------------------------
145 // Global -> Local
146
148 // Returns the envelope Global -> Local G4AffineTransform
149
150 // Local -> Global
152 // Returns the envelope Local -> Global G4AffineTransform
153
154 //-----------------
155 // Private members
156 //-----------------
157private:
158
159 // Current G4Track pointer
160 const G4Track* fTrack;
161
162 //------------------------------------------------
163 // Records the Affine/InverseAffine transformation
164 // of the envelope.
165 //------------------------------------------------
166 void FRecordsAffineTransformation(const G4Navigator*);
167 G4bool fAffineTransformationDefined;
168 G4Envelope* fEnvelope;
169 G4bool fIsUnique;
170 G4LogicalVolume* fEnvelopeLogicalVolume;
171 G4VPhysicalVolume* fEnvelopePhysicalVolume;
172 G4VSolid* fEnvelopeSolid;
173 G4ThreeVector fLocalTrackPosition,
174 fLocalTrackMomentum,
175 fLocalTrackDirection,
176 fLocalTrackPolarization;
177 G4AffineTransform fAffineTransformation,
178 fInverseAffineTransformation;
179};
180
181
182// -----------------
183// -- Inline methods
184// -----------------
185
187{
188 return fEnvelope;
189}
190
192{
193 return fEnvelopeLogicalVolume;
194}
195
197{
198 return fEnvelopePhysicalVolume;
199}
200
202{
203 return fEnvelopeSolid;
204}
205
207{
208 return fTrack;
209}
210
212{
213 return fLocalTrackPosition;
214}
215
217{
218 return fLocalTrackMomentum;
219}
220
222{
223 return fLocalTrackDirection;
224}
225
227{
228 return fLocalTrackPolarization;
229}
230
232{
233 return &fAffineTransformation;
234}
235
237{
238 return &fInverseAffineTransformation;
239}
240
242{
243 // tests if particle are on the boundary and leaving.
244 return GetEnvelopeSolid()->
245 DistanceToOut(GetPrimaryTrackLocalPosition(),
247}
248
249#endif
G4Region G4Envelope
Definition: G4FastTrack.hh:55
bool G4bool
Definition: G4Types.hh:86
G4ThreeVector GetPrimaryTrackLocalPosition() const
Definition: G4FastTrack.hh:211
G4Envelope * GetEnvelope() const
Definition: G4FastTrack.hh:186
const G4Track * GetPrimaryTrack() const
Definition: G4FastTrack.hh:206
G4ThreeVector GetPrimaryTrackLocalPolarization() const
Definition: G4FastTrack.hh:226
const G4AffineTransform * GetInverseAffineTransformation() const
Definition: G4FastTrack.hh:236
G4VPhysicalVolume * GetEnvelopePhysicalVolume() const
Definition: G4FastTrack.hh:196
G4ThreeVector GetPrimaryTrackLocalDirection() const
Definition: G4FastTrack.hh:221
G4ThreeVector GetPrimaryTrackLocalMomentum() const
Definition: G4FastTrack.hh:216
G4VSolid * GetEnvelopeSolid() const
Definition: G4FastTrack.hh:201
G4bool OnTheBoundaryButExiting() const
Definition: G4FastTrack.hh:241
G4LogicalVolume * GetEnvelopeLogicalVolume() const
Definition: G4FastTrack.hh:191
void SetCurrentTrack(const G4Track &, const G4Navigator *a=0)
Definition: G4FastTrack.cc:70
const G4AffineTransform * GetAffineTransformation() const
Definition: G4FastTrack.hh:231