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

#include <G4Navigator.hh>

+ Inheritance diagram for G4Navigator:

Public Member Functions

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

Protected Member Functions

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

Protected Attributes

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

Friends

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

Detailed Description

Definition at line 69 of file G4Navigator.hh.

Constructor & Destructor Documentation

◆ G4Navigator() [1/2]

G4Navigator::G4Navigator ( )

Definition at line 53 of file G4Navigator.cc.

54{
56 // Initialises also all
57 // - exit / entry flags
58 // - flags & variables for exit normals
59 // - zero step counters
60 // - blocked volume
61
62 if( fVerbose > 2 )
63 {
64 G4cout << " G4Navigator parameters: Action Threshold (No Zero Steps) = "
65 << fActionThreshold_NoZeroSteps
66 << " Abandon Threshold (No Zero Steps) = "
67 << fAbandonThreshold_NoZeroSteps << G4endl;
68 }
72
73 fregularNav.SetNormalNavigation( &fnormalNav );
74
75 fStepEndPoint = G4ThreeVector( kInfinity, kInfinity, kInfinity );
76 fLastStepEndPointLocal = G4ThreeVector( kInfinity, kInfinity, kInfinity );
77
78 fpVoxelSafety = new G4VoxelSafety();
79}
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:371
G4ThreeVector fStepEndPoint
Definition: G4Navigator.hh:382
G4int fVerbose
Definition: G4Navigator.hh:389
G4ThreeVector fLastStepEndPointLocal
Definition: G4Navigator.hh:385
G4double fSqTol
Definition: G4Navigator.hh:371
G4double kCarTolerance
Definition: G4Navigator.hh:371
void ResetStackAndState()
void SetNormalNavigation(G4NormalNavigation *fnormnav)
T sqr(const T &x)
Definition: templates.hh:128

◆ G4Navigator() [2/2]

G4Navigator::G4Navigator ( const G4Navigator )
delete

◆ ~G4Navigator()

G4Navigator::~G4Navigator ( )
virtual

Definition at line 85 of file G4Navigator.cc.

86{
87 delete fpVoxelSafety;
88 delete fpExternalNav;
89}

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

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

Referenced by G4SafetyHelper::CheckNextStep().

◆ CheckOverlapsIterative()

G4bool G4Navigator::CheckOverlapsIterative ( G4VPhysicalVolume vol)
protected

Definition at line 2313 of file G4Navigator.cc.

2314{
2315 // Check and report overlaps
2316 //
2317 G4bool foundOverlap = false;
2318 G4int nPoints = 300000, ntrials = 9, numOverlaps = 5;
2319 G4double trialLength = 1.0 * CLHEP::centimeter;
2320 while ( ntrials-- > 0 && !foundOverlap )
2321 {
2322 if ( fVerbose > 1 )
2323 {
2324 G4cout << " ** Running overlap checks in volume "
2325 << vol->GetName()
2326 << " with length = " << trialLength << G4endl;
2327 }
2328 foundOverlap = vol->CheckOverlaps(nPoints, trialLength,
2329 fVerbose, numOverlaps);
2330 trialLength *= 0.1;
2331 if ( trialLength <= 1.0e-5 ) { numOverlaps= 1;}
2332 }
2333 return foundOverlap;
2334}
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual G4bool CheckOverlaps(G4int res=1000, G4double tol=0., G4bool verbose=true, G4int errMax=1)
const G4String & GetName() const

Referenced by ComputeStep().

◆ Clone()

◆ ComputeLocalAxis()

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

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

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

Referenced by 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 764 of file G4Navigator.cc.

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

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

◆ CreateGRSSolid()

G4GRSSolid * G4Navigator::CreateGRSSolid ( ) const
inline

◆ CreateGRSVolume()

G4GRSVolume * G4Navigator::CreateGRSVolume ( ) const
inline

◆ CreateTouchableHistory() [1/2]

◆ CreateTouchableHistory() [2/2]

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

◆ CreateTouchableHistoryHandle()

