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

#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)
 
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)
 
G4VExternalNavigationGetExternalNavigation () const
 
void SetExternalNavigation (G4VExternalNavigation *externalNav)
 
G4NavigatorClone () const
 
G4ThreeVector GetLastStepEndPoint () const
 
void InformLastStep (G4double lastStep, G4bool entersDaughtVol, G4bool exitsMotherVol)
 
void SetVoxelNavigation (G4VoxelNavigation *voxelNav)
 

Protected Member Functions

void SetSavedState ()
 
void RestoreSavedState ()
 
virtual void ResetState ()
 
G4ThreeVector ComputeLocalPoint (const G4ThreeVector &rGlobPoint) const
 
G4ThreeVector ComputeLocalAxis (const G4ThreeVector &pVec) const
 
EVolume VolumeType (const G4VPhysicalVolume *pVol) const
 
EVolume CharacteriseDaughters (const G4LogicalVolume *pLog) const
 
G4int GetDaughtersRegularStructureId (const G4LogicalVolume *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

Definition at line 71 of file G4Navigator.hh.

Constructor & Destructor Documentation

◆ G4Navigator() [1/2]

G4Navigator::G4Navigator ( )

Definition at line 53 of file G4Navigator.cc.

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

◆ G4Navigator() [2/2]

G4Navigator::G4Navigator ( const G4Navigator )
delete

◆ ~G4Navigator()

G4Navigator::~G4Navigator ( )
virtual

Definition at line 88 of file G4Navigator.cc.

89{
90 delete fpVoxelSafety;
91 delete fpExternalNav;
92#ifdef ALTERNATIVE_VOXEL_NAV
93 delete fpvoxelNav;
94#endif
95}

Member Function Documentation

◆ Activate()

void G4Navigator::Activate ( G4bool  flag)
inline

◆ CharacteriseDaughters()

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

◆ CheckMode()

void G4Navigator::CheckMode ( G4bool  mode)
inline

◆ CheckNextStep()

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

Definition at line 1231 of file G4Navigator.cc.

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

Referenced by G4SafetyHelper::CheckNextStep().

◆ CheckOverlapsIterative()

G4bool G4Navigator::CheckOverlapsIterative ( G4VPhysicalVolume vol)
protected

Definition at line 2103 of file G4Navigator.cc.

2104{
2105 // Check and report overlaps
2106 //
2107 G4bool foundOverlap = false;
2108 G4int nPoints = 300000, ntrials = 9, numOverlaps = 5;
2109 G4double trialLength = 1.0 * CLHEP::centimeter;
2110 while ( ntrials-- > 0 && !foundOverlap )
2111 {
2112 if ( fVerbose > 1 )
2113 {
2114 G4cout << " ** Running overlap checks in volume "
2115 << vol->GetName()
2116 << " with length = " << trialLength << G4endl;
2117 }
2118 foundOverlap = vol->CheckOverlaps(nPoints, trialLength,
2119 fVerbose, numOverlaps);
2120 trialLength *= 0.1;
2121 if ( trialLength <= 1.0e-5 ) { numOverlaps= 1;}
2122 }
2123 return foundOverlap;
2124}
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()

◆ ComputeLocalAxis()

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

Referenced by ComputeStep().

◆ ComputeLocalPoint()

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

◆ ComputeSafety()

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

Reimplemented in G4MultiNavigator, and G4ErrorPropagationNavigator.

Definition at line 1784 of file G4Navigator.cc.

1787{
1788#ifdef G4DEBUG_NAVIGATION
1789 G4int oldcoutPrec = G4cout.precision(8);
1790 if( fVerbose > 0 )
1791 {
1792 G4cout << "*** G4Navigator::ComputeSafety: ***" << G4endl
1793 << " Called at point: " << pGlobalpoint << G4endl;
1794
1795 G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
1796 G4cout << " Volume = " << motherPhysical->GetName()
1797 << " - Maximum length = " << pMaxLength << G4endl;
1798 if( fVerbose >= 4 )
1799 {
1800 G4cout << " ----- Upon entering Compute Safety:" << G4endl;
1801 PrintState();
1802 }
1803 }
1804#endif
1805
1806 G4double distEndpointSq = (pGlobalpoint-fStepEndPoint).mag2();
1807 G4bool stayedOnEndpoint = distEndpointSq < sqr(kCarTolerance);
1808 G4bool endpointOnSurface = fEnteredDaughter || fExitedMother;
1809
1810 if( endpointOnSurface && stayedOnEndpoint )
1811 {
1812#ifdef G4DEBUG_NAVIGATION
1813 if( fVerbose >= 2 )
1814 {
1815 G4cout << " G4Navigator::ComputeSafety() finds that point - "
1816 << pGlobalpoint << " - is on surface " << G4endl;
1817 if( fEnteredDaughter ) { G4cout << " entered new daughter volume"; }
1818 if( fExitedMother ) { G4cout << " and exited previous volume."; }
1819 G4cout << G4endl;
1820 G4cout << " EndPoint was = " << fStepEndPoint << G4endl;
1821 G4cout << " ---- Exiting ComputeSafety " << G4endl;
1822 PrintState();
1823 G4cout << " Returned value of Safety is zero " << G4endl;
1824 G4cout.precision(oldcoutPrec);
1825 }
1826#endif
1827 return 0.0;
1828 }
1829
1830 G4double newSafety = 0.0;
1831
1832 if (keepState) { SetSavedState(); }
1833
1834 // Pseudo-relocate to this point (updates voxel information only)
1835 //
1836 LocateGlobalPointWithinVolume( pGlobalpoint );
1837 // --->> DANGER: Side effects on sub-navigator voxel information <<---
1838 // Could be replaced again by 'granular' calls to sub-navigator
1839 // locates (similar side-effects, but faster.
1840 // Solutions:
1841 // 1) Re-locate (to where?)
1842 // 2) Insure that the methods using (G4ComputeStep?)
1843 // does a relocation (if information is disturbed only ?)
1844
1845#ifdef G4DEBUG_NAVIGATION
1846 if( fVerbose >= 2 )
1847 {
1848 G4cout << " G4Navigator::ComputeSafety() relocates-in-volume to point: "
1849 << pGlobalpoint << G4endl;
1850 }
1851#endif
1852 G4VPhysicalVolume* motherPhysical = fHistory.GetTopVolume();
1853 G4LogicalVolume* motherLogical = motherPhysical->GetLogicalVolume();
1854 G4SmartVoxelHeader* pVoxelHeader = motherLogical->GetVoxelHeader();
1855 G4ThreeVector localPoint = ComputeLocalPoint(pGlobalpoint);
1856
1858 {
1859 switch(CharacteriseDaughters(motherLogical))
1860 {
1861 case kNormal:
1862 if ( pVoxelHeader )
1863 {
1864 newSafety = fpVoxelSafety->ComputeSafety(localPoint,
1865 *motherPhysical, pMaxLength);
1866 // = VoxelNav().ComputeSafety(localPoint,fHistory,pMaxLength); // - Old method
1867 }
1868 else
1869 {
1870 newSafety=fnormalNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1871 }
1872 break;
1873 case kParameterised:
1874 if( GetDaughtersRegularStructureId(motherLogical) != 1 )
1875 {
1876 newSafety=fparamNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1877 }
1878 else // Regular structure
1879 {
1880 newSafety=fregularNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1881 }
1882 break;
1883 case kReplica:
1884 G4Exception("G4Navigator::ComputeSafety()", "GeomNav0001",
1885 FatalException, "Not applicable for replicated volumes.");
1886 break;
1887 case kExternal:
1888 newSafety = fpExternalNav->ComputeSafety(localPoint, fHistory,
1889 pMaxLength);
1890 break;
1891 }
1892 }
1893 else
1894 {
1895 newSafety = freplicaNav.ComputeSafety(pGlobalpoint, localPoint,
1896 fHistory, pMaxLength);
1897 }
1898
1899 if (keepState)
1900 {
1902 // This now overwrites the values of the Safety 'sphere' (correction)
1903 }
1904
1905 // Remember last safety origin & value
1906 //
1907 // We overwrite the Safety 'sphere' - keeping old behaviour
1908 fPreviousSftOrigin = pGlobalpoint;
1909 fPreviousSafety = newSafety;
1910
1911#ifdef G4DEBUG_NAVIGATION
1912 if( fVerbose > 1 )
1913 {
1914 G4cout << " ---- Exiting ComputeSafety " << G4endl;
1915 if( fVerbose > 2 ) { PrintState(); }
1916 G4cout << " Returned value of Safety = " << newSafety << G4endl;
1917 }
1918 G4cout.precision(oldcoutPrec);
1919#endif
1920
1921 return newSafety;
1922}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
G4SmartVoxelHeader * GetVoxelHeader() const
EVolume GetTopVolumeType() const
G4VPhysicalVolume * GetTopVolume() const
G4bool fExitedMother
Definition: G4Navigator.hh:404
G4bool fEnteredDaughter
Definition: G4Navigator.hh:398
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
Definition: G4Navigator.cc:600
EVolume CharacteriseDaughters(const G4LogicalVolume *pLog) const
G4ThreeVector ComputeLocalPoint(const G4ThreeVector &rGlobPoint) const
void PrintState() const
G4int GetDaughtersRegularStructureId(const G4LogicalVolume *pLv) const
G4NavigationHistory fHistory
Definition: G4Navigator.hh:384
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)
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4NavigationHistory &history, const G4double pMaxLength=DBL_MAX)=0
G4LogicalVolume * GetLogicalVolume() const
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 G4ImportanceProcess::AlongStepGetPhysicalInteractionLength(), G4WeightCutOffProcess::AlongStepGetPhysicalInteractionLength(), G4WeightWindowProcess::AlongStepGetPhysicalInteractionLength(), G4ParallelWorldProcess::AlongStepGetPhysicalInteractionLength(), G4ParallelWorldScoringProcess::AlongStepGetPhysicalInteractionLength(), G4Transportation::AlongStepGetPhysicalInteractionLength(), G4FastSimulationManagerProcess::AlongStepGetPhysicalInteractionLength(), G4ErrorPropagationNavigator::ComputeSafety(), G4SafetyHelper::ComputeSafety(), and G4PathFinder::DoNextCurvedStep().

