Geant4-11
Data Structures | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | Friends
G4Navigator Class Reference

#include <G4Navigator.hh>

Inheritance diagram for G4Navigator:
G4ErrorPropagationNavigator G4MultiNavigator

Data Structures

struct  G4SaveNavigatorState
 

Public Member Functions

void Activate (G4bool flag)
 
void CheckMode (G4bool mode)
 
G4double CheckNextStep (const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
 
G4NavigatorClone () const
 
virtual G4double ComputeSafety (const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
 
virtual G4double ComputeStep (const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
 
G4GRSSolidCreateGRSSolid () const
 
G4GRSVolumeCreateGRSVolume () const
 
G4TouchableHistoryCreateTouchableHistory () const
 
G4TouchableHistoryCreateTouchableHistory (const G4NavigationHistory *) const
 
virtual G4TouchableHistoryHandle CreateTouchableHistoryHandle () const
 
void EnableBestSafety (G4bool value=false)
 
G4bool EnteredDaughterVolume () const
 
G4bool ExitedMotherVolume () const
 
 G4Navigator ()
 
 G4Navigator (const G4Navigator &)=delete
 
G4ThreeVector GetCurrentLocalCoordinate () const
 
G4VExternalNavigationGetExternalNavigation () const
 
virtual G4ThreeVector GetGlobalExitNormal (const G4ThreeVector &point, G4bool *valid)
 
const G4AffineTransformGetGlobalToLocalTransform () const
 
G4ThreeVector GetLastStepEndPoint () const
 
virtual G4ThreeVector GetLocalExitNormal (G4bool *valid)
 
virtual G4ThreeVector GetLocalExitNormalAndCheck (const G4ThreeVector &point, G4bool *valid)
 
const G4AffineTransform GetLocalToGlobalTransform () const
 
G4AffineTransform GetMotherToDaughterTransform (G4VPhysicalVolume *dVolume, G4int dReplicaNo, EVolume dVolumeType)
 
G4int GetVerboseLevel () const
 
G4VPhysicalVolumeGetWorldVolume () const
 
void InformLastStep (G4double lastStep, G4bool entersDaughtVol, G4bool exitsMotherVol)
 
G4bool IsActive () const
 
G4bool IsCheckModeActive () const
 
virtual G4VPhysicalVolumeLocateGlobalPointAndSetup (const G4ThreeVector &point, const G4ThreeVector *direction=nullptr, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=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 LocateGlobalPointAndUpdateTouchableHandle (const G4ThreeVector &position, const G4ThreeVector &direction, G4TouchableHandle &oldTouchableToUpdate, const G4bool RelativeSearch=true)
 
virtual void LocateGlobalPointWithinVolume (const G4ThreeVector &position)
 
G4RotationMatrix NetRotation () const
 
G4ThreeVector NetTranslation () const
 
G4Navigatoroperator= (const G4Navigator &)=delete
 
void PrintState () const
 
virtual G4VPhysicalVolumeResetHierarchyAndLocate (const G4ThreeVector &point, const G4ThreeVector &direction, const G4TouchableHistory &h)
 
void ResetStackAndState ()
 
void SetExternalNavigation (G4VExternalNavigation *externalNav)
 
void SetGeometricallyLimitedStep ()
 
void SetPushVerbosity (G4bool mode)
 
void SetVerboseLevel (G4int level)
 
void SetVoxelNavigation (G4VoxelNavigation *voxelNav)
 
void SetWorldVolume (G4VPhysicalVolume *pWorld)
 
G4int SeverityOfZeroStepping (G4int *noZeroSteps) const
 
virtual ~G4Navigator ()
 

Protected Member Functions

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

Protected Attributes

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

Private Member Functions

void ComputeStepLog (const G4ThreeVector &pGlobalpoint, G4double moveLenSq) const
 
G4VoxelNavigationGetVoxelNavigator ()
 

Private Attributes

G4int fAbandonThreshold_NoZeroSteps = 25
 
G4int fActionThreshold_NoZeroSteps = 10
 
G4bool fActive = false
 
G4VPhysicalVolumefBlockedPhysicalVolume
 
G4int fBlockedReplicaNo
 
G4bool fCalculatedExitNormal
 
G4bool fChangedGrandMotherRefFrame
 
G4bool fCheck = false
 
G4bool fEntering
 
G4bool fExiting
 
G4ThreeVector fExitNormal
 
G4ThreeVector fExitNormalGlobalFrame
 
G4ThreeVector fGrandMotherExitNormal
 
G4ThreeVector fLastLocatedPointLocal
 
G4VPhysicalVolumefLastMotherPhys = nullptr
 
G4bool fLastStepWasZero
 
G4bool fLastTriedStepComputation = false
 
G4bool fLocatedOnEdge
 
G4bool fLocatedOutsideWorld
 
G4NormalNavigation fnormalNav
 
G4int fNumberZeroSteps
 
G4ParameterisedNavigation fparamNav
 
G4VExternalNavigationfpExternalNav = nullptr
 
G4double fPreviousSafety
 
G4ThreeVector fPreviousSftOrigin
 
G4bool fPushed = false
 
G4VoxelNavigationfpvoxelNav
 
G4VoxelSafetyfpVoxelSafety
 
G4RegularNavigation fregularNav
 
G4ReplicaNavigation freplicaNav
 
struct G4Navigator::G4SaveNavigatorState fSaveState
 
G4VPhysicalVolumefTopPhysical = nullptr
 
G4bool fValidExitNormal
 
G4bool fWarnPush = true
 

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) = "
66 << " Abandon Threshold (No Zero Steps) = "
68 }
72
74
77
79#ifdef ALTERNATIVE_VOXEL_NAV
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
G4int fActionThreshold_NoZeroSteps
Definition: G4Navigator.hh:451
G4VoxelSafety * fpVoxelSafety
Definition: G4Navigator.hh:533
G4ThreeVector fLastStepEndPointLocal
Definition: G4Navigator.hh:391
G4double fSqTol
Definition: G4Navigator.hh:377
G4RegularNavigation fregularNav
Definition: G4Navigator.hh:531
G4double kCarTolerance
Definition: G4Navigator.hh:377
G4int fAbandonThreshold_NoZeroSteps
Definition: G4Navigator.hh:453
G4NormalNavigation fnormalNav
Definition: G4Navigator.hh:523
void ResetStackAndState()
G4VoxelNavigation * fpvoxelNav
Definition: G4Navigator.hh:525
void SetNormalNavigation(G4NormalNavigation *fnormnav)
static const G4double kInfinity
Definition: geomdefs.hh:41
T sqr(const T &x)
Definition: templates.hh:128

References fAbandonThreshold_NoZeroSteps, fActionThreshold_NoZeroSteps, fLastStepEndPointLocal, fMinStep, fnormalNav, fpvoxelNav, fpVoxelSafety, fregularNav, fSqTol, fStepEndPoint, fVerbose, G4cout, G4endl, G4GeometryTolerance::GetInstance(), G4GeometryTolerance::GetSurfaceTolerance(), kCarTolerance, kInfinity, ResetStackAndState(), G4RegularNavigation::SetNormalNavigation(), and sqr().

◆ 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}
G4VExternalNavigation * fpExternalNav
Definition: G4Navigator.hh:532

References fpExternalNav, fpvoxelNav, and fpVoxelSafety.

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 1255 of file G4Navigator.cc.

1259{
1260 G4double step;
1261
1262 // Save the state, for this parasitic call
1263 //
1264 SetSavedState();
1265
1266 step = ComputeStep ( pGlobalpoint,
1267 pDirection,
1268 pCurrentProposedStepLength,
1269 pNewSafety );
1270
1271 // It is a parasitic call, so attempt to restore the key parts of the state
1272 //
1274 // NOTE: the state of the current subnavigator is NOT restored.
1275 // ***> TODO: restore subnavigator state
1276 // if( last_located) Need Position of last location
1277 // if( last_computed step) Need Endposition of last step
1278
1279 return step;
1280}
double G4double
Definition: G4Types.hh:83
void RestoreSavedState()
Definition: G4Navigator.cc:716
void SetSavedState()
Definition: G4Navigator.cc:682
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
Definition: G4Navigator.cc:770

References ComputeStep(), RestoreSavedState(), and SetSavedState().

Referenced by G4SafetyHelper::CheckNextStep().

◆ CheckOverlapsIterative()

G4bool G4Navigator::CheckOverlapsIterative ( G4VPhysicalVolume vol)
protected

Definition at line 2130 of file G4Navigator.cc.

2131{
2132 // Check and report overlaps
2133 //
2134 G4bool foundOverlap = false;
2135 G4int nPoints = 300000, ntrials = 9, numOverlaps = 5;
2136 G4double trialLength = 1.0 * CLHEP::centimeter;
2137 while ( ntrials-- > 0 && !foundOverlap )
2138 {
2139 if ( fVerbose > 1 )
2140 {
2141 G4cout << " ** Running overlap checks in volume "
2142 << vol->GetName()
2143 << " with length = " << trialLength << G4endl;
2144 }
2145 foundOverlap = vol->CheckOverlaps(nPoints, trialLength,
2146 fVerbose, numOverlaps);
2147 trialLength *= 0.1;
2148 if ( trialLength <= 1.0e-5 ) { numOverlaps= 1;}
2149 }
2150 return foundOverlap;
2151}
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
static constexpr double centimeter
Definition: SystemOfUnits.h:67

References CLHEP::centimeter, G4VPhysicalVolume::CheckOverlaps(), fVerbose, G4cout, G4endl, and G4VPhysicalVolume::GetName().

Referenced by ComputeStep().

◆ Clone()

G4Navigator * G4Navigator::Clone ( ) const
inline

◆ 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 1813 of file G4Navigator.cc.

1816{
1817 G4double newSafety = 0.0;
1818
1819#ifdef G4DEBUG_NAVIGATION
1820 G4int oldcoutPrec = G4cout.precision(8);
1821 if( fVerbose > 0 )
1822 {
1823 G4cout << "*** G4Navigator::ComputeSafety: ***" << G4endl
1824 << " Called at point: " << pGlobalpoint << G4endl;
1825
1826 G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
1827 G4cout << " Volume = " << motherPhysical->GetName()
1828 << " - Maximum length = " << pMaxLength << G4endl;
1829 if( fVerbose >= 4 )
1830 {
1831 G4cout << " ----- Upon entering Compute Safety:" << G4endl;
1832 PrintState();
1833 }
1834 }
1835#endif
1836
1837 G4double distEndpointSq = (pGlobalpoint-fStepEndPoint).mag2();
1838 G4bool stayedOnEndpoint = distEndpointSq < sqr(kCarTolerance);
1839 G4bool endpointOnSurface = fEnteredDaughter || fExitedMother;
1840
1841 if( endpointOnSurface && stayedOnEndpoint )
1842 {
1843#ifdef G4DEBUG_NAVIGATION
1844 if( fVerbose >= 2 )
1845 {
1846 G4cout << " G4Navigator::ComputeSafety() finds that point - "
1847 << pGlobalpoint << " - is on surface " << G4endl;
1848 if( fEnteredDaughter ) { G4cout << " entered new daughter volume"; }
1849 if( fExitedMother ) { G4cout << " and exited previous volume."; }
1850 G4cout << G4endl;
1851 G4cout << " EndPoint was = " << fStepEndPoint << G4endl;
1852 }
1853#endif
1854 newSafety = 0.0;
1855 }
1856 else // if( !(endpointOnSurface && stayedOnEndpoint) )
1857 {
1858 if (keepState) { SetSavedState(); }
1859
1860 // Pseudo-relocate to this point (updates voxel information only)
1861 //
1862 LocateGlobalPointWithinVolume( pGlobalpoint );
1863 // --->> DANGER: Side effects on sub-navigator voxel information <<---
1864 // Could be replaced again by 'granular' calls to sub-navigator
1865 // locates (similar side-effects, but faster.
1866 // Solutions:
1867 // 1) Re-locate (to where?)
1868 // 2) Insure that the methods using (G4ComputeStep?)
1869 // does a relocation (if information is disturbed only ?)
1870
1871#ifdef G4DEBUG_NAVIGATION
1872 if( fVerbose >= 2 )
1873 {
1874 G4cout << " G4Navigator::ComputeSafety() relocates-in-volume to point: "
1875 << pGlobalpoint << G4endl;
1876 }
1877#endif
1878 G4VPhysicalVolume* motherPhysical = fHistory.GetTopVolume();
1879 G4LogicalVolume* motherLogical = motherPhysical->GetLogicalVolume();
1880 G4SmartVoxelHeader* pVoxelHeader = motherLogical->GetVoxelHeader();
1881 G4ThreeVector localPoint = ComputeLocalPoint(pGlobalpoint);
1882
1884 {
1885 switch(CharacteriseDaughters(motherLogical))
1886 {
1887 case kNormal:
1888 if ( pVoxelHeader )
1889 {
1890 newSafety = fpVoxelSafety->ComputeSafety(localPoint,
1891 *motherPhysical, pMaxLength);
1892 // = VoxelNav().ComputeSafety(localPoint,fHistory,pMaxLength); // - Old method
1893 }
1894 else
1895 {
1896 newSafety=fnormalNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1897 }
1898 break;
1899 case kParameterised:
1900 if( GetDaughtersRegularStructureId(motherLogical) != 1 )
1901 {
1902 newSafety=fparamNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1903 }
1904 else // Regular structure
1905 {
1906 newSafety=fregularNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1907 }
1908 break;
1909 case kReplica:
1910 G4Exception("G4Navigator::ComputeSafety()", "GeomNav0001",
1911 FatalException, "Not applicable for replicated volumes.");
1912 break;
1913 case kExternal:
1914 newSafety = fpExternalNav->ComputeSafety(localPoint, fHistory,
1915 pMaxLength);
1916 break;
1917 }
1918 }
1919 else
1920 {
1921 newSafety = freplicaNav.ComputeSafety(pGlobalpoint, localPoint,
1922 fHistory, pMaxLength);
1923 }
1924
1925 if (keepState)
1926 {
1928 // This now overwrites the values of the Safety 'sphere' (correction)
1929 }
1930
1931 // Remember last safety origin & value
1932 //
1933 // We overwrite the Safety 'sphere' - keeping old behaviour
1934 fPreviousSftOrigin = pGlobalpoint;
1935 fPreviousSafety = newSafety;
1936 }
1937
1938#ifdef G4DEBUG_NAVIGATION
1939 if( fVerbose > 1 )
1940 {
1941 G4cout << " ---- Exiting ComputeSafety " << G4endl;
1942 if( fVerbose > 2 ) { PrintState(); }
1943 G4cout << " Returned value of Safety = " << newSafety << G4endl;
1944 }
1945 G4cout.precision(oldcoutPrec);
1946#endif
1947
1948 return newSafety;
1949}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
G4SmartVoxelHeader * GetVoxelHeader() const
EVolume GetTopVolumeType() const
G4VPhysicalVolume * GetTopVolume() const
G4double fPreviousSafety
Definition: G4Navigator.hh:429
G4bool fExitedMother
Definition: G4Navigator.hh:404
G4bool fEnteredDaughter
Definition: G4Navigator.hh:398
G4ParameterisedNavigation fparamNav
Definition: G4Navigator.hh:529
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
Definition: G4Navigator.cc:614
G4ReplicaNavigation freplicaNav
Definition: G4Navigator.hh:530
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
G4ThreeVector fPreviousSftOrigin
Definition: G4Navigator.hh:428
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

References CharacteriseDaughters(), ComputeLocalPoint(), G4NormalNavigation::ComputeSafety(), G4VExternalNavigation::ComputeSafety(), G4ReplicaNavigation::ComputeSafety(), G4ParameterisedNavigation::ComputeSafety(), G4RegularNavigation::ComputeSafety(), G4VoxelSafety::ComputeSafety(), FatalException, fEnteredDaughter, fExitedMother, fHistory, fnormalNav, fparamNav, fpExternalNav, fPreviousSafety, fPreviousSftOrigin, fpVoxelSafety, fregularNav, freplicaNav, fStepEndPoint, fVerbose, G4cout, G4endl, G4Exception(), GetDaughtersRegularStructureId(), G4VPhysicalVolume::GetLogicalVolume(), G4VPhysicalVolume::GetName(), G4NavigationHistory::GetTopVolume(), G4NavigationHistory::GetTopVolumeType(), G4LogicalVolume::GetVoxelHeader(), kCarTolerance, kExternal, kNormal, kParameterised, kReplica, LocateGlobalPointWithinVolume(), PrintState(), RestoreSavedState(), SetSavedState(), and sqr().

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 770 of file G4Navigator.cc.

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

References CharacteriseDaughters(), CheckOverlapsIterative(), ComputeLocalAxis(), ComputeLocalPoint(), G4ParameterisedNavigation::ComputeStep(), G4RegularNavigation::ComputeStep(), G4VoxelNavigation::ComputeStep(), G4ReplicaNavigation::ComputeStep(), G4NormalNavigation::ComputeStep(), G4VExternalNavigation::ComputeStep(), ComputeStepLog(), EventMustBeAborted, fAbandonThreshold_NoZeroSteps, fActionThreshold_NoZeroSteps, FatalException, fBlockedPhysicalVolume, fBlockedReplicaNo, fCalculatedExitNormal, fChangedGrandMotherRefFrame, fEnteredDaughter, fEntering, fExitedMother, fExiting, fExitNormal, fExitNormalGlobalFrame, fGrandMotherExitNormal, fHistory, fLastLocatedPointLocal, fLastMotherPhys, fLastStepEndPointLocal, fLastStepWasZero, fLastTriedStepComputation, fLocatedOnEdge, fMinStep, fnormalNav, fNumberZeroSteps, fparamNav, fpExternalNav, fPreviousSafety, fPreviousSftOrigin, fPushed, fregularNav, freplicaNav, fSqTol, fStepEndPoint, fValidExitNormal, fVerbose, fWarnPush, G4cout, G4endl, G4Exception(), G4ThreadLocal, G4VPhysicalVolume::GetCopyNo(), GetDaughtersRegularStructureId(), G4NavigationHistory::GetDepth(), G4VPhysicalVolume::GetLogicalVolume(), G4VPhysicalVolume::GetName(), G4VPhysicalVolume::GetRegularStructureId(), G4VPhysicalVolume::GetRotation(), G4LogicalVolume::GetSolid(), G4NavigationHistory::GetTopVolume(), G4NavigationHistory::GetTopVolumeType(), G4NavigationHistory::GetTransform(), G4LogicalVolume::GetVoxelHeader(), GetVoxelNavigator(), G4AffineTransform::InverseTransformAxis(), JustWarning, kCarTolerance, kExternal, kInfinity, kNormal, kParameterised, kReplica, LocateGlobalPointAndSetup(), LocateGlobalPointWithinVolume(), G4INCL::Math::min(), PrintState(), and G4VSolid::SurfaceNormal().

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

◆ ComputeStepLog()

void G4Navigator::ComputeStepLog ( const G4ThreeVector pGlobalpoint,
G4double  moveLenSq 
) const
private

Definition at line 2027 of file G4Navigator.cc.

2029{
2030 // The following checks only make sense if the move is larger
2031 // than the tolerance.
2032
2033 const G4double fAccuracyForWarning = kCarTolerance,
2034 fAccuracyForException = 1000*kCarTolerance;
2035
2036 G4ThreeVector OriginalGlobalpoint = fHistory.GetTopTransform().
2037 InverseTransformPoint(fLastLocatedPointLocal);
2038
2039 G4double shiftOriginSafSq = (fPreviousSftOrigin-pGlobalpoint).mag2();
2040
2041 // Check that the starting point of this step is
2042 // within the isotropic safety sphere of the last point
2043 // to a accuracy/precision given by fAccuracyForWarning.
2044 // If so give warning.
2045 // If it fails by more than fAccuracyForException exit with error.
2046 //
2047 if( shiftOriginSafSq >= sqr(fPreviousSafety) )
2048 {
2049 G4double shiftOrigin = std::sqrt(shiftOriginSafSq);
2050 G4double diffShiftSaf = shiftOrigin - fPreviousSafety;
2051
2052 if( diffShiftSaf > fAccuracyForWarning )
2053 {
2054 G4int oldcoutPrec = G4cout.precision(8);
2055 G4int oldcerrPrec = G4cerr.precision(10);
2056 std::ostringstream message, suggestion;
2057 message << "Accuracy error or slightly inaccurate position shift."
2058 << G4endl
2059 << " The Step's starting point has moved "
2060 << std::sqrt(moveLenSq)/mm << " mm " << G4endl
2061 << " since the last call to a Locate method." << G4endl
2062 << " This has resulted in moving "
2063 << shiftOrigin/mm << " mm "
2064 << " from the last point at which the safety "
2065 << " was calculated " << G4endl
2066 << " which is more than the computed safety= "
2067 << fPreviousSafety/mm << " mm at that point." << G4endl
2068 << " This difference is "
2069 << diffShiftSaf/mm << " mm." << G4endl
2070 << " The tolerated accuracy is "
2071 << fAccuracyForException/mm << " mm.";
2072
2073 suggestion << " ";
2074 static G4ThreadLocal G4int warnNow = 0;
2075 if( ((++warnNow % 100) == 1) )
2076 {
2077 message << G4endl
2078 << " This problem can be due to either " << G4endl
2079 << " - a process that has proposed a displacement"
2080 << " larger than the current safety , or" << G4endl
2081 << " - inaccuracy in the computation of the safety";
2082 suggestion << "We suggest that you " << G4endl
2083 << " - find i) what particle is being tracked, and "
2084 << " ii) through what part of your geometry " << G4endl
2085 << " for example by re-running this event with "
2086 << G4endl
2087 << " /tracking/verbose 1 " << G4endl
2088 << " - check which processes you declare for"
2089 << " this particle (and look at non-standard ones)"
2090 << G4endl
2091 << " - in case, create a detailed logfile"
2092 << " of this event using:" << G4endl
2093 << " /tracking/verbose 6 ";
2094 }
2095 G4Exception("G4Navigator::ComputeStep()",
2096 "GeomNav1002", JustWarning,
2097 message, G4String(suggestion.str()));
2098 G4cout.precision(oldcoutPrec);
2099 G4cerr.precision(oldcerrPrec);
2100 }
2101#ifdef G4DEBUG_NAVIGATION
2102 else
2103 {
2104 G4cerr << "WARNING - G4Navigator::ComputeStep()" << G4endl
2105 << " The Step's starting point has moved "
2106 << std::sqrt(moveLenSq) << "," << G4endl
2107 << " which has taken it to the limit of"
2108 << " the current safety. " << G4endl;
2109 }
2110#endif
2111 }
2112 G4double safetyPlus = fPreviousSafety + fAccuracyForException;
2113 if ( shiftOriginSafSq > sqr(safetyPlus) )
2114 {
2115 std::ostringstream message;
2116 message << "May lead to a crash or unreliable results." << G4endl
2117 << " Position has shifted considerably without"
2118 << " notifying the navigator !" << G4endl
2119 << " Tolerated safety: " << safetyPlus << G4endl
2120 << " Computed shift : " << shiftOriginSafSq;
2121 G4Exception("G4Navigator::ComputeStep()", "GeomNav1002",
2122 JustWarning, message);
2123 }
2124}
static constexpr double mm
Definition: G4SIunits.hh:95
G4GLOB_DLL std::ostream G4cerr
const G4AffineTransform & GetTopTransform() const

References fHistory, fLastLocatedPointLocal, fPreviousSafety, fPreviousSftOrigin, G4cerr, G4cout, G4endl, G4Exception(), G4ThreadLocal, G4NavigationHistory::GetTopTransform(), JustWarning, kCarTolerance, mm, and sqr().

Referenced by ComputeStep().

◆ CreateGRSSolid()

G4GRSSolid * G4Navigator::CreateGRSSolid ( ) const
inline

◆ CreateGRSVolume()

G4GRSVolume * G4Navigator::CreateGRSVolume ( ) const
inline

◆ CreateTouchableHistory() [1/2]

G4TouchableHistory * G4Navigator::CreateTouchableHistory ( ) const
inline

◆ CreateTouchableHistory() [2/2]

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

◆ CreateTouchableHistoryHandle()

G4TouchableHistoryHandle G4Navigator::CreateTouchableHistoryHandle ( ) const
virtual

Reimplemented in G4MultiNavigator.

Definition at line 1955 of file G4Navigator.cc.

1956{
1958}
G4ReferenceCountedHandle< G4TouchableHistory > G4TouchableHistoryHandle
G4TouchableHistory * CreateTouchableHistory() const

References CreateTouchableHistory().

Referenced by G4FastTrack::FRecordsAffineTransformation(), and 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()

G4VExternalNavigation * G4Navigator::GetExternalNavigation ( ) const
inline

◆ GetGlobalExitNormal()

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

Reimplemented in G4MultiNavigator, and G4ErrorPropagationNavigator.

Definition at line 1642 of file G4Navigator.cc.

1644{
1645 G4bool validNormal;
1646 G4ThreeVector localNormal, globalNormal;
1647
1648 G4bool usingStored = fCalculatedExitNormal && (
1649 ( fLastTriedStepComputation && fExiting ) // Just calculated it
1650 || // No locate in between
1652 && (IntersectPointGlobal-fStepEndPoint).mag2() < 10.0*fSqTol ) );
1653 // Calculated it 'just' before & then called locate
1654 // but it did not move position
1655
1656 if( usingStored )
1657 {
1658 // This was computed in last call to ComputeStep
1659 // and only if it arrived at boundary
1660 //
1661 globalNormal = fExitNormalGlobalFrame;
1662 G4double normMag2 = globalNormal.mag2();
1663 if( std::fabs ( normMag2 - 1.0 ) < perThousand ) // was perMillion
1664 {
1665 *pNormalCalculated = true; // ComputeStep always computes it if Exiting
1666 // (fExiting==true)
1667 }
1668 else
1669 {
1670 G4ExceptionDescription message;
1671 message.precision(10);
1672 message << " WARNING> Expected normal-global-frame to be valid, "
1673 << " i.e. a unit vector!" << G4endl
1674 << " - but |normal| = " << std::sqrt(normMag2)
1675 << " - and |normal|^2 = " << normMag2 << G4endl
1676 << " which differs from 1.0 by " << normMag2 - 1.0 << G4endl
1677 << " n = " << fExitNormalGlobalFrame << G4endl
1678 << " Global point: " << IntersectPointGlobal << G4endl
1679 << " Volume: " << fHistory.GetTopVolume()->GetName() << G4endl;
1680#ifdef G4VERBOSE
1682 if ( candLog )
1683 {
1684 message << " Solid: " << candLog->GetSolid()->GetName()
1685 << ", Type: " << candLog->GetSolid()->GetEntityType() << G4endl
1686 << *candLog->GetSolid() << G4endl;
1687 }
1688#endif
1689 message << "============================================================"
1690 << G4endl;
1691 G4int oldVerbose = fVerbose;
1692 fVerbose = 4;
1693 message << " State of Navigator: " << G4endl;
1694 message << *this << G4endl;
1695 fVerbose = oldVerbose;
1696 message << "============================================================"
1697 << G4endl;
1698
1699 G4Exception("G4Navigator::GetGlobalExitNormal()",
1700 "GeomNav0003",JustWarning, message,
1701 "Value obtained from stored global-normal is not a unit vector.");
1702
1703 // (Re)Compute it now -- as either it was not computed, or it is wrong.
1704 //
1705 localNormal = GetLocalExitNormalAndCheck(IntersectPointGlobal,
1706 &validNormal);
1707 *pNormalCalculated = fCalculatedExitNormal;
1708 globalNormal = fHistory.GetTopTransform()
1709 .InverseTransformAxis(localNormal);
1710 }
1711 }
1712 else
1713 {
1714 localNormal = GetLocalExitNormalAndCheck(IntersectPointGlobal,&validNormal);
1715 *pNormalCalculated = fCalculatedExitNormal;
1716
1717#ifdef G4DEBUG_NAVIGATION
1718 usingStored = false;
1719
1720 if( (!validNormal) && !fCalculatedExitNormal )
1721 {
1723 edN << " Calculated = " << fCalculatedExitNormal << G4endl;
1724 edN << " Entering= " << fEntering << G4endl;
1725 G4int oldVerbose = this->GetVerboseLevel();
1726 this->SetVerboseLevel(4);
1727 edN << " State of Navigator: " << G4endl;
1728 edN << *this << G4endl;
1729 this->SetVerboseLevel( oldVerbose );
1730
1731 G4Exception("G4Navigator::GetGlobalExitNormal()",
1732 "GeomNav0003", JustWarning, edN,
1733 "LocalExitNormalAndCheck() did not calculate Normal.");
1734 }
1735#endif
1736
1737 G4double localMag2 = localNormal.mag2();
1738 if( validNormal && (std::fabs(localMag2-1.0)) > kToleranceNormalCheck )
1739 {
1741 edN.precision(10);
1742 edN << "G4Navigator::GetGlobalExitNormal: "
1743 << " Using Local Normal - from call to GetLocalExitNormalAndCheck. "
1744 << G4endl
1745 << " Local Exit Normal : " << " || = " << std::sqrt(localMag2)
1746 << " vec = " << localNormal << G4endl
1747 << " Global Exit Normal : " << " || = " << globalNormal.mag()
1748 << " vec = " << globalNormal << G4endl
1749 << " Global point: " << IntersectPointGlobal << G4endl;
1750 edN << " Calculated It = " << fCalculatedExitNormal << G4endl
1751 << " Volume: " << fHistory.GetTopVolume()->GetName() << G4endl;
1752#ifdef G4VERBOSE
1754 if ( candLog )
1755 {
1756 edN << " Solid: " << candLog->GetSolid()->GetName()
1757 << ", Type: " << candLog->GetSolid()->GetEntityType() << G4endl
1758 << *candLog->GetSolid();
1759 }
1760#endif
1761 G4Exception("G4Navigator::GetGlobalExitNormal()",
1762 "GeomNav0003",JustWarning, edN,
1763 "Value obtained from new local *solid* is incorrect.");
1764 localNormal = localNormal.unit(); // Should we correct it ??
1765 }
1766 globalNormal = fHistory.GetTopTransform()
1767 .InverseTransformAxis(localNormal);
1768 }
1769
1770#ifdef G4DEBUG_NAVIGATION
1771 if( usingStored )
1772 {
1773 G4ThreeVector globalNormAgn;
1774
1775 localNormal = GetLocalExitNormalAndCheck(IntersectPointGlobal,&validNormal);
1776
1777 globalNormAgn = fHistory.GetTopTransform()
1778 .InverseTransformAxis(localNormal);
1779
1780 // Check the value computed against fExitNormalGlobalFrame
1781 G4ThreeVector diffNorm = globalNormAgn - fExitNormalGlobalFrame;
1782 if( diffNorm.mag2() > kToleranceNormalCheck )
1783 {
1785 edDfn << "Found difference in normals in case of exiting mother "
1786 << "- when Get is called after ComputingStep " << G4endl;
1787 edDfn << " Magnitude of diff = " << diffNorm.mag() << G4endl;
1788 edDfn << " Normal stored (Global) = " << fExitNormalGlobalFrame
1789 << G4endl;
1790 edDfn << " Global Computed from Local = " << globalNormAgn << G4endl;
1791 G4Exception("G4Navigator::GetGlobalExitNormal()", "GeomNav0003",
1792 JustWarning, edDfn);
1793 }
1794 }
1795#endif
1796
1797 // Synchronise stored global exit normal as possibly re-computed here
1798 //
1799 fExitNormalGlobalFrame = globalNormal;
1800
1801 return globalNormal;
1802}
static const G4double kToleranceNormalCheck
Definition: G4Navigator.cc:47
static constexpr double perThousand
Definition: G4SIunits.hh:326
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