G4TouchableHistoryHandle G4Navigator::CreateTouchableHistoryHandle ( ) const
virtual

Reimplemented in G4MultiNavigator.

Definition at line 2138 of file G4Navigator.cc.

2139{
2141}
G4ReferenceCountedHandle< G4TouchableHistory > G4TouchableHistoryHandle
G4TouchableHistory * CreateTouchableHistory() const

Referenced by G4VIntersectionLocator::LocateGlobalPointWithinVolumeAndCheck().

◆ EnableBestSafety()

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

◆ EnteredDaughterVolume()

G4bool G4Navigator::EnteredDaughterVolume ( ) const
inline

◆ ExitedMotherVolume()

G4bool G4Navigator::ExitedMotherVolume ( ) const
inline

◆ GetCurrentLocalCoordinate()

G4ThreeVector G4Navigator::GetCurrentLocalCoordinate ( ) const
inline

◆ GetDaughtersRegularStructureId()

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

◆ GetExternalNavigation()

◆ GetGlobalExitNormal()

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

Reimplemented in G4MultiNavigator, and G4ErrorPropagationNavigator.

Definition at line 1636 of file G4Navigator.cc.

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

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

◆ GetGlobalToLocalTransform()

const G4AffineTransform & G4Navigator::GetGlobalToLocalTransform ( ) const
inline

◆ GetLocalExitNormal()

G4ThreeVector G4Navigator::GetLocalExitNormal ( G4bool valid)
virtual

Reimplemented in G4MultiNavigator.

Definition at line 1383 of file G4Navigator.cc.

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

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

◆ GetLocalExitNormalAndCheck()

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

Reimplemented in G4MultiNavigator.

Definition at line 1600 of file G4Navigator.cc.

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

Referenced by GetGlobalExitNormal().

◆ GetLocalToGlobalTransform()

const G4AffineTransform G4Navigator::GetLocalToGlobalTransform ( ) const
inline

◆ GetMotherToDaughterTransform()

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

Definition at line 1549 of file G4Navigator.cc.

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

Referenced by GetLocalExitNormal().

◆ GetVerboseLevel()

G4int G4Navigator::GetVerboseLevel ( ) const
inline

Referenced by GetGlobalExitNormal().

◆ GetWorldVolume()

◆ IsActive()

G4bool G4Navigator::IsActive ( ) const
inline

◆ IsCheckModeActive()

G4bool G4Navigator::IsCheckModeActive ( ) const
inline

◆ LocateGlobalPointAndSetup()

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

Reimplemented in G4MultiNavigator.

Definition at line 126 of file G4Navigator.cc.

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

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

◆ LocateGlobalPointAndUpdateTouchable() [1/2]

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

◆ LocateGlobalPointAndUpdateTouchable() [2/2]

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

◆ LocateGlobalPointAndUpdateTouchableHandle()

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

◆ LocateGlobalPointWithinVolume()

void G4Navigator::LocateGlobalPointWithinVolume ( const G4ThreeVector position)
virtual

Reimplemented in G4MultiNavigator.

Definition at line 608 of file G4Navigator.cc.

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

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

