Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ITNavigator1 Class Reference

#include <G4ITNavigator1.hh>

Collaboration diagram for G4ITNavigator1:

Public Member Functions

 G4ITNavigator1 ()
 
virtual ~G4ITNavigator1 ()
 
G4ITNavigatorState_Lock1GetNavigatorState ()
 
void SetNavigatorState (G4ITNavigatorState_Lock1 *)
 
void NewNavigatorState ()
 
virtual G4double ComputeStep (const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
 
G4double CheckNextStep (const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
 
virtual G4VPhysicalVolumeResetHierarchyAndLocate (const G4ThreeVector &point, const G4ThreeVector &direction, const G4TouchableHistory &h)
 
virtual G4VPhysicalVolumeLocateGlobalPointAndSetup (const G4ThreeVector &point, const G4ThreeVector *direction=0, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
 
virtual void LocateGlobalPointWithinVolume (const G4ThreeVector &position)
 
void LocateGlobalPointAndUpdateTouchableHandle (const G4ThreeVector &position, const G4ThreeVector &direction, G4TouchableHandle &oldTouchableToUpdate, const G4bool RelativeSearch=true)
 
void LocateGlobalPointAndUpdateTouchable (const G4ThreeVector &position, const G4ThreeVector &direction, G4VTouchable *touchableToUpdate, const G4bool RelativeSearch=true)
 
void LocateGlobalPointAndUpdateTouchable (const G4ThreeVector &position, G4VTouchable *touchableToUpdate, const G4bool RelativeSearch=true)
 
void SetGeometricallyLimitedStep ()
 
virtual G4double ComputeSafety (const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
 
G4VPhysicalVolumeGetWorldVolume () const
 
void SetWorldVolume (G4VPhysicalVolume *pWorld)
 
G4GRSVolumeCreateGRSVolume () const
 
G4GRSSolidCreateGRSSolid () const
 
G4TouchableHistoryCreateTouchableHistory () const
 
G4TouchableHistoryCreateTouchableHistory (const G4NavigationHistory *) const
 
virtual G4TouchableHistoryHandle CreateTouchableHistoryHandle () const
 
virtual G4ThreeVector GetLocalExitNormal (G4bool *valid)
 
virtual G4ThreeVector GetLocalExitNormalAndCheck (const G4ThreeVector &point, G4bool *valid)
 
virtual G4ThreeVector GetGlobalExitNormal (const G4ThreeVector &point, G4bool *valid)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int level)
 
G4bool IsActive () const
 
void Activate (G4bool flag)
 
G4bool EnteredDaughterVolume () const
 
G4bool ExitedMotherVolume () const
 
void CheckMode (G4bool mode)
 
G4bool IsCheckModeActive () const
 
void SetPushVerbosity (G4bool mode)
 
void PrintState () const
 
const G4AffineTransformGetGlobalToLocalTransform () const
 
const G4AffineTransform GetLocalToGlobalTransform () const
 
G4AffineTransform GetMotherToDaughterTransform (G4VPhysicalVolume *dVolume, G4int dReplicaNo, EVolume dVolumeType)
 
void ResetStackAndState ()
 
G4int SeverityOfZeroStepping (G4int *noZeroSteps) const
 
void SetSavedState ()
 
void RestoreSavedState ()
 
G4ThreeVector GetCurrentLocalCoordinate () const
 
G4ThreeVector NetTranslation () const
 
G4RotationMatrix NetRotation () const
 
void EnableBestSafety (G4bool value=false)
 
virtual void ResetState ()
 

Static Public Attributes

static const G4int fMaxNav = 8
 

Protected Member Functions

G4ThreeVector ComputeLocalPoint (const G4ThreeVector &rGlobPoint) const
 
G4ThreeVector ComputeLocalAxis (const G4ThreeVector &pVec) const
 
EVolume VolumeType (const G4VPhysicalVolume *pVol) const
 
EVolume CharacteriseDaughters (const G4LogicalVolume *pLog) const
 
G4int GetDaughtersRegularStructureId (const G4LogicalVolume *pLog) const
 
virtual void SetupHierarchy ()
 

Protected Attributes

G4double kCarTolerance
 
G4NavigationHistory fHistory
 
G4bool fEnteredDaughter
 
G4bool fExitedMother
 
G4bool fWasLimitedByGeometry
 
G4ThreeVector fStepEndPoint
 
G4ThreeVector fLastStepEndPointLocal
 
G4int fVerbose
 

Friends

std::ostream & operator<< (std::ostream &os, const G4ITNavigator1 &n)
 

Detailed Description

Definition at line 91 of file G4ITNavigator1.hh.

Constructor & Destructor Documentation

G4ITNavigator1::G4ITNavigator1 ( )

Definition at line 60 of file G4ITNavigator1.cc.

61  : fWasLimitedByGeometry(false), fVerbose(0),
62  fTopPhysical(0), fCheck(false), fPushed(false), fWarnPush(true)
63 {
64  fActive= false;
65  fLastTriedStepComputation= false;
66 
68  // Initialises also all
69  // - exit / entry flags
70  // - flags & variables for exit normals
71  // - zero step counters
72  // - blocked volume
73 
74  fActionThreshold_NoZeroSteps = 10;
75  fAbandonThreshold_NoZeroSteps = 25;
76 
78  fregularNav.SetNormalNavigation( &fnormalNav );
79 
82 
83  fpVoxelSafety= new G4VoxelSafety();
84  fpSaveState = 0;
85 
86  // this->SetVerboseLevel(3);
87  // this->CheckMode(true);
88 }
static const G4double kInfinity
Definition: geomdefs.hh:42
CLHEP::Hep3Vector G4ThreeVector
G4double GetSurfaceTolerance() const
G4ThreeVector fLastStepEndPointLocal
void ResetStackAndState()
void SetNormalNavigation(G4NormalNavigation *fnormnav)
G4ThreeVector fStepEndPoint
G4double kCarTolerance
static G4GeometryTolerance * GetInstance()
G4bool fWasLimitedByGeometry

Here is the call graph for this function:

G4ITNavigator1::~G4ITNavigator1 ( )
virtual

Definition at line 124 of file G4ITNavigator1.cc.

125 { delete fpVoxelSafety; }

Member Function Documentation

void G4ITNavigator1::Activate ( G4bool  flag)
inline
EVolume G4ITNavigator1::CharacteriseDaughters ( const G4LogicalVolume pLog) const
inlineprotected

Here is the caller graph for this function:

void G4ITNavigator1::CheckMode ( G4bool  mode)
inline
G4double G4ITNavigator1::CheckNextStep ( const G4ThreeVector pGlobalPoint,
const G4ThreeVector pDirection,
const G4double  pCurrentProposedStepLength,
G4double pNewSafety 
)

Definition at line 1230 of file G4ITNavigator1.cc.

1234 {
1235  G4double step;
1236 
1237  // Save the state, for this parasitic call
1238  //
1239  SetSavedState();
1240 
1241  step = ComputeStep ( pGlobalpoint,
1242  pDirection,
1243  pCurrentProposedStepLength,
1244  pNewSafety );
1245 
1246  // If a parasitic call, then attempt to restore the key parts of the state
1247  //
1248  RestoreSavedState();
1249 
1250  return step;
1251 }
void RestoreSavedState()
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4ThreeVector G4ITNavigator1::ComputeLocalAxis ( const G4ThreeVector pVec) const
inlineprotected

Here is the caller graph for this function:

G4ThreeVector G4ITNavigator1::ComputeLocalPoint ( const G4ThreeVector rGlobPoint) const
inlineprotected

Here is the caller graph for this function:

G4double G4ITNavigator1::ComputeSafety ( const G4ThreeVector globalpoint,
const G4double  pProposedMaxLength = DBL_MAX,
const G4bool  keepState = true 
)
virtual

Definition at line 1692 of file G4ITNavigator1.cc.

1695 {
1696  G4double newSafety = 0.0;
1697 
1698 #ifdef G4DEBUG_NAVIGATION
1699  G4int oldcoutPrec = G4cout.precision(8);
1700  if( fVerbose > 0 )
1701  {
1702  G4cout << "*** G4ITNavigator1::ComputeSafety: ***" << G4endl
1703  << " Called at point: " << pGlobalpoint << G4endl;
1704 
1705  G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
1706  G4cout << " Volume = " << motherPhysical->GetName()
1707  << " - Maximum length = " << pMaxLength << G4endl;
1708  if( fVerbose >= 4 )
1709  {
1710  G4cout << " ----- Upon entering Compute Safety:" << G4endl;
1711  PrintState();
1712  }
1713  }
1714 #endif
1715 
1716  if (keepState) { SetSavedState(); }
1717 
1718  G4double distEndpointSq = (pGlobalpoint-fStepEndPoint).mag2();
1719  G4bool stayedOnEndpoint = distEndpointSq < kCarTolerance*kCarTolerance;
1720  G4bool endpointOnSurface = fEnteredDaughter || fExitedMother;
1721 
1722  if( !(endpointOnSurface && stayedOnEndpoint) )
1723  {
1724  // Pseudo-relocate to this point (updates voxel information only)
1725  //
1726  LocateGlobalPointWithinVolume( pGlobalpoint );
1727  // --->> DANGER: Side effects on sub-navigator voxel information <<---
1728  // Could be replaced again by 'granular' calls to sub-navigator
1729  // locates (similar side-effects, but faster.
1730  // Solutions:
1731  // 1) Re-locate (to where?)
1732  // 2) Insure that the methods using (G4ComputeStep?)
1733  // does a relocation (if information is disturbed only ?)
1734 
1735 #ifdef G4DEBUG_NAVIGATION
1736  if( fVerbose >= 2 )
1737  {
1738  G4cout << " G4ITNavigator1::ComputeSafety() relocates-in-volume to point: "
1739  << pGlobalpoint << G4endl;
1740  }
1741 #endif
1742  G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
1743  G4LogicalVolume *motherLogical = motherPhysical->GetLogicalVolume();
1744  G4SmartVoxelHeader* pVoxelHeader = motherLogical->GetVoxelHeader();
1745  G4ThreeVector localPoint = ComputeLocalPoint(pGlobalpoint);
1746 
1748  {
1749  switch(CharacteriseDaughters(motherLogical))
1750  {
1751  case kNormal:
1752  if ( pVoxelHeader )
1753  {
1754 #ifdef G4NEW_SAFETY
1755  G4double safetyTwo = fpVoxelSafety->ComputeSafety(localPoint,
1756  *motherPhysical, pMaxLength);
1757  newSafety= safetyTwo; // Faster and best available
1758 #else
1759  G4double safetyOldVoxel;
1760  safetyOldVoxel =
1761  fvoxelNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1762  newSafety= safetyOldVoxel;
1763 #endif
1764  }
1765  else
1766  {
1767  newSafety=fnormalNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1768  }
1769  break;
1770  case kParameterised:
1771  if( GetDaughtersRegularStructureId(motherLogical) != 1 )
1772  {
1773  newSafety=fparamNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1774  }
1775  else // Regular structure
1776  {
1777  newSafety=fregularNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1778  }
1779  break;
1780  case kReplica:
1781  G4Exception("G4ITNavigator1::ComputeSafety()", "GeomNav0001",
1782  FatalException, "Not applicable for replicated volumes.");
1783  break;
1784  }
1785  }
1786  else
1787  {
1788  newSafety = freplicaNav.ComputeSafety(pGlobalpoint, localPoint,
1789  fHistory, pMaxLength);
1790  }
1791  }
1792  else // if( endpointOnSurface && stayedOnEndpoint )
1793  {
1794 #ifdef G4DEBUG_NAVIGATION
1795  if( fVerbose >= 2 )
1796  {
1797  G4cout << " G4ITNavigator1::ComputeSafety() finds that point - "
1798  << pGlobalpoint << " - is on surface " << G4endl;
1799  if( fEnteredDaughter ) { G4cout << " entered new daughter volume"; }
1800  if( fExitedMother ) { G4cout << " and exited previous volume."; }
1801  G4cout << G4endl;
1802  G4cout << " EndPoint was = " << fStepEndPoint << G4endl;
1803  }
1804 #endif
1805  newSafety = 0.0;
1806  }
1807 
1808  // Remember last safety origin & value
1809  //
1810  fPreviousSftOrigin = pGlobalpoint;
1811  fPreviousSafety = newSafety;
1812 
1813  if (keepState) { RestoreSavedState(); }
1814 
1815 #ifdef G4DEBUG_NAVIGATION
1816  if( fVerbose > 1 )
1817  {
1818  G4cout << " ---- Exiting ComputeSafety " << G4endl;
1819  if( fVerbose > 2 ) { PrintState(); }
1820  G4cout << " Returned value of Safety = " << newSafety << G4endl;
1821  }
1822  G4cout.precision(oldcoutPrec);
1823 #endif
1824 
1825  return newSafety;
1826 }
G4SmartVoxelHeader * GetVoxelHeader() const
G4VPhysicalVolume * GetTopVolume() const
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
G4int GetDaughtersRegularStructureId(const G4LogicalVolume *pLog) const
G4ThreeVector ComputeLocalPoint(const G4ThreeVector &rGlobPoint) const
G4double ComputeSafety(const G4ThreeVector &localPoint, const G4VPhysicalVolume &currentPhysical, G4double maxLength=DBL_MAX)
void RestoreSavedState()
int G4int
Definition: G4Types.hh:78
G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4NavigationHistory &history, const G4double pMaxLength=DBL_MAX)
G4double ComputeSafety(const G4ThreeVector &localPoint, const G4NavigationHistory &history, const G4double pProposedMaxLength=DBL_MAX)
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
bool G4bool
Definition: G4Types.hh:79
EVolume GetTopVolumeType() const
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4NavigationHistory &history, const G4double pMaxLength=DBL_MAX)
G4ThreeVector fStepEndPoint
G4NavigationHistory fHistory
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4LogicalVolume * GetLogicalVolume() const
EVolume CharacteriseDaughters(const G4LogicalVolume *pLog) const
G4double ComputeSafety(const G4ThreeVector &localPoint, const G4NavigationHistory &history, const G4double pProposedMaxLength=DBL_MAX)
#define G4endl
Definition: G4ios.hh:61
G4double ComputeSafety(const G4ThreeVector &globalPoint, const G4ThreeVector &localPoint, G4NavigationHistory &history, const G4double pProposedMaxLength=DBL_MAX)
double G4double
Definition: G4Types.hh:76
G4double kCarTolerance
void PrintState() const

Here is the call graph for this function:

G4double G4ITNavigator1::ComputeStep ( const G4ThreeVector pGlobalPoint,
const G4ThreeVector pDirection,
const G4double  pCurrentProposedStepLength,
G4double pNewSafety 
)
virtual

Definition at line 787 of file G4ITNavigator1.cc.

791 {
792  G4ThreeVector localDirection = ComputeLocalAxis(pDirection);
793  G4double Step = kInfinity;
794  G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
795  G4LogicalVolume *motherLogical = motherPhysical->GetLogicalVolume();
796 
797  // All state relating to exiting normals must be reset
798  //
799  fExitNormalGlobalFrame= G4ThreeVector( 0., 0., 0.);
800  // Reset value - to erase its memory
801  fChangedGrandMotherRefFrame= false;
802  // Reset - used for local exit normal
803  fGrandMotherExitNormal= G4ThreeVector( 0., 0., 0.);
804  fCalculatedExitNormal = false;
805  // Reset for new step
806 
807  static G4ThreadLocal G4int sNavCScalls=0;
808  sNavCScalls++;
809 
810  fLastTriedStepComputation= true;
811 
812 #ifdef G4VERBOSE
813  if( fVerbose > 0 )
814  {
815  G4cout << "*** G4ITNavigator1::ComputeStep: ***" << G4endl;
816  G4cout << " Volume = " << motherPhysical->GetName()
817  << " - Proposed step length = " << pCurrentProposedStepLength
818  << G4endl;
819 #ifdef G4DEBUG_NAVIGATION
820  if( fVerbose >= 2 )
821  {
822  G4cout << " Called with the arguments: " << G4endl
823  << " Globalpoint = " << std::setw(25) << pGlobalpoint << G4endl
824  << " Direction = " << std::setw(25) << pDirection << G4endl;
825  if( fVerbose >= 4 )
826  {
827  G4cout << " ---- Upon entering : State" << G4endl;
828  PrintState();
829  }
830  }
831 #endif
832  }
833 #endif
834 
835  G4ThreeVector newLocalPoint = ComputeLocalPoint(pGlobalpoint);
836  if( newLocalPoint != fLastLocatedPointLocal )
837  {
838  // Check whether the relocation is within safety
839  //
840  G4ThreeVector oldLocalPoint = fLastLocatedPointLocal;
841  G4double moveLenSq = (newLocalPoint-oldLocalPoint).mag2();
842 
843  if ( moveLenSq >= kCarTolerance*kCarTolerance )
844  {
845 #ifdef G4VERBOSE
846  ComputeStepLog(pGlobalpoint, moveLenSq);
847 #endif
848  // Relocate the point within the same volume
849  //
850  LocateGlobalPointWithinVolume( pGlobalpoint );
851  fLastTriedStepComputation= true; // Ensure that this is set again !!
852  }
853  }
855  {
856  switch( CharacteriseDaughters(motherLogical) )
857  {
858  case kNormal:
859  if ( motherLogical->GetVoxelHeader() )
860  {
861  Step = fvoxelNav.ComputeStep(fLastLocatedPointLocal,
862  localDirection,
863  pCurrentProposedStepLength,
864  pNewSafety,
865  fHistory,
866  fValidExitNormal,
867  fExitNormal,
868  fExiting,
869  fEntering,
870  &fBlockedPhysicalVolume,
871  fBlockedReplicaNo);
872 
873  }
874  else
875  {
876  if( motherPhysical->GetRegularStructureId() == 0 )
877  {
878  Step = fnormalNav.ComputeStep(fLastLocatedPointLocal,
879  localDirection,
880  pCurrentProposedStepLength,
881  pNewSafety,
882  fHistory,
883  fValidExitNormal,
884  fExitNormal,
885  fExiting,
886  fEntering,
887  &fBlockedPhysicalVolume,
888  fBlockedReplicaNo);
889  }
890  else // Regular (non-voxelised) structure
891  {
892  LocateGlobalPointAndSetup( pGlobalpoint, &pDirection, true, true );
893  fLastTriedStepComputation= true; // Ensure that this is set again !!
894  //
895  // if physical process limits the step, the voxel will not be the
896  // one given by ComputeStepSkippingEqualMaterials() and the local
897  // point will be wrongly calculated.
898 
899  // There is a problem: when msc limits the step and the point is
900  // assigned wrongly to phantom in previous step (while it is out
901  // of the container volume). Then LocateGlobalPointAndSetup() has
902  // reset the history topvolume to world.
903  //
905  {
906  G4Exception("G4ITNavigator1::ComputeStep()",
907  "GeomNav1001", JustWarning,
908  "Point is relocated in voxels, while it should be outside!");
909  Step = fnormalNav.ComputeStep(fLastLocatedPointLocal,
910  localDirection,
911  pCurrentProposedStepLength,
912  pNewSafety,
913  fHistory,
914  fValidExitNormal,
915  fExitNormal,
916  fExiting,
917  fEntering,
918  &fBlockedPhysicalVolume,
919  fBlockedReplicaNo);
920  }
921  else
922  {
923  Step = fregularNav.
924  ComputeStepSkippingEqualMaterials(fLastLocatedPointLocal,
925  localDirection,
926  pCurrentProposedStepLength,
927  pNewSafety,
928  fHistory,
929  fValidExitNormal,
930  fExitNormal,
931  fExiting,
932  fEntering,
933  &fBlockedPhysicalVolume,
934  fBlockedReplicaNo,
935  motherPhysical);
936  }
937  }
938  }
939  break;
940  case kParameterised:
941  if( GetDaughtersRegularStructureId(motherLogical) != 1 )
942  {
943  Step = fparamNav.ComputeStep(fLastLocatedPointLocal,
944  localDirection,
945  pCurrentProposedStepLength,
946  pNewSafety,
947  fHistory,
948  fValidExitNormal,
949  fExitNormal,
950  fExiting,
951  fEntering,
952  &fBlockedPhysicalVolume,
953  fBlockedReplicaNo);
954  }
955  else // Regular structure
956  {
957  Step = fregularNav.ComputeStep(fLastLocatedPointLocal,
958  localDirection,
959  pCurrentProposedStepLength,
960  pNewSafety,
961  fHistory,
962  fValidExitNormal,
963  fExitNormal,
964  fExiting,
965  fEntering,
966  &fBlockedPhysicalVolume,
967  fBlockedReplicaNo);
968  }
969  break;
970  case kReplica:
971  G4Exception("G4ITNavigator1::ComputeStep()", "GeomNav0001",
972  FatalException, "Not applicable for replicated volumes.");
973  break;
974  }
975  }
976  else
977  {
978  // In the case of a replica, it must handle the exiting
979  // edge/corner problem by itself
980  //
981  G4bool exitingReplica = fExitedMother;
982  G4bool calculatedExitNormal= false;
983 
984  Step = freplicaNav.ComputeStep(pGlobalpoint,
985  pDirection,
986  fLastLocatedPointLocal,
987  localDirection,
988  pCurrentProposedStepLength,
989  pNewSafety,
990  fHistory,
991  fValidExitNormal,
992  calculatedExitNormal,
993  fExitNormal,
994  exitingReplica,
995  fEntering,
996  &fBlockedPhysicalVolume,
997  fBlockedReplicaNo);
998  fExiting= exitingReplica;
999  fCalculatedExitNormal= calculatedExitNormal;
1000  }
1001 
1002  // Remember last safety origin & value.
1003  //
1004  fPreviousSftOrigin = pGlobalpoint;
1005  fPreviousSafety = pNewSafety;
1006 
1007  // Count zero steps - one can occur due to changing momentum at a boundary
1008  // - one, two (or a few) can occur at common edges between
1009  // volumes
1010  // - more than two is likely a problem in the geometry
1011  // description or the Navigation
1012 
1013  // Rule of thumb: likely at an Edge if two consecutive steps are zero,
1014  // because at least two candidate volumes must have been
1015  // checked
1016  //
1017  fLocatedOnEdge = fLastStepWasZero && (Step==0.0);
1018  fLastStepWasZero = (Step==0.0);
1019  if (fPushed) { fPushed = fLastStepWasZero; }
1020 
1021  // Handle large number of consecutive zero steps
1022  //
1023  if ( fLastStepWasZero )
1024  {
1025  fNumberZeroSteps++;
1026 #ifdef G4DEBUG_NAVIGATION
1027  if( fNumberZeroSteps > 1 )
1028  {
1029  G4cout << "G4ITNavigator1::ComputeStep(): another zero step, # "
1030  << fNumberZeroSteps
1031  << " at " << pGlobalpoint
1032  << " in volume " << motherPhysical->GetName()
1033  << " nav-comp-step calls # " << sNavCScalls
1034  << G4endl;
1035  }
1036 #endif
1037  if( fNumberZeroSteps > fActionThreshold_NoZeroSteps-1 )
1038  {
1039  // Act to recover this stuck track. Pushing it along direction
1040  //
1041  Step += 100*kCarTolerance;
1042 #ifdef G4VERBOSE
1043  if ((!fPushed) && (fWarnPush))
1044  {
1045  std::ostringstream message;
1046  message << "Track stuck or not moving." << G4endl
1047  << " Track stuck, not moving for "
1048  << fNumberZeroSteps << " steps" << G4endl
1049  << " in volume -" << motherPhysical->GetName()
1050  << "- at point " << pGlobalpoint << G4endl
1051  << " direction: " << pDirection << "." << G4endl
1052  << " Potential geometry or navigation problem !"
1053  << G4endl
1054  << " Trying pushing it of " << Step << " mm ...";
1055  G4Exception("G4ITNavigator1::ComputeStep()", "GeomNav1002",
1056  JustWarning, message, "Potential overlap in geometry!");
1057  }
1058 #endif
1059  fPushed = true;
1060  }
1061  if( fNumberZeroSteps > fAbandonThreshold_NoZeroSteps-1 )
1062  {
1063  // Must kill this stuck track
1064  //
1065  std::ostringstream message;
1066  message << "Stuck Track: potential geometry or navigation problem."
1067  << G4endl
1068  << " Track stuck, not moving for "
1069  << fNumberZeroSteps << " steps" << G4endl
1070  << " in volume -" << motherPhysical->GetName()
1071  << "- at point " << pGlobalpoint << G4endl
1072  << " direction: " << pDirection << ".";
1073  motherPhysical->CheckOverlaps(5000, false);
1074  G4Exception("G4ITNavigator1::ComputeStep()", "GeomNav0003",
1075  EventMustBeAborted, message);
1076  }
1077  }
1078  else
1079  {
1080  if (!fPushed) fNumberZeroSteps = 0;
1081  }
1082 
1083  fEnteredDaughter = fEntering; // I expect to enter a volume in this Step
1084  fExitedMother = fExiting;
1085 
1086  fStepEndPoint = pGlobalpoint + Step * pDirection;
1087  fLastStepEndPointLocal = fLastLocatedPointLocal + Step * localDirection;
1088 
1089  if( fExiting )
1090  {
1091 #ifdef G4DEBUG_NAVIGATION
1092  if( fVerbose > 2 )
1093  {
1094  G4cout << " At G4Nav CompStep End - if(exiting) - fExiting= " << fExiting
1095  << " fValidExitNormal = " << fValidExitNormal << G4endl;
1096  G4cout << " fExitNormal= " << fExitNormal << G4endl;
1097  }
1098 #endif
1099 
1100  if(fValidExitNormal || fCalculatedExitNormal)
1101  {
1103  {
1104  // Convention: fExitNormal is in the 'grand-mother' coordinate system
1105  //
1106  fGrandMotherExitNormal= fExitNormal;
1107  fCalculatedExitNormal= true;
1108  }
1109  else
1110  {
1111  fGrandMotherExitNormal = fExitNormal;
1112  }
1113  }
1114  else
1115  {
1116  // We must calculate the normal anyway (in order to have it if requested)
1117  //
1118  G4ThreeVector finalLocalPoint =
1119  fLastLocatedPointLocal + localDirection*Step;
1120 
1122  {
1123  // Find normal in the 'mother' coordinate system
1124  //
1125  G4ThreeVector exitNormalMotherFrame=
1126  motherLogical->GetSolid()->SurfaceNormal(finalLocalPoint);
1127 
1128  // Transform it to the 'grand-mother' coordinate system
1129  //
1130  const G4RotationMatrix* mRot = motherPhysical->GetRotation();
1131  if( mRot )
1132  {
1133  fChangedGrandMotherRefFrame= true;
1134  fGrandMotherExitNormal = (*mRot).inverse() * exitNormalMotherFrame;
1135  }
1136  else
1137  {
1138  fGrandMotherExitNormal = exitNormalMotherFrame;
1139  }
1140 
1141  // Do not set fValidExitNormal -- this signifies
1142  // that the solid is convex!
1143  //
1144  fCalculatedExitNormal= true;
1145  }
1146  else
1147  {
1148  fCalculatedExitNormal = false;
1149  //
1150  // Nothing can be done at this stage currently - to solve this
1151  // Replica Navigation must have calculated the normal for this case
1152  // already.
1153  // Cases: mother is not convex, and exit is at previous replica level
1154 
1155 #ifdef G4DEBUG_NAVIGATION
1157 
1158  desc << "Problem in ComputeStep: Replica Navigation did not provide"
1159  << " valid exit Normal. " << G4endl;
1160  desc << " Do not know how calculate it in this case." << G4endl;
1161  desc << " Location = " << finalLocalPoint << G4endl;
1162  desc << " Volume name = " << motherPhysical->GetName()
1163  << " copy/replica No = " << motherPhysical->GetCopyNo() << G4endl;
1164  G4Exception("G4ITNavigator1::ComputeStep()", "GeomNav0003",
1165  JustWarning, desc, "Normal not available for exiting.");
1166 #endif
1167  }
1168  }
1169 
1170  // Now transform it to the global reference frame !!
1171  //
1172  if( fValidExitNormal || fCalculatedExitNormal )
1173  {
1174  G4int depth= fHistory.GetDepth();
1175  if( depth > 0 )
1176  {
1177  G4AffineTransform GrandMotherToGlobalTransf =
1178  fHistory.GetTransform(depth-1).Inverse();
1179  fExitNormalGlobalFrame =
1180  GrandMotherToGlobalTransf.TransformAxis( fGrandMotherExitNormal );
1181  }
1182  else
1183  {
1184  fExitNormalGlobalFrame= fGrandMotherExitNormal;
1185  }
1186  }
1187  else
1188  {
1189  fExitNormalGlobalFrame= G4ThreeVector( 0., 0., 0.);
1190  }
1191  }
1192  fStepEndPoint= pGlobalpoint+Step*pDirection;
1193 
1194  if( (Step == pCurrentProposedStepLength) && (!fExiting) && (!fEntering) )
1195  {
1196  // This if Step is not really limited by the geometry.
1197  // The Navigator is obliged to return "infinity"
1198  //
1199  Step = kInfinity;
1200  }
1201 
1202 #ifdef G4VERBOSE
1203  if( fVerbose > 1 )
1204  {
1205  if( fVerbose >= 4 )
1206  {
1207  G4cout << " ----- Upon exiting :" << G4endl;
1208  PrintState();
1209  }
1210  G4cout << " Returned step= " << Step;
1211  if( fVerbose > 5 ) G4cout << G4endl;
1212  if( Step == kInfinity )
1213  {
1214  G4cout << " Requested step= " << pCurrentProposedStepLength ;
1215  if( fVerbose > 5) G4cout << G4endl;
1216  }
1217  G4cout << " Safety = " << pNewSafety << G4endl;
1218  }
1219 #endif
1220 
1221  return Step;
1222 }
G4SmartVoxelHeader * GetVoxelHeader() const
G4VPhysicalVolume * GetTopVolume() const
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
G4double ComputeStep(const G4ThreeVector &globalPoint, const G4ThreeVector &globalDirection, const G4double currentProposedStepLength, G4double &newSafety, G4NavigationHistory &history, G4bool &validExitNormal, G4ThreeVector &exitNormal, G4bool &exiting, G4bool &entering, G4VPhysicalVolume *(*pBlockedPhysical), G4int &blockedReplicaNo)
static const G4double kInfinity
Definition: geomdefs.hh:42
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4int GetDaughtersRegularStructureId(const G4LogicalVolume *pLog) const
CLHEP::Hep3Vector G4ThreeVector
G4AffineTransform Inverse() const
virtual G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=0, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
G4int GetDepth() const
G4VSolid * GetSolid() const
G4ThreeVector ComputeLocalPoint(const G4ThreeVector &rGlobPoint) const
const G4RotationMatrix * GetRotation() const
G4ThreeVector fLastStepEndPointLocal
G4double ComputeStep(const G4ThreeVector &globalPoint, const G4ThreeVector &globalDirection, const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4double currentProposedStepLength, G4double &newSafety, G4NavigationHistory &history, G4bool &validExitNormal, G4bool &calculatedExitNormal, G4ThreeVector &exitNormal, G4bool &exiting, G4bool &entering, G4VPhysicalVolume *(*pBlockedPhysical), G4int &blockedReplicaNo)
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4ThreeVector ComputeLocalAxis(const G4ThreeVector &pVec) const
const G4String & GetName() const
G4double ComputeStep(const G4ThreeVector &globalPoint, const G4ThreeVector &globalDirection, const G4double currentProposedStepLength, G4double &newSafety, G4NavigationHistory &history, G4bool &validExitNormal, G4ThreeVector &exitNormal, G4bool &exiting, G4bool &entering, G4VPhysicalVolume *(*pBlockedPhysical), G4int &blockedReplicaNo)
bool G4bool
Definition: G4Types.hh:79
EVolume GetTopVolumeType() const
virtual G4int GetRegularStructureId() const =0
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
G4ThreeVector fStepEndPoint
G4NavigationHistory fHistory
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double ComputeStep(const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4double currentProposedStepLength, G4double &newSafety, G4NavigationHistory &history, G4bool &validExitNormal, G4ThreeVector &exitNormal, G4bool &exiting, G4bool &entering, G4VPhysicalVolume *(*pBlockedPhysical), G4int &blockedReplicaNo)
G4LogicalVolume * GetLogicalVolume() const
const G4AffineTransform & GetTransform(G4int n) const
EVolume CharacteriseDaughters(const G4LogicalVolume *pLog) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
virtual G4int GetCopyNo() const =0
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
virtual G4double ComputeStep(const G4ThreeVector &globalPoint, const G4ThreeVector &globalDirection, const G4double currentProposedStepLength, G4double &newSafety, G4NavigationHistory &history, G4bool &validExitNormal, G4ThreeVector &exitNormal, G4bool &exiting, G4bool &entering, G4VPhysicalVolume *(*pBlockedPhysical), G4int &blockedReplicaNo)
G4double kCarTolerance
virtual G4bool CheckOverlaps(G4int res=1000, G4double tol=0., G4bool verbose=true, G4int errMax=1)
void PrintState() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4GRSSolid* G4ITNavigator1::CreateGRSSolid ( ) const
inline
G4GRSVolume* G4ITNavigator1::CreateGRSVolume ( ) const
inline
G4TouchableHistory* G4ITNavigator1::CreateTouchableHistory ( ) const
inline

Here is the caller graph for this function:

G4TouchableHistory* G4ITNavigator1::CreateTouchableHistory ( const G4NavigationHistory ) const
inline
G4TouchableHistoryHandle G4ITNavigator1::CreateTouchableHistoryHandle ( ) const
virtual

Definition at line 1832 of file G4ITNavigator1.cc.

1833 {
1835 }
G4TouchableHistory * CreateTouchableHistory() const
G4ReferenceCountedHandle< G4TouchableHistory > G4TouchableHistoryHandle

Here is the call graph for this function:

void G4ITNavigator1::EnableBestSafety ( G4bool  value = false)
inline
G4bool G4ITNavigator1::EnteredDaughterVolume ( ) const
inline

Here is the caller graph for this function:

G4bool G4ITNavigator1::ExitedMotherVolume ( ) const
inline
G4ThreeVector G4ITNavigator1::GetCurrentLocalCoordinate ( ) const
inline
G4int G4ITNavigator1::GetDaughtersRegularStructureId ( const G4LogicalVolume pLog) const
inlineprotected

Here is the caller graph for this function:

G4ThreeVector G4ITNavigator1::GetGlobalExitNormal ( const G4ThreeVector point,
G4bool valid 
)
virtual

Definition at line 1590 of file G4ITNavigator1.cc.

1592 {
1593  G4bool validNormal;
1594  G4ThreeVector localNormal, globalNormal;
1595 
1596  if( fLastTriedStepComputation && fExiting )
1597  {
1598  // This was computed in ComputeStep -- and only on arrival at boundary
1599  //
1600  globalNormal = fExitNormalGlobalFrame;
1601  *pNormalCalculated = true; // ComputeStep always computes it if Exiting
1602  // (fExiting==true)
1603  }
1604  else
1605  {
1606  localNormal = GetLocalExitNormalAndCheck(IntersectPointGlobal,&validNormal);
1607  *pNormalCalculated = fCalculatedExitNormal;
1608 
1609 #ifdef G4DEBUG_NAVIGATION
1610  if( (!validNormal) && !fCalculatedExitNormal)
1611  {
1613  edN << " Calculated = " << fCalculatedExitNormal << G4endl;
1614  edN << " Entering= " << fEntering << G4endl;
1615  G4int oldVerbose= this->GetVerboseLevel();
1616  this->SetVerboseLevel(4);
1617  edN << " State of Navigator: " << G4endl;
1618  edN << *this << G4endl;
1619  this->SetVerboseLevel( oldVerbose );
1620 
1621  G4Exception("G4ITNavigator1::GetGlobalExitNormal()",
1622  "GeomNav0003", JustWarning, edN,
1623  "LocalExitNormalAndCheck() did not calculate Normal.");
1624  }
1625 #endif
1626 
1627  G4double localMag2= localNormal.mag2();
1628  if( validNormal && (std::fabs(localMag2-1.0)) > CLHEP::perMillion )
1629  {
1631 
1632  edN << "G4ITNavigator1::GetGlobalExitNormal: "
1633  << " Using Local Normal - from call to GetLocalExitNormalAndCheck. "
1634  << G4endl
1635  << " Local Exit Normal = " << localNormal << " || = "
1636  << std::sqrt(localMag2) << G4endl
1637  << " Global Exit Normal = " << globalNormal << " || = "
1638  << globalNormal.mag() << G4endl;
1639  edN << " Calculated It = " << fCalculatedExitNormal << G4endl;
1640 
1641  G4Exception("G4ITNavigator1::GetGlobalExitNormal()",
1642  "GeomNav0003",JustWarning, edN,
1643  "Value obtained from new local *solid* is incorrect.");
1644  localNormal = localNormal.unit(); // Should we correct it ??
1645  }
1646  G4AffineTransform localToGlobal = GetLocalToGlobalTransform();
1647  globalNormal = localToGlobal.TransformAxis( localNormal );
1648  }
1649 
1650 #ifdef G4DEBUG_NAVIGATION
1651  // Temporary extra checks
1652  if( fLastTriedStepComputation && fExiting)
1653  {
1654  localNormal = GetLocalExitNormalAndCheck( IntersectPointGlobal, &validNormal);
1655  *pNormalCalculated = fCalculatedExitNormal;
1656 
1657  G4AffineTransform localToGlobal = GetLocalToGlobalTransform();
1658  globalNormal = localToGlobal.TransformAxis( localNormal );
1659 
1660  // Check the value computed against fExitNormalGlobalFrame
1661  G4ThreeVector diffNorm = globalNormal - fExitNormalGlobalFrame;
1662  if( diffNorm.mag2() > perMillion*CLHEP::perMillion)
1663  {
1664  G4ExceptionDescription edDfn;
1665  edDfn << "Found difference in normals in case of exiting mother "
1666  << "- when Get is called after ComputingStep " << G4endl;
1667  edDfn << " Magnitude of diff = " << diffNorm.mag() << G4endl;
1668  edDfn << " Normal stored (Global) = " << fExitNormalGlobalFrame
1669  << G4endl;
1670  edDfn << " Global Computed from Local = " << globalNormal << G4endl;
1671  G4Exception("G4ITNavigator1::GetGlobalExitNormal()", "GeomNav0003",
1672  JustWarning, edDfn);
1673  }
1674  }
1675 #endif
1676 
1677  return globalNormal;
1678 }
static constexpr double perMillion
Definition: G4SIunits.hh:334
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void SetVerboseLevel(G4int level)
int G4int
Definition: G4Types.hh:78
G4int GetVerboseLevel() const
bool G4bool
Definition: G4Types.hh:79
virtual G4ThreeVector GetLocalExitNormalAndCheck(const G4ThreeVector &point, G4bool *valid)
static constexpr double perMillion
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4AffineTransform GetLocalToGlobalTransform() const
Hep3Vector unit() const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
double mag2() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
double mag() const

Here is the call graph for this function:

const G4AffineTransform& G4ITNavigator1::GetGlobalToLocalTransform ( ) const
inline

Here is the caller graph for this function:

G4ThreeVector G4ITNavigator1::GetLocalExitNormal ( G4bool valid)
virtual

Definition at line 1348 of file G4ITNavigator1.cc.

1349 {
1350  G4ThreeVector ExitNormal(0.,0.,0.);
1351  G4VSolid *currentSolid=0;
1352  G4LogicalVolume *candidateLogical;
1353  if ( fLastTriedStepComputation )
1354  {
1355  // use fLastLocatedPointLocal and next candidate volume
1356  //
1357  G4ThreeVector nextSolidExitNormal(0.,0.,0.);
1358 
1359  if( fEntering && (fBlockedPhysicalVolume!=0) )
1360  {
1361  candidateLogical= fBlockedPhysicalVolume->GetLogicalVolume();
1362  if( candidateLogical )
1363  {
1364  // fLastStepEndPointLocal is in the coordinates of the mother
1365  // we need it in the daughter's coordinate system.
1366 
1367  // The following code should also work in case of Replica
1368  {
1369  // First transform fLastLocatedPointLocal to the new daughter
1370  // coordinates
1371  //
1372  G4AffineTransform MotherToDaughterTransform=
1373  GetMotherToDaughterTransform( fBlockedPhysicalVolume,
1374  fBlockedReplicaNo,
1375  VolumeType(fBlockedPhysicalVolume) );
1376  G4ThreeVector daughterPointOwnLocal=
1377  MotherToDaughterTransform.TransformPoint( fLastStepEndPointLocal );
1378 
1379  // OK if it is a parameterised volume
1380  //
1381  EInside inSideIt;
1382  G4bool onSurface;
1383  G4double safety= -1.0;
1384  currentSolid= candidateLogical->GetSolid();
1385  inSideIt = currentSolid->Inside(daughterPointOwnLocal);
1386  onSurface = (inSideIt == kSurface);
1387  if( ! onSurface )
1388  {
1389  if( inSideIt == kOutside )
1390  {
1391  safety = (currentSolid->DistanceToIn(daughterPointOwnLocal));
1392  onSurface = safety < 100.0 * kCarTolerance;
1393  }
1394  else if (inSideIt == kInside )
1395  {
1396  safety = (currentSolid->DistanceToOut(daughterPointOwnLocal));
1397  onSurface = safety < 100.0 * kCarTolerance;
1398  }
1399  }
1400 
1401  if( onSurface )
1402  {
1403  nextSolidExitNormal =
1404  currentSolid->SurfaceNormal(daughterPointOwnLocal);
1405 
1406  // Entering the solid ==> opposite
1407  //
1408  ExitNormal = -nextSolidExitNormal;
1409  fCalculatedExitNormal= true;
1410  }
1411  else
1412  {
1413 #ifdef G4VERBOSE
1414  if(( fVerbose == 1 ) && ( fCheck ))
1415  {
1416  std::ostringstream message;
1417  message << "Point not on surface ! " << G4endl
1418  << " Point = "
1419  << daughterPointOwnLocal << G4endl
1420  << " Physical volume = "
1421  << fBlockedPhysicalVolume->GetName() << G4endl
1422  << " Logical volume = "
1423  << candidateLogical->GetName() << G4endl
1424  << " Solid = " << currentSolid->GetName()
1425  << " Type = "
1426  << currentSolid->GetEntityType() << G4endl
1427  << *currentSolid << G4endl;
1428  if( inSideIt == kOutside )
1429  {
1430  message << "Point is Outside. " << G4endl
1431  << " Safety (from outside) = " << safety << G4endl;
1432  }
1433  else // if( inSideIt == kInside )
1434  {
1435  message << "Point is Inside. " << G4endl
1436  << " Safety (from inside) = " << safety << G4endl;
1437  }
1438  G4Exception("G4ITNavigator1::GetLocalExitNormal()", "GeomNav1001",
1439  JustWarning, message);
1440  }
1441 #endif
1442  }
1443  *valid = onSurface; // was =true;
1444  }
1445  }
1446  }
1447  else if ( fExiting )
1448  {
1449  ExitNormal = fGrandMotherExitNormal;
1450  *valid = true;
1451  fCalculatedExitNormal= true; // Should be true already
1452  }
1453  else // i.e. ( fBlockedPhysicalVolume == 0 )
1454  {
1455  *valid = false;
1456  G4Exception("G4ITNavigator1::GetLocalExitNormal()",
1457  "GeomNav0003", JustWarning,
1458  "Incorrect call to GetLocalSurfaceNormal." );
1459  }
1460  }
1461  else // ( ! fLastTriedStepComputation ) ie. last call was to Locate
1462  {
1463  if ( EnteredDaughterVolume() )
1464  {
1465  G4VSolid* daughterSolid =fHistory.GetTopVolume()->GetLogicalVolume()
1466  ->GetSolid();
1467  ExitNormal= -(daughterSolid->SurfaceNormal(fLastLocatedPointLocal));
1468  if( std::fabs(ExitNormal.mag2()-1.0 ) > CLHEP::perMillion )
1469  {
1471  desc << " Parameters of solid: " << *daughterSolid
1472  << " Point for surface = " << fLastLocatedPointLocal << std::endl;
1473  G4Exception("G4ITNavigator1::GetLocalExitNormal()",
1474  "GeomNav0003", FatalException, desc,
1475  "Surface Normal returned by Solid is not a Unit Vector." );
1476  }
1477  fCalculatedExitNormal= true;
1478  *valid = true;
1479  }
1480  else
1481  {
1482  if( fExitedMother )
1483  {
1484  ExitNormal = fGrandMotherExitNormal;
1485  *valid = true;
1486  fCalculatedExitNormal= true;
1487  }
1488  else // We are not at a boundary. ExitNormal remains (0,0,0)
1489  {
1490  *valid = false;
1491  fCalculatedExitNormal= false;
1492  G4ExceptionDescription message;
1493  message << "Function called when *NOT* at a Boundary." << G4endl;
1494  G4Exception("G4ITNavigator1::GetLocalExitNormal()",
1495  "GeomNav0003", JustWarning, message);
1496  }
1497  }
1498  }
1499  return ExitNormal;
1500 }
G4String GetName() const
G4VPhysicalVolume * GetTopVolume() const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4VSolid * GetSolid() const
EVolume VolumeType(const G4VPhysicalVolume *pVol) const
G4ThreeVector fLastStepEndPointLocal
virtual G4GeometryType GetEntityType() const =0
const G4String & GetName() const
virtual EInside Inside(const G4ThreeVector &p) const =0
bool G4bool
Definition: G4Types.hh:79
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
G4NavigationHistory fHistory
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
static constexpr double perMillion
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4LogicalVolume * GetLogicalVolume() const
EInside
Definition: geomdefs.hh:58
G4bool EnteredDaughterVolume() const
#define G4endl
Definition: G4ios.hh:61
const G4String & GetName() const
double G4double
Definition: G4Types.hh:76
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
G4AffineTransform GetMotherToDaughterTransform(G4VPhysicalVolume *dVolume, G4int dReplicaNo, EVolume dVolumeType)
G4double kCarTolerance

Here is the call graph for this function:

Here is the caller graph for this function:

G4ThreeVector G4ITNavigator1::GetLocalExitNormalAndCheck ( const G4ThreeVector point,
G4bool valid 
)
virtual

Definition at line 1555 of file G4ITNavigator1.cc.

1562 {
1563 #ifdef G4DEBUG_NAVIGATION
1564  // Check Current point against expected 'local' value
1565  //
1566  if ( fLastTriedStepComputation )
1567  {
1568  G4ThreeVector ExpectedBoundaryPointLocal;
1569 
1570  const G4AffineTransform& GlobalToLocal= GetGlobalToLocalTransform();
1571  ExpectedBoundaryPointLocal =
1572  GlobalToLocal.TransformPoint( ExpectedBoundaryPointGlobal );
1573 
1574  // Add here: Comparison against expected position,
1575  // i.e. the endpoint of ComputeStep
1576  }
1577 #endif
1578 
1579  return GetLocalExitNormal( pValid);
1580 }
virtual G4ThreeVector GetLocalExitNormal(G4bool *valid)
const G4AffineTransform & GetGlobalToLocalTransform() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const

Here is the call graph for this function:

Here is the caller graph for this function:

const G4AffineTransform G4ITNavigator1::GetLocalToGlobalTransform ( ) const
inline

Here is the caller graph for this function:

G4AffineTransform G4ITNavigator1::GetMotherToDaughterTransform ( G4VPhysicalVolume dVolume,
G4int  dReplicaNo,
EVolume  dVolumeType 
)

Definition at line 1509 of file G4ITNavigator1.cc.

1512 {
1513  switch (enteringVolumeType)
1514  {
1515  case kNormal: // Nothing is needed to prepare the transformation
1516  break; // It is stored already in the physical volume (placement)
1517  case kReplica: // Sets the transform in the Replica - tbc
1518  G4Exception("G4ITNavigator1::GetMotherToDaughterTransform()",
1519  "GeomNav0001", FatalException,
1520  "Method NOT Implemented yet for replica volumes.");
1521  break;
1522  case kParameterised:
1523  if( pEnteringPhysVol->GetRegularStructureId() == 0 )
1524  {
1525  G4VPVParameterisation *pParam =
1526  pEnteringPhysVol->GetParameterisation();
1527  G4VSolid* pSolid =
1528  pParam->ComputeSolid(enteringReplicaNo, pEnteringPhysVol);
1529  pSolid->ComputeDimensions(pParam, enteringReplicaNo, pEnteringPhysVol);
1530 
1531  // Sets the transform in the Parameterisation
1532  //
1533  pParam->ComputeTransformation(enteringReplicaNo, pEnteringPhysVol);
1534 
1535  // Set the correct solid and material in Logical Volume
1536  //
1537  G4LogicalVolume* pLogical = pEnteringPhysVol->GetLogicalVolume();
1538  pLogical->SetSolid( pSolid );
1539  }
1540  break;
1541  }
1542  return G4AffineTransform(pEnteringPhysVol->GetRotation(),
1543  pEnteringPhysVol->GetTranslation()).Invert();
1544 }
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:138
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
void SetSolid(G4VSolid *pSolid)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual void ComputeTransformation(const G4int, G4VPhysicalVolume *) const =0

Here is the call graph for this function:

Here is the caller graph for this function:

G4ITNavigatorState_Lock1 * G4ITNavigator1::GetNavigatorState ( )

Definition at line 656 of file G4ITNavigator1.cc.

657 {
658  SetSavedState();
659  return fpSaveState;
660 }

Here is the call graph for this function:

G4int G4ITNavigator1::GetVerboseLevel ( ) const
inline

Here is the caller graph for this function:

G4VPhysicalVolume* G4ITNavigator1::GetWorldVolume ( ) const
inline
G4bool G4ITNavigator1::IsActive ( ) const
inline
G4bool G4ITNavigator1::IsCheckModeActive ( ) const
inline
G4VPhysicalVolume * G4ITNavigator1::LocateGlobalPointAndSetup ( const G4ThreeVector point,
const G4ThreeVector direction = 0,
const G4bool  pRelativeSearch = true,
const G4bool  ignoreDirection = true 
)
virtual

Definition at line 160 of file G4ITNavigator1.cc.

164 {
165  G4bool notKnownContained=true, noResult;
166  G4VPhysicalVolume *targetPhysical;
167  G4LogicalVolume *targetLogical;
168  G4VSolid *targetSolid=0;
169  G4ThreeVector localPoint, globalDirection;
170  EInside insideCode;
171 
172  G4bool considerDirection = (!ignoreDirection) || fLocatedOnEdge;
173 
174  fLastTriedStepComputation= false;
175  fChangedGrandMotherRefFrame= false; // For local exit normal
176 
177  if( considerDirection && pGlobalDirection != 0 )
178  {
179  globalDirection=*pGlobalDirection;
180  }
181 
182 
183 #ifdef G4VERBOSE
184  if( fVerbose > 2 )
185  {
186  G4int oldcoutPrec = G4cout.precision(8);
187  G4cout << "*** G4ITNavigator1::LocateGlobalPointAndSetup: ***" << G4endl;
188  G4cout << " Called with arguments: " << G4endl
189  << " Globalpoint = " << globalPoint << G4endl
190  << " RelativeSearch = " << relativeSearch << G4endl;
191  if( fVerbose == 4 )
192  {
193  G4cout << " ----- Upon entering:" << G4endl;
194  PrintState();
195  }
196  G4cout.precision(oldcoutPrec);
197  }
198 #endif
199 
200  if ( !relativeSearch )
201  {
203  }
204  else
205  {
206  if ( fWasLimitedByGeometry )
207  {
208  fWasLimitedByGeometry = false;
209  fEnteredDaughter = fEntering; // Remember
210  fExitedMother = fExiting; // Remember
211  if ( fExiting )
212  {
213  if ( fHistory.GetDepth() )
214  {
215  fBlockedPhysicalVolume = fHistory.GetTopVolume();
216  fBlockedReplicaNo = fHistory.GetTopReplicaNo();
218  }
219  else
220  {
221  fLastLocatedPointLocal = localPoint;
222  fLocatedOutsideWorld = true;
223  return 0; // Have exited world volume
224  }
225  // A fix for the case where a volume is "entered" at an edge
226  // and a coincident surface exists outside it.
227  // - This stops it from exiting further volumes and cycling
228  // - However ReplicaNavigator treats this case itself
229  //
230  if ( fLocatedOnEdge && (VolumeType(fBlockedPhysicalVolume)!=kReplica ))
231  {
232  fExiting= false;
233  }
234  }
235  else
236  if ( fEntering )
237  {
238  switch (VolumeType(fBlockedPhysicalVolume))
239  {
240  case kNormal:
241  fHistory.NewLevel(fBlockedPhysicalVolume, kNormal,
242  fBlockedPhysicalVolume->GetCopyNo());
243  break;
244  case kReplica:
245  freplicaNav.ComputeTransformation(fBlockedReplicaNo,
246  fBlockedPhysicalVolume);
247  fHistory.NewLevel(fBlockedPhysicalVolume, kReplica,
248  fBlockedReplicaNo);
249  fBlockedPhysicalVolume->SetCopyNo(fBlockedReplicaNo);
250  break;
251  case kParameterised:
252  if( fBlockedPhysicalVolume->GetRegularStructureId() == 0 )
253  {
254  G4VSolid *pSolid;
255  G4VPVParameterisation *pParam;
256  G4TouchableHistory parentTouchable( fHistory );
257  pParam = fBlockedPhysicalVolume->GetParameterisation();
258  pSolid = pParam->ComputeSolid(fBlockedReplicaNo,
259  fBlockedPhysicalVolume);
260  pSolid->ComputeDimensions(pParam, fBlockedReplicaNo,
261  fBlockedPhysicalVolume);
262  pParam->ComputeTransformation(fBlockedReplicaNo,
263  fBlockedPhysicalVolume);
264  fHistory.NewLevel(fBlockedPhysicalVolume, kParameterised,
265  fBlockedReplicaNo);
266  fBlockedPhysicalVolume->SetCopyNo(fBlockedReplicaNo);
267  //
268  // Set the correct solid and material in Logical Volume
269  //
270  G4LogicalVolume *pLogical;
271  pLogical = fBlockedPhysicalVolume->GetLogicalVolume();
272  pLogical->SetSolid( pSolid );
273  pLogical->UpdateMaterial(pParam ->
274  ComputeMaterial(fBlockedReplicaNo,
275  fBlockedPhysicalVolume,
276  &parentTouchable));
277  }
278  break;
279  }
280  fEntering = false;
281  fBlockedPhysicalVolume = 0;
282  localPoint = fHistory.GetTopTransform().TransformPoint(globalPoint);
283  notKnownContained = false;
284  }
285  }
286  else
287  {
288  fBlockedPhysicalVolume = 0;
289  fEntering = false;
290  fEnteredDaughter = false; // Full Step was not taken, did not enter
291  fExiting = false;
292  fExitedMother = false; // Full Step was not taken, did not exit
293  }
294  }
295  //
296  // Search from top of history up through geometry until
297  // containing volume found:
298  // If on
299  // o OUTSIDE - Back up level, not/no longer exiting volumes
300  // o SURFACE and EXITING - Back up level, setting new blocking no.s
301  // else
302  // o containing volume found
303  //
304  G4int noLevelsExited=0 ;
305 
306  while (notKnownContained)
307  {
309  {
310  targetSolid = fHistory.GetTopVolume()->GetLogicalVolume()->GetSolid();
311  localPoint = fHistory.GetTopTransform().TransformPoint(globalPoint);
312  insideCode = targetSolid->Inside(localPoint);
313 #ifdef G4VERBOSE
314  if(( fVerbose == 1 ) && ( fCheck ))
315  {
316  G4String solidResponse = "-kInside-";
317  if (insideCode == kOutside)
318  solidResponse = "-kOutside-";
319  else if (insideCode == kSurface)
320  solidResponse = "-kSurface-";
321  G4cout << "*** G4ITNavigator1::LocateGlobalPointAndSetup(): ***" << G4endl
322  << " Invoked Inside() for solid: " << targetSolid->GetName()
323  << ". Solid replied: " << solidResponse << G4endl
324  << " For local point p: " << localPoint << G4endl;
325  }
326 #endif
327  }
328  else
329  {
330  insideCode = freplicaNav.BackLocate(fHistory, globalPoint, localPoint,
331  fExiting, notKnownContained);
332  // !CARE! if notKnownContained returns false then the point is within
333  // the containing placement volume of the replica(s). If insidecode
334  // will result in the history being backed up one level, then the
335  // local point returned is the point in the system of this new level
336  }
337 
338 
339  if ( insideCode==kOutside )
340  {
341  noLevelsExited++;
342  if ( fHistory.GetDepth() )
343  {
344  fBlockedPhysicalVolume = fHistory.GetTopVolume();
345  fBlockedReplicaNo = fHistory.GetTopReplicaNo();
347  fExiting = false;
348 
349  if( noLevelsExited > 1 )
350  {
351  // The first transformation was done by the sub-navigator
352  //
353  const G4RotationMatrix* mRot = fBlockedPhysicalVolume->GetRotation();
354  if( mRot )
355  {
356  fGrandMotherExitNormal *= (*mRot).inverse();
357  fChangedGrandMotherRefFrame= true;
358  }
359  }
360  }
361  else
362  {
363  fLastLocatedPointLocal = localPoint;
364  fLocatedOutsideWorld = true;
365  // No extra transformation for ExitNormal - is in frame of Top Volume
366  return 0; // Have exited world volume
367  }
368  }
369  else
370  if ( insideCode==kSurface )
371  {
372  G4bool isExiting = fExiting;
373  if( (!fExiting)&&considerDirection )
374  {
375  // Figure out whether we are exiting this level's volume
376  // by using the direction
377  //
378  G4bool directionExiting = false;
379  G4ThreeVector localDirection =
380  fHistory.GetTopTransform().TransformAxis(globalDirection);
381 
382  // Make sure localPoint in correct reference frame
383  // ( Was it already correct ? How ? )
384  //
385  localPoint= fHistory.GetTopTransform().TransformPoint(globalPoint);
387  {
388  G4ThreeVector normal = targetSolid->SurfaceNormal(localPoint);
389  directionExiting = normal.dot(localDirection) > 0.0;
390  isExiting = isExiting || directionExiting;
391  }
392  }
393  if( isExiting )
394  {
395  noLevelsExited++;
396  if ( fHistory.GetDepth() )
397  {
398  fBlockedPhysicalVolume = fHistory.GetTopVolume();
399  fBlockedReplicaNo = fHistory.GetTopReplicaNo();
401  //
402  // Still on surface but exited volume not necessarily convex
403  //
404  fValidExitNormal = false;
405 
406  if( noLevelsExited > 1 )
407  {
408  // The first transformation was done by the sub-navigator
409  //
410  const G4RotationMatrix* mRot=fBlockedPhysicalVolume->GetRotation();
411  if( mRot )
412  {
413  fGrandMotherExitNormal *= (*mRot).inverse();
414  fChangedGrandMotherRefFrame= true;
415  }
416  }
417  }
418  else
419  {
420  fLastLocatedPointLocal = localPoint;
421  fLocatedOutsideWorld = true;
422  // No extra transformation for ExitNormal, is in frame of Top Vol
423  return 0; // Have exited world volume
424  }
425  }
426  else
427  {
428  notKnownContained=false;
429  }
430  }
431  else
432  {
433  notKnownContained=false;
434  }
435  } // END while (notKnownContained)
436  //
437  // Search downwards until deepest containing volume found,
438  // blocking fBlockedPhysicalVolume/BlockedReplicaNum
439  //
440  // 3 Cases:
441  //
442  // o Parameterised daughters
443  // =>Must be one G4PVParameterised daughter & voxels
444  // o Positioned daughters & voxels
445  // o Positioned daughters & no voxels
446 
447  noResult = true; // noResult should be renamed to
448  // something like enteredLevel, as that is its meaning.
449  do
450  {
451  // Determine `type' of current mother volume
452  //
453  targetPhysical = fHistory.GetTopVolume();
454  if (!targetPhysical) { break; }
455  targetLogical = targetPhysical->GetLogicalVolume();
456  switch( CharacteriseDaughters(targetLogical) )
457  {
458  case kNormal:
459  if ( targetLogical->GetVoxelHeader() ) // use optimised navigation
460  {
461  noResult = fvoxelNav.LevelLocate(fHistory,
462  fBlockedPhysicalVolume,
463  fBlockedReplicaNo,
464  globalPoint,
465  pGlobalDirection,
466  considerDirection,
467  localPoint);
468  }
469  else // do not use optimised navigation
470  {
471  noResult = fnormalNav.LevelLocate(fHistory,
472  fBlockedPhysicalVolume,
473  fBlockedReplicaNo,
474  globalPoint,
475  pGlobalDirection,
476  considerDirection,
477  localPoint);
478  }
479  break;
480  case kReplica:
481  noResult = freplicaNav.LevelLocate(fHistory,
482  fBlockedPhysicalVolume,
483  fBlockedReplicaNo,
484  globalPoint,
485  pGlobalDirection,
486  considerDirection,
487  localPoint);
488  break;
489  case kParameterised:
490  if( GetDaughtersRegularStructureId(targetLogical) != 1 )
491  {
492  noResult = fparamNav.LevelLocate(fHistory,
493  fBlockedPhysicalVolume,
494  fBlockedReplicaNo,
495  globalPoint,
496  pGlobalDirection,
497  considerDirection,
498  localPoint);
499  }
500  else // Regular structure
501  {
502  noResult = fregularNav.LevelLocate(fHistory,
503  fBlockedPhysicalVolume,
504  fBlockedReplicaNo,
505  globalPoint,
506  pGlobalDirection,
507  considerDirection,
508  localPoint);
509  }
510  break;
511  }
512 
513  // LevelLocate returns true if it finds a daughter volume
514  // in which globalPoint is inside (or on the surface).
515 
516  if ( noResult )
517  {
518  // Entering a daughter after ascending
519  //
520  // The blocked volume is no longer valid - it was for another level
521  //
522  fBlockedPhysicalVolume = 0;
523  fBlockedReplicaNo = -1;
524 
525  // fEntering should be false -- else blockedVolume is assumed good.
526  // fEnteredDaughter is used for ExitNormal
527  //
528  fEntering = false;
529  fEnteredDaughter = true;
530 
531  if( fExitedMother )
532  {
533  G4VPhysicalVolume* enteredPhysical = fHistory.GetTopVolume();
534  const G4RotationMatrix* mRot = enteredPhysical->GetRotation();
535  if( mRot )
536  {
537  fGrandMotherExitNormal *= (*mRot).inverse();
538  }
539  }
540 
541 #ifdef G4DEBUG_NAVIGATION
542  if( fVerbose > 2 )
543  {
544  G4VPhysicalVolume* enteredPhysical = fHistory.GetTopVolume();
545  G4cout << "*** G4ITNavigator1::LocateGlobalPointAndSetup() ***" << G4endl;
546  G4cout << " Entering volume: " << enteredPhysical->GetName()
547  << G4endl;
548  }
549 #endif
550  }
551  } while (noResult);
552 
553  fLastLocatedPointLocal = localPoint;
554 
555 #ifdef G4VERBOSE
556  if( fVerbose >= 4 )
557  {
558  G4int oldcoutPrec = G4cout.precision(8);
559  G4String curPhysVol_Name("None");
560  if (targetPhysical) { curPhysVol_Name = targetPhysical->GetName(); }
561  G4cout << " Return value = new volume = " << curPhysVol_Name << G4endl;
562  G4cout << " ----- Upon exiting:" << G4endl;
563  PrintState();
564  if( fVerbose == 5 )
565  {
566  G4cout << "Upon exiting LocateGlobalPointAndSetup():" << G4endl;
567  G4cout << " History = " << G4endl << fHistory << G4endl << G4endl;
568  }
569  G4cout.precision(oldcoutPrec);
570  }
571 #endif
572 
573  fLocatedOutsideWorld= false;
574 
575  return targetPhysical;
576 }
G4String GetName() const
G4SmartVoxelHeader * GetVoxelHeader() const
G4VPhysicalVolume * GetTopVolume() const
virtual G4bool LevelLocate(G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint)
void UpdateMaterial(G4Material *pMaterial)
G4int GetDaughtersRegularStructureId(const G4LogicalVolume *pLog) const
double dot(const Hep3Vector &) const
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:138
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
G4int GetDepth() const
G4VSolid * GetSolid() const
EVolume VolumeType(const G4VPhysicalVolume *pVol) const
void SetSolid(G4VSolid *pSolid)
const G4RotationMatrix * GetRotation() const
G4bool LevelLocate(G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint)
void ComputeTransformation(const G4int replicaNo, G4VPhysicalVolume *pVol, G4ThreeVector &point) const
int G4int
Definition: G4Types.hh:78
void ResetStackAndState()
G4int GetTopReplicaNo() const
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
virtual EInside Inside(const G4ThreeVector &p) const =0
bool G4bool
Definition: G4Types.hh:79
EVolume GetTopVolumeType() const
virtual G4int GetRegularStructureId() const =0
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
virtual G4VPVParameterisation * GetParameterisation() const =0
virtual void SetCopyNo(G4int CopyNo)=0
G4NavigationHistory fHistory
G4bool LevelLocate(G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint)
void NewLevel(G4VPhysicalVolume *pNewMother, EVolume vType=kNormal, G4int nReplica=-1)
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4LogicalVolume * GetLogicalVolume() const
EInside
Definition: geomdefs.hh:58
EVolume CharacteriseDaughters(const G4LogicalVolume *pLog) const
virtual void ComputeTransformation(const G4int, G4VPhysicalVolume *) const =0
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
virtual G4int GetCopyNo() const =0
#define G4endl
Definition: G4ios.hh:61
const G4AffineTransform & GetTopTransform() const
G4bool LevelLocate(G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint)
G4bool LevelLocate(G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint)
EInside BackLocate(G4NavigationHistory &history, const G4ThreeVector &globalPoint, G4ThreeVector &localPoint, const G4bool &exiting, G4bool &notKnownInside) const
void PrintState() const
G4bool fWasLimitedByGeometry

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ITNavigator1::LocateGlobalPointAndUpdateTouchable ( const G4ThreeVector position,
const G4ThreeVector direction,
G4VTouchable touchableToUpdate,
const G4bool  RelativeSearch = true 
)
inline
void G4ITNavigator1::LocateGlobalPointAndUpdateTouchable ( const G4ThreeVector position,
G4VTouchable touchableToUpdate,
const G4bool  RelativeSearch = true 
)
inline
void G4ITNavigator1::LocateGlobalPointAndUpdateTouchableHandle ( const G4ThreeVector position,
const G4ThreeVector direction,
G4TouchableHandle oldTouchableToUpdate,
const G4bool  RelativeSearch = true 
)
inline
void G4ITNavigator1::LocateGlobalPointWithinVolume ( const G4ThreeVector position)
virtual

Definition at line 592 of file G4ITNavigator1.cc.

593 {
594  fLastLocatedPointLocal = ComputeLocalPoint(pGlobalpoint);
595  fLastTriedStepComputation= false;
596  fChangedGrandMotherRefFrame= false; // Frame for Exit Normal
597 
598 #ifdef G4DEBUG_NAVIGATION
599  if( fVerbose > 2 )
600  {
601  G4cout << "Entering LocateGlobalWithinVolume(): History = " << G4endl;
602  G4cout << fHistory << G4endl;
603  }
604 #endif
605 
606  // For the case of Voxel (or Parameterised) volume the respective
607  // Navigator must be messaged to update its voxel information etc
608 
609  // Update the state of the Sub Navigators
610  // - in particular any voxel information they store/cache
611  //
612  G4VPhysicalVolume* motherPhysical = fHistory.GetTopVolume();
613  G4LogicalVolume* motherLogical = motherPhysical->GetLogicalVolume();
614  G4SmartVoxelHeader* pVoxelHeader = motherLogical->GetVoxelHeader();
615 
617  {
618  switch( CharacteriseDaughters(motherLogical) )
619  {
620  case kNormal:
621  if ( pVoxelHeader )
622  {
623  fvoxelNav.VoxelLocate( pVoxelHeader, fLastLocatedPointLocal );
624  }
625  break;
626  case kParameterised:
627  if( GetDaughtersRegularStructureId(motherLogical) != 1 )
628  {
629  // Resets state & returns voxel node
630  //
631  fparamNav.ParamVoxelLocate( pVoxelHeader, fLastLocatedPointLocal );
632  }
633  break;
634  case kReplica:
635  G4Exception("G4ITNavigator1::LocateGlobalPointWithinVolume()",
636  "GeomNav0001", FatalException,
637  "Not applicable for replicated volumes.");
638  break;
639  }
640  }
641 
642  // Reset the state variables
643  // - which would have been affected
644  // by the 'equivalent' call to LocateGlobalPointAndSetup
645  // - who's values have been invalidated by the 'move'.
646  //
647  fBlockedPhysicalVolume = 0;
648  fBlockedReplicaNo = -1;
649  fEntering = false;
650  fEnteredDaughter = false; // Boundary not encountered, did not enter
651  fExiting = false;
652  fExitedMother = false; // Boundary not encountered, did not exit
653 }
G4SmartVoxelHeader * GetVoxelHeader() const
G4VPhysicalVolume * GetTopVolume() const
G4int GetDaughtersRegularStructureId(const G4LogicalVolume *pLog) const
G4ThreeVector ComputeLocalPoint(const G4ThreeVector &rGlobPoint) const
G4SmartVoxelNode * ParamVoxelLocate(G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint)
G4GLOB_DLL std::ostream G4cout
EVolume GetTopVolumeType() const
G4NavigationHistory fHistory
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4LogicalVolume * GetLogicalVolume() const
EVolume CharacteriseDaughters(const G4LogicalVolume *pLog) const
G4SmartVoxelNode * VoxelLocate(G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint)
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:

G4RotationMatrix G4ITNavigator1::NetRotation ( ) const
inline
G4ThreeVector G4ITNavigator1::NetTranslation ( ) const
inline
void G4ITNavigator1::NewNavigatorState ( )

Definition at line 668 of file G4ITNavigator1.cc.

669 {
670  fpSaveState = new G4SaveNavigatorState();
671  ResetState();
672 }
virtual void ResetState()

Here is the call graph for this function:

void G4ITNavigator1::PrintState ( ) const

Definition at line 1841 of file G4ITNavigator1.cc.

1842 {
1843  G4int oldcoutPrec = G4cout.precision(4);
1844  if( fVerbose == 4 )
1845  {
1846  G4cout << "The current state of G4ITNavigator1 is: " << G4endl;
1847  G4cout << " ValidExitNormal= " << fValidExitNormal << G4endl
1848  << " ExitNormal = " << fExitNormal << G4endl
1849  << " Exiting = " << fExiting << G4endl
1850  << " Entering = " << fEntering << G4endl
1851  << " BlockedPhysicalVolume= " ;
1852  if (fBlockedPhysicalVolume==0)
1853  G4cout << "None";
1854  else
1855  G4cout << fBlockedPhysicalVolume->GetName();
1856  G4cout << G4endl
1857  << " BlockedReplicaNo = " << fBlockedReplicaNo << G4endl
1858  << " LastStepWasZero = " << fLastStepWasZero << G4endl
1859  << G4endl;
1860  }
1861  if( ( 1 < fVerbose) && (fVerbose < 4) )
1862  {
1863  G4cout << G4endl; // Make sure to line up
1864  G4cout << std::setw(30) << " ExitNormal " << " "
1865  << std::setw( 5) << " Valid " << " "
1866  << std::setw( 9) << " Exiting " << " "
1867  << std::setw( 9) << " Entering" << " "
1868  << std::setw(15) << " Blocked:Volume " << " "
1869  << std::setw( 9) << " ReplicaNo" << " "
1870  << std::setw( 8) << " LastStepZero " << " "
1871  << G4endl;
1872  G4cout << "( " << std::setw(7) << fExitNormal.x()
1873  << ", " << std::setw(7) << fExitNormal.y()
1874  << ", " << std::setw(7) << fExitNormal.z() << " ) "
1875  << std::setw( 5) << fValidExitNormal << " "
1876  << std::setw( 9) << fExiting << " "
1877  << std::setw( 9) << fEntering << " ";
1878  if ( fBlockedPhysicalVolume==0 )
1879  {
1880  G4cout << std::setw(15) << "None";
1881  }
1882  else
1883  {
1884  G4cout << std::setw(15)<< fBlockedPhysicalVolume->GetName();
1885  }
1886  G4cout << std::setw( 9) << fBlockedReplicaNo << " "
1887  << std::setw( 8) << fLastStepWasZero << " "
1888  << G4endl;
1889  }
1890  if( fVerbose > 2 )
1891  {
1892  G4cout.precision(8);
1893  G4cout << " Current Localpoint = " << fLastLocatedPointLocal << G4endl;
1894  G4cout << " PreviousSftOrigin = " << fPreviousSftOrigin << G4endl;
1895  G4cout << " PreviousSafety = " << fPreviousSafety << G4endl;
1896  }
1897  G4cout.precision(oldcoutPrec);
1898 }
double x() const
int G4int
Definition: G4Types.hh:78
double z() const
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
double y() const
#define G4endl
Definition: G4ios.hh:61

Here is the caller graph for this function:

G4VPhysicalVolume * G4ITNavigator1::ResetHierarchyAndLocate ( const G4ThreeVector point,
const G4ThreeVector direction,
const G4TouchableHistory h 
)
virtual

Definition at line 132 of file G4ITNavigator1.cc.

135 {
136  ResetState();
137  fHistory = *h.GetHistory();
138  SetupHierarchy();
139  fLastTriedStepComputation= false; // Redundant, but best
140  return LocateGlobalPointAndSetup(p, &direction, true, false);
141 }
virtual void ResetState()
virtual G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=0, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
const char * p
Definition: xmltok.h:285
virtual void SetupHierarchy()
G4NavigationHistory fHistory
const G4NavigationHistory * GetHistory() const

Here is the call graph for this function:

void G4ITNavigator1::ResetStackAndState ( )
inline

Here is the caller graph for this function:

void G4ITNavigator1::ResetState ( )
virtual

Definition at line 1259 of file G4ITNavigator1.cc.

1260 {
1261  fWasLimitedByGeometry = false;
1262  fEntering = false;
1263  fExiting = false;
1264  fLocatedOnEdge = false;
1265  fLastStepWasZero = false;
1266  fEnteredDaughter = false;
1267  fExitedMother = false;
1268  fPushed = false;
1269 
1270  fValidExitNormal = false;
1271  fChangedGrandMotherRefFrame= false;
1272  fCalculatedExitNormal = false;
1273 
1274  fExitNormal = G4ThreeVector(0,0,0);
1275  fGrandMotherExitNormal = G4ThreeVector(0,0,0);
1276  fExitNormalGlobalFrame = G4ThreeVector(0,0,0);
1277 
1278  fPreviousSftOrigin = G4ThreeVector(0,0,0);
1279  fPreviousSafety = 0.0;
1280 
1281  fNumberZeroSteps = 0;
1282 
1283  fBlockedPhysicalVolume = 0;
1284  fBlockedReplicaNo = -1;
1285 
1286  fLastLocatedPointLocal = G4ThreeVector( kInfinity, -kInfinity, 0.0 );
1287  fLocatedOutsideWorld = false;
1288 }
static const G4double kInfinity
Definition: geomdefs.hh:42
CLHEP::Hep3Vector G4ThreeVector
G4bool fWasLimitedByGeometry

Here is the caller graph for this function:

void G4ITNavigator1::RestoreSavedState ( )

Definition at line 727 of file G4ITNavigator1.cc.

728 {
729  fExitNormal = fpSaveState->sExitNormal;
730  fValidExitNormal = fpSaveState->sValidExitNormal;
731  fExiting = fpSaveState->sExiting;
732  fEntering = fpSaveState->sEntering;
733 
734  fBlockedPhysicalVolume = fpSaveState->spBlockedPhysicalVolume;
735  fBlockedReplicaNo = fpSaveState->sBlockedReplicaNo,
736 
737  fLastStepWasZero = fpSaveState->sLastStepWasZero;
738 
739  // !>
740  fPreviousSftOrigin = fpSaveState->sPreviousSftOrigin ;
741  fPreviousSafety = fpSaveState->sPreviousSafety ;
742  fNumberZeroSteps = fpSaveState->sNumberZeroSteps ;
743  fLocatedOnEdge = fpSaveState->sLocatedOnEdge ;
744  fWasLimitedByGeometry = fpSaveState->sWasLimitedByGeometry;
745  fPushed = fpSaveState->sPushed;
746  fNumberZeroSteps = fpSaveState->sNumberZeroSteps;
747  fEnteredDaughter= fpSaveState->sEnteredDaughter ;
748  fExitedMother = fpSaveState->sExitedMother ;
749 
750  fLastLocatedPointLocal = fpSaveState->sLastLocatedPointLocal ;
751  fLocatedOutsideWorld = fpSaveState->sLocatedOutsideWorld;
752  // <!
753 }
G4bool fWasLimitedByGeometry

Here is the caller graph for this function:

void G4ITNavigator1::SetGeometricallyLimitedStep ( )
inline
void G4ITNavigator1::SetNavigatorState ( G4ITNavigatorState_Lock1 navState)

Definition at line 662 of file G4ITNavigator1.cc.

663 {
664  fpSaveState = (G4SaveNavigatorState*) navState;
665  if(navState) RestoreSavedState();
666 }
void RestoreSavedState()

Here is the call graph for this function:

void G4ITNavigator1::SetPushVerbosity ( G4bool  mode)
inline
void G4ITNavigator1::SetSavedState ( )

Definition at line 683 of file G4ITNavigator1.cc.

684 {
685  // !>
686  // This check can be avoid if instead, at every first step of a track,
687  // the IT tracking uses NewNavigatorSate
688  // The normal tracking would just call once NewNavigatorState() before tracking
689 
690 // if(fpSaveState == 0)
691 // fpSaveState = new G4SaveNavigatorState;
692  // <!
693 
694  // fSaveExitNormal = fExitNormal;
695  fpSaveState->sExitNormal = fExitNormal;
696  fpSaveState->sValidExitNormal = fValidExitNormal;
697  fpSaveState->sExiting = fExiting;
698  fpSaveState->sEntering = fEntering;
699 
700  fpSaveState->spBlockedPhysicalVolume = fBlockedPhysicalVolume;
701  fpSaveState->sBlockedReplicaNo = fBlockedReplicaNo,
702 
703  fpSaveState->sLastStepWasZero = fLastStepWasZero;
704 
705  // !>
706  fpSaveState->sPreviousSftOrigin = fPreviousSftOrigin;
707  fpSaveState->sPreviousSafety = fPreviousSafety;
708  fpSaveState->sNumberZeroSteps = fNumberZeroSteps;
709  fpSaveState->sLocatedOnEdge = fLocatedOnEdge;
710  fpSaveState->sWasLimitedByGeometry= fWasLimitedByGeometry;
711  fpSaveState->sPushed=fPushed;
712  fpSaveState->sNumberZeroSteps=fNumberZeroSteps;
713  fpSaveState->sEnteredDaughter = fEnteredDaughter;
714  fpSaveState->sExitedMother = fExitedMother;
715 
716  fpSaveState->sLastLocatedPointLocal = fLastLocatedPointLocal;
717  fpSaveState->sLocatedOutsideWorld = fLocatedOutsideWorld;
718  // <!
719 }
G4bool fWasLimitedByGeometry

Here is the caller graph for this function:

void G4ITNavigator1::SetupHierarchy ( )
protectedvirtual

Definition at line 1298 of file G4ITNavigator1.cc.

1299 {
1300  G4int i;
1301  const G4int cdepth = fHistory.GetDepth();
1302  G4VPhysicalVolume *current;
1303  G4VSolid *pSolid;
1304  G4VPVParameterisation *pParam;
1305 
1306  for ( i=1; i<=cdepth; i++ )
1307  {
1308  current = fHistory.GetVolume(i);
1309  switch ( fHistory.GetVolumeType(i) )
1310  {
1311  case kNormal:
1312  break;
1313  case kReplica:
1314  freplicaNav.ComputeTransformation(fHistory.GetReplicaNo(i), current);
1315  break;
1316  case kParameterised:
1317  G4int replicaNo;
1318  pParam = current->GetParameterisation();
1319  replicaNo = fHistory.GetReplicaNo(i);
1320  pSolid = pParam->ComputeSolid(replicaNo, current);
1321 
1322  // Set up dimensions & transform in solid/physical volume
1323  //
1324  pSolid->ComputeDimensions(pParam, replicaNo, current);
1325  pParam->ComputeTransformation(replicaNo, current);
1326 
1327  G4TouchableHistory touchable( fHistory );
1328  touchable.MoveUpHistory(); // move up to the parent level
1329 
1330  // Set up the correct solid and material in Logical Volume
1331  //
1332  G4LogicalVolume *pLogical = current->GetLogicalVolume();
1333  pLogical->SetSolid( pSolid );
1334  pLogical->UpdateMaterial( pParam ->
1335  ComputeMaterial(replicaNo, current, &touchable) );
1336  break;
1337  }
1338  }
1339 }
void UpdateMaterial(G4Material *pMaterial)
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:138
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
G4int GetDepth() const
void SetSolid(G4VSolid *pSolid)
EVolume GetVolumeType(G4int n) const
void ComputeTransformation(const G4int replicaNo, G4VPhysicalVolume *pVol, G4ThreeVector &point) const
int G4int
Definition: G4Types.hh:78
virtual G4VPVParameterisation * GetParameterisation() const =0
G4int GetReplicaNo(G4int n) const
G4NavigationHistory fHistory
G4LogicalVolume * GetLogicalVolume() const
virtual void ComputeTransformation(const G4int, G4VPhysicalVolume *) const =0
G4VPhysicalVolume * GetVolume(G4int n) const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ITNavigator1::SetVerboseLevel ( G4int  level)
inline

Here is the caller graph for this function:

void G4ITNavigator1::SetWorldVolume ( G4VPhysicalVolume pWorld)
inline
G4int G4ITNavigator1::SeverityOfZeroStepping ( G4int noZeroSteps) const
inline
EVolume G4ITNavigator1::VolumeType ( const G4VPhysicalVolume pVol) const
inlineprotected

Here is the caller graph for this function:

Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  os,
const G4ITNavigator1 n 
)
friend

