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