2148{
2149 G4int oldcoutPrec = G4cout.precision(4);
2150 if( fVerbose >= 4 )
2151 {
2152 G4cout << "The current state of G4Navigator is: " << G4endl;
2153 G4cout << " ValidExitNormal= " << fValidExitNormal // << G4endl
2154 << " ExitNormal = " << fExitNormal // << G4endl
2155 << " Exiting = " << fExiting // << G4endl
2156 << " Entering = " << fEntering // << G4endl
2157 << " BlockedPhysicalVolume= " ;
2158 if (fBlockedPhysicalVolume==0)
2159 {
2160 G4cout << "None";
2161 }
2162 else
2163 {
2164 G4cout << fBlockedPhysicalVolume->GetName();
2165 }
2166 G4cout << G4endl
2167 << " BlockedReplicaNo = " << fBlockedReplicaNo // << G4endl
2168 << " LastStepWasZero = " << fLastStepWasZero // << G4endl
2169 << G4endl;
2170 }
2171 if( ( 1 < fVerbose) && (fVerbose < 4) )
2172 {
2173 G4cout << G4endl; // Make sure to line up
2174 G4cout << std::setw(30) << " ExitNormal " << " "
2175 << std::setw( 5) << " Valid " << " "
2176 << std::setw( 9) << " Exiting " << " "
2177 << std::setw( 9) << " Entering" << " "
2178 << std::setw(15) << " Blocked:Volume " << " "
2179 << std::setw( 9) << " ReplicaNo" << " "
2180 << std::setw( 8) << " LastStepZero " << " "
2181 << G4endl;
2182 G4cout << "( " << std::setw(7) << fExitNormal.x()
2183 << ", " << std::setw(7) << fExitNormal.y()
2184 << ", " << std::setw(7) << fExitNormal.z() << " ) "
2185 << std::setw( 5) << fValidExitNormal << " "
2186 << std::setw( 9) << fExiting << " "
2187 << std::setw( 9) << fEntering << " ";
2188 if ( fBlockedPhysicalVolume == nullptr )
2189 { G4cout << std::setw(15) << "None"; }
2190 else
2191 { G4cout << std::setw(15)<< fBlockedPhysicalVolume->GetName(); }
2192 G4cout << std::setw( 9) << fBlockedReplicaNo << " "
2193 << std::setw( 8) << fLastStepWasZero << " "
2194 << G4endl;
2195 }
2196 if( fVerbose > 2 )
2197 {
2198 G4cout.precision(8);
2199 G4cout << " Current Localpoint = " << fLastLocatedPointLocal << G4endl;
2200 G4cout << " PreviousSftOrigin = " << fPreviousSftOrigin << G4endl;
2201 G4cout << " PreviousSafety = " << fPreviousSafety << G4endl;
2202 }
2203 G4cout.precision(oldcoutPrec);
2204}
#define fEntering
double z() const
double x() const
double y() const

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

◆ RecheckDistanceToCurrentBoundary()

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

Definition at line 1973 of file G4Navigator.cc.

