Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ITNavigator2 Class Reference

#include <G4ITNavigator2.hh>

Classes

struct  G4NavigatorState
 
struct  G4SaveNavigatorState
 

Public Member Functions

 G4ITNavigator2 ()
 
virtual ~G4ITNavigator2 ()
 
G4ITNavigatorState_Lock2GetNavigatorState ()
 
void SetNavigatorState (G4ITNavigatorState_Lock2 *)
 
void NewNavigatorState ()
 
void NewNavigatorState (const G4TouchableHistory &h)
 
void ResetNavigatorState ()
 
G4VPhysicalVolumeNewNavigatorStateAndLocate (const G4ThreeVector &p, const G4ThreeVector &direction)
 
void CheckNavigatorState () const
 
std::shared_ptr< G4ITNavigatorState_Lock2GetSnapshotOfState ()
 
void ResetFromSnapshot (std::shared_ptr< G4ITNavigatorState_Lock2 >)
 
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)
 
EInside InsideCurrentVolume (const G4ThreeVector &globalPoint) const
 
void GetRandomInCurrentVolume (G4ThreeVector &rndmPoint) const
 
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)
 
virtual G4bool RecheckDistanceToCurrentBoundary (const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double CurrentProposedStepLength, G4double *prDistance, G4double *prNewSafety=0) const
 
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
 
G4ThreeVector GetCurrentLocalCoordinate () const
 
G4ThreeVector NetTranslation () const
 
G4RotationMatrix NetRotation () const
 
void EnableBestSafety (G4bool value=false)
 

Public Attributes

G4NavigatorStatefpNavigatorState
 
G4VPhysicalVolumefTopPhysical
 
G4bool fCheck
 
G4bool fWarnPush
 
G4NormalNavigation fnormalNav
 
G4VoxelNavigation fvoxelNav
 
G4ParameterisedNavigation fparamNav
 
G4ReplicaNavigation freplicaNav
 
G4RegularNavigation fregularNav
 
G4VoxelSafetyfpVoxelSafety
 

Static Public Attributes

static const G4int fMaxNav = 8
 

Protected Member Functions

G4ThreeVector ComputeLocalPoint (const G4ThreeVector &rGlobPoint) const
 
G4ThreeVector ComputeLocalAxis (const G4ThreeVector &pVec) const
 
virtual void ResetState ()
 
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
 
G4int fVerbose
 

Friends

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

Detailed Description

Definition at line 101 of file G4ITNavigator2.hh.

Constructor & Destructor Documentation

◆ G4ITNavigator2()

G4ITNavigator2::G4ITNavigator2 ( )

Definition at line 99 of file G4ITNavigator2.cc.

100 : fVerbose(0),
101 fTopPhysical(0), fCheck(false),
102 fWarnPush(true)
103{
104 fActive= false;
105 fActionThreshold_NoZeroSteps = 1000;
106 fAbandonThreshold_NoZeroSteps = 2500;
108 fregularNav.SetNormalNavigation( &fnormalNav );
109
111// fSaveState = 0;
113
114 // this->SetVerboseLevel(3);
115 // this->CheckMode(true);
116}
const G4double kCarTolerance
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4RegularNavigation fregularNav
G4VoxelSafety * fpVoxelSafety
G4NavigatorState * fpNavigatorState
G4VPhysicalVolume * fTopPhysical
void SetNormalNavigation(G4NormalNavigation *fnormnav)

◆ ~G4ITNavigator2()

G4ITNavigator2::~G4ITNavigator2 ( )
virtual

Definition at line 122 of file G4ITNavigator2.cc.

123{ delete fpVoxelSafety; }

Member Function Documentation

◆ Activate()

void G4ITNavigator2::Activate ( G4bool  flag)
inline

◆ CharacteriseDaughters()

EVolume G4ITNavigator2::CharacteriseDaughters ( const G4LogicalVolume pLog) const
inlineprotected

◆ CheckMode()

void G4ITNavigator2::CheckMode ( G4bool  mode)
inline

◆ CheckNavigatorState()

void G4ITNavigator2::CheckNavigatorState ( ) const

Definition at line 700 of file G4ITNavigator2.cc.

701{
702 if(fpNavigatorState == 0)
703 {
704 G4ExceptionDescription exceptionDescription;
705 exceptionDescription << "The navigator state is NULL. ";
706 exceptionDescription << "Either NewNavigatorStateAndLocate was not called ";
707 exceptionDescription << "or the provided navigator state was already NULL.";
708
709 G4Exception("G4ITNavigator::CheckNavigatorStateIsValid",
710 "NavigatorStateNotValid",FatalException,exceptionDescription);
711 return;
712 }
713}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40

◆ CheckNextStep()

G4double G4ITNavigator2::CheckNextStep ( const G4ThreeVector pGlobalPoint,
const G4ThreeVector pDirection,
const G4double  pCurrentProposedStepLength,
G4double pNewSafety 
)

Definition at line 1363 of file G4ITNavigator2.cc.

1367{
1369 G4double step;
1370
1371 // Save the state, for this parasitic call
1372 //
1373 //SetSavedState();
1374// G4SaveNavigatorState savedState(fpNavigatorState);
1375 G4NavigatorState savedState(*fpNavigatorState);
1376
1377 step = ComputeStep ( pGlobalpoint,
1378 pDirection,
1379 pCurrentProposedStepLength,
1380 pNewSafety );
1381
1382 // If a parasitic call, then attempt to restore the key parts of the state
1383 //
1384 *fpNavigatorState = savedState;
1385 //RestoreSavedState();
1386 // NOTE: the state of the current subnavigator is NOT restored.
1387 // ***> TODO: restore subnavigator state
1388 // if( last_located) Need Position of last location
1389 // if( last_computed step) Need Endposition of last step
1390
1391
1392 return step;
1393}
#define CheckNavigatorStateIsValid()
double G4double
Definition: G4Types.hh:83
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)

Referenced by ResetFromSnapshot().

◆ ComputeLocalAxis()

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

◆ ComputeLocalPoint()

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

◆ ComputeSafety()

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

Definition at line 1912 of file G4ITNavigator2.cc.