References fCalculatedExitNormal, fEntering, fExiting, fExitNormalGlobalFrame, fHistory, fLastTriedStepComputation, fSqTol, fStepEndPoint, fVerbose, G4endl, G4Exception(), G4VSolid::GetEntityType(), GetLocalExitNormalAndCheck(), G4VPhysicalVolume::GetLogicalVolume(), G4VPhysicalVolume::GetName(), G4VSolid::GetName(), G4LogicalVolume::GetSolid(), G4NavigationHistory::GetTopTransform(), G4NavigationHistory::GetTopVolume(), GetVerboseLevel(), G4AffineTransform::InverseTransformAxis(), JustWarning, kToleranceNormalCheck, CLHEP::Hep3Vector::mag(), CLHEP::Hep3Vector::mag2(), perThousand, SetVerboseLevel(), and CLHEP::Hep3Vector::unit().

Referenced by G4MultiNavigator::GetGlobalExitNormal(), G4ErrorPropagationNavigator::GetGlobalExitNormal(), G4VIntersectionLocator::GetLastSurfaceNormal(), 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;}

References fStepEndPoint.

◆ GetLocalExitNormal()

G4ThreeVector G4Navigator::GetLocalExitNormal ( G4bool valid)
virtual

Reimplemented in G4MultiNavigator.