◆ ComputeStep()

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

Reimplemented in G4ErrorPropagationNavigator, and G4MultiNavigator.

Definition at line 756 of file G4Navigator.cc.

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

Referenced by G4Transportation::AlongStepGetPhysicalInteractionLength(), CheckNextStep(), G4AdjointPrimaryGenerator::ComputeAccumulatedDepthVectorAlongBackRay(), and G4ErrorPropagationNavigator::ComputeStep().

◆ CreateGRSSolid()

G4GRSSolid * G4Navigator::CreateGRSSolid ( ) const
inline

◆ CreateGRSVolume()

G4GRSVolume * G4Navigator::CreateGRSVolume ( ) const
inline

◆ CreateTouchableHistory() [1/2]

◆ CreateTouchableHistory() [2/2]

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

◆ CreateTouchableHistoryHandle()

G4TouchableHistoryHandle G4Navigator::CreateTouchableHistoryHandle ( ) const
virtual

Reimplemented in G4MultiNavigator.

Definition at line 1928 of file G4Navigator.cc.

1929{
1931}
G4ReferenceCountedHandle< G4TouchableHistory > G4TouchableHistoryHandle
G4TouchableHistory * CreateTouchableHistory() const

Referenced by G4VIntersectionLocator::LocateGlobalPointWithinVolumeAndCheck().