1979{
1980 G4ThreeVector localPosition = ComputeLocalPoint(aDisplacedGlobalPoint);
1981 G4ThreeVector localDirection = ComputeLocalAxis(aNewDirection);
1982 // G4double Step = kInfinity;
1983
1984 G4bool validExitNormal;
1985 G4ThreeVector exitNormal;
1986 // Check against mother solid
1987 G4VPhysicalVolume* motherPhysical = fHistory.GetTopVolume();
1988 G4LogicalVolume* motherLogical = motherPhysical->GetLogicalVolume();
1989
1990#ifdef CHECK_ORDER_OF_METHODS
1991 if( ! fLastTriedStepComputation )
1992 {
1993 G4Exception("G4Navigator::RecheckDistanceToCurrentBoundary()",
1994 "GeomNav0001", FatalException,
1995 "Method must be called after ComputeStep(), before call to LocateMethod.");
1996 }
1997#endif
1998
1999 EInside locatedDaug; // = kUndefined;
2000 G4double daughterStep = DBL_MAX;
2001 G4double daughterSafety = DBL_MAX;
2002
2003 if( fEnteredDaughter )
2004 {
2005 if( motherLogical->CharacteriseDaughters() == kReplica ) { return false; }
2006
2007 // Track arrived at boundary of a daughter volume at
2008 // the last call of ComputeStep().
2009 // In case the proposed displaced point is inside this daughter,
2010 // it must backtrack at least to the entry point.
2011 // NOTE: No check is made against other daughter volumes. It is
2012 // assumed that the proposed displacement is small enough that
2013 // this is not needed.
2014
2015 // Must check boundary of current daughter
2016 //
2017 G4VPhysicalVolume* candPhysical = fBlockedPhysicalVolume;
2018 G4LogicalVolume* candLogical = candPhysical->GetLogicalVolume();
2019 G4VSolid* candSolid = candLogical->GetSolid();
2020
2021 G4AffineTransform nextLevelTrf(candPhysical->GetRotation(),
2022 candPhysical->GetTranslation());
2023
2024 G4ThreeVector dgPosition = nextLevelTrf.TransformPoint(localPosition);
2025 G4ThreeVector dgDirection = nextLevelTrf.TransformAxis(localDirection);
2026 locatedDaug = candSolid->Inside(dgPosition);
2027
2028 if( locatedDaug == kInside )
2029 {
2030 // Reverse direction - and find first exit. ( Is it valid?)
2031 // Must backtrack
2032 G4double distanceBackOut =
2033 candSolid->DistanceToOut(dgPosition,
2034 - dgDirection, // Reverse direction
2035 true, &validExitNormal, &exitNormal);
2036 daughterStep= -distanceBackOut;
2037 // No check is made whether the particle could have arrived at
2038 // at this location without encountering another volume or
2039 // a different psurface of the current volume
2040 if( prNewSafety )
2041 {
2042 daughterSafety = candSolid->DistanceToOut(dgPosition);
2043 }
2044 }
2045 else
2046 {
2047 if( locatedDaug == kOutside )
2048 {
2049 // See whether it still intersects it
2050 //
2051 daughterStep = candSolid->DistanceToIn(dgPosition, dgDirection);
2052 if( prNewSafety )
2053 {
2054 daughterSafety = candSolid->DistanceToIn(dgPosition);
2055 }
2056 }
2057 else
2058 {
2059 // The point remains on the surface of candidate solid
2060 //
2061 daughterStep = 0.0;
2062 daughterSafety = 0.0;
2063 }
2064 }
2065
2066 // If trial point is in daughter (or on its surface) we have the
2067 // answer, the rest is not relevant
2068 //
2069 if( locatedDaug != kOutside )
2070 {
2071 *prDistance = daughterStep;
2072 if( prNewSafety ) { *prNewSafety= daughterSafety; }
2073 return true;
2074 }
2075 // If ever extended, so that some type of mother cut daughter,
2076 // this would change
2077 }
2078
2079 G4VSolid* motherSolid = motherLogical->GetSolid();
2080
2081 G4double motherStep = DBL_MAX, motherSafety = DBL_MAX;
2082
2083 // Check distance to boundary of mother
2084 //
2086 {
2087 EInside locatedMoth = motherSolid->Inside(localPosition);
2088
2089 if( locatedMoth == kInside )
2090 {
2091 motherSafety = motherSolid->DistanceToOut(localPosition);
2092 if( ProposedMove >= motherSafety )
2093 {
2094 motherStep = motherSolid->DistanceToOut(localPosition,
2095 localDirection,
2096 true, &validExitNormal, &exitNormal);
2097 }
2098 else
2099 {
2100 motherStep= ProposedMove;
2101 }
2102 }
2103 else if( locatedMoth == kOutside)
2104 {
2105 motherSafety= motherSolid->DistanceToIn(localPosition);
2106 if( ProposedMove >= motherSafety )
2107 {
2108 motherStep = -motherSolid->DistanceToIn(localPosition,
2109 -localDirection);
2110 }
2111 }
2112 else
2113 {
2114 motherSafety = 0.0;
2115 *prDistance = 0.0; // On surface - no move // motherStep;
2116 if( prNewSafety ) { *prNewSafety= motherSafety; }
2117 return false;
2118 }
2119 }
2120 else
2121 {
2122 return false;
2123 }
2124
2125 *prDistance = std::min( motherStep, daughterStep );
2126 if( prNewSafety )
2127 {
2128 *prNewSafety = std::min( motherSafety, daughterSafety );
2129 }
2130 return true;
2131}
EVolume CharacteriseDaughters() const
const G4ThreeVector GetTranslation() const
#define DBL_MAX
Definition: templates.hh:62

Referenced by G4PathFinder::RecheckDistanceToCurrentBoundary(), and G4SafetyHelper::RecheckDistanceToCurrentBoundary().

◆ ResetHierarchyAndLocate()

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

Reimplemented in G4MultiNavigator.

