Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4FastTrack.cc
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.cc
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#include "G4ios.hh"
42#include "G4FastTrack.hh"
45
46// -----------
47// Constructor
48// -----------
49//
51 : fTrack ( nullptr ),
52 fAffineTransformationDefined( false ),
53 fEnvelope ( anEnvelope ),
54 fIsUnique ( IsUnique ),
55 fEnvelopeLogicalVolume ( nullptr ),
56 fEnvelopePhysicalVolume ( nullptr ),
57 fEnvelopeSolid ( nullptr )
58{}
59
60// -----------
61// Destructor:
62// -----------
64{}
65
66//------------------------------------------------------------
67// The parameterised simulation manager uses the SetCurrentTrack
68// method to setup the current G4FastTrack object
69//------------------------------------------------------------
71 const G4Navigator* theNavigator)
72{
73
74 // -- Register track pointer (used everywhere):
75 fTrack = &track;
76
77 //-----------------------------------------------------
78 // First time the track enters the volume or if the
79 // Logical Volume was placed n-Times in the geometry :
80 //
81 // Records the Rotation+Translation for the Envelope !
82 // When the particle is inside or on the boundary, the
83 // NavigationHistory IS UP TO DATE.
84 //------------------------------------------------------
85 if (!fAffineTransformationDefined || !fIsUnique) FRecordsAffineTransformation(theNavigator);
86
87 //-------------------------------------------
88 // Records local position/momentum/direction
89 // of the Track.
90 // They are accessible to the user through a
91 // set of Get functions and should be useful
92 // to decide to trigger or not.
93 //-------------------------------------------
94 // -- local position:
95 fLocalTrackPosition = fAffineTransformation.TransformPoint(fTrack->GetPosition());
96 // -- local momentum:
97 fLocalTrackMomentum = fAffineTransformation.TransformAxis(fTrack->GetMomentum());
98 // -- local direction:
99 fLocalTrackDirection = fLocalTrackMomentum.unit();
100 // -- local polarization:
101 fLocalTrackPolarization = fAffineTransformation.TransformAxis(fTrack->GetPolarization());
102}
103
104//------------------------------------
105//
106// 3D transformation of the envelope
107// This is Done only one time.
108//
109//------------------------------------
110void
111G4FastTrack::FRecordsAffineTransformation(const G4Navigator* theNavigator)
112{
113
114 //--------------------------------------------------------
115 // Get the touchable history which represents the current
116 // volume hierachy the particle is in.
117 // Note that TouchableHistory allocated by the Navigator
118 // must be deleted by G4FastTrack.
119 //--------------------------------------------------------
120 const G4Navigator* NavigatorToUse;
121 if(theNavigator != nullptr ) NavigatorToUse = theNavigator;
123
124 G4TouchableHistoryHandle history = NavigatorToUse->CreateTouchableHistoryHandle();
125
126 //-----------------------------------------------------
127 // Run accross the hierarchy to find the physical volume
128 // associated with the envelope
129 //-----------------------------------------------------
130 G4int depth = (G4int)history->GetHistory()->GetDepth();
131 G4int idepth;
132 G4bool Done = false;
133 for (idepth = 0; idepth <= depth; ++idepth)
134 {
135 G4VPhysicalVolume* currPV = history->GetHistory()->GetVolume(idepth);
136 G4LogicalVolume* currLV = currPV->GetLogicalVolume();
137 if ( (currLV->GetRegion() == fEnvelope) && (currLV->IsRootRegion()) )
138 {
139 fEnvelopePhysicalVolume = currPV;
140 fEnvelopeLogicalVolume = currLV;
141 fEnvelopeSolid = currLV->GetSolid();
142 Done = true;
143 break;
144 }
145 }
146 //---------------------------------------------
147 //-- Verification: should be removed in future:
148 //---------------------------------------------
149 if ( Done == false )
150 {
152 ed << "Can't find transformation for `" << fEnvelopePhysicalVolume->GetName() << "'" << G4endl;
153 G4Exception("G4FastTrack::FRecordsAffineTransformation()",
154 "FastSim011",
155 JustWarning, ed);
156 }
157 else
158 {
159 //-------------------------------------------------------
160 // Records the transformation and inverse transformation:
161 //-------------------------------------------------------
162 fAffineTransformation = history->GetHistory()->GetTransform(idepth);
163 fInverseAffineTransformation = fAffineTransformation.Inverse();
164
165 fAffineTransformationDefined = true;
166 }
167}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
Hep3Vector unit() const
G4AffineTransform Inverse() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
G4FastTrack(G4Envelope *anEnvelope, G4bool IsUnique)
Definition: G4FastTrack.cc:50
void SetCurrentTrack(const G4Track &, const G4Navigator *a=0)
Definition: G4FastTrack.cc:70
G4VSolid * GetSolid() const
G4bool IsRootRegion() const
G4Region * GetRegion() const
virtual G4TouchableHistoryHandle CreateTouchableHistoryHandle() const
const G4ThreeVector & GetPosition() const
G4ThreeVector GetMomentum() const
const G4ThreeVector & GetPolarization() const
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const