◆ EnableBestSafety()

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

◆ EnteredDaughterVolume()

G4bool G4Navigator::EnteredDaughterVolume ( ) const
inline

◆ ExitedMotherVolume()

G4bool G4Navigator::ExitedMotherVolume ( ) const
inline

◆ GetCurrentLocalCoordinate()

G4ThreeVector G4Navigator::GetCurrentLocalCoordinate ( ) const
inline

◆ GetDaughtersRegularStructureId()

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

◆ GetExternalNavigation()

◆ GetGlobalExitNormal()

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

Reimplemented in G4MultiNavigator, and G4ErrorPropagationNavigator.

Definition at line 1613 of file G4Navigator.cc.

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

◆ GetGlobalToLocalTransform()

const G4AffineTransform & G4Navigator::GetGlobalToLocalTransform ( ) const
inline

◆ GetLastStepEndPoint()

G4ThreeVector G4Navigator::GetLastStepEndPoint ( ) const
inline

Definition at line 304 of file G4Navigator.hh.

304{ return fStepEndPoint;}

◆ GetLocalExitNormal()

G4ThreeVector G4Navigator::GetLocalExitNormal ( G4bool valid)
virtual

Reimplemented in G4MultiNavigator.

Definition at line 1360 of file G4Navigator.cc.

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

◆ GetLocalExitNormalAndCheck()

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

Reimplemented in G4MultiNavigator.

Definition at line 1577 of file G4Navigator.cc.

1584{
1585#ifdef G4DEBUG_NAVIGATION
1586 // Check Current point against expected 'local' value
1587 //
1588 if ( fLastTriedStepComputation )
1589 {
1590 G4ThreeVector ExpectedBoundaryPointLocal;
1591
1592 const G4AffineTransform& GlobalToLocal = GetGlobalToLocalTransform();
1593 ExpectedBoundaryPointLocal =
1594 GlobalToLocal.TransformPoint( ExpectedBoundaryPointGlobal );
1595
1596 // Add here: Comparison against expected position,
1597 // i.e. the endpoint of ComputeStep
1598 }
1599#endif
1600
1601 return GetLocalExitNormal( pValid );
1602}
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 
)