1915{
1917 G4double newSafety = 0.0;
1918
1919#ifdef G4DEBUG_NAVIGATION
1920 G4int oldcoutPrec = G4cout.precision(8);
1921 if( fVerbose > 0 )
1922 {
1923 G4cout << "*** G4ITNavigator2::ComputeSafety: ***" << G4endl
1924 << " Called at point: " << pGlobalpoint << G4endl;
1925
1926 G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
1927 G4cout << " Volume = " << motherPhysical->GetName()
1928 << " - Maximum length = " << pMaxLength << G4endl;
1929 if( fVerbose >= 4 )
1930 {
1931 G4cout << " ----- Upon entering Compute Safety:" << G4endl;
1932 PrintState();
1933 }
1934 }
1935#endif
1936
1937 G4SaveNavigatorState* savedState(0);
1938
1939 G4double distEndpointSq = (pGlobalpoint-fStepEndPoint).mag2();
1940 G4bool stayedOnEndpoint = distEndpointSq < kCarTolerance*kCarTolerance;
1941 G4bool endpointOnSurface = fEnteredDaughter || fExitedMother;
1942
1943 if( endpointOnSurface && stayedOnEndpoint )
1944 {
1945#ifdef G4DEBUG_NAVIGATION
1946 if( fVerbose >= 2 )
1947 {
1948 G4cout << " G4Navigator::ComputeSafety() finds that point - "
1949 << pGlobalpoint << " - is on surface " << G4endl;
1950 if( fEnteredDaughter ) { G4cout << " entered new daughter volume"; }
1951 if( fExitedMother ) { G4cout << " and exited previous volume."; }
1952 G4cout << G4endl;
1953 G4cout << " EndPoint was = " << fStepEndPoint << G4endl;
1954 }
1955#endif
1956 newSafety = 0.0;
1957 // return newSafety;
1958 }
1959 else // if( !(endpointOnSurface && stayedOnEndpoint) )
1960 {
1961
1962 if (keepState)
1963 {
1964// SetSavedState();
1965 savedState = new G4SaveNavigatorState(fpNavigatorState);
1966 }
1967
1968 // Pseudo-relocate to this point (updates voxel information only)
1969 //
1970 LocateGlobalPointWithinVolume( pGlobalpoint );
1971 // --->> DANGER: Side effects on sub-navigator voxel information <<---
1972 // Could be replaced again by 'granular' calls to sub-navigator
1973 // locates (similar side-effects, but faster.
1974 // Solutions:
1975 // 1) Re-locate (to where?)
1976 // 2) Insure that the methods using (G4ComputeStep?)
1977 // does a relocation (if information is disturbed only ?)
1978
1979#ifdef G4DEBUG_NAVIGATION
1980 if( fVerbose >= 2 )
1981 {
1982 G4cout << " G4ITNavigator2::ComputeSafety() relocates-in-volume to point: "
1983 << pGlobalpoint << G4endl;
1984 }
1985#endif
1986 G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
1987 G4LogicalVolume *motherLogical = motherPhysical->GetLogicalVolume();
1988 G4SmartVoxelHeader* pVoxelHeader = motherLogical->GetVoxelHeader();
1989 G4ThreeVector localPoint = ComputeLocalPoint(pGlobalpoint);
1990
1991 if ( fHistory.GetTopVolumeType()!=kReplica )
1992 {
1993 switch(CharacteriseDaughters(motherLogical))
1994 {
1995 case kNormal:
1996 if ( pVoxelHeader )
1997 {
1998#ifdef G4NEW_SAFETY
1999 G4double safetyTwo = fpVoxelSafety->ComputeSafety(localPoint,
2000 *motherPhysical, pMaxLength);
2001 newSafety= safetyTwo; // Faster and best available
2002#else
2003 G4double safetyOldVoxel;
2004 LocateGlobalPointWithinVolume(pGlobalpoint);
2005 safetyOldVoxel =
2006 fvoxelNav.ComputeSafety(localPoint,fHistory,pMaxLength);
2007 newSafety= safetyOldVoxel;
2008#endif
2009 }
2010 else
2011 {
2012 newSafety=fnormalNav.ComputeSafety(localPoint,fHistory,pMaxLength);
2013 }
2014 break;
2015 case kParameterised:
2016 if( GetDaughtersRegularStructureId(motherLogical) != 1 )
2017 {
2018 newSafety=fparamNav.ComputeSafety(localPoint,fHistory,pMaxLength);
2019 }
2020 else // Regular structure
2021 {
2022 newSafety=fregularNav.ComputeSafety(localPoint,fHistory,pMaxLength);
2023 }
2024 break;
2025 case kReplica:
2026 G4Exception("G4ITNavigator2::ComputeSafety()", "GeomNav0001",
2027 FatalException, "Not applicable for replicated volumes.");
2028 break;
2029 case kExternal:
2030 G4Exception("G4ITNavigator2::ComputeSafety()",
2031 "GeomNav0001", FatalException,
2032 "Not applicable for external volumes.");
2033 break;
2034 }
2035 }
2036 else
2037 {
2038 newSafety = freplicaNav.ComputeSafety(pGlobalpoint, localPoint,
2039 fHistory, pMaxLength);
2040 }
2041
2042 if (keepState)
2043 {
2044 *fpNavigatorState = *savedState;
2045 delete savedState;
2046 // RestoreSavedState();
2047 // This now overwrites the values of the Safety 'sphere' (correction)
2048 }
2049
2050 // Remember last safety origin & value
2051 //
2052 // We overwrite the Safety 'sphere' - keeping old behaviour
2053 fPreviousSftOrigin = pGlobalpoint;
2054 fPreviousSafety = newSafety;
2055 }
2056
2057#ifdef G4DEBUG_NAVIGATION
2058 if( fVerbose > 1 )
2059 {
2060 G4cout << " ---- Exiting ComputeSafety " << G4endl;
2061 if( fVerbose > 2 ) { PrintState(); }
2062 G4cout << " Returned value of Safety = " << newSafety << G4endl;
2063 }
2064 G4cout.precision(oldcoutPrec);
2065#endif
2066
2067 return newSafety;
2068}
#define fStepEndPoint
#define fHistory
#define fPreviousSftOrigin
#define fEnteredDaughter
#define fExitedMother
#define fPreviousSafety
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ReplicaNavigation freplicaNav
G4ThreeVector ComputeLocalPoint(const G4ThreeVector &rGlobPoint) const
G4VoxelNavigation fvoxelNav
G4NormalNavigation fnormalNav
G4int GetDaughtersRegularStructureId(const G4LogicalVolume *pLog) const
EVolume CharacteriseDaughters(const G4LogicalVolume *pLog) const
G4ParameterisedNavigation fparamNav
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
void PrintState() const
G4SmartVoxelHeader * GetVoxelHeader() const
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)
G4double ComputeSafety(const G4ThreeVector &localPoint, const G4NavigationHistory &history, const G4double pProposedMaxLength=DBL_MAX)
G4double ComputeSafety(const G4ThreeVector &globalPoint, const G4ThreeVector &localPoint, G4NavigationHistory &history, const G4double pProposedMaxLength=DBL_MAX)
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4NavigationHistory &history, const G4double pMaxLength=DBL_MAX)
G4double ComputeSafety(const G4ThreeVector &localPoint, const G4VPhysicalVolume &currentPhysical, G4double maxLength=DBL_MAX)
@ kNormal
Definition: geomdefs.hh:84
@ kParameterised
Definition: geomdefs.hh:86
@ kExternal
Definition: geomdefs.hh:87
@ kReplica
Definition: geomdefs.hh:85

Referenced by SetGeometricallyLimitedStep().

◆ ComputeStep()

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

Definition at line 916 of file G4ITNavigator2.cc.