Definition at line 1389 of file G4Navigator.cc.

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

References G4VSolid::DistanceToIn(), G4VSolid::DistanceToOut(), EnteredDaughterVolume(), FatalException, fBlockedPhysicalVolume, fBlockedReplicaNo, fCalculatedExitNormal, fCheck, fEntering, fExitedMother, fExiting, fGrandMotherExitNormal, fHistory, fLastLocatedPointLocal, fLastStepEndPointLocal, fLastTriedStepComputation, fVerbose, G4endl, G4Exception(), G4VSolid::GetEntityType(), G4VPhysicalVolume::GetLogicalVolume(), GetMotherToDaughterTransform(), G4LogicalVolume::GetName(), G4VPhysicalVolume::GetName(), G4VSolid::GetName(), G4LogicalVolume::GetSolid(), G4NavigationHistory::GetTopVolume(), G4VSolid::Inside(), G4AffineTransform::InverseTransformAxis(), JustWarning, kCarTolerance, kInside, kOutside, kSurface, kToleranceNormalCheck, CLHEP::Hep3Vector::mag2(), G4VSolid::SurfaceNormal(), G4AffineTransform::TransformPoint(), and VolumeType().

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 1606 of file G4Navigator.cc.

1613{
1614#ifdef G4DEBUG_NAVIGATION
1615 // Check Current point against expected 'local' value
1616 //
1618 {
1619 G4ThreeVector ExpectedBoundaryPointLocal;
1620
1621 const G4AffineTransform& GlobalToLocal = GetGlobalToLocalTransform();
1622 ExpectedBoundaryPointLocal =
1623 GlobalToLocal.TransformPoint( ExpectedBoundaryPointGlobal );
1624
1625 // Add here: Comparison against expected position,
1626 // i.e. the endpoint of ComputeStep
1627 }
1628#endif
1629
1630 return GetLocalExitNormal( pValid );
1631}
virtual G4ThreeVector GetLocalExitNormal(G4bool *valid)
const G4AffineTransform & GetGlobalToLocalTransform() const

