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

G4Navigator is a class for use by the tracking management, able to obtain/calculate dynamic tracking time information such as the distance to the next volume, or to find the physical volume containing a given point in the world reference system. The navigator maintains a transformation history and other information to optimise the tracking time performance. More...

#include <G4Navigator.hh>

Inheritance diagram for G4Navigator:

Public Member Functions

 G4Navigator ()
 G4Navigator (const G4Navigator &)=delete
G4Navigatoroperator= (const G4Navigator &)=delete
virtual ~G4Navigator ()
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=nullptr, 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)
G4TouchableHistoryCreateTouchableHistory () const
G4TouchableHistoryCreateTouchableHistory (const G4NavigationHistory *h) const
virtual G4TouchableHandle 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)
G4VExternalNavigationGetExternalNavigation () const
void SetExternalNavigation (G4VExternalNavigation *externalNav)
G4VoxelNavigationGetVoxelNavigator ()
void SetVoxelNavigation (G4VoxelNavigation *voxelNav)
G4NavigatorClone () const
G4ThreeVector GetLastStepEndPoint () const
void InformLastStep (G4double lastStep, G4bool entersDaughtVol, G4bool exitsMotherVol)

Protected Member Functions

void SetSavedState ()
void RestoreSavedState ()
virtual void ResetState ()
G4ThreeVector ComputeLocalPoint (const G4ThreeVector &rGlobP) const
G4ThreeVector ComputeLocalAxis (const G4ThreeVector &pVec) const
EVolume VolumeType (const G4VPhysicalVolume *pVol) const
EVolume CharacteriseDaughters (const G4LogicalVolume *pLog) const
G4int GetDaughtersRegularStructureId (const G4LogicalVolume *pLv) const
virtual void SetupHierarchy ()
G4bool CheckOverlapsIterative (G4VPhysicalVolume *vol)

Protected Attributes

G4double kCarTolerance
G4double fMinStep
G4double fSqTol
G4NavigationHistory fHistory
G4ThreeVector fStepEndPoint
G4ThreeVector fLastStepEndPointLocal
G4int fVerbose = 0
G4bool fEnteredDaughter
G4bool fExitedMother
G4bool fWasLimitedByGeometry = false

Friends

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

Detailed Description

G4Navigator is a class for use by the tracking management, able to obtain/calculate dynamic tracking time information such as the distance to the next volume, or to find the physical volume containing a given point in the world reference system. The navigator maintains a transformation history and other information to optimise the tracking time performance.

Definition at line 75 of file G4Navigator.hh.

Constructor & Destructor Documentation

◆ G4Navigator() [1/2]

G4Navigator::G4Navigator ( )

Constructor - initialisers and setup.

Definition at line 54 of file G4Navigator.cc.

55{
57 // Initialises also all
58 // - exit / entry flags
59 // - flags & variables for exit normals
60 // - zero step counters
61 // - blocked volume
62
63 if( fVerbose > 2 )
64 {
65 G4cout << " G4Navigator parameters: Action Threshold (No Zero Steps) = "
66 << fActionThreshold_NoZeroSteps
67 << " Abandon Threshold (No Zero Steps) = "
68 << fAbandonThreshold_NoZeroSteps << G4endl;
69 }
73
74 fregularNav.SetNormalNavigation( &fnormalNav );
75
76 fStepEndPoint = G4ThreeVector( kInfinity, kInfinity, kInfinity );
77 fLastStepEndPointLocal = G4ThreeVector( kInfinity, kInfinity, kInfinity );
78
79 fpVoxelSafety = new G4VoxelSafety();
80 fpvoxelNav = new G4VoxelNavigation();
81 fpSafetyCalculator = new G4SafetyCalculator( *this, fHistory );
82 fpSafetyCalculator->SetExternalNavigation(fpExternalNav);
83}
CLHEP::Hep3Vector G4ThreeVector
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4double fMinStep
G4ThreeVector fStepEndPoint
G4ThreeVector fLastStepEndPointLocal
G4double fSqTol
G4double kCarTolerance
void ResetStackAndState()
G4NavigationHistory fHistory
T sqr(const T &x)
Definition templates.hh:128

Referenced by Clone(), G4MultiNavigator::G4MultiNavigator(), G4Navigator(), G4MultiNavigator::GetNavigator(), G4MultiNavigator::operator<<, operator=(), G4MultiNavigator::PrintLimited(), G4MultiNavigator::ResetHierarchyAndLocate(), and G4MultiNavigator::~G4MultiNavigator().

◆ G4Navigator() [2/2]

G4Navigator::G4Navigator ( const G4Navigator & )
delete

Copy constructor & assignment operator not allowed.

◆ ~G4Navigator()

G4Navigator::~G4Navigator ( )
virtual

Destructor.

Definition at line 89 of file G4Navigator.cc.

90{
91 delete fpVoxelSafety;
92 delete fpExternalNav;
93 delete fpvoxelNav;
94 delete fpSafetyCalculator;
95}

Member Function Documentation

◆ Activate()

void G4Navigator::Activate ( G4bool flag)
inline

Activate/inactivate the navigator.

Referenced by G4TransportationManager::ActivateNavigator().

◆ CharacteriseDaughters()

EVolume G4Navigator::CharacteriseDaughters ( const G4LogicalVolume * pLog) const
inlineprotected

Characterises the daughters of given logical volume.

Referenced by ComputeStep(), LocateGlobalPointAndSetup(), and LocateGlobalPointWithinVolume().

◆ CheckMode()

void G4Navigator::CheckMode ( G4bool mode)
inline

Run navigation in "check-mode", therefore using additional verifications and more strict correctness conditions.

Note
Is effective only with G4VERBOSE set.

Referenced by G4VIntersectionLocator::LocateGlobalPointWithinVolumeAndCheck().

◆ CheckNextStep()

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

Same as ComputeStep() above, but does not affect/modify the state of the Navigator.

Definition at line 1228 of file G4Navigator.cc.

1232{
1233 G4double step;
1234
1235 // Save the state, for this parasitic call
1236 //
1237 SetSavedState();
1238
1239 step = ComputeStep ( pGlobalpoint,
1240 pDirection,
1241 pCurrentProposedStepLength,
1242 pNewSafety );
1243
1244 // It is a parasitic call, so attempt to restore the key parts of the state
1245 //
1247 // NOTE: the state of the current subnavigator is NOT restored.
1248 // ***> TODO: restore subnavigator state
1249 // if( last_located) Need Position of last location
1250 // if( last_computed step) Need Endposition of last step
1251
1252 return step;
1253}
double G4double
Definition G4Types.hh:83
void RestoreSavedState()
void SetSavedState()
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)

◆ CheckOverlapsIterative()

G4bool G4Navigator::CheckOverlapsIterative ( G4VPhysicalVolume * vol)
protected

Utility method to trigger overlaps check on a volume with reported overlaps ordered by relevance. Used in ComputeStep() when loopings with zero step are detected.

Definition at line 1986 of file G4Navigator.cc.