920{
922 G4ThreeVector localDirection = ComputeLocalAxis(pDirection);
923 G4double Step = kInfinity;
924 G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
925 G4LogicalVolume *motherLogical = motherPhysical->GetLogicalVolume();
926
927 // All state relating to exiting normals must be reset
928 //
930 // Reset value - to erase its memory
932 // Reset - used for local exit normal
934 fCalculatedExitNormal = false;
935 // Reset for new step
936
938
939#ifdef G4VERBOSE
940 if( fVerbose > 0 )
941 {
942 G4cout << "*** G4ITNavigator2::ComputeStep: ***" << G4endl;
943 G4cout << " Volume = " << motherPhysical->GetName()
944 << " - Proposed step length = " << pCurrentProposedStepLength
945 << G4endl;
946#ifdef G4DEBUG_NAVIGATION
947 if( fVerbose >= 2 )
948 {
949 G4cout << " Called with the arguments: " << G4endl
950 << " Globalpoint = " << std::setw(25) << pGlobalpoint << G4endl
951 << " Direction = " << std::setw(25) << pDirection << G4endl;
952 if( fVerbose >= 4 )
953 {
954 G4cout << " ---- Upon entering : State" << G4endl;
955 PrintState();
956 }
957 }
958#endif
959 }
960#endif
961
962 G4ThreeVector newLocalPoint = ComputeLocalPoint(pGlobalpoint);
963 if( newLocalPoint != fLastLocatedPointLocal )
964 {
965 // Check whether the relocation is within safety
966 //
968 G4double moveLenSq = (newLocalPoint-oldLocalPoint).mag2();
969
970 if ( moveLenSq >= kCarTolerance*kCarTolerance )
971 {
972#ifdef G4VERBOSE
973 ComputeStepLog(pGlobalpoint, moveLenSq);
974#endif
975 // Relocate the point within the same volume
976 //
977 LocateGlobalPointWithinVolume( pGlobalpoint );
978 fLastTriedStepComputation= true; // Ensure that this is set again !!
979 }
980 }
981 if ( fHistory.GetTopVolumeType()!=kReplica )
982 {
983 switch( CharacteriseDaughters(motherLogical) )
984 {
985 case kNormal:
986 if ( motherLogical->GetVoxelHeader() )
987 {
988 LocateGlobalPointWithinVolume(pGlobalpoint);
990 localDirection,
991 pCurrentProposedStepLength,
992 pNewSafety,
993 fHistory,
996 fExiting,
997 fEntering,
1000
1001 }
1002 else
1003 {
1004 if( motherPhysical->GetRegularStructureId() == 0 )
1005 {
1007 localDirection,
1008 pCurrentProposedStepLength,
1009 pNewSafety,
1010 fHistory,
1013 fExiting,
1014 fEntering,
1017 }
1018 else // Regular (non-voxelised) structure
1019 {
1020 LocateGlobalPointAndSetup( pGlobalpoint, &pDirection, true, true );
1021 fLastTriedStepComputation= true; // Ensure that this is set again !!
1022 //
1023 // if physical process limits the step, the voxel will not be the
1024 // one given by ComputeStepSkippingEqualMaterials() and the local
1025 // point will be wrongly calculated.
1026
1027 // There is a problem: when msc limits the step and the point is
1028 // assigned wrongly to phantom in previous step (while it is out
1029 // of the container volume). Then LocateGlobalPointAndSetup() has
1030 // reset the history topvolume to world.
1031 //
1032 if(fHistory.GetTopVolume()->GetRegularStructureId() == 0 )
1033 {
1034 G4Exception("G4ITNavigator2::ComputeStep()",
1035 "GeomNav1001", JustWarning,
1036 "Point is relocated in voxels, while it should be outside!");
1038 localDirection,
1039 pCurrentProposedStepLength,
1040 pNewSafety,
1041 fHistory,
1044 fExiting,
1045 fEntering,
1048 }
1049 else
1050 {
1051 Step = fregularNav.
1052 ComputeStepSkippingEqualMaterials(fLastLocatedPointLocal,
1053 localDirection,
1054 pCurrentProposedStepLength,
1055 pNewSafety,
1056 fHistory,
1059 fExiting,
1060 fEntering,
1063 motherPhysical);
1064 }
1065 }
1066 }
1067 break;
1068 case kParameterised:
1069 if( GetDaughtersRegularStructureId(motherLogical) != 1 )
1070 {
1072 localDirection,
1073 pCurrentProposedStepLength,
1074 pNewSafety,
1075 fHistory,
1078 fExiting,
1079 fEntering,
1082 }
1083 else // Regular structure
1084 {
1086 localDirection,
1087 pCurrentProposedStepLength,
1088 pNewSafety,
1089 fHistory,
1092 fExiting,
1093 fEntering,
1096 }
1097 break;
1098 case kReplica:
1099 G4Exception("G4ITNavigator2::ComputeStep()", "GeomNav0001",
1100 FatalException, "Not applicable for replicated volumes.");
1101 break;
1102 case kExternal:
1103 G4Exception("G4ITNavigator2::ComputeStep()", "GeomNav0001",
1104 FatalException, "Not applicable for external volumes.");
1105 break;
1106 }
1107 }
1108 else
1109 {
1110 // In the case of a replica, it must handle the exiting
1111 // edge/corner problem by itself
1112 //
1113 G4bool exitingReplica = fExitedMother;
1114 G4bool calculatedExitNormal;
1115 Step = freplicaNav.ComputeStep(pGlobalpoint,
1116 pDirection,
1118 localDirection,
1119 pCurrentProposedStepLength,
1120 pNewSafety,
1121 fHistory,
1123 calculatedExitNormal,
1125 exitingReplica,
1126 fEntering,
1129 fExiting= exitingReplica;
1130 fCalculatedExitNormal= calculatedExitNormal;
1131 }
1132
1133
1134// G4cout << " !!!! Step = " << Step << G4endl;
1135
1136 // Remember last safety origin & value.
1137 //
1138 fPreviousSftOrigin = pGlobalpoint;
1139 fPreviousSafety = pNewSafety;
1140
1141 // Count zero steps - one can occur due to changing momentum at a boundary
1142 // - one, two (or a few) can occur at common edges between
1143 // volumes
1144 // - more than two is likely a problem in the geometry
1145 // description or the Navigation
1146
1147 // Rule of thumb: likely at an Edge if two consecutive steps are zero,
1148 // because at least two candidate volumes must have been
1149 // checked
1150 //
1151 fLocatedOnEdge = fLastStepWasZero && (Step==0.0);
1152 fLastStepWasZero = (Step==0.0);
1153 if (fPushed) { fPushed = fLastStepWasZero; }
1154
1155 // Handle large number of consecutive zero steps
1156 //
1157 if ( fLastStepWasZero )
1158 {
1160#ifdef G4DEBUG_NAVIGATION
1161 if( fNumberZeroSteps > 1 )
1162 {
1163 G4cout << "G4ITNavigator2::ComputeStep(): another zero step, # "
1165 << " at " << pGlobalpoint
1166 << " in volume " << motherPhysical->GetName()
1167 << G4endl;
1168 }
1169#endif
1170 if( fNumberZeroSteps > fActionThreshold_NoZeroSteps-1 )
1171 {
1172 // Act to recover this stuck track. Pushing it along direction
1173 //
1174 Step += 100*kCarTolerance;
1175#ifdef G4VERBOSE
1176 if ((!fPushed) && (fWarnPush))
1177 {
1178 std::ostringstream message;
1179 message << "Track stuck or not moving." << G4endl
1180 << " Track stuck, not moving for "
1181 << fNumberZeroSteps << " steps" << G4endl
1182 << " in volume -" << motherPhysical->GetName()
1183 << "- at point " << pGlobalpoint << G4endl
1184 << " direction: " << pDirection << "." << G4endl
1185 << " Potential geometry or navigation problem !"
1186 << G4endl
1187 << " Trying pushing it of " << Step << " mm ...";
1188 G4Exception("G4ITNavigator2::ComputeStep()", "GeomNav1002",
1189 JustWarning, message, "Potential overlap in geometry!");
1190 }
1191#endif
1192 fPushed = true;
1193 }
1194 if( fNumberZeroSteps > fAbandonThreshold_NoZeroSteps-1 )
1195 {
1196 // Must kill this stuck track
1197 //
1198 std::ostringstream message;
1199 message << "Stuck Track: potential geometry or navigation problem."
1200 << G4endl
1201 << " Track stuck, not moving for "
1202 << fNumberZeroSteps << " steps" << G4endl
1203 << " in volume -" << motherPhysical->GetName()
1204 << "- at point " << pGlobalpoint << G4endl
1205 << " direction: " << pDirection << ".";
1206 motherPhysical->CheckOverlaps(5000, false);
1207 G4Exception("G4ITNavigator2::ComputeStep()", "GeomNav0003",
1208 EventMustBeAborted, message);
1209 }
1210 }
1211 else
1212 {
1213 if (!fPushed) fNumberZeroSteps = 0;
1214 }
1215
1216 fEnteredDaughter = fEntering; // I expect to enter a volume in this Step
1218
1219 fStepEndPoint = pGlobalpoint
1220 + std::min(Step,pCurrentProposedStepLength) * pDirection;
1221 fLastStepEndPointLocal = fLastLocatedPointLocal + Step * localDirection;
1222
1223 if( fExiting )
1224 {
1225#ifdef G4DEBUG_NAVIGATION
1226 if( fVerbose > 2 )
1227 {
1228 G4cout << " At G4Nav CompStep End - if(exiting) - fExiting= " << fExiting
1229 << " fValidExitNormal = " << fValidExitNormal << G4endl;
1230 G4cout << " fExitNormal= " << fExitNormal << G4endl;
1231 }
1232#endif
1233
1235 {
1236 if ( fHistory.GetTopVolumeType()!=kReplica )
1237 {
1238 // Convention: fExitNormal is in the 'grand-mother' coordinate system
1239 //
1242 }
1243 else
1244 {
1246 }
1247 }
1248 else
1249 {
1250 // We must calculate the normal anyway (in order to have it if requested)
1251 //
1252 G4ThreeVector finalLocalPoint =
1253 fLastLocatedPointLocal + localDirection*Step;
1254
1255 if ( fHistory.GetTopVolumeType()!=kReplica )
1256 {
1257 // Find normal in the 'mother' coordinate system
1258 //
1259 G4ThreeVector exitNormalMotherFrame=
1260 motherLogical->GetSolid()->SurfaceNormal(finalLocalPoint);
1261
1262 // Transform it to the 'grand-mother' coordinate system
1263 //
1264 const G4RotationMatrix* mRot = motherPhysical->GetRotation();
1265 if( mRot )
1266 {
1268 fGrandMotherExitNormal = (*mRot).inverse() * exitNormalMotherFrame;
1269 }
1270 else
1271 {
1272 fGrandMotherExitNormal = exitNormalMotherFrame;
1273 }
1274
1275 // Do not set fValidExitNormal -- this signifies
1276 // that the solid is convex!
1277 //
1279 }
1280 else
1281 {
1282 fCalculatedExitNormal = false;
1283 //
1284 // Nothing can be done at this stage currently - to solve this
1285 // Replica Navigation must have calculated the normal for this case
1286 // already.
1287 // Cases: mother is not convex, and exit is at previous replica level
1288
1289#ifdef G4DEBUG_NAVIGATION
1291
1292 desc << "Problem in ComputeStep: Replica Navigation did not provide"
1293 << " valid exit Normal. " << G4endl;
1294 desc << " Do not know how calculate it in this case." << G4endl;
1295 desc << " Location = " << finalLocalPoint << G4endl;
1296 desc << " Volume name = " << motherPhysical->GetName()
1297 << " copy/replica No = " << motherPhysical->GetCopyNo() << G4endl;
1298 G4Exception("G4ITNavigator2::ComputeStep()", "GeomNav0003",
1299 JustWarning, desc, "Normal not available for exiting.");
1300#endif
1301 }
1302 }
1303
1304 // Now transform it to the global reference frame !!
1305 //
1307 {
1308 G4int depth= (G4int)fHistory.GetDepth();
1309 if( depth > 0 )
1310 {
1311 G4AffineTransform GrandMotherToGlobalTransf =
1312 fHistory.GetTransform(depth-1).Inverse();
1314 GrandMotherToGlobalTransf.TransformAxis( fGrandMotherExitNormal );
1315 }
1316 else
1317 {
1319 }
1320 }
1321 else
1322 {
1324 }
1325 }
1326
1327 if( (Step == pCurrentProposedStepLength) && (!fExiting) && (!fEntering) )
1328 {
1329 // This if Step is not really limited by the geometry.
1330 // The Navigator is obliged to return "infinity"
1331 //
1332 Step = kInfinity;
1333 }
1334
1335#ifdef G4VERBOSE
1336 if( fVerbose > 1 )
1337 {
1338 if( fVerbose >= 4 )
1339 {
1340 G4cout << " ----- Upon exiting :" << G4endl;
1341 PrintState();
1342 }
1343 G4cout << " Returned step= " << Step;
1344 if( fVerbose > 5 ) G4cout << G4endl;
1345 if( Step == kInfinity )
1346 {
1347 G4cout << " Requested step= " << pCurrentProposedStepLength ;
1348 if( fVerbose > 5) G4cout << G4endl;
1349 }
1350 G4cout << " Safety = " << pNewSafety << G4endl;
1351 }
1352#endif
1353
1354 return Step;
1355}
@ JustWarning
@ EventMustBeAborted
#define fLastStepEndPointLocal
#define fExiting
#define fExitNormalGlobalFrame
#define fBlockedReplicaNo
#define fLocatedOnEdge
#define fGrandMotherExitNormal
#define fLastLocatedPointLocal
#define fValidExitNormal
#define fChangedGrandMotherRefFrame
#define fBlockedPhysicalVolume
#define fCalculatedExitNormal
#define fNumberZeroSteps
#define fPushed
#define fExitNormal
#define fLastStepWasZero
#define fEntering
#define fLastTriedStepComputation
CLHEP::Hep3Vector G4ThreeVector
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
virtual G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=0, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
G4ThreeVector ComputeLocalAxis(const G4ThreeVector &pVec) const
G4VSolid * GetSolid() const
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)
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 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 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)
const G4RotationMatrix * GetRotation() const
virtual G4bool CheckOverlaps(G4int res=1000, G4double tol=0., G4bool verbose=true, G4int errMax=1)
virtual G4int GetCopyNo() const =0
virtual G4int GetRegularStructureId() const =0
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
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)

Referenced by CheckNextStep(), and ResetFromSnapshot().

◆ CreateGRSSolid()