Definition at line 1526 of file G4Navigator.cc.

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

◆ GetVerboseLevel()

G4int G4Navigator::GetVerboseLevel ( ) const
inline

Referenced by GetGlobalExitNormal().

◆ GetWorldVolume()

◆ InformLastStep()

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

Definition at line 2213 of file G4Navigator.cc.

2214{
2215 G4bool zeroStep = ( lastStep == 0.0 );
2216 fLocatedOnEdge = fLastStepWasZero && zeroStep;
2217 fLastStepWasZero = zeroStep;
2218
2219 fExiting = exitsMotherVol;
2220 fEntering = entersDaughtVol;
2221}

◆ IsActive()

G4bool G4Navigator::IsActive ( ) const
inline

◆ IsCheckModeActive()

G4bool G4Navigator::IsCheckModeActive ( ) const
inline

◆ LocateGlobalPointAndSetup()

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

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 = 0;
141 G4ThreeVector localPoint, globalDirection;
142 EInside insideCode;
143
144 G4bool considerDirection = pGlobalDirection && ((!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() )
189 {
190 fBlockedPhysicalVolume = fHistory.GetTopVolume();
191 fBlockedReplicaNo = fHistory.GetTopReplicaNo();
193 }
194 else
195 {
196 fLastLocatedPointLocal = localPoint;
197 fLocatedOutsideWorld = true;
198 fBlockedPhysicalVolume = 0; // 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 solidResponse = "-kOutside-";
310 else if (insideCode == kSurface)
311 solidResponse = "-kSurface-";
312 G4cout << "*** G4Navigator::LocateGlobalPointAndSetup(): ***" << G4endl
313 << " Invoked Inside() for solid: " << targetSolid->GetName()
314 << ". Solid replied: " << solidResponse << G4endl
315 << " For local point p: " << localPoint << G4endl;
316 }
317#endif
318 }
319 else
320 {
321 if( topVolumeType == kReplica )
322 {
323 insideCode = freplicaNav.BackLocate(fHistory, globalPoint, localPoint,
324 fExiting, notKnownContained);
325 // !CARE! if notKnownContained returns false then the point is within
326 // the containing placement volume of the replica(s). If insidecode
327 // will result in the history being backed up one level, then the
328 // local point returned is the point in the system of this new level
329 }
330 else
331 {
332 targetSolid = fHistory.GetTopVolume()->GetLogicalVolume()->GetSolid();
333 localPoint = fHistory.GetTopTransform().TransformPoint(globalPoint);
334 G4ThreeVector localDirection =
335 fHistory.GetTopTransform().TransformAxis(globalDirection);
336 insideCode = fpExternalNav->Inside(targetSolid, localPoint, localDirection);
337 }
338 }
339
340 // Point is inside current volume, break out of the loop
341 if ( insideCode == kInside )
342 break;
343
344 // Point is outside current volume, move up a level in the hierarchy
345 if ( insideCode == kOutside )
346 {
347 ++noLevelsExited;
348
349 // Exiting world volume
350 if ( fHistory.GetDepth() == 0 )
351 {
352 fLocatedOutsideWorld = true;
353 fLastLocatedPointLocal = localPoint;
354 return nullptr;
355 }
356
357 fBlockedPhysicalVolume = fHistory.GetTopVolume();
358 fBlockedReplicaNo = fHistory.GetTopReplicaNo();
360 fExiting = false;
361
362 if( noLevelsExited > 1 )
363 {
364 // The first transformation was done by the sub-navigator
365 //
366 if(const auto *mRot = fBlockedPhysicalVolume->GetRotation())
367 {
368 fGrandMotherExitNormal *= (*mRot).inverse();
369 fChangedGrandMotherRefFrame = true;
370 }
371 }
372 continue;
373 }
374
375 // Point is on the surface of a volume
376 G4bool isExiting = fExiting;
377 if( (!fExiting) && considerDirection )
378 {
379 // Figure out whether we are exiting this level's volume
380 // by using the direction
381 //
382 G4bool directionExiting = false;
383 G4ThreeVector localDirection =
384 fHistory.GetTopTransform().TransformAxis(globalDirection);
385
386 // Make sure localPoint in correct reference frame
387 // ( Was it already correct ? How ? )
388 //
389 localPoint= fHistory.GetTopTransform().TransformPoint(globalPoint);
391 {
392 G4ThreeVector normal = targetSolid->SurfaceNormal(localPoint);
393 directionExiting = normal.dot(localDirection) > 0.0;
394 isExiting = isExiting || directionExiting;
395 }
396 }
397
398 // Point is on a surface, but no longer exiting, break out of the loop
399 if ( !isExiting )
400 break;
401
402 ++noLevelsExited;
403
404 // Point is on the outer surface, leaving world volume
405 if ( fHistory.GetDepth() == 0 )
406 {
407 fLocatedOutsideWorld = true;
408 fLastLocatedPointLocal = localPoint;
409 return nullptr;
410 }
411
412 // Point is still on a surface, but exited a volume not necessarily convex
413 fValidExitNormal = false;
414 fBlockedPhysicalVolume = fHistory.GetTopVolume();
415 fBlockedReplicaNo = fHistory.GetTopReplicaNo();
417
418 if( noLevelsExited > 1 )
419 {
420 // The first transformation was done by the sub-navigator
421 //
422 const G4RotationMatrix* mRot =
423 fBlockedPhysicalVolume->GetRotation();
424 if( mRot )
425 {
426 fGrandMotherExitNormal *= (*mRot).inverse();
427 fChangedGrandMotherRefFrame = true;
428 }
429 }
430 } // END while (notKnownContained)
431 //
432 // Search downwards until deepest containing volume found,
433 // blocking fBlockedPhysicalVolume/BlockedReplicaNum
434 //
435 // 3 Cases:
436 //
437 // o Parameterised daughters
438 // =>Must be one G4PVParameterised daughter & voxels
439 // o Positioned daughters & voxels
440 // o Positioned daughters & no voxels
441
442 noResult = true; // noResult should be renamed to
443 // something like enteredLevel, as that is its meaning.
444 do
445 {
446 // Determine `type' of current mother volume
447 //
448 targetPhysical = fHistory.GetTopVolume();
449 if (!targetPhysical) { break; }
450 targetLogical = targetPhysical->GetLogicalVolume();
451 switch( CharacteriseDaughters(targetLogical) )
452 {
453 case kNormal:
454 if ( targetLogical->GetVoxelHeader() ) // use optimised navigation
455 {
456 noResult = GetVoxelNavigator().LevelLocate(fHistory,
457 fBlockedPhysicalVolume,
458 fBlockedReplicaNo,
459 globalPoint,
460 pGlobalDirection,
461 considerDirection,
462 localPoint);
463 }
464 else // do not use optimised navigation
465 {
466 noResult = fnormalNav.LevelLocate(fHistory,
467 fBlockedPhysicalVolume,
468 fBlockedReplicaNo,
469 globalPoint,
470 pGlobalDirection,
471 considerDirection,
472 localPoint);
473 }
474 break;
475 case kReplica:
476 noResult = freplicaNav.LevelLocate(fHistory,
477 fBlockedPhysicalVolume,
478 fBlockedReplicaNo,
479 globalPoint,
480 pGlobalDirection,
481 considerDirection,
482 localPoint);
483 break;
484 case kParameterised:
485 if( GetDaughtersRegularStructureId(targetLogical) != 1 )
486 {
487 noResult = fparamNav.LevelLocate(fHistory,
488 fBlockedPhysicalVolume,
489 fBlockedReplicaNo,
490 globalPoint,
491 pGlobalDirection,
492 considerDirection,
493 localPoint);
494 }
495 else // Regular structure
496 {
497 noResult = fregularNav.LevelLocate(fHistory,
498 fBlockedPhysicalVolume,
499 fBlockedReplicaNo,
500 globalPoint,
501 pGlobalDirection,
502 considerDirection,
503 localPoint);
504 }
505 break;
506 case kExternal:
507 noResult = fpExternalNav->LevelLocate(fHistory,
508 fBlockedPhysicalVolume,
509 fBlockedReplicaNo,
510 globalPoint,
511 pGlobalDirection,
512 considerDirection,
513 localPoint);
514 break;
515 }
516
517 // LevelLocate returns true if it finds a daughter volume
518 // in which globalPoint is inside (or on the surface).
519
520 if ( noResult )
521 {
522 // Entering a daughter after ascending
523 //
524 // The blocked volume is no longer valid - it was for another level
525 //
526 fBlockedPhysicalVolume = nullptr;
527 fBlockedReplicaNo = -1;
528
529 // fEntering should be false -- else blockedVolume is assumed good.
530 // fEnteredDaughter is used for ExitNormal
531 //
532 fEntering = false;
533 fEnteredDaughter = true;
534
535 if( fExitedMother )
536 {
537 G4VPhysicalVolume* enteredPhysical = fHistory.GetTopVolume();
538 const G4RotationMatrix* mRot = enteredPhysical->GetRotation();
539 if( mRot )
540 {
541 // Go deeper, i.e. move 'down' in the hierarchy
542 // Apply direct rotation, not inverse
543 //
544 fGrandMotherExitNormal *= (*mRot);
545 fChangedGrandMotherRefFrame= true;
546 }
547 }
548
549#ifdef G4DEBUG_NAVIGATION
550 if( fVerbose > 2 )
551 {
552 G4VPhysicalVolume* enteredPhysical = fHistory.GetTopVolume();
553 G4cout << "*** G4Navigator::LocateGlobalPointAndSetup() ***" << G4endl;
554 G4cout << " Entering volume: " << enteredPhysical->GetName()
555 << G4endl;
556 }
557#endif
558 }
559 } while (noResult); // Loop checking, 07.10.2016, J.Apostolakis
560
561 fLastLocatedPointLocal = localPoint;
562
563#ifdef G4VERBOSE
564 if( fVerbose >= 4 )
565 {
566 G4long oldcoutPrec = G4cout.precision(8);
567 G4String curPhysVol_Name("None");
568 if (targetPhysical) { curPhysVol_Name = targetPhysical->GetName(); }
569 G4cout << " Return value = new volume = " << curPhysVol_Name << G4endl;
570 G4cout << " ----- Upon exiting:" << G4endl;
571 PrintState();
572 if( fVerbose >= 5 )
573 {
574 G4cout << "Upon exiting LocateGlobalPointAndSetup():" << G4endl;
575 G4cout << " History = " << G4endl << fHistory << G4endl << G4endl;
576 }
577 G4cout.precision(oldcoutPrec);
578 }
579#endif
580
581 fLocatedOutsideWorld = false;
582
583 return targetPhysical;
584}
long G4long
Definition: G4Types.hh:87
double dot(const Hep3Vector &) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
void UpdateMaterial(G4Material *pMaterial)
void NewLevel(G4VPhysicalVolume *pNewMother, EVolume vType=kNormal, G4int nReplica=-1)
G4int GetTopReplicaNo() const
G4bool fWasLimitedByGeometry
Definition: G4Navigator.hh:408
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 EInside Inside(const G4VSolid *solid, const G4ThreeVector &position, const G4ThreeVector &direction)
virtual G4bool LevelLocate(G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint)=0
virtual void SetCopyNo(G4int CopyNo)=0
virtual G4VPVParameterisation * GetParameterisation() const =0
virtual G4bool LevelLocate(G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint)
EVolume
Definition: geomdefs.hh:83