1987{
1988 // Check and report overlaps
1989 //
1990 G4bool foundOverlap = false;
1991 G4int nPoints = 300000, ntrials = 9, numOverlaps = 5;
1992 G4double trialLength = 1.0 * CLHEP::centimeter;
1993 while ( ntrials-- > 0 && !foundOverlap )
1994 {
1995 if ( fVerbose > 1 )
1996 {
1997 G4cout << " ** Running overlap checks in volume "
1998 << vol->GetName()
1999 << " with length = " << trialLength << G4endl;
2000 }
2001 foundOverlap = vol->CheckOverlaps(nPoints, trialLength,
2002 fVerbose != 0, numOverlaps);
2003 trialLength *= 0.1;
2004 if ( trialLength <= 1.0e-5 ) { numOverlaps= 1;}
2005 }
2006 return foundOverlap;
2007}
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
virtual G4bool CheckOverlaps(G4int res=1000, G4double tol=0., G4bool verbose=true, G4int errMax=1)
const G4String & GetName() const

Referenced by ComputeStep().

◆ Clone()

G4Navigator * G4Navigator::Clone ( ) const
inline

Cloning feature for use in MT applications to clone the navigator, including external sub-navigator.

Note
Client has responsibility for ownership of the returned allocated pointer.
Returns
A pointer to the cloned navigator object.

◆ ComputeLocalAxis()

G4ThreeVector G4Navigator::ComputeLocalAxis ( const G4ThreeVector & pVec) const
inlineprotected

Returns the local direction of the specified vector in the reference system of the volume that was found by LocalGlobalPointAndSetup(). The Local Coordinates of point in world coordinate system.

Referenced by ComputeStep().

◆ ComputeLocalPoint()

G4ThreeVector G4Navigator::ComputeLocalPoint ( const G4ThreeVector & rGlobP) const
inlineprotected

Returns position vector in local coordinate system, given a position vector in world coordinate system.

Referenced by ComputeStep(), and LocateGlobalPointWithinVolume().

◆ ComputeSafety()

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

Calculates the isotropic distance to the nearest boundary from the specified point in the global coordinate system.

Note
The geometry must be closed.
Parameters
[in]globalpointThe point in global coordinates system. The point must be within the current volume.
[in]pProposedMaxLengthThe proposed maximum length is used to avoid volume safety calculations.
[in]keepStateFlag to instruct keeping the state (default true) to ensure minimum side effects from the call.
Returns
Length from current point to closest boundary surface. The value returned is usually an underestimate.

Reimplemented in G4ErrorPropagationNavigator, and G4MultiNavigator.

Definition at line 1780 of file G4Navigator.cc.

1783{
1784 G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
1785 G4double safety = 0.0;
1786
1787 G4double distEndpointSq = (pGlobalpoint-fStepEndPoint).mag2();
1788 G4bool stayedOnEndpoint = distEndpointSq < sqr(kCarTolerance);
1789 G4bool endpointOnSurface = fEnteredDaughter || fExitedMother;
1790
1791 G4bool onSurface = endpointOnSurface && stayedOnEndpoint;
1792 if( ! onSurface )
1793 {
1794 safety= fpSafetyCalculator->SafetyInCurrentVolume(pGlobalpoint, motherPhysical, pMaxLength);
1795 // offload to G4SafetyCalculator - avoids need to save / reload state
1796
1797 // Remember last safety origin & value
1798 //
1799 fPreviousSftOrigin = pGlobalpoint;
1800 fPreviousSafety = safety;
1801 // We overwrite the Safety 'sphere' - keeping old behaviour
1802 }
1803
1804 return safety;
1805}
G4bool fExitedMother
G4bool fEnteredDaughter

Referenced by G4ErrorPropagationNavigator::ComputeSafety().

◆ ComputeStep()

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

Calculates the distance to the next boundary intersected along the specified NORMALISED vector direction and from the specified point in the global coordinate system. LocateGlobalPointAndSetup() or LocateGlobalPointWithinVolume() must have been called with the same global point prior to this call. The isotropic distance to the nearest boundary is also calculated (usually an underestimate). The current proposed Step length is used to avoid intersection calculations: if it can be determined that the nearest boundary is >pCurrentProposedStepLength away, kInfinity is returned together with the computed isotropic safety distance.

Note
Geometry must be closed.
Parameters
[in]pGlobalPointThe point in global coordinates system.
[in]pDirectionThe normalised vector direction.
[in]pCurrentProposedStepLengthCurrent proposed step length.
[in,out]newSafetyNew safety.
Returns
Length from current point to next boundary surface along pDirection.

Reimplemented in G4ErrorPropagationNavigator, and G4MultiNavigator.

Definition at line 751 of file G4Navigator.cc.

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

Referenced by CheckNextStep(), and G4ErrorPropagationNavigator::ComputeStep().

◆ CreateTouchableHistory() [1/2]

G4TouchableHistory * G4Navigator::CreateTouchableHistory ( ) const
inline

Touchable creation method.

Note
Caller has deletion responsibility.
Returns
A pointer to the allocated touchable history.

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

◆ CreateTouchableHistory() [2/2]

G4TouchableHistory * G4Navigator::CreateTouchableHistory ( const G4NavigationHistory * h) const
inline

Touchable creation method, given a history.

Note
Caller has deletion responsibility.
Parameters
[in]hPointer to a navigation history to copy from.
Returns
A pointer to the allocated touchable history.

◆ CreateTouchableHistoryHandle()

G4TouchableHandle G4Navigator::CreateTouchableHistoryHandle ( ) const
virtual

Returns a reference counted handle to a touchable history.

Reimplemented in G4MultiNavigator.

Definition at line 1811 of file G4Navigator.cc.

1812{
1813 return { CreateTouchableHistory() };
1814}
G4TouchableHistory * CreateTouchableHistory() const

Referenced by G4VIntersectionLocator::LocateGlobalPointWithinVolumeAndCheck().

◆ EnableBestSafety()

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

Enables best-possible evaluation of isotropic safety.

◆ EnteredDaughterVolume()

G4bool G4Navigator::EnteredDaughterVolume ( ) const
inline

The purpose of this function is to inform the caller if the track is entering a daughter volume while exiting from the current volume.

Note
It is not guaranteed to work if SetGeometricallyLimitedStep() was not called when it should have been called.
Returns
True only in case when the Step has caused the track to arrive at a boundary of a daughter. False, in all other cases.

Referenced by GetLocalExitNormal().

◆ ExitedMotherVolume()

G4bool G4Navigator::ExitedMotherVolume ( ) const
inline

Verify if the step has exited the mother volume.

◆ GetCurrentLocalCoordinate()

G4ThreeVector G4Navigator::GetCurrentLocalCoordinate ( ) const
inline

Returns the local coordinate of the point in the reference system of its containing volume that was found by LocalGlobalPointAndSetup(). The local coordinate of the last located track.

◆ GetDaughtersRegularStructureId()

G4int G4Navigator::GetDaughtersRegularStructureId ( const G4LogicalVolume * pLv) const
inlineprotected

Gets regular structure ID of first daughter.

Referenced by ComputeStep(), and LocateGlobalPointAndSetup().

◆ GetExternalNavigation()

G4VExternalNavigation * G4Navigator::GetExternalNavigation ( ) const
inline

Accessor & modifier for custom external navigation.

◆ GetGlobalExitNormal()

G4ThreeVector G4Navigator::GetGlobalExitNormal ( const G4ThreeVector & point,
G4bool * valid )
virtual

Obtains the Normal vector to a surface (in global coordinates) pointing out of previous volume and into current volume The method takes full care about how to calculate the normal, but if the surfaces are not convex it will return valid=false.