G4GRSSolid * G4ITNavigator2::CreateGRSSolid ( ) const
inline

◆ CreateGRSVolume()

G4GRSVolume * G4ITNavigator2::CreateGRSVolume ( ) const
inline

◆ CreateTouchableHistory() [1/2]

G4TouchableHistory * G4ITNavigator2::CreateTouchableHistory ( ) const
inline

◆ CreateTouchableHistory() [2/2]

G4TouchableHistory * G4ITNavigator2::CreateTouchableHistory ( const G4NavigationHistory ) const
inline

◆ CreateTouchableHistoryHandle()

G4TouchableHistoryHandle G4ITNavigator2::CreateTouchableHistoryHandle ( ) const
virtual

Definition at line 2254 of file G4ITNavigator2.cc.

2255{
2258}
G4ReferenceCountedHandle< G4TouchableHistory > G4TouchableHistoryHandle
G4TouchableHistory * CreateTouchableHistory() const

Referenced by CreateTouchableHistory().

◆ EnableBestSafety()

void G4ITNavigator2::EnableBestSafety ( G4bool  value = false)
inline

◆ EnteredDaughterVolume()

G4bool G4ITNavigator2::EnteredDaughterVolume ( ) const
inline

Referenced by GetLocalExitNormal().

◆ ExitedMotherVolume()

G4bool G4ITNavigator2::ExitedMotherVolume ( ) const
inline

◆ GetCurrentLocalCoordinate()

G4ThreeVector G4ITNavigator2::GetCurrentLocalCoordinate ( ) const
inline

◆ GetDaughtersRegularStructureId()

G4int G4ITNavigator2::GetDaughtersRegularStructureId ( const G4LogicalVolume pLog) const
inlineprotected

◆ GetGlobalExitNormal()

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

Definition at line 1759 of file G4ITNavigator2.cc.

1761{
1763 G4bool validNormal;
1764 G4ThreeVector localNormal, globalNormal;
1765 G4bool usingStored;
1766
1767 usingStored=
1769 ( ( fLastTriedStepComputation && fExiting) // Just calculated it
1770 || // No locate in between
1772 && (IntersectPointGlobal-fStepEndPoint).mag2()
1774 ) // Calculated it 'just' before & then called locate
1775 // but it did not move position
1776 );
1777
1778 if( usingStored )
1779 {
1780 // This was computed in ComputeStep -- and only on arrival at boundary
1781 //
1782 globalNormal = fExitNormalGlobalFrame;
1783 G4double normMag2 = globalNormal.mag2();
1784 if( std::fabs ( normMag2 - 1.0 ) < perMillion ) // Value is good
1785 {
1786 *pNormalCalculated = true; // ComputeStep always computes it if Exiting
1787 // (fExiting==true)
1788 }
1789 else
1790 {
1791 G4ExceptionDescription message;
1792
1793 message << " ERROR> Expected normal-global-frame to valid (unit vector) "
1794 << " - but |normal| = " << std::sqrt(normMag2)
1795 << " - and |normal|^ = " << normMag2
1796 << " which differs from 1.0 by " << normMag2 - 1.0 << G4endl
1797 << " n = " << fExitNormalGlobalFrame << G4endl;
1798 message << "============================================================"
1799 << G4endl;
1800 G4int oldVerbose = fVerbose;
1801 fVerbose=4;
1802 message << " State of Navigator: " << G4endl;
1803 message << *this << G4endl;
1804 fVerbose = oldVerbose;
1805 message << "============================================================"
1806 << G4endl;
1807
1808 G4Exception("G4ITNavigator2::GetGlobalExitNormal()",
1809 "GeomNav0003",JustWarning, message,
1810 "Value obtained from stored global-normal is not a unit vector.");
1811
1812 // (Re)Compute it now -- as either it was not computed, or it is wrong.
1813 //
1814 localNormal = GetLocalExitNormalAndCheck(IntersectPointGlobal,
1815 &validNormal);
1816 *pNormalCalculated = fCalculatedExitNormal;
1817
1819 globalNormal = localToGlobal.TransformAxis( localNormal );
1820 }
1821 }
1822 else
1823 {
1824 localNormal = GetLocalExitNormalAndCheck(IntersectPointGlobal,&validNormal);
1825 *pNormalCalculated = fCalculatedExitNormal;
1826
1827#ifdef G4DEBUG_NAVIGATION
1828 usingStored= false;
1829
1830 if( (!validNormal) && !fCalculatedExitNormal)
1831 {
1833 edN << " Calculated = " << fCalculatedExitNormal << G4endl;
1834 edN << " Entering= " << fEntering << G4endl;
1835 G4int oldVerbose= this->GetVerboseLevel();
1836 this->SetVerboseLevel(4);
1837 edN << " State of Navigator: " << G4endl;
1838 edN << *this << G4endl;
1839 this->SetVerboseLevel( oldVerbose );
1840
1841 G4Exception("G4ITNavigator2::GetGlobalExitNormal()",
1842 "GeomNav0003", JustWarning, edN,
1843 "LocalExitNormalAndCheck() did not calculate Normal.");
1844 }
1845#endif
1846
1847 G4double localMag2= localNormal.mag2();
1848 if( validNormal && (std::fabs(localMag2-1.0)) > CLHEP::perMillion )
1849 {
1851
1852 edN << "G4ITNavigator2::GetGlobalExitNormal: "
1853 << " Using Local Normal - from call to GetLocalExitNormalAndCheck. "
1854 << G4endl
1855 << " Local Exit Normal : " << " || = " << std::sqrt(localMag2)
1856 << " vec = " << localNormal << G4endl
1857 << " Global Exit Normal : " << " || = " << globalNormal.mag()
1858 << " vec = " << globalNormal << G4endl;
1859 edN << " Calculated It = " << fCalculatedExitNormal << G4endl;
1860
1861 G4Exception("G4ITNavigator2::GetGlobalExitNormal()",
1862 "GeomNav0003",JustWarning, edN,
1863 "Value obtained from new local *solid* is incorrect.");
1864 localNormal = localNormal.unit(); // Should we correct it ??
1865 }
1867 globalNormal = localToGlobal.TransformAxis( localNormal );
1868 }
1869
1870#ifdef G4DEBUG_NAVIGATION
1871 if( usingStored )
1872 {
1873 G4ThreeVector globalNormAgn;
1874
1875 localNormal= GetLocalExitNormalAndCheck(IntersectPointGlobal, &validNormal);
1876
1878 globalNormAgn = localToGlobal.TransformAxis( localNormal );
1879
1880 // Check the value computed against fExitNormalGlobalFrame
1881 G4ThreeVector diffNorm = globalNormAgn - fExitNormalGlobalFrame;
1882 if( diffNorm.mag2() > perMillion*CLHEP::perMillion)
1883 {
1885 edDfn << "Found difference in normals in case of exiting mother "
1886 << "- when Get is called after ComputingStep " << G4endl;
1887 edDfn << " Magnitude of diff = " << diffNorm.mag() << G4endl;
1888 edDfn << " Normal stored (Global) = " << fExitNormalGlobalFrame
1889 << G4endl;
1890 edDfn << " Global Computed from Local = " << globalNormAgn << G4endl;
1891 G4Exception("G4ITNavigator::GetGlobalExitNormal()", "GeomNav0003",
1892 JustWarning, edDfn);
1893 }
1894 }
1895#endif
1896
1897 return globalNormal;
1898}
Hep3Vector unit() const
double mag2() const
double mag() const
void SetVerboseLevel(G4int level)
virtual G4ThreeVector GetLocalExitNormalAndCheck(const G4ThreeVector &point, G4bool *valid)
const G4AffineTransform GetLocalToGlobalTransform() const
G4int GetVerboseLevel() const

Referenced by CreateTouchableHistory().

◆ GetGlobalToLocalTransform()

const G4AffineTransform & G4ITNavigator2::GetGlobalToLocalTransform ( ) const
inline

◆ GetLocalExitNormal()

G4ThreeVector G4ITNavigator2::GetLocalExitNormal ( G4bool valid)
virtual

Definition at line 1509 of file G4ITNavigator2.cc.