Referenced by G4AdjointPrimaryGenerator::ComputeAccumulatedDepthVectorAlongBackRay(), ComputeStep(), G4TheRayTracer::CreateBitMap(), G4VFieldModel::DescribeYourselfTo(), G4TrajectoryDrawByOriginVolume::Draw(), G4TrajectoryOriginVolumeFilter::Evaluate(), G4SafetyHelper::Locate(), G4VIntersectionLocator::LocateGlobalPointWithinVolumeAndCheck(), ResetHierarchyAndLocate(), and G4SteppingManager::SetInitialStep().

◆ LocateGlobalPointAndUpdateTouchable() [1/2]

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

◆ LocateGlobalPointAndUpdateTouchable() [2/2]

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

◆ LocateGlobalPointAndUpdateTouchableHandle()

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

◆ LocateGlobalPointWithinVolume()

void G4Navigator::LocateGlobalPointWithinVolume ( const G4ThreeVector position)
virtual

Reimplemented in G4MultiNavigator.

Definition at line 600 of file G4Navigator.cc.

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

Referenced by G4VIntersectionLocator::AdjustmentOfFoundIntersection(), G4TransportationWithMsc::AlongStepGetPhysicalInteractionLength(), ComputeSafety(), ComputeStep(), G4PropagatorInField::ComputeStep(), G4BrentLocator::EstimateIntersectionPoint(), G4MultiLevelLocator::EstimateIntersectionPoint(), G4SimpleLocator::EstimateIntersectionPoint(), G4VIntersectionLocator::LocateGlobalPointWithinVolumeAndCheck(), G4Transportation::PostStepDoIt(), and G4SafetyHelper::ReLocateWithinVolume().