Note
Can only be called if the Navigator's last Step has crossed a volume geometrical boundary.
Normals are not available for replica volumes (i.e. valid=false).
Parameters
[in]pointPoint in global coordinates system to compare to.
[in,out]validFlag indicating if normal is valid.
Returns
A Exit Surface Normal vector and validity too.

Reimplemented in G4ErrorPropagationNavigator, and G4MultiNavigator.

Definition at line 1609 of file G4Navigator.cc.

1611{
1612 G4bool validNormal;
1613 G4ThreeVector localNormal, globalNormal;
1614
1615 G4bool usingStored = fCalculatedExitNormal && (
1616 ( fLastTriedStepComputation && fExiting ) // Just calculated it
1617 || // No locate in between
1618 ( !fLastTriedStepComputation
1619 && (IntersectPointGlobal-fStepEndPoint).mag2() < 10.0*fSqTol ) );
1620 // Calculated it 'just' before & then called locate
1621 // but it did not move position
1622
1623 if( usingStored )
1624 {
1625 // This was computed in last call to ComputeStep
1626 // and only if it arrived at boundary
1627 //
1628 globalNormal = fExitNormalGlobalFrame;
1629 G4double normMag2 = globalNormal.mag2();
1630 if( std::fabs ( normMag2 - 1.0 ) < perThousand ) // was perMillion
1631 {
1632 *pNormalCalculated = true; // ComputeStep always computes it if Exiting
1633 // (fExiting==true)
1634 }
1635 else
1636 {
1637 G4ExceptionDescription message;
1638 message.precision(10);
1639 message << " WARNING> Expected normal-global-frame to be valid, "
1640 << " i.e. a unit vector!" << G4endl
1641 << " - but |normal| = " << std::sqrt(normMag2)
1642 << " - and |normal|^2 = " << normMag2 << G4endl
1643 << " which differs from 1.0 by " << normMag2 - 1.0 << G4endl
1644 << " n = " << fExitNormalGlobalFrame << G4endl
1645 << " Global point: " << IntersectPointGlobal << G4endl
1646 << " Volume: " << fHistory.GetTopVolume()->GetName() << G4endl;
1647#ifdef G4VERBOSE
1648 G4LogicalVolume* candLog = fHistory.GetTopVolume()->GetLogicalVolume();
1649 if ( candLog != nullptr )
1650 {
1651 message << " Solid: " << candLog->GetSolid()->GetName()
1652 << ", Type: " << candLog->GetSolid()->GetEntityType() << G4endl
1653 << *candLog->GetSolid() << G4endl;
1654 }
1655#endif
1656 message << "============================================================"
1657 << G4endl;
1658 G4int oldVerbose = fVerbose;
1659 fVerbose = 4;
1660 message << " State of Navigator: " << G4endl;
1661 message << *this << G4endl;
1662 fVerbose = oldVerbose;
1663 message << "============================================================"
1664 << G4endl;
1665
1666 G4Exception("G4Navigator::GetGlobalExitNormal()",
1667 "GeomNav0003",JustWarning, message,
1668 "Value obtained from stored global-normal is not a unit vector.");
1669
1670 // (Re)Compute it now -- as either it was not computed, or it is wrong.
1671 //
1672 localNormal = GetLocalExitNormalAndCheck(IntersectPointGlobal,
1673 &validNormal);
1674 *pNormalCalculated = fCalculatedExitNormal;
1675 globalNormal = fHistory.GetTopTransform()
1676 .InverseTransformAxis(localNormal);
1677 }
1678 }
1679 else
1680 {
1681 localNormal = GetLocalExitNormalAndCheck(IntersectPointGlobal,&validNormal);
1682 *pNormalCalculated = fCalculatedExitNormal;
1683
1684#ifdef G4DEBUG_NAVIGATION
1685 usingStored = false;
1686
1687 if( (!validNormal) && !fCalculatedExitNormal )
1688 {
1690 edN << " Calculated = " << fCalculatedExitNormal << G4endl;
1691 edN << " Entering= " << fEntering << G4endl;
1692 G4int oldVerbose = this->GetVerboseLevel();
1693 this->SetVerboseLevel(4);
1694 edN << " State of Navigator: " << G4endl;
1695 edN << *this << G4endl;
1696 this->SetVerboseLevel( oldVerbose );
1697
1698 G4Exception("G4Navigator::GetGlobalExitNormal()",
1699 "GeomNav0003", JustWarning, edN,
1700 "LocalExitNormalAndCheck() did not calculate Normal.");
1701 }
1702#endif
1703
1704 G4double localMag2 = localNormal.mag2();
1705 if( validNormal && (std::fabs(localMag2-1.0)) > kToleranceNormalCheck )
1706 {
1708 edN.precision(10);
1709 edN << "G4Navigator::GetGlobalExitNormal: "
1710 << " Using Local Normal - from call to GetLocalExitNormalAndCheck. "
1711 << G4endl
1712 << " Local Exit Normal : " << " || = " << std::sqrt(localMag2)
1713 << " vec = " << localNormal << G4endl
1714 << " Global Exit Normal : " << " || = " << globalNormal.mag()
1715 << " vec = " << globalNormal << G4endl
1716 << " Global point: " << IntersectPointGlobal << G4endl;
1717 edN << " Calculated It = " << fCalculatedExitNormal << G4endl
1718 << " Volume: " << fHistory.GetTopVolume()->GetName() << G4endl;
1719#ifdef G4VERBOSE
1720 G4LogicalVolume* candLog = fHistory.GetTopVolume()->GetLogicalVolume();
1721 if ( candLog != nullptr )
1722 {
1723 edN << " Solid: " << candLog->GetSolid()->GetName()
1724 << ", Type: " << candLog->GetSolid()->GetEntityType() << G4endl
1725 << *candLog->GetSolid();
1726 }
1727#endif
1728 G4Exception("G4Navigator::GetGlobalExitNormal()",
1729 "GeomNav0003",JustWarning, edN,
1730 "Value obtained from new local *solid* is incorrect.");
1731 localNormal = localNormal.unit(); // Should we correct it ??
1732 }
1733 globalNormal = fHistory.GetTopTransform()
1734 .InverseTransformAxis(localNormal);
1735 }
1736
1737#ifdef G4DEBUG_NAVIGATION
1738 if( usingStored )
1739 {
1740 G4ThreeVector globalNormAgn;
1741
1742 localNormal = GetLocalExitNormalAndCheck(IntersectPointGlobal,&validNormal);
1743
1744 globalNormAgn = fHistory.GetTopTransform()
1745 .InverseTransformAxis(localNormal);
1746
1747 // Check the value computed against fExitNormalGlobalFrame
1748 G4ThreeVector diffNorm = globalNormAgn - fExitNormalGlobalFrame;
1749 if( diffNorm.mag2() > kToleranceNormalCheck )
1750 {
1752 edDfn << "Found difference in normals in case of exiting mother "
1753 << "- when Get is called after ComputingStep " << G4endl;
1754 edDfn << " Magnitude of diff = " << diffNorm.mag() << G4endl;
1755 edDfn << " Normal stored (Global) = " << fExitNormalGlobalFrame
1756 << G4endl;
1757 edDfn << " Global Computed from Local = " << globalNormAgn << G4endl;
1758 G4Exception("G4Navigator::GetGlobalExitNormal()", "GeomNav0003",
1759 JustWarning, edDfn);
1760 }
1761 }
1762#endif
1763
1764 // Synchronise stored global exit normal as possibly re-computed here
1765 //
1766 fExitNormalGlobalFrame = globalNormal;
1767
1768 return globalNormal;
1769}
Hep3Vector unit() const
double mag2() const
double mag() const
void SetVerboseLevel(G4int level)
G4int GetVerboseLevel() const
virtual G4ThreeVector GetLocalExitNormalAndCheck(const G4ThreeVector &point, G4bool *valid)
G4String GetName() const
virtual G4GeometryType GetEntityType() const =0