1510{
1512 G4ThreeVector ExitNormal(0.,0.,0.);
1513 G4VSolid *currentSolid=0;
1514 G4LogicalVolume *candidateLogical;
1516 {
1517 // use fLastLocatedPointLocal and next candidate volume
1518 //
1519 G4ThreeVector nextSolidExitNormal(0.,0.,0.);
1520
1521 if( fEntering && (fBlockedPhysicalVolume!=0) )
1522 {
1523 candidateLogical= fBlockedPhysicalVolume->GetLogicalVolume();
1524 if( candidateLogical )
1525 {
1526 // fLastStepEndPointLocal is in the coordinates of the mother
1527 // we need it in the daughter's coordinate system.
1528
1529 // The following code should also work in case of Replica
1530 {
1531 // First transform fLastLocatedPointLocal to the new daughter
1532 // coordinates
1533 //
1534 G4AffineTransform MotherToDaughterTransform=
1538 G4ThreeVector daughterPointOwnLocal=
1539 MotherToDaughterTransform.TransformPoint( fLastStepEndPointLocal );
1540
1541 // OK if it is a parameterised volume
1542 //
1543 EInside inSideIt;
1544 G4bool onSurface;
1545 G4double safety= -1.0;
1546 currentSolid= candidateLogical->GetSolid();
1547 inSideIt = currentSolid->Inside(daughterPointOwnLocal);
1548 onSurface = (inSideIt == kSurface);
1549 if( ! onSurface )
1550 {
1551 if( inSideIt == kOutside )
1552 {
1553 safety = (currentSolid->DistanceToIn(daughterPointOwnLocal));
1554 onSurface = safety < 100.0 * kCarTolerance;
1555 }
1556 else if (inSideIt == kInside )
1557 {
1558 safety = (currentSolid->DistanceToOut(daughterPointOwnLocal));
1559 onSurface = safety < 100.0 * kCarTolerance;
1560 }
1561 }
1562
1563 if( onSurface )
1564 {
1565 nextSolidExitNormal =
1566 currentSolid->SurfaceNormal(daughterPointOwnLocal);
1567
1568 // Entering the solid ==> opposite
1569 //
1570 ExitNormal = -nextSolidExitNormal;
1572 }
1573 else
1574 {
1575#ifdef G4VERBOSE
1576 if(( fVerbose == 1 ) && ( fCheck ))
1577 {
1578 std::ostringstream message;
1579 message << "Point not on surface ! " << G4endl
1580 << " Point = "
1581 << daughterPointOwnLocal << G4endl
1582 << " Physical volume = "
1583 << fBlockedPhysicalVolume->GetName() << G4endl
1584 << " Logical volume = "
1585 << candidateLogical->GetName() << G4endl
1586 << " Solid = " << currentSolid->GetName()
1587 << " Type = "
1588 << currentSolid->GetEntityType() << G4endl
1589 << *currentSolid << G4endl;
1590 if( inSideIt == kOutside )
1591 {
1592 message << "Point is Outside. " << G4endl
1593 << " Safety (from outside) = " << safety << G4endl;
1594 }
1595 else // if( inSideIt == kInside )
1596 {
1597 message << "Point is Inside. " << G4endl
1598 << " Safety (from inside) = " << safety << G4endl;
1599 }
1600 G4Exception("G4ITNavigator2::GetLocalExitNormal()", "GeomNav1001",
1601 JustWarning, message);
1602 }
1603#endif
1604 }
1605 *valid = onSurface; // was =true;
1606 }
1607 }
1608 }
1609 else if ( fExiting )
1610 {
1611 ExitNormal = fGrandMotherExitNormal;
1612 *valid = true;
1613 fCalculatedExitNormal= true; // Should be true already
1614 }
1615 else // i.e. ( fBlockedPhysicalVolume == 0 )
1616 {
1617 *valid = false;
1618 G4Exception("G4ITNavigator2::GetLocalExitNormal()",
1619 "GeomNav0003", JustWarning,
1620 "Incorrect call to GetLocalSurfaceNormal." );
1621 }
1622 }
1623 else // ( ! fLastTriedStepComputation ) ie. last call was to Locate
1624 {
1625 if ( EnteredDaughterVolume() )
1626 {
1627 G4VSolid* daughterSolid =fHistory.GetTopVolume()->GetLogicalVolume()
1628 ->GetSolid();
1629 ExitNormal= -(daughterSolid->SurfaceNormal(fLastLocatedPointLocal));
1630 if( std::fabs(ExitNormal.mag2()-1.0 ) > CLHEP::perMillion )
1631 {
1633 desc << " Parameters of solid: " << *daughterSolid
1634 << " Point for surface = " << fLastLocatedPointLocal << std::endl;
1635 G4Exception("G4ITNavigator2::GetLocalExitNormal()",
1636 "GeomNav0003", FatalException, desc,
1637 "Surface Normal returned by Solid is not a Unit Vector." );
1638 }
1640 *valid = true;
1641 }
1642 else
1643 {
1644 if( fExitedMother )
1645 {
1646 ExitNormal = fGrandMotherExitNormal;
1647 *valid = true;
1649 }
1650 else // We are not at a boundary. ExitNormal remains (0,0,0)
1651 {
1652 *valid = false;
1653 fCalculatedExitNormal= false;
1654 G4ExceptionDescription message;
1655 message << "Function called when *NOT* at a Boundary." << G4endl;
1656 G4Exception("G4ITNavigator2::GetLocalExitNormal()",
1657 "GeomNav0003", JustWarning, message);
1658 }
1659 }
1660 }
1661 return ExitNormal;
1662}
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4AffineTransform GetMotherToDaughterTransform(G4VPhysicalVolume *dVolume, G4int dReplicaNo, EVolume dVolumeType)
EVolume VolumeType(const G4VPhysicalVolume *pVol) const
G4bool EnteredDaughterVolume() const
const G4String & GetName() const
G4String GetName() const
virtual EInside Inside(const G4ThreeVector &p) const =0
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
virtual G4GeometryType GetEntityType() const =0
EInside
Definition: geomdefs.hh:67
@ kInside
Definition: geomdefs.hh:70
@ kOutside
Definition: geomdefs.hh:68
@ kSurface
Definition: geomdefs.hh:69

Referenced by CreateTouchableHistory(), and GetLocalExitNormalAndCheck().

◆ GetLocalExitNormalAndCheck()

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

Definition at line 1722 of file G4ITNavigator2.cc.

1730{
1732#ifdef G4DEBUG_NAVIGATION
1733 // Check Current point against expected 'local' value
1734 //
1736 {
1737 G4ThreeVector ExpectedBoundaryPointLocal;
1738
1739 const G4AffineTransform& GlobalToLocal= GetGlobalToLocalTransform();
1740 ExpectedBoundaryPointLocal =
1741 GlobalToLocal.TransformPoint( ExpectedBoundaryPointGlobal );
1742
1743 // Add here: Comparison against expected position,
1744 // i.e. the endpoint of ComputeStep
1745 }
1746#endif
1747
1748 return GetLocalExitNormal( pValid);
1749}
const G4AffineTransform & GetGlobalToLocalTransform() const
virtual G4ThreeVector GetLocalExitNormal(G4bool *valid)

Referenced by CreateTouchableHistory(), and GetGlobalExitNormal().

◆ GetLocalToGlobalTransform()

const G4AffineTransform G4ITNavigator2::GetLocalToGlobalTransform ( ) const
inline

◆ GetMotherToDaughterTransform()

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

Definition at line 1671 of file G4ITNavigator2.cc.

1674{
1676 switch (enteringVolumeType)
1677 {
1678 case kNormal: // Nothing is needed to prepare the transformation
1679 break; // It is stored already in the physical volume (placement)
1680 case kReplica: // Sets the transform in the Replica - tbc
1681 G4Exception("G4ITNavigator2::GetMotherToDaughterTransform()",
1682 "GeomNav0001", FatalException,
1683 "Method NOT Implemented yet for replica volumes.");
1684 break;
1685 case kParameterised:
1686 if( pEnteringPhysVol->GetRegularStructureId() == 0 )
1687 {
1688 G4VPVParameterisation *pParam =
1689 pEnteringPhysVol->GetParameterisation();
1690 G4VSolid* pSolid =
1691 pParam->ComputeSolid(enteringReplicaNo, pEnteringPhysVol);
1692 pSolid->ComputeDimensions(pParam, enteringReplicaNo, pEnteringPhysVol);
1693
1694 // Sets the transform in the Parameterisation
1695 //
1696 pParam->ComputeTransformation(enteringReplicaNo, pEnteringPhysVol);
1697
1698 // Set the correct solid and material in Logical Volume
1699 //
1700 G4LogicalVolume* pLogical = pEnteringPhysVol->GetLogicalVolume();
1701 pLogical->SetSolid( pSolid );
1702 }
1703 break;
1704 case kExternal:
1705 G4Exception("G4ITNavigator2::GetMotherToDaughterTransform()",
1706 "GeomNav0001", FatalException,
1707 "Not applicable for external volumes.");
1708 break;
1709 }
1710 return G4AffineTransform(pEnteringPhysVol->GetRotation(),
1711 pEnteringPhysVol->GetTranslation()).Invert();
1712}
G4AffineTransform & Invert()
void SetSolid(G4VSolid *pSolid)
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
virtual void ComputeTransformation(const G4int, G4VPhysicalVolume *) const =0
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:137

Referenced by GetLocalExitNormal(), and GetLocalToGlobalTransform().

◆ GetNavigatorState()

G4ITNavigatorState_Lock2 * G4ITNavigator2::GetNavigatorState ( )

Definition at line 715 of file G4ITNavigator2.cc.

716{
717 return fpNavigatorState;
718}

◆ GetRandomInCurrentVolume()

void G4ITNavigator2::GetRandomInCurrentVolume ( G4ThreeVector rndmPoint) const

Definition at line 2522 of file G4ITNavigator2.cc.

2523{
2524 const G4AffineTransform& local2Global = GetLocalToGlobalTransform();
2525 G4VSolid* solid = fHistory.GetTopVolume()->GetLogicalVolume()->GetSolid();
2526
2527 G4VoxelLimits voxelLimits; // Defaults to "infinite" limits.
2528 G4double vmin, vmax;
2529 G4AffineTransform dummy;
2530 std::vector<std::vector<G4double> > fExtend;
2531
2532 solid->CalculateExtent(kXAxis,voxelLimits,dummy,vmin,vmax);
2533 fExtend[kXAxis][BoundingBox::kMin] = vmin;
2534 fExtend[kXAxis][BoundingBox::kMax] = vmax;
2535
2536 solid->CalculateExtent(kYAxis,voxelLimits,dummy,vmin,vmax);
2537 fExtend[kYAxis][BoundingBox::kMin] = vmin;
2538 fExtend[kYAxis][BoundingBox::kMax] = vmax;
2539
2540 solid->CalculateExtent(kZAxis,voxelLimits,dummy,vmin,vmax);
2541 fExtend[kZAxis][BoundingBox::kMin] = vmin;
2542 fExtend[kZAxis][BoundingBox::kMax] = vmax;
2543
2544 G4ThreeVector rndmPos;
2545
2546 while(1)
2547 {
2548 for(G4int i = 0 ; i < 3 ; ++i)
2549 {
2550 G4double min = fExtend[i][BoundingBox::kMin];
2551 G4double max = fExtend[i][BoundingBox::kMax];
2552 rndmPos[i] = G4UniformRand()*(max-min)+min;
2553 }
2554
2555 if(solid->Inside(rndmPos) == kInside) break;
2556 }
2557
2558 _rndmPoint = local2Global.TransformPoint(rndmPos);
2559}
#define G4UniformRand()
Definition: Randomize.hh:52
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
@ kYAxis
Definition: geomdefs.hh:56
@ kXAxis
Definition: geomdefs.hh:55
@ kZAxis
Definition: geomdefs.hh:57
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