Definition at line 96 of file G4Navigator.cc.

99{
100 ResetState();
101 fHistory = *h.GetHistory();
103 fLastTriedStepComputation = false; // Redundant, but best
104 return LocateGlobalPointAndSetup(p, &direction, true, false);
105}
virtual void SetupHierarchy()
virtual void ResetState()
const G4NavigationHistory * GetHistory() const

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

◆ ResetStackAndState()

void G4Navigator::ResetStackAndState ( )
inline

◆ ResetState()

void G4Navigator::ResetState ( )
protectedvirtual

Reimplemented in G4MultiNavigator.

Definition at line 1282 of file G4Navigator.cc.

1283{
1284 fWasLimitedByGeometry = false;
1285 fEntering = false;
1286 fExiting = false;
1287 fLocatedOnEdge = false;
1288 fLastStepWasZero = false;
1289 fEnteredDaughter = false;
1290 fExitedMother = false;
1291 fPushed = false;
1292
1293 fValidExitNormal = false;
1294 fChangedGrandMotherRefFrame = false;
1295 fCalculatedExitNormal = false;
1296
1297 fExitNormal = G4ThreeVector(0,0,0);
1298 fGrandMotherExitNormal = G4ThreeVector(0,0,0);
1299 fExitNormalGlobalFrame = G4ThreeVector(0,0,0);
1300
1301 fPreviousSftOrigin = G4ThreeVector(0,0,0);
1302 fPreviousSafety = 0.0;
1303
1304 fNumberZeroSteps = 0;
1305
1306 fBlockedPhysicalVolume = nullptr;
1307 fBlockedReplicaNo = -1;
1308
1309 fLastLocatedPointLocal = G4ThreeVector( kInfinity, -kInfinity, 0.0 );
1310 fLocatedOutsideWorld = false;
1311
1312 fLastMotherPhys = nullptr;
1313}

Referenced by ResetHierarchyAndLocate().

◆ RestoreSavedState()

void G4Navigator::RestoreSavedState ( )
protected

Definition at line 710 of file G4Navigator.cc.

711{
712 fExitNormal = fSaveState.sExitNormal;
713 fValidExitNormal = fSaveState.sValidExitNormal;
714 fExiting = fSaveState.sExiting;
715 fEntering = fSaveState.sEntering;
716
717 fBlockedPhysicalVolume = fSaveState.spBlockedPhysicalVolume;
718 fBlockedReplicaNo = fSaveState.sBlockedReplicaNo;
719
720 fLastStepWasZero = fSaveState.sLastStepWasZero;
721
722 fLocatedOutsideWorld = fSaveState.sLocatedOutsideWorld;
723 fLastLocatedPointLocal = fSaveState.sLastLocatedPointLocal;
724 fEnteredDaughter = fSaveState.sEnteredDaughter;
725 fExitedMother = fSaveState.sExitedMother;
726 fWasLimitedByGeometry = fSaveState.sWasLimitedByGeometry;
727
728 // The 'expected' behaviour is to restore these too (fix 2014.05.26)
729 fPreviousSftOrigin = fSaveState.sPreviousSftOrigin;
730 fPreviousSafety = fSaveState.sPreviousSafety;
731}

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

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

Referenced by CheckNextStep(), and ComputeSafety().

◆ SetupHierarchy()

void G4Navigator::SetupHierarchy ( )
protectedvirtual

Reimplemented in G4MultiNavigator.

Definition at line 1323 of file G4Navigator.cc.

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

Referenced by ResetHierarchyAndLocate().

◆ SetVerboseLevel()

void G4Navigator::SetVerboseLevel ( G4int  level)
inline

Referenced by GetGlobalExitNormal().

◆ SetWorldVolume()

◆ SeverityOfZeroStepping()

G4int G4Navigator::SeverityOfZeroStepping ( G4int noZeroSteps) const
inline

◆ VolumeType()

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

Friends And Related Function Documentation

◆ operator<<

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

Definition at line 2340 of file G4Navigator.cc.