Referenced by G4ErrorPropagationNavigator::GetGlobalExitNormal(), G4XrayReflection::GetMeanFreePath(), and G4MicroElecSurface::PostStepDoIt().

◆ GetGlobalToLocalTransform()

const G4AffineTransform & G4Navigator::GetGlobalToLocalTransform ( ) const
inline

Obtains the transformations Global/Local (and inverse).

Note
Clients of these methods must copy the data if they need to keep it.

Referenced by GetLocalExitNormalAndCheck(), and G4VIntersectionLocator::LocateGlobalPointWithinVolumeAndCheck().

◆ GetLastStepEndPoint()

G4ThreeVector G4Navigator::GetLastStepEndPoint ( ) const
inline

Gets endpoint of last step.

◆ GetLocalExitNormal()

G4ThreeVector G4Navigator::GetLocalExitNormal ( G4bool * valid)
virtual

Obtains the Normal vector to a surface (in local coordinates) pointing out of previous volume and into current volume Convention: the local normal is in the coordinate system of the final volume. The method takes full care about how to calculate this normal, but if the surfaces are not convex it will return valid=false.

Note
Can only be called if the Navigator's last Step has crossed a volume geometrical boundary.
Normals are not available for replica volumes (i.e. valid=false).
Parameters
[in,out]validFlag indicating if normal is valid.
Returns
A Exit Surface Normal vector and validity too.

Reimplemented in G4MultiNavigator.

Definition at line 1357 of file G4Navigator.cc.

1358{
1359 G4ThreeVector ExitNormal(0.,0.,0.);
1360 G4VSolid* currentSolid = nullptr;
1361 G4LogicalVolume* candidateLogical;
1362
1363 if ( fLastTriedStepComputation )
1364 {
1365 // use fLastLocatedPointLocal and next candidate volume
1366 //
1367 G4ThreeVector nextSolidExitNormal(0.,0.,0.);
1368
1369 if( fEntering && (fBlockedPhysicalVolume!=nullptr) )
1370 {
1371 candidateLogical = fBlockedPhysicalVolume->GetLogicalVolume();
1372 if( candidateLogical != nullptr )
1373 {
1374 // fLastStepEndPointLocal is in the coordinates of the mother
1375 // we need it in the daughter's coordinate system.
1376
1377 // The following code should also work in case of Replica
1378 {
1379 // First transform fLastLocatedPointLocal to the new daughter
1380 // coordinates
1381 //
1382 G4AffineTransform MotherToDaughterTransform=
1383 GetMotherToDaughterTransform( fBlockedPhysicalVolume,
1384 fBlockedReplicaNo,
1385 VolumeType(fBlockedPhysicalVolume) );
1386 G4ThreeVector daughterPointOwnLocal =
1387 MotherToDaughterTransform.TransformPoint( fLastStepEndPointLocal );
1388
1389 // OK if it is a parameterised volume
1390 //
1391 EInside inSideIt;
1392 G4bool onSurface;
1393 G4double safety = -1.0;
1394 currentSolid = candidateLogical->GetSolid();
1395 inSideIt = currentSolid->Inside(daughterPointOwnLocal);
1396 onSurface = (inSideIt == kSurface);
1397 if( !onSurface )
1398 {
1399 if( inSideIt == kOutside )
1400 {
1401 safety = (currentSolid->DistanceToIn(daughterPointOwnLocal));
1402 onSurface = safety < 100.0 * kCarTolerance;
1403 }
1404 else if (inSideIt == kInside )
1405 {
1406 safety = (currentSolid->DistanceToOut(daughterPointOwnLocal));
1407 onSurface = safety < 100.0 * kCarTolerance;
1408 }
1409 }
1410
1411 if( onSurface )
1412 {
1413 nextSolidExitNormal =
1414 currentSolid->SurfaceNormal(daughterPointOwnLocal);
1415
1416 // Entering the solid ==> opposite
1417 //
1418 // First flip ( ExitNormal = -nextSolidExitNormal; )
1419 // and then rotate the the normal to the frame of the mother (current volume)
1420 ExitNormal = MotherToDaughterTransform
1421 .InverseTransformAxis( -nextSolidExitNormal );
1422 fCalculatedExitNormal = true;
1423 }
1424 else
1425 {
1426#ifdef G4VERBOSE
1427 if(( fVerbose == 1 ) && ( fCheck ))
1428 {
1429 std::ostringstream message;
1430 message << "Point not on surface ! " << G4endl
1431 << " Point = "
1432 << daughterPointOwnLocal << G4endl
1433 << " Physical volume = "
1434 << fBlockedPhysicalVolume->GetName() << G4endl
1435 << " Logical volume = "
1436 << candidateLogical->GetName() << G4endl
1437 << " Solid = " << currentSolid->GetName()
1438 << " Type = "
1439 << currentSolid->GetEntityType() << G4endl
1440 << *currentSolid << G4endl;
1441 if( inSideIt == kOutside )
1442 {
1443 message << "Point is Outside. " << G4endl
1444 << " Safety (from outside) = " << safety << G4endl;
1445 }
1446 else // if( inSideIt == kInside )
1447 {
1448 message << "Point is Inside. " << G4endl
1449 << " Safety (from inside) = " << safety << G4endl;
1450 }
1451 G4Exception("G4Navigator::GetLocalExitNormal()", "GeomNav1001",
1452 JustWarning, message);
1453 }
1454#endif
1455 }
1456 *valid = onSurface; // was =true;
1457 }
1458 }
1459 }
1460 else if ( fExiting )
1461 {
1462 ExitNormal = fGrandMotherExitNormal;
1463 *valid = true;
1464 fCalculatedExitNormal = true; // Should be true already
1465 }
1466 else // i.e. ( fBlockedPhysicalVolume == 0 )
1467 {
1468 *valid = false;
1469 G4Exception("G4Navigator::GetLocalExitNormal()",
1470 "GeomNav0003", JustWarning,
1471 "Incorrect call to GetLocalSurfaceNormal." );
1472 }
1473 }
1474 else // ( ! fLastTriedStepComputation ) i.e. last call was to Locate
1475 {
1476 if ( EnteredDaughterVolume() )
1477 {
1478 G4VSolid* daughterSolid = fHistory.GetTopVolume()->GetLogicalVolume()
1479 ->GetSolid();
1480 ExitNormal = -(daughterSolid->SurfaceNormal(fLastLocatedPointLocal));
1481 if( std::fabs(ExitNormal.mag2()-1.0 ) > kToleranceNormalCheck )
1482 {
1484 desc << " Parameters of solid: " << *daughterSolid
1485 << " Point for surface = " << fLastLocatedPointLocal << std::endl;
1486 G4Exception("G4Navigator::GetLocalExitNormal()",
1487 "GeomNav0003", FatalException, desc,
1488 "Surface Normal returned by Solid is not a Unit Vector." );
1489 }
1490 fCalculatedExitNormal = true;
1491 *valid = true;
1492 }
1493 else
1494 {
1495 if( fExitedMother )
1496 {
1497 ExitNormal = fGrandMotherExitNormal;
1498 *valid = true;
1499 fCalculatedExitNormal = true;
1500 }
1501 else // We are not at a boundary. ExitNormal remains (0,0,0)
1502 {
1503 *valid = false;
1504 fCalculatedExitNormal = false;
1505 G4ExceptionDescription message;
1506 message << "Function called when *NOT* at a Boundary." << G4endl;
1507 message << "Exit Normal not calculated." << G4endl;
1508 G4Exception("G4Navigator::GetLocalExitNormal()",
1509 "GeomNav0003", JustWarning, message);
1510 }
1511 }
1512 }
1513 return ExitNormal;
1514}
EVolume VolumeType() const override
G4ThreeVector InverseTransformAxis(const G4ThreeVector &axis) const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
const G4String & GetName() const
G4AffineTransform GetMotherToDaughterTransform(G4VPhysicalVolume *dVolume, G4int dReplicaNo, EVolume dVolumeType)
G4bool EnteredDaughterVolume() 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
EInside
Definition geomdefs.hh:67
@ kInside
Definition geomdefs.hh:70
@ kOutside
Definition geomdefs.hh:68
@ kSurface
Definition geomdefs.hh:69