◆ NetRotation()

G4RotationMatrix G4Navigator::NetRotation ( ) const
inline

◆ NetTranslation()

G4ThreeVector G4Navigator::NetTranslation ( ) const
inline

◆ operator=()

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

◆ PrintState()

void G4Navigator::PrintState ( ) const

Definition at line 1937 of file G4Navigator.cc.

1938{
1939 G4long oldcoutPrec = G4cout.precision(4);
1940 if( fVerbose >= 4 )
1941 {
1942 G4cout << "The current state of G4Navigator is: " << G4endl;
1943 G4cout << " ValidExitNormal= " << fValidExitNormal // << G4endl
1944 << " ExitNormal = " << fExitNormal // << G4endl
1945 << " Exiting = " << fExiting // << G4endl
1946 << " Entering = " << fEntering // << G4endl
1947 << " BlockedPhysicalVolume= " ;
1948 if (fBlockedPhysicalVolume==0)
1949 {
1950 G4cout << "None";
1951 }
1952 else
1953 {
1954 G4cout << fBlockedPhysicalVolume->GetName();
1955 }
1956 G4cout << G4endl
1957 << " BlockedReplicaNo = " << fBlockedReplicaNo // << G4endl
1958 << " LastStepWasZero = " << fLastStepWasZero // << G4endl
1959 << G4endl;
1960 }
1961 if( ( 1 < fVerbose) && (fVerbose < 4) )
1962 {
1963 G4cout << G4endl; // Make sure to line up
1964 G4cout << std::setw(30) << " ExitNormal " << " "
1965 << std::setw( 5) << " Valid " << " "
1966 << std::setw( 9) << " Exiting " << " "
1967 << std::setw( 9) << " Entering" << " "
1968 << std::setw(15) << " Blocked:Volume " << " "
1969 << std::setw( 9) << " ReplicaNo" << " "
1970 << std::setw( 8) << " LastStepZero " << " "
1971 << G4endl;
1972 G4cout << "( " << std::setw(7) << fExitNormal.x()
1973 << ", " << std::setw(7) << fExitNormal.y()
1974 << ", " << std::setw(7) << fExitNormal.z() << " ) "
1975 << std::setw( 5) << fValidExitNormal << " "
1976 << std::setw( 9) << fExiting << " "
1977 << std::setw( 9) << fEntering << " ";
1978 if ( fBlockedPhysicalVolume == nullptr )
1979 { G4cout << std::setw(15) << "None"; }
1980 else
1981 { G4cout << std::setw(15)<< fBlockedPhysicalVolume->GetName(); }
1982 G4cout << std::setw( 9) << fBlockedReplicaNo << " "
1983 << std::setw( 8) << fLastStepWasZero << " "
1984 << G4endl;
1985 }
1986 if( fVerbose > 2 )
1987 {
1988 G4cout.precision(8);
1989 G4cout << " Current Localpoint = " << fLastLocatedPointLocal << G4endl;
1990 G4cout << " PreviousSftOrigin = " << fPreviousSftOrigin << G4endl;
1991 G4cout << " PreviousSafety = " << fPreviousSafety << G4endl;
1992 }
1993 G4cout.precision(oldcoutPrec);
1994}
#define fEntering
double z() const
double x() const
double y() const

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