Referenced by ResetFromSnapshot().

◆ GetSnapshotOfState()

std::shared_ptr< G4ITNavigatorState_Lock2 > G4ITNavigator2::GetSnapshotOfState ( )
inline

◆ GetVerboseLevel()

G4int G4ITNavigator2::GetVerboseLevel ( ) const
inline

Referenced by GetGlobalExitNormal().

◆ GetWorldVolume()

G4VPhysicalVolume * G4ITNavigator2::GetWorldVolume ( ) const
inline

◆ InsideCurrentVolume()

EInside G4ITNavigator2::InsideCurrentVolume ( const G4ThreeVector globalPoint) const

Definition at line 2510 of file G4ITNavigator2.cc.

2511{
2512 const G4AffineTransform& transform = GetGlobalToLocalTransform();
2513 G4ThreeVector localPoint(transform.TransformPoint(globalPoint));
2514
2515 G4VSolid* solid = fHistory.GetTopVolume()->GetLogicalVolume()->GetSolid();
2516
2517 return solid->Inside(localPoint);
2518}

Referenced by ResetFromSnapshot().

◆ IsActive()

G4bool G4ITNavigator2::IsActive ( ) const
inline

◆ IsCheckModeActive()

G4bool G4ITNavigator2::IsCheckModeActive ( ) const
inline

◆ LocateGlobalPointAndSetup()

G4VPhysicalVolume * G4ITNavigator2::LocateGlobalPointAndSetup ( const G4ThreeVector point,
const G4ThreeVector direction = 0,
const G4bool  pRelativeSearch = true,
const G4bool  ignoreDirection = true 
)
virtual

Definition at line 158 of file G4ITNavigator2.cc.

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

Referenced by ComputeStep(), NewNavigatorStateAndLocate(), ResetFromSnapshot(), and ResetHierarchyAndLocate().

◆ LocateGlobalPointAndUpdateTouchable() [1/2]

void G4ITNavigator2::LocateGlobalPointAndUpdateTouchable ( const G4ThreeVector position,
const G4ThreeVector direction,
G4VTouchable touchableToUpdate,
const G4bool  RelativeSearch = true 
)
inline

◆ LocateGlobalPointAndUpdateTouchable() [2/2]

void G4ITNavigator2::LocateGlobalPointAndUpdateTouchable ( const G4ThreeVector position,
G4VTouchable touchableToUpdate,
const G4bool  RelativeSearch = true 
)
inline

◆ LocateGlobalPointAndUpdateTouchableHandle()

void G4ITNavigator2::LocateGlobalPointAndUpdateTouchableHandle ( const G4ThreeVector position,
const G4ThreeVector direction,
G4TouchableHandle oldTouchableToUpdate,
const G4bool  RelativeSearch = true 
)
inline

◆ LocateGlobalPointWithinVolume()

void G4ITNavigator2::LocateGlobalPointWithinVolume ( const G4ThreeVector position)
virtual

Definition at line 622 of file G4ITNavigator2.cc.

623{
625
626#ifdef G4DEBUG_NAVIGATION
627 // Check: Either step was not limited by a boundary
628 // or else the full step is no longer being taken
629 assert( !fWasLimitedByGeometry );
630#endif
631
634 fChangedGrandMotherRefFrame= false; // Frame for Exit Normal
635
636#ifdef G4DEBUG_NAVIGATION
637 if( fVerbose > 2 )
638 {
639 G4cout << "Entering LocateGlobalWithinVolume(): History = " << G4endl;
640 G4cout << fHistory << G4endl;
641 }
642#endif
643
644 // For the case of Voxel (or Parameterised) volume the respective
645 // Navigator must be messaged to update its voxel information etc
646
647 // Update the state of the Sub Navigators
648 // - in particular any voxel information they store/cache
649 //
650 G4VPhysicalVolume* motherPhysical = fHistory.GetTopVolume();
651 G4LogicalVolume* motherLogical = motherPhysical->GetLogicalVolume();
652 G4SmartVoxelHeader* pVoxelHeader = motherLogical->GetVoxelHeader();
653
654 if ( fHistory.GetTopVolumeType()!=kReplica )
655 {
656 switch( CharacteriseDaughters(motherLogical) )
657 {
658 case kNormal:
659 if ( pVoxelHeader )
660 {
662 }
663 break;
664 case kParameterised:
665 if( GetDaughtersRegularStructureId(motherLogical) != 1 )
666 {
667 // Resets state & returns voxel node
668 //
670 }
671 break;
672 case kReplica:
673 G4Exception("G4ITNavigator2::LocateGlobalPointWithinVolume()",
674 "GeomNav0001", FatalException,
675 "Not applicable for replicated volumes.");
676 break;
677 case kExternal:
678 G4Exception("G4ITNavigator2::LocateGlobalPointWithinVolume()",
679 "GeomNav0001", FatalException,
680 "Not applicable for external volumes.");
681 break;
682 }
683 }
684
685 // Reset the state variables
686 // - which would have been affected
687 // by the 'equivalent' call to LocateGlobalPointAndSetup
688 // - who's values have been invalidated by the 'move'.
689 //
692 fEntering = false;
693 fEnteredDaughter = false; // Boundary not encountered, did not enter
694 fExiting = false;
695 fExitedMother = false; // Boundary not encountered, did not exit
696}
G4SmartVoxelNode * ParamVoxelLocate(G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint)
G4SmartVoxelNode * VoxelLocate(G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint)

Referenced by ComputeSafety(), ComputeStep(), and ResetFromSnapshot().

◆ NetRotation()

G4RotationMatrix G4ITNavigator2::NetRotation ( ) const
inline

◆ NetTranslation()

G4ThreeVector G4ITNavigator2::NetTranslation ( ) const
inline

◆ NewNavigatorState() [1/2]

void G4ITNavigator2::NewNavigatorState ( )

Definition at line 733 of file G4ITNavigator2.cc.

734{
735 fpNavigatorState = new G4NavigatorState();
736 if(fTopPhysical == 0)
737 {
738 G4ExceptionDescription exceptionDescription;
739 exceptionDescription << "No World Volume";
740
741 G4Exception("G4ITNavigator::NewNavigatorState",
742 "NoWorldVolume",FatalException,exceptionDescription);
743 return;
744 }
745
746 fHistory.SetFirstEntry(fTopPhysical );
748}
virtual void SetupHierarchy()

◆ NewNavigatorState() [2/2]

void G4ITNavigator2::NewNavigatorState ( const G4TouchableHistory h)

Definition at line 750 of file G4ITNavigator2.cc.

751{
752 fpNavigatorState = new G4NavigatorState();
753 if(fTopPhysical == 0)
754 {
755 G4ExceptionDescription exceptionDescription;
756 exceptionDescription << "No World Volume";
757
758 G4Exception("G4ITNavigator::NewNavigatorState",
759 "NoWorldVolume",FatalException,exceptionDescription);
760 return;
761 }
762
763 fHistory = *h.GetHistory();
764 fLastTriedStepComputation= false; // Redundant, but best
766}
const G4NavigationHistory * GetHistory() const

◆ NewNavigatorStateAndLocate()

G4VPhysicalVolume * G4ITNavigator2::NewNavigatorStateAndLocate ( const G4ThreeVector p,
const G4ThreeVector direction 
)

Definition at line 768 of file G4ITNavigator2.cc.

770{
771 fpNavigatorState = new G4NavigatorState();
772
773 if(fTopPhysical == 0)
774 {
775 G4ExceptionDescription exceptionDescription;
776 exceptionDescription << "No World Volume";
777
778 G4Exception("G4ITNavigator::NewNavigatorStateAndLocate",
779 "NoWorldVolume",FatalException,exceptionDescription);
780 return 0;
781 }
782
783 fHistory.SetFirstEntry(fTopPhysical );
785 return LocateGlobalPointAndSetup(p, &direction, false, false);
786}

◆ PrintState()

void G4ITNavigator2::PrintState ( ) const

Definition at line 2264 of file G4ITNavigator2.cc.