Referenced by G4RayTrajectory::AppendStep(), GetLocalExitNormalAndCheck(), and G4VTransitionRadiation::PostStepDoIt().

◆ GetLocalExitNormalAndCheck()

G4ThreeVector G4Navigator::GetLocalExitNormalAndCheck ( const G4ThreeVector & point,
G4bool * valid )
virtual

Obtains the Normal vector to a surface (in local coordinates) pointing out of previous volume and into current volume, and checks the current point against expected 'local' value. Convention: the local normal is in the coordinate system of the final volume. The method takes full care about how to calculate this normal, but if the surfaces are not convex it will return valid=false.

Note
Can only be called if the Navigator's last Step has crossed a volume geometrical boundary.
Normals are not available for replica volumes (i.e. valid=false).
Parameters
[in]pointPoint in global coordinates system to compare to.
[in,out]validFlag indicating if normal is valid.
Returns
A Exit Surface Normal vector and validity too.

Reimplemented in G4MultiNavigator.

Definition at line 1574 of file G4Navigator.cc.

1581{
1582#ifdef G4DEBUG_NAVIGATION
1583 // Check Current point against expected 'local' value
1584 //
1585 if ( fLastTriedStepComputation )
1586 {
1587 G4ThreeVector ExpectedBoundaryPointLocal;
1588
1589 const G4AffineTransform& GlobalToLocal = GetGlobalToLocalTransform();
1590 ExpectedBoundaryPointLocal =
1591 GlobalToLocal.TransformPoint( ExpectedBoundaryPointGlobal );
1592
1593 // Add here: Comparison against expected position,
1594 // i.e. the endpoint of ComputeStep
1595 }
1596#endif
1597
1598 return GetLocalExitNormal( pValid );
1599}
virtual G4ThreeVector GetLocalExitNormal(G4bool *valid)
const G4AffineTransform & GetGlobalToLocalTransform() const

Referenced by GetGlobalExitNormal().

◆ GetLocalToGlobalTransform()

const G4AffineTransform G4Navigator::GetLocalToGlobalTransform ( ) const
inline

◆ GetMotherToDaughterTransform()

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

Obtains mother to daughter transformation.

Definition at line 1523 of file G4Navigator.cc.

1526{
1527 switch (enteringVolumeType)
1528 {
1529 case kNormal: // Nothing is needed to prepare the transformation
1530 break; // It is stored already in the physical volume (placement)
1531 case kReplica: // Sets the transform in the Replica - tbc
1532 G4Exception("G4Navigator::GetMotherToDaughterTransform()",
1533 "GeomNav0001", FatalException,
1534 "Method NOT Implemented yet for replica volumes.");
1535 break;
1536 case kParameterised:
1537 if( pEnteringPhysVol->GetRegularStructureId() == 0 )
1538 {
1539 G4VPVParameterisation *pParam =
1540 pEnteringPhysVol->GetParameterisation();
1541 G4VSolid* pSolid =
1542 pParam->ComputeSolid(enteringReplicaNo, pEnteringPhysVol);
1543 pSolid->ComputeDimensions(pParam, enteringReplicaNo, pEnteringPhysVol);
1544
1545 // Sets the transform in the Parameterisation
1546 //
1547 pParam->ComputeTransformation(enteringReplicaNo, pEnteringPhysVol);
1548
1549 // Set the correct solid and material in Logical Volume
1550 //
1551 G4LogicalVolume* pLogical = pEnteringPhysVol->GetLogicalVolume();
1552 pLogical->SetSolid( pSolid );
1553 }
1554 break;
1555 case kExternal:
1556 // Expect that nothing is needed to prepare the transformation.
1557 // It is stored already in the physical volume (placement)
1558 break;
1559 }
1560 return G4AffineTransform(pEnteringPhysVol->GetRotation(),
1561 pEnteringPhysVol->GetTranslation()).Invert();
1562}
void SetSolid(G4VSolid *pSolid)
virtual void ComputeTransformation(const G4int no, G4VPhysicalVolume *pv) const =0
virtual G4VSolid * ComputeSolid(const G4int no, G4VPhysicalVolume *pv)
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition G4VSolid.cc:136

Referenced by GetLocalExitNormal().

◆ GetVerboseLevel()

G4int G4Navigator::GetVerboseLevel ( ) const
inline

Verbosity control.

Note
If level>0 && G4VERBOSE, printout can occur.

Referenced by GetGlobalExitNormal().

◆ GetVoxelNavigator()

G4VoxelNavigation & G4Navigator::GetVoxelNavigator ( )
inline

Gets/sets alternative navigator for voxel volumes.

Referenced by ComputeStep(), LocateGlobalPointAndSetup(), and LocateGlobalPointWithinVolume().

◆ GetWorldVolume()

◆ InformLastStep()

void G4Navigator::InformLastStep ( G4double lastStep,
G4bool entersDaughtVol,
G4bool exitsMotherVol )

Derived navigators which rely on LocateGlobalPointAndSetup() need to inform size of step, to maintain logic about arriving on boundary for challenging cases. Required in order to cope with multiple trials at boundaries => Locate with use direction rather than simple, fast logic.

Definition at line 2098 of file G4Navigator.cc.

2100{
2101 G4bool zeroStep = ( lastStep == 0.0 );
2102 fLocatedOnEdge = fLastStepWasZero && zeroStep;
2103 fLastStepWasZero = zeroStep;
2104
2105 fExiting = exitsMotherVol;
2106 fEntering = entersDaughtVol;
2107}

◆ IsActive()

G4bool G4Navigator::IsActive ( ) const
inline

Verify if the navigator is active.

◆ IsCheckModeActive()

G4bool G4Navigator::IsCheckModeActive ( ) const
inline

Set/unset verbosity for pushed tracks (default is true).

Referenced by G4VIntersectionLocator::LocateGlobalPointWithinVolumeAndCheck().

◆ LocateGlobalPointAndSetup()

G4VPhysicalVolume * G4Navigator::LocateGlobalPointAndSetup ( const G4ThreeVector & point,
const G4ThreeVector * direction = nullptr,
const G4bool pRelativeSearch = true,
const G4bool ignoreDirection = true )
virtual