References fLastTriedStepComputation, GetGlobalToLocalTransform(), GetLocalExitNormal(), and G4AffineTransform::TransformPoint().

Referenced by GetGlobalExitNormal().

◆ GetLocalToGlobalTransform()

const G4AffineTransform G4Navigator::GetLocalToGlobalTransform ( ) const
inline

◆ GetMotherToDaughterTransform()

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

Definition at line 1555 of file G4Navigator.cc.

1558{
1559 switch (enteringVolumeType)
1560 {
1561 case kNormal: // Nothing is needed to prepare the transformation
1562 break; // It is stored already in the physical volume (placement)
1563 case kReplica: // Sets the transform in the Replica - tbc
1564 G4Exception("G4Navigator::GetMotherToDaughterTransform()",
1565 "GeomNav0001", FatalException,
1566 "Method NOT Implemented yet for replica volumes.");
1567 break;
1568 case kParameterised:
1569 if( pEnteringPhysVol->GetRegularStructureId() == 0 )
1570 {
1571 G4VPVParameterisation *pParam =
1572 pEnteringPhysVol->GetParameterisation();
1573 G4VSolid* pSolid =
1574 pParam->ComputeSolid(enteringReplicaNo, pEnteringPhysVol);
1575 pSolid->ComputeDimensions(pParam, enteringReplicaNo, pEnteringPhysVol);
1576
1577 // Sets the transform in the Parameterisation
1578 //
1579 pParam->ComputeTransformation(enteringReplicaNo, pEnteringPhysVol);
1580
1581 // Set the correct solid and material in Logical Volume
1582 //
1583 G4LogicalVolume* pLogical = pEnteringPhysVol->GetLogicalVolume();
1584 pLogical->SetSolid( pSolid );
1585 }
1586 break;
1587 case kExternal:
1588 // Expect that nothing is needed to prepare the transformation.
1589 // It is stored already in the physical volume (placement)
1590 break;
1591 }
1592 return G4AffineTransform(pEnteringPhysVol->GetRotation(),
1593 pEnteringPhysVol->GetTranslation()).Invert();
1594}
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