2265{
2267 G4long oldcoutPrec = G4cout.precision(4);
2268 if( fVerbose >= 4 )
2269 {
2270 G4cout << "The current state of G4Navigator is: " << G4endl;
2271 G4cout << " ValidExitNormal= " << fValidExitNormal // << G4endl
2272 << " ExitNormal = " << fExitNormal // << G4endl
2273 << " Exiting = " << fExiting // << G4endl
2274 << " Entering = " << fEntering // << G4endl
2275 << " BlockedPhysicalVolume= " ;
2277 G4cout << "None";
2278 else
2279 G4cout << fBlockedPhysicalVolume->GetName();
2280 G4cout << G4endl
2281 << " BlockedReplicaNo = " << fBlockedReplicaNo // << G4endl
2282 << " LastStepWasZero = " << fLastStepWasZero // << G4endl
2283 << G4endl;
2284 }
2285 if( ( 1 < fVerbose) && (fVerbose < 4) )
2286 {
2287 G4cout << G4endl; // Make sure to line up
2288 G4cout << std::setw(30) << " ExitNormal " << " "
2289 << std::setw( 5) << " Valid " << " "
2290 << std::setw( 9) << " Exiting " << " "
2291 << std::setw( 9) << " Entering" << " "
2292 << std::setw(15) << " Blocked:Volume " << " "
2293 << std::setw( 9) << " ReplicaNo" << " "
2294 << std::setw( 8) << " LastStepZero " << " "
2295 << G4endl;
2296 G4cout << "( " << std::setw(7) << fExitNormal.x()
2297 << ", " << std::setw(7) << fExitNormal.y()
2298 << ", " << std::setw(7) << fExitNormal.z() << " ) "
2299 << std::setw( 5) << fValidExitNormal << " "
2300 << std::setw( 9) << fExiting << " "
2301 << std::setw( 9) << fEntering << " ";
2302 if ( fBlockedPhysicalVolume==0 )
2303 {
2304 G4cout << std::setw(15) << "None";
2305 }
2306 else
2307 {
2308 G4cout << std::setw(15)<< fBlockedPhysicalVolume->GetName();
2309 }
2310 G4cout << std::setw( 9) << fBlockedReplicaNo << " "
2311 << std::setw( 8) << fLastStepWasZero << " "
2312 << G4endl;
2313 }
2314 if( fVerbose > 2 )
2315 {
2316 G4cout.precision(8);
2317 G4cout << " Current Localpoint = " << fLastLocatedPointLocal << G4endl;
2318 G4cout << " PreviousSftOrigin = " << fPreviousSftOrigin << G4endl;
2319 G4cout << " PreviousSafety = " << fPreviousSafety << G4endl;
2320 }
2321 G4cout.precision(oldcoutPrec);
2322}

Referenced by ComputeSafety(), ComputeStep(), LocateGlobalPointAndSetup(), and SetPushVerbosity().

◆ RecheckDistanceToCurrentBoundary()

G4bool G4ITNavigator2::RecheckDistanceToCurrentBoundary ( const G4ThreeVector pGlobalPoint,
const G4ThreeVector pDirection,
const G4double  CurrentProposedStepLength,
G4double prDistance,
G4double prNewSafety = 0 
) const
virtual

Definition at line 2088 of file G4ITNavigator2.cc.

2094{
2095 G4ThreeVector localPosition = ComputeLocalPoint(aDisplacedGlobalPoint);
2096 G4ThreeVector localDirection = ComputeLocalAxis(aNewDirection);
2097 // G4double Step = kInfinity;
2098
2099 G4bool validExitNormal;
2100 G4ThreeVector exitNormal;
2101 // Check against mother solid
2102 G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
2103 G4LogicalVolume *motherLogical = motherPhysical->GetLogicalVolume();
2104
2105#ifdef CHECK_ORDER_OF_METHODS
2107 {
2108 G4Exception("G4Navigator::RecheckDistanceToCurrentBoundary()",
2109 "GeomNav0001", FatalException,
2110 "Method must be called after ComputeStep(), before call to LocateMethod.");
2111 }
2112#endif
2113
2114 EInside locatedDaug; // = kUndefined;
2115 G4double daughterStep= DBL_MAX;
2116 G4double daughterSafety= DBL_MAX;
2117
2118 if( fEnteredDaughter )
2119 {
2120 if( motherLogical->CharacteriseDaughters() ==kReplica ) { return false; }
2121
2122 // Track arrived at boundary of a daughter volume at
2123 // the last call of ComputeStep().
2124 // In case the proposed displaced point is inside this daughter,
2125 // it must backtrack at least to the entry point.
2126 // NOTE: No check is made against other daughter volumes. It is
2127 // assumed that the proposed displacement is small enough that
2128 // this is not needed.
2129
2130 // Must check boundary of current daughter
2131 //
2133 G4LogicalVolume *candLogical= candPhysical->GetLogicalVolume();
2134 G4VSolid *candSolid= candLogical->GetSolid();
2135
2136 G4AffineTransform nextLevelTrf(candPhysical->GetRotation(),
2137 candPhysical->GetTranslation());
2138
2139 G4ThreeVector dgPosition= nextLevelTrf.TransformPoint(localPosition);
2140 G4ThreeVector dgDirection= nextLevelTrf.TransformAxis(localDirection);
2141 locatedDaug = candSolid->Inside(dgPosition);
2142
2143 if( locatedDaug == kInside )
2144 {
2145 // Reverse direction - and find first exit. ( Is it valid?)
2146 // Must backtrack
2147 G4double distanceBackOut =
2148 candSolid->DistanceToOut(dgPosition,
2149 - dgDirection, // Reverse direction
2150 true, &validExitNormal, &exitNormal);
2151 daughterStep= - distanceBackOut;
2152 // No check is made whether the particle could have arrived at
2153 // at this location without encountering another volume or
2154 // a different psurface of the current volume
2155 if( prNewSafety )
2156 {
2157 daughterSafety= candSolid->DistanceToOut(dgPosition);
2158 }
2159 }
2160 else
2161 {
2162 if( locatedDaug == kOutside )
2163 {
2164 // See whether it still intersects it
2165 //
2166 daughterStep= candSolid->DistanceToIn(dgPosition,
2167 dgDirection);
2168 if( prNewSafety )
2169 {
2170 daughterSafety= candSolid->DistanceToIn(dgPosition);
2171 }
2172 }
2173 else
2174 {
2175 // The point remains on the surface of candidate solid
2176 //
2177 daughterStep= 0.0;
2178 daughterSafety= 0.0;
2179 }
2180 }
2181
2182 // If trial point is in daughter (or on its surface) we have the
2183 // answer, the rest is not relevant
2184 //
2185 if( locatedDaug != kOutside )
2186 {
2187 *prDistance= daughterStep;
2188 if( prNewSafety ) { *prNewSafety= daughterSafety; }
2189 return true;
2190 }
2191 // If ever extended, so that some type of mother cut daughter,
2192 // this would change
2193 }
2194
2195 G4VSolid *motherSolid= motherLogical->GetSolid();
2196
2197 G4double motherStep= DBL_MAX, motherSafety= DBL_MAX;
2198
2199 // Check distance to boundary of mother
2200 //
2201 if ( fHistory.GetTopVolumeType()!=kReplica )
2202 {
2203 EInside locatedMoth = motherSolid->Inside(localPosition);
2204
2205 if( locatedMoth == kInside )
2206 {
2207 motherSafety= motherSolid->DistanceToOut(localPosition);
2208 if( ProposedMove >= motherSafety )
2209 {
2210 motherStep= motherSolid->DistanceToOut(localPosition,
2211 localDirection,
2212 true, &validExitNormal, &exitNormal);
2213 }
2214 else
2215 {
2216 motherStep= ProposedMove;
2217 }
2218 }
2219 else if( locatedMoth == kOutside)
2220 {
2221 motherSafety= motherSolid->DistanceToIn(localPosition);
2222 if( ProposedMove >= motherSafety )
2223 {
2224 motherStep= - motherSolid->DistanceToIn(localPosition,
2225 -localDirection);
2226 }
2227 }
2228 else
2229 {
2230 motherSafety= 0.0;
2231 *prDistance= 0.0; // On surface - no move // motherStep;
2232 if( prNewSafety ) { *prNewSafety= motherSafety; }
2233 return false;
2234 }
2235 }
2236 else
2237 {
2238 return false;
2239 }
2240
2241 *prDistance= std::min( motherStep, daughterStep );
2242 if( prNewSafety )
2243 {
2244 *prNewSafety= std::min( motherSafety, daughterSafety );
2245 }
2246 return true;
2247}
EVolume CharacteriseDaughters() const
const G4ThreeVector GetTranslation() const
#define DBL_MAX
Definition: templates.hh:62

Referenced by SetGeometricallyLimitedStep().

◆ ResetFromSnapshot()

void G4ITNavigator2::ResetFromSnapshot ( std::shared_ptr< G4ITNavigatorState_Lock2 )
inline

◆ ResetHierarchyAndLocate()

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

Definition at line 130 of file G4ITNavigator2.cc.

133{
134// ResetState();
135 fHistory = *h.GetHistory();
137 fLastTriedStepComputation= false; // Redundant, but best
138 return LocateGlobalPointAndSetup(p, &direction, true, false);
139}

Referenced by ResetFromSnapshot().

◆ ResetNavigatorState()

void G4ITNavigator2::ResetNavigatorState ( )

Definition at line 728 of file G4ITNavigator2.cc.

729{
731}

◆ ResetStackAndState()

void G4ITNavigator2::ResetStackAndState ( )
inline

◆ ResetState()

void G4ITNavigator2::ResetState ( )
protectedvirtual

Definition at line 1401 of file G4ITNavigator2.cc.