Searches the geometrical hierarchy for the volumes deepest in hierarchy containing the point in the global coordinate space. Two main cases are: i) If pRelativeSearch=false it makes use of no previous/state information. Returns the physical volume containing the point, with all previous mothers correctly set up. ii) If pRelativeSearch is set to true, the search begin is the geometrical hierarchy at the location of the last located point, or the endpoint of previous Step if SetGeometricallyLimitedStep() has been called immediately before. The direction is used (to check if a volume is entered) if either

  • the argument ignoreDirection is false, or
  • the Navigator has determined that it is on an edge shared by two or more volumes (this is state information).
    Note
    In order to call this the geometry MUST be closed.
    Parameters
    [in]pointThe point in global coordinates system.
    [in]directionThe normalised vector direction.
    [in]pRelativeSearchFlag to specify where search starts from.
    [in]ignoreDirectionFlag to specify if to use direction or not.
    Returns
    The pointer to the physical volume where point is located.

Reimplemented in G4MultiNavigator.

Definition at line 132 of file G4Navigator.cc.

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

Referenced by ComputeStep(), G4TheRayTracer::CreateBitMap(), G4VFieldModel::DescribeYourselfTo(), G4TrajectoryDrawByOriginVolume::Draw(), G4TrajectoryOriginVolumeFilter::Evaluate(), G4VIntersectionLocator::LocateGlobalPointWithinVolumeAndCheck(), and ResetHierarchyAndLocate().

◆ LocateGlobalPointAndUpdateTouchable() [1/2]

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

Same as the method above LocateGlobalPointAndUpdateTouchableHandle(), except that a pointer to G4VTouchable is used for updating the touchable.

◆ LocateGlobalPointAndUpdateTouchable() [2/2]

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

Same as the method above LocateGlobalPointAndUpdateTouchable(), except that direction is not specified.

◆ LocateGlobalPointAndUpdateTouchableHandle()

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

It first searches the geometrical hierarchy like the above method LocateGlobalPointAndSetup(), then it uses the volume found and its navigation history to update the touchable handle.

Parameters
[in]positionThe point in global coordinates system.
[in]directionThe normalised vector direction.
[in,out]oldTouchableToUpdateTouchable handle to update.
[in]RelativeSearchFlag to specify where search starts from.

◆ LocateGlobalPointWithinVolume()

void G4Navigator::LocateGlobalPointWithinVolume ( const G4ThreeVector & position)
virtual

Notifies the Navigator that a track has moved to the new Global point 'position', that is known to be within the current safety. No check is performed to ensure that it is within the volume. This method can be called instead of LocateGlobalPointAndSetup() ONLY if the caller is certain that the new global point (position) is inside the same volume as the previous position. Usually this can be guaranteed only if the point is within safety.

Parameters
[in]positionThe position point in global coordinates system.

Reimplemented in G4MultiNavigator.

Definition at line 602 of file G4Navigator.cc.

603{
604#ifdef G4DEBUG_NAVIGATION
605 assert( !fWasLimitedByGeometry );
606 // Check: Either step was not limited by a boundary or
607 // else the full step is no longer being taken
608#endif
609
610 fLastLocatedPointLocal = ComputeLocalPoint(pGlobalpoint);
611 fLastTriedStepComputation = false;
612 fChangedGrandMotherRefFrame = false; // Frame for Exit Normal
613
614 // For the case of Voxel (or Parameterised) volume the respective
615 // Navigator must be messaged to update its voxel information etc
616
617 // Update the state of the Sub Navigators
618 // - in particular any voxel information they store/cache
619 //
620 G4VPhysicalVolume* motherPhysical = fHistory.GetTopVolume();
621 G4LogicalVolume* motherLogical = motherPhysical->GetLogicalVolume();
622
623 switch( CharacteriseDaughters(motherLogical) )
624 {
625 case kNormal:
626 GetVoxelNavigator().RelocateWithinVolume( motherPhysical, fLastLocatedPointLocal );
627 break;
628 case kParameterised:
629 fparamNav.RelocateWithinVolume( motherPhysical, fLastLocatedPointLocal );
630 break;
631 case kReplica:
632 // Nothing to do
633 break;
634 case kExternal:
635 fpExternalNav->RelocateWithinVolume( motherPhysical,
636 fLastLocatedPointLocal );
637 break;
638 }
639
640 // Reset the state variables
641 // - which would have been affected
642 // by the 'equivalent' call to LocateGlobalPointAndSetup
643 // - who's values have been invalidated by the 'move'.
644 //
645 fBlockedPhysicalVolume = nullptr;
646 fBlockedReplicaNo = -1;
647 fEntering = false;
648 fEnteredDaughter = false; // Boundary not encountered, did not enter
649 fExiting = false;
650 fExitedMother = false; // Boundary not encountered, did not exit
651}
void RelocateWithinVolume(G4VPhysicalVolume *motherPhysical, const G4ThreeVector &localPoint) override

Referenced by G4VIntersectionLocator::AdjustmentOfFoundIntersection(), ComputeStep(), G4BrentLocator::EstimateIntersectionPoint(), G4MultiLevelLocator::EstimateIntersectionPoint(), G4SimpleLocator::EstimateIntersectionPoint(), and G4VIntersectionLocator::LocateGlobalPointWithinVolumeAndCheck().

◆ NetRotation()

G4RotationMatrix G4Navigator::NetRotation ( ) const
inline

◆ NetTranslation()

G4ThreeVector G4Navigator::NetTranslation ( ) const
inline

Computes and returns the local->global translation/rotation of current volume.

◆ operator=()

G4Navigator & G4Navigator::operator= ( const G4Navigator & )
delete

◆ PrintState()

void G4Navigator::PrintState ( ) const

Prints the internal state of the Navigator (for debugging). The level of detail is according to the verbosity.

Definition at line 1820 of file G4Navigator.cc.

1821{
1822 G4long oldcoutPrec = G4cout.precision(4);
1823 if( fVerbose >= 4 )
1824 {
1825 G4cout << "The current state of G4Navigator is: " << G4endl;
1826 G4cout << " ValidExitNormal= " << fValidExitNormal // << G4endl
1827 << " ExitNormal = " << fExitNormal // << G4endl
1828 << " Exiting = " << fExiting // << G4endl
1829 << " Entering = " << fEntering // << G4endl
1830 << " BlockedPhysicalVolume= " ;
1831 if (fBlockedPhysicalVolume==nullptr)
1832 {
1833 G4cout << "None";
1834 }
1835 else
1836 {
1837 G4cout << fBlockedPhysicalVolume->GetName();
1838 }
1839 G4cout << G4endl
1840 << " BlockedReplicaNo = " << fBlockedReplicaNo // << G4endl
1841 << " LastStepWasZero = " << fLastStepWasZero // << G4endl
1842 << G4endl;
1843 }
1844 if( ( 1 < fVerbose) && (fVerbose < 4) )
1845 {
1846 G4cout << G4endl; // Make sure to line up
1847 G4cout << std::setw(30) << " ExitNormal " << " "
1848 << std::setw( 5) << " Valid " << " "
1849 << std::setw( 9) << " Exiting " << " "
1850 << std::setw( 9) << " Entering" << " "
1851 << std::setw(15) << " Blocked:Volume " << " "
1852 << std::setw( 9) << " ReplicaNo" << " "
1853 << std::setw( 8) << " LastStepZero " << " "
1854 << G4endl;
1855 G4cout << "( " << std::setw(7) << fExitNormal.x()
1856 << ", " << std::setw(7) << fExitNormal.y()
1857 << ", " << std::setw(7) << fExitNormal.z() << " ) "
1858 << std::setw( 5) << fValidExitNormal << " "
1859 << std::setw( 9) << fExiting << " "
1860 << std::setw( 9) << fEntering << " ";
1861 if ( fBlockedPhysicalVolume == nullptr )
1862 { G4cout << std::setw(15) << "None"; }
1863 else
1864 { G4cout << std::setw(15)<< fBlockedPhysicalVolume->GetName(); }
1865 G4cout << std::setw( 9) << fBlockedReplicaNo << " "
1866 << std::setw( 8) << fLastStepWasZero << " "
1867 << G4endl;
1868 }
1869 if( fVerbose > 2 )
1870 {
1871 G4cout.precision(8);
1872 G4cout << " Current Localpoint = " << fLastLocatedPointLocal << G4endl;
1873 G4cout << " PreviousSftOrigin = " << fPreviousSftOrigin << G4endl;
1874 G4cout << " PreviousSafety = " << fPreviousSafety << G4endl;
1875 }
1876 G4cout.precision(oldcoutPrec);
1877}

