G4TouchableHistory Class Reference

#include <G4TouchableHistory.hh>

Inheritance diagram for G4TouchableHistory:

G4VTouchable

Public Member Functions

 G4TouchableHistory (const G4NavigationHistory &history)
 G4TouchableHistory ()
 ~G4TouchableHistory ()
G4VPhysicalVolumeGetVolume (G4int depth=0) const
G4VSolidGetSolid (G4int depth=0) const
const G4ThreeVectorGetTranslation (G4int depth=0) const
const G4RotationMatrixGetRotation (G4int depth=0) const
G4int GetReplicaNumber (G4int depth=0) const
G4int GetHistoryDepth () const
G4int MoveUpHistory (G4int num_levels=1)
void UpdateYourself (G4VPhysicalVolume *pPhysVol, const G4NavigationHistory *history=0)
const G4NavigationHistoryGetHistory () const
void * operator new (size_t)
void operator delete (void *aTH)

Detailed Description

Definition at line 51 of file G4TouchableHistory.hh.


Constructor & Destructor Documentation

G4TouchableHistory::G4TouchableHistory ( const G4NavigationHistory history  ) 

Definition at line 47 of file G4TouchableHistory.cc.

References G4NavigationHistory::GetTopTransform(), and G4AffineTransform::Inverse().

00048   : fhistory(history)
00049 {
00050   G4AffineTransform tf(fhistory.GetTopTransform().Inverse());
00051   ftlate = tf.NetTranslation();
00052   frot = tf.NetRotation();
00053 }

G4TouchableHistory::G4TouchableHistory (  ) 

Definition at line 38 of file G4TouchableHistory.cc.

References G4NavigationHistory::SetFirstEntry().

00039   : frot(G4RotationMatrix()),
00040     ftlate(G4ThreeVector(0.,0.,0.)),
00041     fhistory()
00042 {
00043    G4VPhysicalVolume* pPhysVol=0;
00044    fhistory.SetFirstEntry(pPhysVol);
00045 }

G4TouchableHistory::~G4TouchableHistory (  ) 

Definition at line 55 of file G4TouchableHistory.cc.

00056 {
00057 }


Member Function Documentation

const G4NavigationHistory * G4TouchableHistory::GetHistory (  )  const [inline, virtual]

Reimplemented from G4VTouchable.

Definition at line 110 of file G4TouchableHistory.icc.

Referenced by G4PathFinder::CreateTouchableHandle(), G4MultiNavigator::CreateTouchableHistoryHandle(), G4Navigator::ResetHierarchyAndLocate(), and G4ITNavigator::ResetHierarchyAndLocate().

00111 {
00112   return &fhistory;
00113 }

G4int G4TouchableHistory::GetHistoryDepth (  )  const [inline, virtual]

Reimplemented from G4VTouchable.

Definition at line 83 of file G4TouchableHistory.icc.

References G4NavigationHistory::GetDepth().

00084 {
00085   return  fhistory.GetDepth();
00086 }

G4int G4TouchableHistory::GetReplicaNumber ( G4int  depth = 0  )  const [inline, virtual]

Reimplemented from G4VTouchable.

Definition at line 77 of file G4TouchableHistory.icc.

References G4NavigationHistory::GetReplicaNo().

Referenced by G4VPrimitiveScorer::GetIndex().

00078 {
00079   return fhistory.GetReplicaNo(CalculateHistoryIndex(depth));
00080 }

const G4RotationMatrix * G4TouchableHistory::GetRotation ( G4int  depth = 0  )  const [virtual]

Implements G4VTouchable.

Definition at line 79 of file G4TouchableHistory.cc.

References G4NavigationHistory::GetTransform(), and G4AffineTransform::NetRotation().

00080 {
00081   // The value returned will change at the next call
00082   // Copy it if you want to use it!
00083   //
00084   static G4RotationMatrix rotM;
00085 
00086   if(depth==0.0)
00087   {
00088     return &frot;
00089   }
00090   else
00091   {
00092     rotM = fhistory.GetTransform(CalculateHistoryIndex(depth)).NetRotation();
00093     return &rotM;
00094   }
00095 }

G4VSolid * G4TouchableHistory::GetSolid ( G4int  depth = 0  )  const [inline, virtual]

Reimplemented from G4VTouchable.