References G4VSolid::ComputeDimensions(), G4VPVParameterisation::ComputeSolid(), G4VPVParameterisation::ComputeTransformation(), FatalException, G4Exception(), G4VPhysicalVolume::GetLogicalVolume(), G4VPhysicalVolume::GetParameterisation(), G4VPhysicalVolume::GetRegularStructureId(), G4VPhysicalVolume::GetRotation(), G4VPhysicalVolume::GetTranslation(), G4AffineTransform::Invert(), kExternal, kNormal, kParameterised, kReplica, and G4LogicalVolume::SetSolid().

Referenced by GetLocalExitNormal().

◆ GetVerboseLevel()

G4int G4Navigator::GetVerboseLevel ( ) const
inline

◆ GetVoxelNavigator()

G4VoxelNavigation & G4Navigator::GetVoxelNavigator ( )
inlineprivate

◆ GetWorldVolume()

G4VPhysicalVolume * G4Navigator::GetWorldVolume ( ) const
inline

◆ InformLastStep()

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

Definition at line 2240 of file G4Navigator.cc.

2241{
2242 G4bool zeroStep = ( lastStep == 0.0 );
2243 fLocatedOnEdge = fLastStepWasZero && zeroStep;
2244 fLastStepWasZero = zeroStep;
2245
2246 fExiting = exitsMotherVol;
2247 fEntering = entersDaughtVol;
2248}

References fEntering, fExiting, fLastStepWasZero, and fLocatedOnEdge.

◆ IsActive()

G4bool G4Navigator::IsActive ( ) const
inline