Referenced by ComputeStep(), and LocateGlobalPointAndSetup().

◆ ResetHierarchyAndLocate()

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

Resets the geometrical hierarchy and searches for the volumes deepest in the hierarchy containing the point in the global coordinates space. The direction is used to check if a volume is entered. The search begin is the geometrical hierarchy at the location of the last located point, or the endpoint of the previous Step if SetGeometricallyLimitedStep() has been called immediately before.

Note
: In order to call this the geometry MUST be closed.
Parameters
[in]pointThe point in global coordinates system.
[in]directionThe normalised vector direction.
[in]hThe touchable history to be used for initialisation.
Returns
The pointer to the physical volume where point is located.

Reimplemented in G4MultiNavigator.

Definition at line 102 of file G4Navigator.cc.

105{
106 ResetState();
107 fHistory = *h.GetHistory();
109 fLastTriedStepComputation = false; // Redundant, but best
110 return LocateGlobalPointAndSetup(p, &direction, true, false);
111}
virtual void SetupHierarchy()
virtual void ResetState()
const G4NavigationHistory * GetHistory() const

Referenced by G4MultiNavigator::ResetHierarchyAndLocate().

◆ ResetStackAndState()

void G4Navigator::ResetStackAndState ( )
inline

Resets stack and minimum or navigator state machine necessary for reset as needed by LocalGlobalPointAndSetup().

Note
Does not perform clears, resizes, or reset fLastLocatedPointLocal.

Referenced by G4Navigator(), and LocateGlobalPointAndSetup().

◆ ResetState()

void G4Navigator::ResetState ( )
protectedvirtual

Utility method to reset the navigator state machine.

Reimplemented in G4MultiNavigator.

Definition at line 1261 of file G4Navigator.cc.

1262{
1263 fWasLimitedByGeometry = false;
1264 fEntering = false;
1265 fExiting = false;
1266 fLocatedOnEdge = false;
1267 fLastStepWasZero = false;
1268 fEnteredDaughter = false;
1269 fExitedMother = false;
1270 fPushed = false;
1271
1272 fValidExitNormal = false;
1273 fChangedGrandMotherRefFrame = false;
1274 fCalculatedExitNormal = false;
1275
1276 fExitNormal = G4ThreeVector(0,0,0);
1277 fGrandMotherExitNormal = G4ThreeVector(0,0,0);
1278 fExitNormalGlobalFrame = G4ThreeVector(0,0,0);
1279
1280 fPreviousSftOrigin = G4ThreeVector(0,0,0);
1281 fPreviousSafety = 0.0;
1282
1283 fNumberZeroSteps = 0;
1284
1285 fBlockedPhysicalVolume = nullptr;
1286 fBlockedReplicaNo = -1;
1287
1288 fLastLocatedPointLocal = G4ThreeVector( kInfinity, -kInfinity, 0.0 );
1289 fLocatedOutsideWorld = false;
1290
1291 fLastMotherPhys = nullptr;
1292}

Referenced by ResetHierarchyAndLocate().

◆ RestoreSavedState()

void G4Navigator::RestoreSavedState ( )
protected

Copy aspects of the state, to enable a non-state changing call to ComputeStep().

Definition at line 697 of file G4Navigator.cc.

698{
699 fExitNormal = fSaveState.sExitNormal;
700 fValidExitNormal = fSaveState.sValidExitNormal;
701 fExiting = fSaveState.sExiting;
702 fEntering = fSaveState.sEntering;
703
704 fBlockedPhysicalVolume = fSaveState.spBlockedPhysicalVolume;
705 fBlockedReplicaNo = fSaveState.sBlockedReplicaNo;
706
707 fLastStepWasZero = (fSaveState.sLastStepWasZero != 0);
708
709 fLocatedOutsideWorld = fSaveState.sLocatedOutsideWorld;
710 fLastLocatedPointLocal = fSaveState.sLastLocatedPointLocal;
711 fEnteredDaughter = fSaveState.sEnteredDaughter;
712 fExitedMother = fSaveState.sExitedMother;
713 fWasLimitedByGeometry = fSaveState.sWasLimitedByGeometry;
714
715 // The 'expected' behaviour is to restore these too (fix 2014.05.26)
716 fPreviousSftOrigin = fSaveState.sPreviousSftOrigin;
717 fPreviousSafety = fSaveState.sPreviousSafety;
718}

Referenced by CheckNextStep().

◆ SetExternalNavigation()

void G4Navigator::SetExternalNavigation ( G4VExternalNavigation * externalNav)

Definition at line 2113 of file G4Navigator.cc.

2114{
2115 fpExternalNav = externalNav;
2116 fpSafetyCalculator->SetExternalNavigation(externalNav);
2117}

◆ SetGeometricallyLimitedStep()

void G4Navigator::SetGeometricallyLimitedStep ( )
inline

Informs the navigator that the previous Step calculated by the geometry was taken in its entirety.

◆ SetPushVerbosity()

void G4Navigator::SetPushVerbosity ( G4bool mode)
inline

◆ SetSavedState()

void G4Navigator::SetSavedState ( )
protected

Saves the state: fValidExitNormal, fExitNormal, fExiting, fEntering, fBlockedPhysicalVolume, fBlockedReplicaNo, fLastStepWasZero, fLastLocatedPointLocal, fLocatedOutsideWorld, fEnteredDaughter, fExitedMother, fPreviousSftOrigin, fPreviousSafety.

Definition at line 663 of file G4Navigator.cc.

664{
665 // Note: the state of dependent objects is not currently saved.
666 // ( This means that the full state is changed by calls between
667 // SetSavedState() and RestoreSavedState();
668
669 fSaveState.sExitNormal = fExitNormal;
670 fSaveState.sValidExitNormal = fValidExitNormal;
671 fSaveState.sExiting = fExiting;
672 fSaveState.sEntering = fEntering;
673
674 fSaveState.spBlockedPhysicalVolume = fBlockedPhysicalVolume;
675 fSaveState.sBlockedReplicaNo = fBlockedReplicaNo;
676
677 fSaveState.sLastStepWasZero = static_cast<G4int>(fLastStepWasZero);
678
679 fSaveState.sLocatedOutsideWorld = fLocatedOutsideWorld;
680 fSaveState.sLastLocatedPointLocal = fLastLocatedPointLocal;
681 fSaveState.sEnteredDaughter = fEnteredDaughter;
682 fSaveState.sExitedMother = fExitedMother;
683 fSaveState.sWasLimitedByGeometry = fWasLimitedByGeometry;
684
685 // Even the safety sphere - if you want to change it do it explicitly!
686 //
687 fSaveState.sPreviousSftOrigin = fPreviousSftOrigin;
688 fSaveState.sPreviousSafety = fPreviousSafety;
689}