◆ ResetHierarchyAndLocate()

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

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(), and G4SteppingManager::SetInitialStep().

◆ ResetStackAndState()

void G4Navigator::ResetStackAndState ( )
inline

◆ ResetState()

void G4Navigator::ResetState ( )
protectedvirtual

Reimplemented in G4MultiNavigator.

Definition at line 1264 of file G4Navigator.cc.

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

Referenced by ResetHierarchyAndLocate().

◆ RestoreSavedState()

void G4Navigator::RestoreSavedState ( )
protected

Definition at line 702 of file G4Navigator.cc.

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

Referenced by CheckNextStep(), and ComputeSafety().

◆ SetExternalNavigation()

void G4Navigator::SetExternalNavigation ( G4VExternalNavigation externalNav)
inline

◆ SetGeometricallyLimitedStep()

void G4Navigator::SetGeometricallyLimitedStep ( )
inline

◆ SetPushVerbosity()

void G4Navigator::SetPushVerbosity ( G4bool  mode)
inline

◆ SetSavedState()

void G4Navigator::SetSavedState ( )
protected

Definition at line 668 of file G4Navigator.cc.

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

Referenced by CheckNextStep(), and ComputeSafety().

◆ SetupHierarchy()

void G4Navigator::SetupHierarchy ( )
protectedvirtual

Reimplemented in G4MultiNavigator.

Definition at line 1305 of file G4Navigator.cc.

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

Referenced by ResetHierarchyAndLocate().

◆ SetVerboseLevel()

void G4Navigator::SetVerboseLevel ( G4int  level)
inline