Referenced by export_G4Navigator().

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

References G4NavigationHistory::BackLevel(), G4ReplicaNavigation::BackLocate(), CharacteriseDaughters(), G4VSolid::ComputeDimensions(), G4VPVParameterisation::ComputeSolid(), G4ReplicaNavigation::ComputeTransformation(), G4VPVParameterisation::ComputeTransformation(), FatalException, fBlockedPhysicalVolume, fBlockedReplicaNo, fChangedGrandMotherRefFrame, fCheck, fEnteredDaughter, fEntering, fExitedMother, fExiting, fGrandMotherExitNormal, fHistory, fLastLocatedPointLocal, fLastTriedStepComputation, fLocatedOnEdge, fLocatedOutsideWorld, fnormalNav, fparamNav, fpExternalNav, fregularNav, freplicaNav, fValidExitNormal, fVerbose, fWasLimitedByGeometry, G4cout, G4endl, G4Exception(), G4VPhysicalVolume::GetCopyNo(), GetDaughtersRegularStructureId(), G4NavigationHistory::GetDepth(), G4VPhysicalVolume::GetLogicalVolume(), G4VPhysicalVolume::GetName(), G4VSolid::GetName(), G4VPhysicalVolume::GetParameterisation(), G4VPhysicalVolume::GetRegularStructureId(), G4VPhysicalVolume::GetRotation(), G4LogicalVolume::GetSolid(), G4NavigationHistory::GetTopReplicaNo(), G4NavigationHistory::GetTopTransform(), G4NavigationHistory::GetTopVolume(), G4NavigationHistory::GetTopVolumeType(), G4LogicalVolume::GetVoxelHeader(), GetVoxelNavigator(), G4VSolid::Inside(), G4VExternalNavigation::Inside(), kExternal, kNormal, kOutside, kParameterised, kReplica, kSurface, G4NormalNavigation::LevelLocate(), G4ParameterisedNavigation::LevelLocate(), G4RegularNavigation::LevelLocate(), G4ReplicaNavigation::LevelLocate(), G4VoxelNavigation::LevelLocate(), G4VExternalNavigation::LevelLocate(), G4NavigationHistory::NewLevel(), CLHEP::normal(), PrintState(), ResetStackAndState(), G4VPhysicalVolume::SetCopyNo(), G4LogicalVolume::SetSolid(), G4VSolid::SurfaceNormal(), G4AffineTransform::TransformAxis(), G4AffineTransform::TransformPoint(), G4LogicalVolume::UpdateMaterial(), and VolumeType().

Referenced by G4AdjointPrimaryGenerator::ComputeAccumulatedDepthVectorAlongBackRay(), ComputeStep(), G4VIntersectionLocator::GetLocalSurfaceNormal(), G4SPSPosDistribution::IsSourceConfined(), 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 614 of file G4Navigator.cc.

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

References CharacteriseDaughters(), ComputeLocalPoint(), fBlockedPhysicalVolume, fBlockedReplicaNo, fChangedGrandMotherRefFrame, fEnteredDaughter, fEntering, fExitedMother, fExiting, fHistory, fLastLocatedPointLocal, fLastTriedStepComputation, fparamNav, fpExternalNav, fWasLimitedByGeometry, GetDaughtersRegularStructureId(), G4VPhysicalVolume::GetLogicalVolume(), G4NavigationHistory::GetTopVolume(), G4LogicalVolume::GetVoxelHeader(), GetVoxelNavigator(), kExternal, kNormal, kParameterised, kReplica, G4ParameterisedNavigation::ParamVoxelLocate(), G4VExternalNavigation::RelocateWithinVolume(), and G4VoxelNavigation::VoxelLocate().

Referenced by G4VIntersectionLocator::AdjustmentOfFoundIntersection(), 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 1964 of file G4Navigator.cc.

1965{
1966 G4int oldcoutPrec = G4cout.precision(4);
1967 if( fVerbose >= 4 )
1968 {
1969 G4cout << "The current state of G4Navigator is: " << G4endl;
1970 G4cout << " ValidExitNormal= " << fValidExitNormal // << G4endl
1971 << " ExitNormal = " << fExitNormal // << G4endl
1972 << " Exiting = " << fExiting // << G4endl
1973 << " Entering = " << fEntering // << G4endl
1974 << " BlockedPhysicalVolume= " ;
1976 {
1977 G4cout << "None";
1978 }
1979 else
1980 {
1982 }
1983 G4cout << G4endl
1984 << " BlockedReplicaNo = " << fBlockedReplicaNo // << G4endl
1985 << " LastStepWasZero = " << fLastStepWasZero // << G4endl
1986 << G4endl;
1987 }
1988 if( ( 1 < fVerbose) && (fVerbose < 4) )
1989 {
1990 G4cout << G4endl; // Make sure to line up
1991 G4cout << std::setw(30) << " ExitNormal " << " "
1992 << std::setw( 5) << " Valid " << " "
1993 << std::setw( 9) << " Exiting " << " "
1994 << std::setw( 9) << " Entering" << " "
1995 << std::setw(15) << " Blocked:Volume " << " "
1996 << std::setw( 9) << " ReplicaNo" << " "
1997 << std::setw( 8) << " LastStepZero " << " "
1998 << G4endl;
1999 G4cout << "( " << std::setw(7) << fExitNormal.x()
2000 << ", " << std::setw(7) << fExitNormal.y()
2001 << ", " << std::setw(7) << fExitNormal.z() << " ) "
2002 << std::setw( 5) << fValidExitNormal << " "
2003 << std::setw( 9) << fExiting << " "
2004 << std::setw( 9) << fEntering << " ";
2005 if ( fBlockedPhysicalVolume == nullptr )
2006 { G4cout << std::setw(15) << "None"; }
2007 else
2008 { G4cout << std::setw(15)<< fBlockedPhysicalVolume->GetName(); }
2009 G4cout << std::setw( 9) << fBlockedReplicaNo << " "
2010 << std::setw( 8) << fLastStepWasZero << " "
2011 << G4endl;
2012 }
2013 if( fVerbose > 2 )
2014 {
2015 G4cout.precision(8);
2016 G4cout << " Current Localpoint = " << fLastLocatedPointLocal << G4endl;
2017 G4cout << " PreviousSftOrigin = " << fPreviousSftOrigin << G4endl;
2018 G4cout << " PreviousSafety = " << fPreviousSafety << G4endl;
2019 }
2020 G4cout.precision(oldcoutPrec);
2021}
double z() const
double x() const
double y() const

