Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups 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 // $Id: G4FastTrack.cc 100945 2016-11-03 11:26:19Z gcosmo $
28 //
29 //---------------------------------------------------------------
30 //
31 // G4FastTrack.cc
32 //
33 // Description:
34 // Keeps the current track information and special features
35 // for Parameterised Simulation Models.
36 //
37 // History:
38 // Oct 97: Verderi && MoraDeFreitas - First Implementation.
39 //
40 //---------------------------------------------------------------
41 
42 #include "G4ios.hh"
43 #include "G4FastTrack.hh"
46 
47 // -----------
48 // Constructor
49 // -----------
50 //
52  : fTrack ( nullptr ),
53  fAffineTransformationDefined( false ),
54  fEnvelope ( anEnvelope ),
55  fIsUnique ( IsUnique ),
56  fEnvelopeLogicalVolume ( nullptr ),
57  fEnvelopePhysicalVolume ( nullptr ),
58  fEnvelopeSolid ( nullptr )
59 {}
60 
61 // -----------
62 // Destructor:
63 // -----------
65 {}
66 
67 //------------------------------------------------------------
68 // The parameterised simulation manager uses the SetCurrentTrack
69 // method to setup the current G4FastTrack object
70 //------------------------------------------------------------
72  const G4Navigator* theNavigator)
73 {
74 
75  // -- Register track pointer (used everywhere):
76  fTrack = &track;
77 
78  //-----------------------------------------------------
79  // First time the track enters the volume or if the
80  // Logical Volume was placed n-Times in the geometry :
81  //
82  // Records the Rotation+Translation for the Envelope !
83  // When the particle is inside or on the boundary, the
84  // NavigationHistory IS UP TO DATE.
85  //------------------------------------------------------
86  if (!fAffineTransformationDefined || !fIsUnique) FRecordsAffineTransformation(theNavigator);
87 
88  //-------------------------------------------
89  // Records local position/momentum/direction
90  // of the Track.
91  // They are accessible to the user through a
92  // set of Get functions and should be useful
93  // to decide to trigger or not.
94  //-------------------------------------------
95  // -- local position:
96  fLocalTrackPosition = fAffineTransformation.TransformPoint(fTrack->GetPosition());
97  // -- local momentum:
98  fLocalTrackMomentum = fAffineTransformation.TransformAxis(fTrack->GetMomentum());
99  // -- local direction:
100  fLocalTrackDirection = fLocalTrackMomentum.unit();
101  // -- local polarization:
102  fLocalTrackPolarization = fAffineTransformation.TransformAxis(fTrack->GetPolarization());
103 }
104 
105 //------------------------------------
106 //
107 // 3D transformation of the envelope
108 // This is Done only one time.
109 //
110 //------------------------------------
111 void
112 G4FastTrack::FRecordsAffineTransformation(const G4Navigator* theNavigator)
113 {
114 
115  //--------------------------------------------------------
116  // Get the touchable history which represents the current
117  // volume hierachy the particle is in.
118  // Note that TouchableHistory allocated by the Navigator
119  // must be deleted by G4FastTrack.
120  //--------------------------------------------------------
121  const G4Navigator* NavigatorToUse;
122  if(theNavigator != 0 ) NavigatorToUse = theNavigator;
124 
125  G4TouchableHistoryHandle history = NavigatorToUse->CreateTouchableHistoryHandle();
126 
127  //-----------------------------------------------------
128  // Run accross the hierarchy to find the physical volume
129  // associated with the envelope
130  //-----------------------------------------------------
131  G4int depth = history->GetHistory()->GetDepth();
132  G4int idepth, Done = 0;
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 = 1;
143  break;
144  }
145  }
146  //---------------------------------------------
147  //-- Verification: should be removed in future:
148  //---------------------------------------------
149  if ( !Done )
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 }
const G4ThreeVector & GetPolarization() const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4AffineTransform Inverse() const
virtual G4TouchableHistoryHandle CreateTouchableHistoryHandle() const
G4int GetDepth() const
const G4ThreeVector & GetPosition() const
G4Navigator * GetNavigatorForTracking() const
G4VSolid * GetSolid() const
G4Region * GetRegion() const
int G4int
Definition: G4Types.hh:78
G4FastTrack(G4Envelope *anEnvelope, G4bool IsUnique)
Definition: G4FastTrack.cc:51
const G4String & GetName() const
bool G4bool
Definition: G4Types.hh:79
G4bool IsRootRegion() const
void SetCurrentTrack(const G4Track &, const G4Navigator *a=0)
Definition: G4FastTrack.cc:71
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4TransportationManager * GetTransportationManager()
G4ThreeVector GetMomentum() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4LogicalVolume * GetLogicalVolume() const
const G4AffineTransform & GetTransform(G4int n) const
Hep3Vector unit() const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
#define G4endl
Definition: G4ios.hh:61
const G4NavigationHistory * GetHistory() const
G4VPhysicalVolume * GetVolume(G4int n) const