2341{
2342 // Old version did only the following:
2343 // os << "Current History: " << G4endl << n.fHistory;
2344 // Old behaviour is recovered for fVerbose = 0
2345
2346 // Adapted from G4Navigator::PrintState() const
2347
2348 G4int oldcoutPrec = os.precision(4);
2349 if( n.fVerbose >= 4 )
2350 {
2351 os << "The current state of G4Navigator is: " << G4endl;
2352 os << " ValidExitNormal= " << n.fValidExitNormal << G4endl
2353 << " ExitNormal = " << n.fExitNormal << G4endl
2354 << " Exiting = " << n.fExiting << G4endl
2355 << " Entering = " << n.fEntering << G4endl
2356 << " BlockedPhysicalVolume= " ;
2357 if (n.fBlockedPhysicalVolume==0)
2358 os << "None";
2359 else
2360 os << n.fBlockedPhysicalVolume->GetName();
2361 os << G4endl
2362 << " BlockedReplicaNo = " << n.fBlockedReplicaNo << G4endl
2363 << " LastStepWasZero = " << n.fLastStepWasZero << G4endl
2364 << G4endl;
2365 }
2366 if( ( 1 < n.fVerbose) && (n.fVerbose < 4) )
2367 {
2368 os << G4endl; // Make sure to line up
2369 os << std::setw(30) << " ExitNormal " << " "
2370 << std::setw( 5) << " Valid " << " "
2371 << std::setw( 9) << " Exiting " << " "
2372 << std::setw( 9) << " Entering" << " "
2373 << std::setw(15) << " Blocked:Volume " << " "
2374 << std::setw( 9) << " ReplicaNo" << " "
2375 << std::setw( 8) << " LastStepZero " << " "
2376 << G4endl;
2377 os << "( " << std::setw(7) << n.fExitNormal.x()
2378 << ", " << std::setw(7) << n.fExitNormal.y()
2379 << ", " << std::setw(7) << n.fExitNormal.z() << " ) "
2380 << std::setw( 5) << n.fValidExitNormal << " "
2381 << std::setw( 9) << n.fExiting << " "
2382 << std::setw( 9) << n.fEntering << " ";
2383 if ( n.fBlockedPhysicalVolume==0 )
2384 { os << std::setw(15) << "None"; }
2385 else
2386 { os << std::setw(15)<< n.fBlockedPhysicalVolume->GetName(); }
2387 os << std::setw( 9) << n.fBlockedReplicaNo << " "
2388 << std::setw( 8) << n.fLastStepWasZero << " "
2389 << G4endl;
2390 }
2391 if( n.fVerbose > 2 )
2392 {
2393 os.precision(8);
2394 os << " Current Localpoint = " << n.fLastLocatedPointLocal << G4endl;
2395 os << " PreviousSftOrigin = " << n.fPreviousSftOrigin << G4endl;
2396 os << " PreviousSafety = " << n.fPreviousSafety << G4endl;
2397 }
2398 if( n.fVerbose > 3 || n.fVerbose == 0 )
2399 {
2400 os << "Current History: " << G4endl << n.fHistory;
2401 }
2402
2403 os.precision(oldcoutPrec);
2404 return os;
2405}

Member Data Documentation

◆ fEnteredDaughter

◆ fExitedMother

◆ fHistory

◆ fLastStepEndPointLocal

G4ThreeVector G4Navigator::fLastStepEndPointLocal
protected

Definition at line 385 of file G4Navigator.hh.

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

◆ fMinStep

G4double G4Navigator::fMinStep
protected

Definition at line 371 of file G4Navigator.hh.

Referenced by ComputeStep(), and G4Navigator().

◆ fSqTol

G4double G4Navigator::fSqTol
protected

Definition at line 371 of file G4Navigator.hh.

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

◆ fStepEndPoint

G4ThreeVector G4Navigator::fStepEndPoint
protected

Definition at line 382 of file G4Navigator.hh.

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

◆ fVerbose

◆ fWasLimitedByGeometry

◆ kCarTolerance

G4double G4Navigator::kCarTolerance
protected

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