Definition at line 2007 of file G4ITNavigator1.cc.

2008 {
2009  // Old version did only the following:
2010  // os << "Current History: " << G4endl << n.fHistory;
2011  // Old behaviour is recovered for fVerbose = 0
2012 
2013  // Adapted from G4ITNavigator1::PrintState() const
2014 
2015  G4int oldcoutPrec = os.precision(4);
2016  if( n.fVerbose >= 4 )
2017  {
2018  os << "The current state of G4ITNavigator1 is: " << G4endl;
2019  os << " ValidExitNormal= " << n.fValidExitNormal << G4endl
2020  << " ExitNormal = " << n.fExitNormal << G4endl
2021  << " Exiting = " << n.fExiting << G4endl
2022  << " Entering = " << n.fEntering << G4endl
2023  << " BlockedPhysicalVolume= " ;
2024  if (n.fBlockedPhysicalVolume==0)
2025  os << "None";
2026  else
2027  os << n.fBlockedPhysicalVolume->GetName();
2028  os << G4endl
2029  << " BlockedReplicaNo = " << n.fBlockedReplicaNo << G4endl
2030  << " LastStepWasZero = " << n.fLastStepWasZero << G4endl
2031  << G4endl;
2032  }
2033  if( ( 1 < n.fVerbose) && (n.fVerbose < 4) )
2034  {
2035  os << G4endl; // Make sure to line up
2036  os << std::setw(30) << " ExitNormal " << " "
2037  << std::setw( 5) << " Valid " << " "
2038  << std::setw( 9) << " Exiting " << " "
2039  << std::setw( 9) << " Entering" << " "
2040  << std::setw(15) << " Blocked:Volume " << " "
2041  << std::setw( 9) << " ReplicaNo" << " "
2042  << std::setw( 8) << " LastStepZero " << " "
2043  << G4endl;
2044  os << "( " << std::setw(7) << n.fExitNormal.x()
2045  << ", " << std::setw(7) << n.fExitNormal.y()
2046  << ", " << std::setw(7) << n.fExitNormal.z() << " ) "
2047  << std::setw( 5) << n.fValidExitNormal << " "
2048  << std::setw( 9) << n.fExiting << " "
2049  << std::setw( 9) << n.fEntering << " ";
2050  if ( n.fBlockedPhysicalVolume==0 )
2051  { os << std::setw(15) << "None"; }
2052  else
2053  { os << std::setw(15)<< n.fBlockedPhysicalVolume->GetName(); }
2054  os << std::setw( 9) << n.fBlockedReplicaNo << " "
2055  << std::setw( 8) << n.fLastStepWasZero << " "
2056  << G4endl;
2057  }
2058  if( n.fVerbose > 2 )
2059  {
2060  os.precision(8);
2061  os << " Current Localpoint = " << n.fLastLocatedPointLocal << G4endl;
2062  os << " PreviousSftOrigin = " << n.fPreviousSftOrigin << G4endl;
2063  os << " PreviousSafety = " << n.fPreviousSafety << G4endl;
2064  }
2065  if( n.fVerbose > 3 || n.fVerbose == 0 )
2066  {
2067  os << "Current History: " << G4endl << n.fHistory;
2068  }
2069 
2070  os.precision(oldcoutPrec);
2071  return os;
2072 }
double x() const
int G4int
Definition: G4Types.hh:78
double z() const
const G4String & GetName() const
G4NavigationHistory fHistory
double y() const
#define G4endl
Definition: G4ios.hh:61

Member Data Documentation

G4bool G4ITNavigator1::fEnteredDaughter
protected

Definition at line 379 of file G4ITNavigator1.hh.

G4bool G4ITNavigator1::fExitedMother
protected

Definition at line 385 of file G4ITNavigator1.hh.

G4NavigationHistory G4ITNavigator1::fHistory
protected

Definition at line 375 of file G4ITNavigator1.hh.

G4ThreeVector G4ITNavigator1::fLastStepEndPointLocal
protected

Definition at line 395 of file G4ITNavigator1.hh.

const G4int G4ITNavigator1::fMaxNav = 8
static

Definition at line 94 of file G4ITNavigator1.hh.

G4ThreeVector G4ITNavigator1::fStepEndPoint
protected

Definition at line 392 of file G4ITNavigator1.hh.

G4int G4ITNavigator1::fVerbose
protected

Definition at line 399 of file G4ITNavigator1.hh.

G4bool G4ITNavigator1::fWasLimitedByGeometry
protected

Definition at line 389 of file G4ITNavigator1.hh.

G4double G4ITNavigator1::kCarTolerance
protected

Definition at line 368 of file G4ITNavigator1.hh.


The documentation for this class was generated from the following files: