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