1402{
1403 fWasLimitedByGeometry = false;
1404 fEntering = false;
1405 fExiting = false;
1406 fLocatedOnEdge = false;
1407 fLastStepWasZero = false;
1408 fEnteredDaughter = false;
1409 fExitedMother = false;
1410 fPushed = false;
1411
1412 fValidExitNormal = false;
1414 fCalculatedExitNormal = false;
1415
1416 fExitNormal = G4ThreeVector(0,0,0);
1419
1421 fPreviousSafety = 0.0;
1422
1423 fNumberZeroSteps = 0;
1424
1426 fBlockedReplicaNo = -1;
1427
1428 fLastLocatedPointLocal = G4ThreeVector( kInfinity, -kInfinity, 0.0 );
1429 fLocatedOutsideWorld = false;
1430}

Referenced by ComputeLocalAxis().

◆ SetGeometricallyLimitedStep()

void G4ITNavigator2::SetGeometricallyLimitedStep ( )
inline

◆ SetNavigatorState()

void G4ITNavigator2::SetNavigatorState ( G4ITNavigatorState_Lock2 navState)

Definition at line 720 of file G4ITNavigator2.cc.

721{
722 fpNavigatorState = (G4NavigatorState*) navState;
723 if(fpNavigatorState) SetupHierarchy(); // ?
724// fpSaveState = (G4SaveNavigatorState*) navState;
725// if(navState) RestoreSavedState();
726}

◆ SetPushVerbosity()

void G4ITNavigator2::SetPushVerbosity ( G4bool  mode)
inline

◆ SetupHierarchy()

void G4ITNavigator2::SetupHierarchy ( )
protectedvirtual

Definition at line 1440 of file G4ITNavigator2.cc.

1441{
1442 G4int i;
1443 const G4int cdepth = (G4int)fHistory.GetDepth();
1444 G4VPhysicalVolume *current;
1445 G4VSolid *pSolid;
1446 G4VPVParameterisation *pParam;
1447
1448 for ( i=1; i<=cdepth; i++ )
1449 {
1450 current = fHistory.GetVolume(i);
1451 switch ( fHistory.GetVolumeType(i) )
1452 {
1453 case kNormal:
1454 break;
1455 case kReplica:
1456 freplicaNav.ComputeTransformation(fHistory.GetReplicaNo(i), current);
1457 break;
1458 case kParameterised:
1459 {
1460 G4int replicaNo;
1461 pParam = current->GetParameterisation();
1462 replicaNo = fHistory.GetReplicaNo(i);
1463 pSolid = pParam->ComputeSolid(replicaNo, current);
1464
1465 // Set up dimensions & transform in solid/physical volume
1466 //
1467 pSolid->ComputeDimensions(pParam, replicaNo, current);
1468 pParam->ComputeTransformation(replicaNo, current);
1469
1470 G4TouchableHistory *pTouchable= 0;
1471 if( pParam->IsNested() )
1472 {
1473 pTouchable= new G4TouchableHistory( fHistory );
1474 pTouchable->MoveUpHistory(); // Move up to the parent level
1475 // Adequate only if Nested at the Branch level (last)
1476 // To extend to other cases:
1477 // pTouchable->MoveUpHistory(cdepth-i-1);
1478 // Move to the parent level of *Current* level
1479 // Could replace this line and constructor with a revised
1480 // c-tor for History(levels to drop)
1481 }
1482 // Set up the correct solid and material in Logical Volume
1483 //
1484 G4LogicalVolume *pLogical = current->GetLogicalVolume();
1485 pLogical->SetSolid( pSolid );
1486 pLogical->UpdateMaterial( pParam ->
1487 ComputeMaterial(replicaNo, current, pTouchable) );
1488 delete pTouchable;
1489 }
1490 break;
1491
1492 case kExternal:
1493 G4Exception("G4ITNavigator2::SetupHierarchy()",
1494 "GeomNav0001", FatalException,
1495 "Not applicable for external volumes.");
1496 break;
1497
1498 }
1499 }
1500}
G4int MoveUpHistory(G4int num_levels=1)
virtual G4bool IsNested() const
virtual G4VPVParameterisation * GetParameterisation() const =0

Referenced by GetDaughtersRegularStructureId(), NewNavigatorState(), NewNavigatorStateAndLocate(), ResetHierarchyAndLocate(), and SetNavigatorState().

◆ SetVerboseLevel()

void G4ITNavigator2::SetVerboseLevel ( G4int  level)
inline

Referenced by GetGlobalExitNormal().

◆ SetWorldVolume()

void G4ITNavigator2::SetWorldVolume ( G4VPhysicalVolume pWorld)
inline

◆ SeverityOfZeroStepping()

G4int G4ITNavigator2::SeverityOfZeroStepping ( G4int noZeroSteps) const
inline

◆ VolumeType()

EVolume G4ITNavigator2::VolumeType ( const G4VPhysicalVolume pVol) const
inlineprotected

Friends And Related Function Documentation

◆ operator<<

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

Definition at line 2432 of file G4ITNavigator2.cc.

2433{
2434 // Old version did only the following:
2435 // os << "Current History: " << G4endl << n.fHistory;
2436 // Old behaviour is recovered for fVerbose = 0
2437
2438 // Adapted from G4ITNavigator2::PrintState() const
2439
2440 G4long oldcoutPrec = os.precision(4);
2441 if( n.fVerbose >= 4 )
2442 {
2443 os << "The current state of G4ITNavigator2 is: " << G4endl;
2444 os << " ValidExitNormal= " << n.fValidExitNormal << G4endl
2445 << " ExitNormal = " << n.fExitNormal << G4endl
2446 << " Exiting = " << n.fExiting << G4endl
2447 << " Entering = " << n.fEntering << G4endl
2448 << " BlockedPhysicalVolume= " ;
2449
2450 if (n.fBlockedPhysicalVolume==0)
2451 {
2452 os << "None";
2453 }
2454 else
2455 {
2456 os << n.fBlockedPhysicalVolume->GetName();
2457 }
2458
2459 os << G4endl
2460 << " BlockedReplicaNo = " << n.fBlockedReplicaNo << G4endl
2461 << " LastStepWasZero = " << n.fLastStepWasZero << G4endl
2462 << G4endl;
2463 }
2464 if( ( 1 < n.fVerbose) && (n.fVerbose < 4) )
2465 {
2466 os << G4endl; // Make sure to line up
2467 os << std::setw(30) << " ExitNormal " << " "
2468 << std::setw( 5) << " Valid " << " "
2469 << std::setw( 9) << " Exiting " << " "
2470 << std::setw( 9) << " Entering" << " "
2471 << std::setw(15) << " Blocked:Volume " << " "
2472 << std::setw( 9) << " ReplicaNo" << " "
2473 << std::setw( 8) << " LastStepZero " << " "
2474 << G4endl;
2475 os << "( " << std::setw(7) << n.fExitNormal.x()
2476 << ", " << std::setw(7) << n.fExitNormal.y()
2477 << ", " << std::setw(7) << n.fExitNormal.z() << " ) "
2478 << std::setw( 5) << n.fValidExitNormal << " "
2479 << std::setw( 9) << n.fExiting << " "
2480 << std::setw( 9) << n.fEntering << " ";
2481
2482 if ( n.fBlockedPhysicalVolume==0 )
2483 { os << std::setw(15) << "None"; }
2484 else
2485 { os << std::setw(15)<< n.fBlockedPhysicalVolume->GetName(); }
2486
2487 os << std::setw( 9) << n.fBlockedReplicaNo << " "
2488 << std::setw( 8) << n.fLastStepWasZero << " "
2489 << G4endl;
2490 }
2491 if( n.fVerbose > 2 )
2492 {
2493 os.precision(8);
2494 os << " Current Localpoint = " << n.fLastLocatedPointLocal << G4endl;
2495 os << " PreviousSftOrigin = " << n.fPreviousSftOrigin << G4endl;
2496 os << " PreviousSafety = " << n.fPreviousSafety << G4endl;
2497 }
2498 if( n.fVerbose > 3 || n.fVerbose == 0 )
2499 {
2500 os << "Current History: " << G4endl << n.fHistory;
2501 }
2502
2503 os.precision(oldcoutPrec);
2504 return os;
2505}

Member Data Documentation

◆ fCheck

G4bool G4ITNavigator2::fCheck

Definition at line 585 of file G4ITNavigator2.hh.

Referenced by GetLocalExitNormal(), and LocateGlobalPointAndSetup().

◆ fMaxNav

const G4int G4ITNavigator2::fMaxNav = 8
static

Definition at line 104 of file G4ITNavigator2.hh.

◆ fnormalNav

G4NormalNavigation G4ITNavigator2::fnormalNav

◆ fparamNav

◆ fpNavigatorState

◆ fpVoxelSafety

G4VoxelSafety* G4ITNavigator2::fpVoxelSafety

Definition at line 598 of file G4ITNavigator2.hh.

Referenced by ComputeSafety(), G4ITNavigator2(), and ~G4ITNavigator2().

◆ fregularNav

G4RegularNavigation G4ITNavigator2::fregularNav

◆ freplicaNav

G4ReplicaNavigation G4ITNavigator2::freplicaNav

◆ fTopPhysical

G4VPhysicalVolume* G4ITNavigator2::fTopPhysical

Definition at line 579 of file G4ITNavigator2.hh.

Referenced by NewNavigatorState(), and NewNavigatorStateAndLocate().

◆ fVerbose

◆ fvoxelNav

G4VoxelNavigation G4ITNavigator2::fvoxelNav

◆ fWarnPush

G4bool G4ITNavigator2::fWarnPush

Definition at line 588 of file G4ITNavigator2.hh.

Referenced by ComputeStep().

◆ kCarTolerance

G4double G4ITNavigator2::kCarTolerance
protected

Definition at line 408 of file G4ITNavigator2.hh.


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