Referenced by GetGlobalExitNormal().

◆ SetVoxelNavigation()

void G4Navigator::SetVoxelNavigation ( G4VoxelNavigation voxelNav)

◆ SetWorldVolume()

◆ SeverityOfZeroStepping()

G4int G4Navigator::SeverityOfZeroStepping ( G4int noZeroSteps) const
inline

◆ VolumeType()

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

Friends And Related Function Documentation

◆ operator<<

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

Definition at line 2130 of file G4Navigator.cc.

2131{
2132 // Old version did only the following:
2133 // os << "Current History: " << G4endl << n.fHistory;
2134 // Old behaviour is recovered for fVerbose = 0
2135
2136 // Adapted from G4Navigator::PrintState() const
2137
2138 G4long oldcoutPrec = os.precision(4);
2139 if( n.fVerbose >= 4 )
2140 {
2141 os << "The current state of G4Navigator is: " << G4endl;
2142 os << " ValidExitNormal= " << n.fValidExitNormal << G4endl
2143 << " ExitNormal = " << n.fExitNormal << G4endl
2144 << " Exiting = " << n.fExiting << G4endl
2145 << " Entering = " << n.fEntering << G4endl
2146 << " BlockedPhysicalVolume= " ;
2147 if (n.fBlockedPhysicalVolume==0)
2148 os << "None";
2149 else
2150 os << n.fBlockedPhysicalVolume->GetName();
2151 os << G4endl
2152 << " BlockedReplicaNo = " << n.fBlockedReplicaNo << G4endl
2153 << " LastStepWasZero = " << n.fLastStepWasZero << G4endl
2154 << G4endl;
2155 }
2156 if( ( 1 < n.fVerbose) && (n.fVerbose < 4) )
2157 {
2158 os << G4endl; // Make sure to line up
2159 os << std::setw(30) << " ExitNormal " << " "
2160 << std::setw( 5) << " Valid " << " "
2161 << std::setw( 9) << " Exiting " << " "
2162 << std::setw( 9) << " Entering" << " "
2163 << std::setw(15) << " Blocked:Volume " << " "
2164 << std::setw( 9) << " ReplicaNo" << " "
2165 << std::setw( 8) << " LastStepZero " << " "
2166 << G4endl;
2167 os << "( " << std::setw(7) << n.fExitNormal.x()
2168 << ", " << std::setw(7) << n.fExitNormal.y()
2169 << ", " << std::setw(7) << n.fExitNormal.z() << " ) "
2170 << std::setw( 5) << n.fValidExitNormal << " "
2171 << std::setw( 9) << n.fExiting << " "
2172 << std::setw( 9) << n.fEntering << " ";
2173 if ( n.fBlockedPhysicalVolume==0 )
2174 { os << std::setw(15) << "None"; }
2175 else
2176 { os << std::setw(15)<< n.fBlockedPhysicalVolume->GetName(); }
2177 os << std::setw( 9) << n.fBlockedReplicaNo << " "
2178 << std::setw( 8) << n.fLastStepWasZero << " "
2179 << G4endl;
2180 }
2181 if( n.fVerbose > 2 )
2182 {
2183 os.precision(8);
2184 os << " Current Localpoint = " << n.fLastLocatedPointLocal << G4endl;
2185 os << " PreviousSftOrigin = " << n.fPreviousSftOrigin << G4endl;
2186 os << " PreviousSafety = " << n.fPreviousSafety << G4endl;
2187 }
2188 if( n.fVerbose > 3 || n.fVerbose == 0 )
2189 {
2190 os << "Current History: " << G4endl << n.fHistory;
2191 }
2192
2193 os.precision(oldcoutPrec);
2194 return os;
2195}

Member Data Documentation

◆ fEnteredDaughter

G4bool G4Navigator::fEnteredDaughter
protected

◆ fExitedMother

◆ fHistory

◆ fLastStepEndPointLocal

G4ThreeVector G4Navigator::fLastStepEndPointLocal
protected

Definition at line 391 of file G4Navigator.hh.

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

◆ fMinStep

G4double G4Navigator::fMinStep
protected

Definition at line 377 of file G4Navigator.hh.

Referenced by ComputeStep(), and G4Navigator().

◆ fSqTol

G4double G4Navigator::fSqTol
protected

Definition at line 377 of file G4Navigator.hh.

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

◆ fStepEndPoint

G4ThreeVector G4Navigator::fStepEndPoint
protected

◆ fVerbose

◆ fWasLimitedByGeometry

◆ kCarTolerance

G4double G4Navigator::kCarTolerance
protected

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