Referenced by CheckNextStep().

◆ SetupHierarchy()

void G4Navigator::SetupHierarchy ( )
protectedvirtual

Renavigates & resets hierarchy described by the current history: Resets volumes and recomputes transforms and/or solids of replicated/parameterised volumes.

Reimplemented in G4MultiNavigator.

Definition at line 1302 of file G4Navigator.cc.

1303{
1304 const auto depth = (G4int)fHistory.GetDepth();
1305 for ( auto i = 1; i <= depth; ++i )
1306 {
1307 switch ( fHistory.GetVolumeType(i) )
1308 {
1309 case kNormal:
1310 case kExternal:
1311 break;
1312 case kReplica:
1313 freplicaNav.ComputeTransformation(fHistory.GetReplicaNo(i), fHistory.GetVolume(i));
1314 break;
1315 case kParameterised:
1316 G4VPhysicalVolume* current = fHistory.GetVolume(i);
1317 G4int replicaNo = fHistory.GetReplicaNo(i);
1318 G4VPVParameterisation* pParam = current->GetParameterisation();
1319 G4VSolid* pSolid = pParam->ComputeSolid(replicaNo, current);
1320
1321 // Set up dimensions & transform in solid/physical volume
1322 //
1323 pSolid->ComputeDimensions(pParam, replicaNo, current);
1324 pParam->ComputeTransformation(replicaNo, current);
1325
1326 G4TouchableHistory* pTouchable = nullptr;
1327 if( pParam->IsNested() )
1328 {
1329 pTouchable= new G4TouchableHistory( fHistory );
1330 pTouchable->MoveUpHistory(); // Move up to the parent level
1331 // Adequate only if Nested at the Branch level (last)
1332 // To extend to other cases:
1333 // pTouchable->MoveUpHistory(cdepth-i-1);
1334 // Move to the parent level of *Current* level
1335 // Could replace this line and constructor with a revised
1336 // c-tor for History(levels to drop)
1337 }
1338 // Set up the correct solid and material in Logical Volume
1339 //
1340 G4LogicalVolume* pLogical = current->GetLogicalVolume();
1341 pLogical->SetSolid( pSolid );
1342 pLogical->UpdateMaterial( pParam ->
1343 ComputeMaterial(replicaNo, current, pTouchable) );
1344 delete pTouchable;
1345 break;
1346 }
1347 }
1348}
G4int MoveUpHistory(G4int num_levels=1)
virtual G4bool IsNested() const
virtual G4VPVParameterisation * GetParameterisation() const =0

Referenced by ResetHierarchyAndLocate().

◆ SetVerboseLevel()

void G4Navigator::SetVerboseLevel ( G4int level)
inline

Referenced by GetGlobalExitNormal().

◆ SetVoxelNavigation()

void G4Navigator::SetVoxelNavigation ( G4VoxelNavigation * voxelNav)

Definition at line 2088 of file G4Navigator.cc.

2089{
2090 delete fpvoxelNav;
2091 fpvoxelNav = voxelNav;
2092}

◆ SetWorldVolume()

void G4Navigator::SetWorldVolume ( G4VPhysicalVolume * pWorld)
inline

Sets the world (topmost) volume. This must be positioned at the origin (0,0,0) and unrotated.

Referenced by G4MultiNavigator::G4MultiNavigator(), G4VisManager::G4VisSubThread(), G4TransportationManager::GetNavigator(), G4TransportationManager::GetNavigator(), and G4MultiNavigator::PrepareNavigators().

◆ SeverityOfZeroStepping()

G4int G4Navigator::SeverityOfZeroStepping ( G4int * noZeroSteps) const
inline

Reports on severity of error and number of zero steps, in case Navigator is stuck and is returning zero steps. Values: 1 (small problem), 5 (correcting), 9 (ready to abandon), 10 (abandoned)

Parameters
[in,out]noZeroStepsReturns the number of zero steps in case pointer is not null.
Returns
The error severity.

◆ VolumeType()

EVolume G4Navigator::VolumeType ( const G4VPhysicalVolume * pVol) const
inlineprotected

Characterises the type of volume - normal/replicated/parameterised.

◆ operator<<

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

Definition at line 2013 of file G4Navigator.cc.

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

Member Data Documentation

◆ fEnteredDaughter

G4bool G4Navigator::fEnteredDaughter
protected

A memory of whether in this Step a daughter volume is entered (set in Compute & Locate). After Compute: it expects to enter a daughter After Locate: it has entered a daughter.

Definition at line 560 of file G4Navigator.hh.

Referenced by ComputeSafety(), ComputeStep(), LocateGlobalPointAndSetup(), LocateGlobalPointWithinVolume(), ResetState(), RestoreSavedState(), and SetSavedState().

◆ fExitedMother

G4bool G4Navigator::fExitedMother
protected

A similar memory whether the Step exited current "mother" volume completely, not entering daughter.

Definition at line 564 of file G4Navigator.hh.

Referenced by ComputeSafety(), ComputeStep(), GetLocalExitNormal(), LocateGlobalPointAndSetup(), LocateGlobalPointWithinVolume(), ResetState(), RestoreSavedState(), and SetSavedState().

◆ fHistory

G4NavigationHistory G4Navigator::fHistory
protected

Transformation and history of the current path through the geometrical hierarchy.

Definition at line 543 of file G4Navigator.hh.

Referenced by ComputeSafety(), ComputeStep(), G4Navigator(), GetGlobalExitNormal(), GetLocalExitNormal(), LocateGlobalPointAndSetup(), LocateGlobalPointWithinVolume(), ResetHierarchyAndLocate(), and SetupHierarchy().

◆ fLastStepEndPointLocal

G4ThreeVector G4Navigator::fLastStepEndPointLocal
protected

Position of the end-point of the last call to ComputeStep() in last local coordinates.

Definition at line 551 of file G4Navigator.hh.

Referenced by ComputeStep(), G4Navigator(), and GetLocalExitNormal().

◆ fMinStep

G4double G4Navigator::fMinStep
protected

Definition at line 536 of file G4Navigator.hh.

Referenced by ComputeStep(), and G4Navigator().

◆ fSqTol

G4double G4Navigator::fSqTol
protected

Definition at line 536 of file G4Navigator.hh.

Referenced by ComputeStep(), G4Navigator(), and GetGlobalExitNormal().

◆ fStepEndPoint

G4ThreeVector G4Navigator::fStepEndPoint
protected

Endpoint of last ComputeStep(). Can be used for optimisation (e.g. when computing safety).

Definition at line 547 of file G4Navigator.hh.

Referenced by ComputeSafety(), ComputeStep(), G4Navigator(), and GetGlobalExitNormal().

◆ fVerbose

◆ fWasLimitedByGeometry

◆ kCarTolerance

G4double G4Navigator::kCarTolerance
protected

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