References fBlockedPhysicalVolume, fBlockedReplicaNo, fEntering, fExiting, fExitNormal, fLastLocatedPointLocal, fLastStepWasZero, fPreviousSafety, fPreviousSftOrigin, fValidExitNormal, fVerbose, G4cout, G4endl, G4VPhysicalVolume::GetName(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by ComputeSafety(), ComputeStep(), export_G4Navigator(), 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

References fHistory, fLastTriedStepComputation, G4TouchableHistory::GetHistory(), LocateGlobalPointAndSetup(), ResetState(), and SetupHierarchy().

Referenced by G4MultiNavigator::ResetHierarchyAndLocate(), and G4SteppingManager::SetInitialStep().

◆ ResetStackAndState()

void G4Navigator::ResetStackAndState ( )
inline

◆ ResetState()

void G4Navigator::ResetState ( )
protectedvirtual

◆ RestoreSavedState()

void G4Navigator::RestoreSavedState ( )
protected

Definition at line 716 of file G4Navigator.cc.

717{
722
725
727
733
734 // The 'expected' behaviour is to restore these too (fix 2014.05.26)
737}
struct G4Navigator::G4SaveNavigatorState fSaveState
G4VPhysicalVolume * spBlockedPhysicalVolume
Definition: G4Navigator.hh:500

References fBlockedPhysicalVolume, fBlockedReplicaNo, fEnteredDaughter, fEntering, fExitedMother, fExiting, fExitNormal, fLastLocatedPointLocal, fLastStepWasZero, fLocatedOutsideWorld, fPreviousSafety, fPreviousSftOrigin, fSaveState, fValidExitNormal, fWasLimitedByGeometry, G4Navigator::G4SaveNavigatorState::sBlockedReplicaNo, G4Navigator::G4SaveNavigatorState::sEnteredDaughter, G4Navigator::G4SaveNavigatorState::sEntering, G4Navigator::G4SaveNavigatorState::sExitedMother, G4Navigator::G4SaveNavigatorState::sExiting, G4Navigator::G4SaveNavigatorState::sExitNormal, G4Navigator::G4SaveNavigatorState::sLastLocatedPointLocal, G4Navigator::G4SaveNavigatorState::sLastStepWasZero, G4Navigator::G4SaveNavigatorState::sLocatedOutsideWorld, G4Navigator::G4SaveNavigatorState::spBlockedPhysicalVolume, G4Navigator::G4SaveNavigatorState::sPreviousSafety, G4Navigator::G4SaveNavigatorState::sPreviousSftOrigin, G4Navigator::G4SaveNavigatorState::sValidExitNormal, and G4Navigator::G4SaveNavigatorState::sWasLimitedByGeometry.

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 682 of file G4Navigator.cc.

683{
684 // Note: the state of dependent objects is not currently saved.
685 // ( This means that the full state is changed by calls between
686 // SetSavedState() and RestoreSavedState();
687
692
695
697
703
704 // Even the safety sphere - if you want to change it do it explicitly!
705 //
708}

References fBlockedPhysicalVolume, fBlockedReplicaNo, fEnteredDaughter, fEntering, fExitedMother, fExiting, fExitNormal, fLastLocatedPointLocal, fLastStepWasZero, fLocatedOutsideWorld, fPreviousSafety, fPreviousSftOrigin, fSaveState, fValidExitNormal, fWasLimitedByGeometry, G4Navigator::G4SaveNavigatorState::sBlockedReplicaNo, G4Navigator::G4SaveNavigatorState::sEnteredDaughter, G4Navigator::G4SaveNavigatorState::sEntering, G4Navigator::G4SaveNavigatorState::sExitedMother, G4Navigator::G4SaveNavigatorState::sExiting, G4Navigator::G4SaveNavigatorState::sExitNormal, G4Navigator::G4SaveNavigatorState::sLastLocatedPointLocal, G4Navigator::G4SaveNavigatorState::sLastStepWasZero, G4Navigator::G4SaveNavigatorState::sLocatedOutsideWorld, G4Navigator::G4SaveNavigatorState::spBlockedPhysicalVolume, G4Navigator::G4SaveNavigatorState::sPreviousSafety, G4Navigator::G4SaveNavigatorState::sPreviousSftOrigin, G4Navigator::G4SaveNavigatorState::sValidExitNormal, and G4Navigator::G4SaveNavigatorState::sWasLimitedByGeometry.

Referenced by CheckNextStep(), and ComputeSafety().

◆ SetupHierarchy()

void G4Navigator::SetupHierarchy ( )
protectedvirtual

Reimplemented in G4MultiNavigator.

Definition at line 1329 of file G4Navigator.cc.

1330{
1331 const G4int cdepth = fHistory.GetDepth();
1332 G4VPhysicalVolume* current;
1333 G4VSolid* pSolid;
1334 G4VPVParameterisation* pParam;
1335
1336 for ( auto i=1; i<=cdepth; ++i )
1337 {
1338 current = fHistory.GetVolume(i);
1339 switch ( fHistory.GetVolumeType(i) )
1340 {
1341 case kNormal:
1342 case kExternal:
1343 break;
1344 case kReplica:
1346 break;
1347 case kParameterised:
1348 G4int replicaNo;
1349 pParam = current->GetParameterisation();
1350 replicaNo = fHistory.GetReplicaNo(i);
1351 pSolid = pParam->ComputeSolid(replicaNo, current);
1352
1353 // Set up dimensions & transform in solid/physical volume
1354 //
1355 pSolid->ComputeDimensions(pParam, replicaNo, current);
1356 pParam->ComputeTransformation(replicaNo, current);
1357
1358 G4TouchableHistory* pTouchable = nullptr;
1359 if( pParam->IsNested() )
1360 {
1361 pTouchable= new G4TouchableHistory( fHistory );
1362 pTouchable->MoveUpHistory(); // Move up to the parent level
1363 // Adequate only if Nested at the Branch level (last)
1364 // To extend to other cases:
1365 // pTouchable->MoveUpHistory(cdepth-i-1);
1366 // Move to the parent level of *Current* level
1367 // Could replace this line and constructor with a revised
1368 // c-tor for History(levels to drop)
1369 }
1370 // Set up the correct solid and material in Logical Volume
1371 //
1372 G4LogicalVolume* pLogical = current->GetLogicalVolume();
1373 pLogical->SetSolid( pSolid );
1374 pLogical->UpdateMaterial( pParam ->
1375 ComputeMaterial(replicaNo, current, pTouchable) );
1376 delete pTouchable;
1377 break;
1378 }
1379 }
1380}
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

References G4VSolid::ComputeDimensions(), G4VPVParameterisation::ComputeSolid(), G4ReplicaNavigation::ComputeTransformation(), G4VPVParameterisation::ComputeTransformation(), fHistory, freplicaNav, G4NavigationHistory::GetDepth(), G4VPhysicalVolume::GetLogicalVolume(), G4VPhysicalVolume::GetParameterisation(), G4NavigationHistory::GetReplicaNo(), G4NavigationHistory::GetVolume(), G4NavigationHistory::GetVolumeType(), G4VPVParameterisation::IsNested(), kExternal, kNormal, kParameterised, kReplica, G4TouchableHistory::MoveUpHistory(), G4LogicalVolume::SetSolid(), and G4LogicalVolume::UpdateMaterial().

Referenced by ResetHierarchyAndLocate().

◆ SetVerboseLevel()

void G4Navigator::SetVerboseLevel ( G4int  level)
inline

◆ SetVoxelNavigation()

void G4Navigator::SetVoxelNavigation ( G4VoxelNavigation voxelNav)

◆ SetWorldVolume()

void G4Navigator::SetWorldVolume ( G4VPhysicalVolume pWorld)
inline

◆ 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 2157 of file G4Navigator.cc.

2158{
2159 // Old version did only the following:
2160 // os << "Current History: " << G4endl << n.fHistory;
2161 // Old behaviour is recovered for fVerbose = 0
2162
2163 // Adapted from G4Navigator::PrintState() const
2164
2165 G4int oldcoutPrec = os.precision(4);
2166 if( n.fVerbose >= 4 )
2167 {
2168 os << "The current state of G4Navigator is: " << G4endl;
2169 os << " ValidExitNormal= " << n.fValidExitNormal << G4endl
2170 << " ExitNormal = " << n.fExitNormal << G4endl
2171 << " Exiting = " << n.fExiting << G4endl
2172 << " Entering = " << n.fEntering << G4endl
2173 << " BlockedPhysicalVolume= " ;
2174 if (n.fBlockedPhysicalVolume==0)
2175 os << "None";
2176 else
2177 os << n.fBlockedPhysicalVolume->GetName();
2178 os << G4endl
2179 << " BlockedReplicaNo = " << n.fBlockedReplicaNo << G4endl
2180 << " LastStepWasZero = " << n.fLastStepWasZero << G4endl
2181 << G4endl;
2182 }
2183 if( ( 1 < n.fVerbose) && (n.fVerbose < 4) )
2184 {
2185 os << G4endl; // Make sure to line up
2186 os << std::setw(30) << " ExitNormal " << " "
2187 << std::setw( 5) << " Valid " << " "
2188 << std::setw( 9) << " Exiting " << " "
2189 << std::setw( 9) << " Entering" << " "
2190 << std::setw(15) << " Blocked:Volume " << " "
2191 << std::setw( 9) << " ReplicaNo" << " "
2192 << std::setw( 8) << " LastStepZero " << " "
2193 << G4endl;
2194 os << "( " << std::setw(7) << n.fExitNormal.x()
2195 << ", " << std::setw(7) << n.fExitNormal.y()
2196 << ", " << std::setw(7) << n.fExitNormal.z() << " ) "
2197 << std::setw( 5) << n.fValidExitNormal << " "
2198 << std::setw( 9) << n.fExiting << " "
2199 << std::setw( 9) << n.fEntering << " ";
2200 if ( n.fBlockedPhysicalVolume==0 )
2201 { os << std::setw(15) << "None"; }
2202 else
2203 { os << std::setw(15)<< n.fBlockedPhysicalVolume->GetName(); }
2204 os << std::setw( 9) << n.fBlockedReplicaNo << " "
2205 << std::setw( 8) << n.fLastStepWasZero << " "
2206 << G4endl;
2207 }
2208 if( n.fVerbose > 2 )
2209 {
2210 os.precision(8);
2211 os << " Current Localpoint = " << n.fLastLocatedPointLocal << G4endl;
2212 os << " PreviousSftOrigin = " << n.fPreviousSftOrigin << G4endl;
2213 os << " PreviousSafety = " << n.fPreviousSafety << G4endl;
2214 }
2215 if( n.fVerbose > 3 || n.fVerbose == 0 )
2216 {
2217 os << "Current History: " << G4endl << n.fHistory;
2218 }
2219
2220 os.precision(oldcoutPrec);
2221 return os;
2222}

Field Documentation

◆ fAbandonThreshold_NoZeroSteps

G4int G4Navigator::fAbandonThreshold_NoZeroSteps = 25
private

Definition at line 453 of file G4Navigator.hh.

Referenced by ComputeStep(), and G4Navigator().

◆ fActionThreshold_NoZeroSteps

G4int G4Navigator::fActionThreshold_NoZeroSteps = 10
private

Definition at line 451 of file G4Navigator.hh.

Referenced by ComputeStep(), and G4Navigator().

◆ fActive

G4bool G4Navigator::fActive = false
private

Definition at line 456 of file G4Navigator.hh.

◆ fBlockedPhysicalVolume

G4VPhysicalVolume* G4Navigator::fBlockedPhysicalVolume
private

◆ fBlockedReplicaNo

G4int G4Navigator::fBlockedReplicaNo
private

◆ fCalculatedExitNormal

G4bool G4Navigator::fCalculatedExitNormal
private

Definition at line 483 of file G4Navigator.hh.

Referenced by ComputeStep(), GetGlobalExitNormal(), GetLocalExitNormal(), and ResetState().

◆ fChangedGrandMotherRefFrame

G4bool G4Navigator::fChangedGrandMotherRefFrame
private

◆ fCheck

G4bool G4Navigator::fCheck = false
private

Definition at line 537 of file G4Navigator.hh.

Referenced by GetLocalExitNormal(), and LocateGlobalPointAndSetup().

◆ fEnteredDaughter

G4bool G4Navigator::fEnteredDaughter
protected

◆ fEntering

G4bool G4Navigator::fEntering
private

◆ fExitedMother

G4bool G4Navigator::fExitedMother
protected

◆ fExiting

G4bool G4Navigator::fExiting
private

◆ fExitNormal

G4ThreeVector G4Navigator::fExitNormal
private

Definition at line 417 of file G4Navigator.hh.

Referenced by ComputeStep(), PrintState(), ResetState(), RestoreSavedState(), and SetSavedState().

◆ fExitNormalGlobalFrame

G4ThreeVector G4Navigator::fExitNormalGlobalFrame
private

Definition at line 425 of file G4Navigator.hh.

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

◆ fGrandMotherExitNormal

G4ThreeVector G4Navigator::fGrandMotherExitNormal
private

◆ fHistory

G4NavigationHistory G4Navigator::fHistory
protected

◆ fLastLocatedPointLocal

G4ThreeVector G4Navigator::fLastLocatedPointLocal
private

◆ fLastMotherPhys

G4VPhysicalVolume* G4Navigator::fLastMotherPhys = nullptr
private

Definition at line 434 of file G4Navigator.hh.

Referenced by ComputeStep(), and ResetState().

◆ fLastStepEndPointLocal

G4ThreeVector G4Navigator::fLastStepEndPointLocal
protected

Definition at line 391 of file G4Navigator.hh.

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

◆ fLastStepWasZero

G4bool G4Navigator::fLastStepWasZero
private

◆ fLastTriedStepComputation

G4bool G4Navigator::fLastTriedStepComputation = false
private

◆ fLocatedOnEdge

G4bool G4Navigator::fLocatedOnEdge
private

Definition at line 477 of file G4Navigator.hh.

Referenced by ComputeStep(), InformLastStep(), LocateGlobalPointAndSetup(), and ResetState().

◆ fLocatedOutsideWorld

G4bool G4Navigator::fLocatedOutsideWorld
private

◆ fMinStep

G4double G4Navigator::fMinStep
protected

Definition at line 377 of file G4Navigator.hh.

Referenced by ComputeStep(), and G4Navigator().

◆ fnormalNav

G4NormalNavigation G4Navigator::fnormalNav
private

Definition at line 523 of file G4Navigator.hh.

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

◆ fNumberZeroSteps

G4int G4Navigator::fNumberZeroSteps
private

Definition at line 449 of file G4Navigator.hh.

Referenced by ComputeStep(), and ResetState().

◆ fparamNav

G4ParameterisedNavigation G4Navigator::fparamNav
private

◆ fpExternalNav

G4VExternalNavigation* G4Navigator::fpExternalNav = nullptr
private

◆ fPreviousSafety

G4double G4Navigator::fPreviousSafety
private

◆ fPreviousSftOrigin

G4ThreeVector G4Navigator::fPreviousSftOrigin
private

◆ fPushed

G4bool G4Navigator::fPushed = false
private

Definition at line 539 of file G4Navigator.hh.

Referenced by ComputeStep(), and ResetState().

◆ fpvoxelNav

G4VoxelNavigation* G4Navigator::fpvoxelNav
private

Definition at line 525 of file G4Navigator.hh.

Referenced by G4Navigator(), and ~G4Navigator().

◆ fpVoxelSafety

G4VoxelSafety* G4Navigator::fpVoxelSafety
private

Definition at line 533 of file G4Navigator.hh.

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

◆ fregularNav

G4RegularNavigation G4Navigator::fregularNav
private

Definition at line 531 of file G4Navigator.hh.

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

◆ freplicaNav

G4ReplicaNavigation G4Navigator::freplicaNav
private

◆ fSaveState

struct G4Navigator::G4SaveNavigatorState G4Navigator::fSaveState
private

Referenced by RestoreSavedState(), and SetSavedState().

◆ fSqTol

G4double G4Navigator::fSqTol
protected

Definition at line 377 of file G4Navigator.hh.

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

◆ fStepEndPoint

G4ThreeVector G4Navigator::fStepEndPoint
protected

◆ fTopPhysical

G4VPhysicalVolume* G4Navigator::fTopPhysical = nullptr
private

Definition at line 517 of file G4Navigator.hh.

◆ fValidExitNormal

G4bool G4Navigator::fValidExitNormal
private

◆ fVerbose

G4int G4Navigator::fVerbose = 0
protected

◆ fWarnPush

G4bool G4Navigator::fWarnPush = true
private

Definition at line 539 of file G4Navigator.hh.

Referenced by ComputeStep().

◆ fWasLimitedByGeometry

G4bool G4Navigator::fWasLimitedByGeometry = false
protected

◆ kCarTolerance

G4double G4Navigator::kCarTolerance
protected

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