Definition at line 70 of file G4TouchableHistory.icc.

References G4VPhysicalVolume::GetLogicalVolume(), G4LogicalVolume::GetSolid(), and G4NavigationHistory::GetVolume().

00071 {
00072   return fhistory.GetVolume(CalculateHistoryIndex(depth))
00073                             ->GetLogicalVolume()->GetSolid();
00074 }

const G4ThreeVector & G4TouchableHistory::GetTranslation ( G4int  depth = 0  )  const [virtual]

Implements G4VTouchable.

Definition at line 60 of file G4TouchableHistory.cc.

References G4NavigationHistory::GetTransform(), and G4AffineTransform::NetTranslation().

00061 {
00062   // The value returned will change at the next call
00063   // Copy it if you want to use it!
00064   //
00065   static G4ThreeVector currTranslation;
00066   if(depth==0.0)
00067   {
00068     return ftlate;
00069   }
00070   else
00071   {
00072     currTranslation =
00073       fhistory.GetTransform(CalculateHistoryIndex(depth)).NetTranslation();
00074     return currTranslation;
00075   }
00076 }

G4VPhysicalVolume * G4TouchableHistory::GetVolume ( G4int  depth = 0  )  const [inline, virtual]

Reimplemented from G4VTouchable.

Definition at line 64 of file G4TouchableHistory.icc.

References G4NavigationHistory::GetVolume().

Referenced by G4VReadOutGeometry::FindROTouchable().

00065 {
00066   return fhistory.GetVolume(CalculateHistoryIndex(depth));
00067 }

G4int G4TouchableHistory::MoveUpHistory ( G4int  num_levels = 1  )  [inline, virtual]

Reimplemented from G4VTouchable.

Definition at line 89 of file G4TouchableHistory.icc.

References G4NavigationHistory::BackLevel(), and G4NavigationHistory::GetDepth().

Referenced by G4Navigator::SetupHierarchy(), and G4ITNavigator::SetupHierarchy().

00090 {
00091   G4int maxLevelsMove = fhistory.GetDepth();
00092   G4int minLevelsMove = 0;              // Cannot redescend today!
00093                                         // Soon it will be possible
00094                                         // by adding a data member here
00095                                         //     fCurrentDepth;
00096   if( num_levels > maxLevelsMove )
00097   {
00098     num_levels = maxLevelsMove;
00099   }
00100   else if( num_levels < minLevelsMove )
00101   {
00102     num_levels = minLevelsMove;
00103   }
00104   fhistory.BackLevel( num_levels ); 
00105 
00106   return num_levels;
00107 }

void G4TouchableHistory::operator delete ( void *  aTH  )  [inline]

Definition at line 125 of file G4TouchableHistory.icc.

References aTouchableHistoryAllocator.

00126 {
00127   aTouchableHistoryAllocator.FreeSingle((G4TouchableHistory *) aTH);
00128 }

void * G4TouchableHistory::operator new ( size_t   )  [inline]

Definition at line 119 of file G4TouchableHistory.icc.

References aTouchableHistoryAllocator.

00120 {
00121   return (void *) aTouchableHistoryAllocator.MallocSingle();
00122 }

void G4TouchableHistory::UpdateYourself ( G4VPhysicalVolume pPhysVol,
const G4NavigationHistory history = 0 
) [inline, virtual]

Reimplemented from G4VTouchable.

Definition at line 40 of file G4TouchableHistory.icc.

References G4NavigationHistory::GetTopTransform(), G4AffineTransform::Inverse(), and G4NavigationHistory::SetFirstEntry().

Referenced by G4PathFinder::CreateTouchableHandle(), and G4MultiNavigator::CreateTouchableHistoryHandle().

00042 {
00043   fhistory = *pHistory;  
00044   G4AffineTransform tf(fhistory.GetTopTransform().Inverse());
00045   if( pPhysVol == 0 )
00046   {
00047     // This means that the track has left the World Volume.
00048     // Since the Navigation History does not already reflect this,
00049     // we must correct this problem here.
00050     //
00051     fhistory.SetFirstEntry(pPhysVol);
00052   }
00053   ftlate = tf.NetTranslation();
00054   frot = tf.NetRotation();
00055 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:32 2013 for Geant4 by  doxygen 1.4.7