Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ITMultiNavigator Class Reference

#include <G4ITMultiNavigator.hh>

Inheritance diagram for G4ITMultiNavigator:
Collaboration diagram for G4ITMultiNavigator:

Public Member Functions

 G4ITMultiNavigator ()
 
 ~G4ITMultiNavigator ()
 
G4double ComputeStep (const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
 
G4double ObtainFinalStep (G4int navigatorId, G4double &pNewSafety, G4double &minStepLast, ELimited &limitedStep)
 
void PrepareNavigators ()
 
void PrepareNewTrack (const G4ThreeVector position, const G4ThreeVector direction)
 
G4VPhysicalVolumeResetHierarchyAndLocate (const G4ThreeVector &point, const G4ThreeVector &direction, const G4TouchableHistory &h)
 
G4VPhysicalVolumeLocateGlobalPointAndSetup (const G4ThreeVector &point, const G4ThreeVector *direction=0, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
 
void LocateGlobalPointWithinVolume (const G4ThreeVector &position)
 
G4double ComputeSafety (const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=false)
 
G4TouchableHistoryHandle CreateTouchableHistoryHandle () const
 
virtual G4ThreeVector GetLocalExitNormal (G4bool *obtained)
 
virtual G4ThreeVector GetLocalExitNormalAndCheck (const G4ThreeVector &CurrentE_Point, G4bool *obtained)
 
virtual G4ThreeVector GetGlobalExitNormal (const G4ThreeVector &CurrentE_Point, G4bool *obtained)
 
G4ITNavigator * GetNavigator (G4int n) const
 
- Public Member Functions inherited from G4TrackStateDependent< G4ITMultiNavigator >
virtual ~G4TrackStateDependent ()
 
virtual void SetTrackState (G4shared_ptr< StateType > state)
 
virtual G4VTrackStateHandle PopTrackState ()
 
virtual G4VTrackStateHandle GetTrackState () const
 
virtual StateTypeHandle GetConcreteTrackState () const
 
virtual void LoadTrackState (G4TrackStateManager &manager)
 
virtual void SaveTrackState (G4TrackStateManager &manager)
 
virtual void NewTrackState ()
 
virtual StateTypeHandle CreateTrackState () const
 
virtual void ResetTrackState ()
 
- Public Member Functions inherited from G4VTrackStateDependent
 G4VTrackStateDependent ()
 
virtual ~G4VTrackStateDependent ()
 

Protected Member Functions

void ResetState ()
 
void SetupHierarchy ()
 
void WhichLimited ()
 
void PrintLimited ()
 
void CheckMassWorld ()
 
- Protected Member Functions inherited from G4TrackStateDependent< G4ITMultiNavigator >
 G4TrackStateDependent ()
 

Friends

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

Additional Inherited Members

- Public Types inherited from G4TrackStateDependent< G4ITMultiNavigator >
typedef G4ITMultiNavigator ClassType
 
typedef G4TrackState
< G4ITMultiNavigator
StateType
 
typedef G4shared_ptr< StateTypeStateTypeHandle
 
- Protected Attributes inherited from G4TrackStateDependent< G4ITMultiNavigator >
StateTypeHandle fpTrackState
 

Detailed Description

Definition at line 144 of file G4ITMultiNavigator.hh.

Constructor & Destructor Documentation

G4ITMultiNavigator::G4ITMultiNavigator ( )

Definition at line 68 of file G4ITMultiNavigator.cc.

69  : G4ITNavigator(), fLastMassWorld(0)
70 {
71  fNoActiveNavigators= 0;
72 
73  for(G4int num=0; num< fMaxNav; ++num )
74  {
75  fpNavigator[num] = 0;
76  }
77 
79 
80  G4ITNavigator* massNav= pTransportManager->GetNavigatorForTracking();
81  if( massNav )
82  {
83  G4VPhysicalVolume* pWorld= massNav->GetWorldVolume();
84  if( pWorld )
85  {
86  SetWorldVolume( pWorld );
87  fLastMassWorld = pWorld;
88  }
89  }
90 }
int G4int
Definition: G4Types.hh:78
static G4ITTransportationManager * GetTransportationManager()
G4ITNavigator * GetNavigatorForTracking() const

Here is the call graph for this function:

G4ITMultiNavigator::~G4ITMultiNavigator ( )

Definition at line 92 of file G4ITMultiNavigator.cc.

93 {
94 }

Member Function Documentation

void G4ITMultiNavigator::CheckMassWorld ( )
protected

Definition at line 641 of file G4ITMultiNavigator.cc.

642 {
643  G4VPhysicalVolume* navTrackWorld=
644  pTransportManager->GetNavigatorForTracking()->GetWorldVolume();
645 
646  if( navTrackWorld != fLastMassWorld )
647  {
648  G4Exception( "G4ITMultiNavigator::CheckMassWorld()",
649  "GeomNav0003", FatalException,
650  "Mass world pointer has been changed." );
651  }
652 }
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4ITNavigator * GetNavigatorForTracking() const

Here is the call graph for this function:

G4double G4ITMultiNavigator::ComputeSafety ( const G4ThreeVector globalpoint,
const G4double  pProposedMaxLength = DBL_MAX,
const G4bool  keepState = false 
)

Definition at line 431 of file G4ITMultiNavigator.cc.

434 {
435  // Recompute safety for the relevant point
436 
437  G4double minSafety = kInfinity, safety = kInfinity;
438 
439  std::vector<G4ITNavigator*>::iterator pNavigatorIter;
440  pNavigatorIter= pTransportManager-> GetActiveNavigatorsIterator();
441 
442  for( G4int num=0; num< fNoActiveNavigators; ++pNavigatorIter,++num )
443  {
444  safety = (*pNavigatorIter)->ComputeSafety( position, maxDistance, state);
445  if( safety < minSafety ) { minSafety = safety; }
446  }
447 
449  fMinSafety_atSafLocation = minSafety;
450 
451 #ifdef G4DEBUG_NAVIGATION
452  if( fVerbose > 1 )
453  {
454  G4cout << " G4ITMultiNavigator::ComputeSafety - returns: "
455  << minSafety << ", at location: " << position << G4endl;
456  }
457 #endif
458  return minSafety;
459 }
#define fSafetyLocation
static const G4double kInfinity
Definition: geomdefs.hh:42
#define fMinSafety_atSafLocation
int G4int
Definition: G4Types.hh:78
#define position
Definition: xmlparse.cc:622
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double G4ITMultiNavigator::ComputeStep ( const G4ThreeVector pGlobalPoint,
const G4ThreeVector pDirection,
const G4double  pCurrentProposedStepLength,
G4double pNewSafety 
)

Definition at line 96 of file G4ITMultiNavigator.cc.

100 {
101  G4double safety= 0.0, step=0.0;
102  G4double minSafety= kInfinity, minStep= kInfinity;
103 
104  fNoLimitingStep= -1;
105  fIdNavLimiting= -1; // Reset for new step
106 
107 #ifdef G4DEBUG_NAVIGATION
108  if( fVerbose > 2 )
109  {
110  G4cout << " G4ITMultiNavigator::ComputeStep : entered " << G4endl;
111  G4cout << " Input position= " << pGlobalPoint
112  << " direction= " << pDirection << G4endl;
113  G4cout << " Requested step= " << proposedStepLength << G4endl;
114  }
115 #endif
116 
117  std::vector<G4ITNavigator*>::iterator pNavigatorIter;
118 
119  pNavigatorIter= pTransportManager-> GetActiveNavigatorsIterator();
120 
121  G4ThreeVector initialPosition = pGlobalPoint;
122  G4ThreeVector initialDirection= pDirection;
123 
124  for( G4int num=0; num< fNoActiveNavigators; ++pNavigatorIter,++num )
125  {
126  safety= kInfinity;
127 
128  step= (*pNavigatorIter)->ComputeStep( initialPosition,
129  initialDirection,
130  proposedStepLength,
131  safety );
132  if( safety < minSafety ){ minSafety = safety; }
133  if( step < minStep ) { minStep= step; }
134 
135  fCurrentStepSize[num] = step;
136  fNewSafety[num]= safety;
137  // This is currently the safety from the last sub-step
138 
139 #ifdef G4DEBUG_NAVIGATION
140  if( fVerbose > 2 )
141  {
142  G4cout << "G4ITMultiNavigator::ComputeStep : Navigator ["
143  << num << "] -- step size " << step
144  << " safety= " << safety << G4endl;
145  }
146 #endif
147  }
148 
149  // Save safety value, related position
150  //
151  fPreStepLocation = initialPosition;
152  fMinSafety_PreStepPt = minSafety;
153  fMinStep = minStep;
154 
155  if( fMinStep == kInfinity )
156  {
157  fTrueMinStep = proposedStepLength; // Use this below for endpoint !!
158  }
159  else
160  {
161  fTrueMinStep = minStep;
162  }
163 
164 #ifdef G4DEBUG_NAVIGATION
165  if( fVerbose > 1 )
166  {
167  G4ThreeVector endPosition = initialPosition+fTrueMinStep*initialDirection;
168 
169  G4int oldPrec = G4cout.precision(8);
170  G4cout << "G4ITMultiNavigator::ComputeStep : "
171  << " initialPosition = " << initialPosition
172  << " and endPosition = " << endPosition << G4endl;
173  G4cout.precision( oldPrec );
174  }
175 #endif
176 
177  pNewSafety = minSafety;
178 
179  this->WhichLimited();
180 
181 #ifdef G4DEBUG_NAVIGATION
182  if( fVerbose > 2 )
183  {
184  G4cout << " G4ITMultiNavigator::ComputeStep : exits returning "
185  << minStep << G4endl;
186  }
187 #endif
188 
189  return minStep; // must return kInfinity if do not limit step
190 }
static const G4double kInfinity
Definition: geomdefs.hh:42
#define fIdNavLimiting
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
#define fTrueMinStep
#define fPreStepLocation
#define fMinStep
#define fCurrentStepSize
#define fMinSafety_PreStepPt
#define G4endl
Definition: G4ios.hh:61
#define fNewSafety
double G4double
Definition: G4Types.hh:76
#define fNoLimitingStep

Here is the call graph for this function:

G4TouchableHistoryHandle G4ITMultiNavigator::CreateTouchableHistoryHandle ( ) const

Definition at line 464 of file G4ITMultiNavigator.cc.

465 {
466  G4Exception( "G4ITMultiNavigator::CreateTouchableHistoryHandle()",
467  "GeomNav0001", FatalException,
468  "Getting a touchable from G4ITMultiNavigator is not defined.");
469 
470  G4TouchableHistory* touchHist;
471  touchHist= fpNavigator[0] -> CreateTouchableHistory();
472 
473  G4VPhysicalVolume* locatedVolume= fLocatedVolume[0];
474  if( locatedVolume == 0 )
475  {
476  // Workaround to ensure that the touchable is fixed !! // TODO: fix
477  //
478  touchHist->UpdateYourself( locatedVolume, touchHist->GetHistory() );
479  }
480 
481  return G4TouchableHistoryHandle(touchHist);
482 }
#define fLocatedVolume
G4ReferenceCountedHandle< G4TouchableHistory > G4TouchableHistoryHandle
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void UpdateYourself(G4VPhysicalVolume *pPhysVol, const G4NavigationHistory *history=0)
const G4NavigationHistory * GetHistory() const

Here is the call graph for this function:

G4ThreeVector G4ITMultiNavigator::GetGlobalExitNormal ( const G4ThreeVector CurrentE_Point,
G4bool obtained 
)
virtual

Definition at line 696 of file G4ITMultiNavigator.cc.

698 {
699  G4ThreeVector normalGlobalCrd(0.0, 0.0, 0.0);
700  G4bool isObtained= false;
701  // These default values will be used if fNoLimitingStep== 0
702  G4int firstNavigatorId= -1;
703  G4bool oneObtained= false;
704 
705  if( fNoLimitingStep==1 )
706  {
707  // Only message the Navigator which limited the step!
708  normalGlobalCrd= fpNavigator[ fIdNavLimiting ]->GetGlobalExitNormal( argPoint, &isObtained);
709  *argpObtained= isObtained;
710  }
711  else
712  {
713  if( fNoLimitingStep > 1 )
714  {
715  std::vector<G4ITNavigator*>::iterator pNavIter=
716  pTransportManager->GetActiveNavigatorsIterator();
717 
718  for ( G4int num=0; num< fNoActiveNavigators ; ++pNavIter,++num )
719  {
720  G4ThreeVector oneNormal;
721  if( fLimitTruth[ num ] ) // Did this geometry limit the step ?
722  {
723  G4ThreeVector newNormal= (*pNavIter)-> GetGlobalExitNormal( argPoint, &oneObtained );
724  if( oneObtained )
725  {
726  // Keep first one - only if it is valid (ie not null)
727  if( !isObtained && (newNormal.mag2() != 0.0) )
728  {
729  normalGlobalCrd= newNormal;
730  isObtained = oneObtained;
731  firstNavigatorId= num;
732  }else{
733  // Check for clash
734  G4double dotNewPrevious= newNormal.dot( normalGlobalCrd );
735  G4double productMagSq= normalGlobalCrd.mag2() * newNormal.mag2();
736  if( productMagSq > 0.0 )
737  {
738  G4double productMag= std::sqrt( productMagSq );
739  dotNewPrevious /= productMag; // Normalise
740  if( dotNewPrevious < (1 - perThousand) )
741  {
742  *argpObtained= false;
743 
744  if( fVerbose > 2 ) // dotNewPrevious <= 0.0 )
745  {
746  std::ostringstream message;
747  message << "Clash of Normal from different Navigators!" << G4endl
748  << " Previous Navigator Id = " << firstNavigatorId << G4endl
749  << " Current Navigator Id = " << num << G4endl;
750  message << " Dot product of 2 normals = " << dotNewPrevious << G4endl;
751  message << " Normal (previous) = " << normalGlobalCrd << G4endl;
752  message << " Normal (current) = " << newNormal << G4endl;
753  G4Exception("G4ITMultiNavigator::GetGlobalExitNormal()", "GeomNav0002",
754  JustWarning, message);
755  }
756  }
757  else
758  {
759  // Close agreement - Do not change
760  }
761  }
762  }
763  }
764  }
765  } // end for over the Navigators
766 
767  // Report if no Normal was obtained
768  if( !oneObtained )
769  {
770  std::ostringstream message;
771  message << "No Normal obtained despite having " << fNoLimitingStep
772  << " candidate Navigators limiting the step!" << G4endl;
773  G4Exception("G4ITMultiNavigator::GetGlobalExitNormal()", "GeomNav0002",
774  JustWarning, message);
775  }
776 
777  } // end if ( fNoLimiting > 1 )
778  } // end else
779 
780  *argpObtained= isObtained;
781  return normalGlobalCrd;
782 }
double dot(const Hep3Vector &) const
#define fIdNavLimiting
int G4int
Definition: G4Types.hh:78
#define fLimitTruth
bool G4bool
Definition: G4Types.hh:79
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual G4ThreeVector GetGlobalExitNormal(const G4ThreeVector &CurrentE_Point, G4bool *obtained)
std::vector< G4ITNavigator * >::iterator GetActiveNavigatorsIterator()
double mag2() const
#define G4endl
Definition: G4ios.hh:61
static constexpr double perThousand
Definition: G4SIunits.hh:333
double G4double
Definition: G4Types.hh:76
#define fNoLimitingStep

Here is the call graph for this function:

G4ThreeVector G4ITMultiNavigator::GetLocalExitNormal ( G4bool obtained)
virtual

Definition at line 787 of file G4ITMultiNavigator.cc.

788 {
789  // If it is the mass navigator, then expect
790  G4ThreeVector normalGlobalCrd(0.0, 0.0, 0.0);
791  G4bool isObtained= false;
792  // These default values will be used if fNoLimitingStep== 0
793 
794  if( fNoLimitingStep==1 )
795  {
796  // Only message the Navigator which limited the step!
797  normalGlobalCrd= fpNavigator[ fIdNavLimiting ]->GetLocalExitNormal( &isObtained);
798  *argpObtained= isObtained;
799 
800  static G4ThreadLocal G4int numberWarnings= 0;
801  G4int noWarningsStart= 10, noModuloWarnings=100;
802  numberWarnings++;
803  if( (numberWarnings < noWarningsStart ) || (numberWarnings%noModuloWarnings==0) )
804  {
805  std::ostringstream message;
806  message << "Cannot obtain normal in local coordinates of two or more coordinate systems." << G4endl;
807  G4Exception("G4ITMultiNavigator::GetGlobalExitNormal()", "GeomNav0002",
808  JustWarning, message);
809  }
810  }
811  else
812  {
813  if( fNoLimitingStep > 1 )
814  {
815  // Does not make sense - cannot obtain *local* normal in several coordinate systems
816  std::ostringstream message;
817  message << "Cannot obtain normal in local coordinates of two or more coordinate systems." << G4endl;
818  G4Exception("G4ITMultiNavigator::GetGlobalExitNormal()", "GeomNav0002",
819  FatalException, message);
820  }
821  }
822 
823  *argpObtained= isObtained;
824  return normalGlobalCrd;
825 }
#define fIdNavLimiting
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
bool G4bool
Definition: G4Types.hh:79
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
#define fNoLimitingStep

Here is the call graph for this function:

Here is the caller graph for this function:

G4ThreeVector G4ITMultiNavigator::GetLocalExitNormalAndCheck ( const G4ThreeVector CurrentE_Point,
G4bool obtained 
)
virtual

Definition at line 831 of file G4ITMultiNavigator.cc.

833 {
834  return G4ITMultiNavigator::GetLocalExitNormal( obtained);
835 }
virtual G4ThreeVector GetLocalExitNormal(G4bool *obtained)

Here is the call graph for this function:

G4ITNavigator* G4ITMultiNavigator::GetNavigator ( G4int  n) const
inline

Definition at line 230 of file G4ITMultiNavigator.hh.

231  {
232  if( (n>fNoActiveNavigators)||(n<0))
233  { n=0;}
234  return fpNavigator[n];
235  }
G4VPhysicalVolume * G4ITMultiNavigator::LocateGlobalPointAndSetup ( const G4ThreeVector point,
const G4ThreeVector direction = 0,
const G4bool  pRelativeSearch = true,
const G4bool  ignoreDirection = true 
)

Definition at line 319 of file G4ITMultiNavigator.cc.

323 {
324  // Locate the point in each geometry
325 
326  G4ThreeVector direction(0.0, 0.0, 0.0);
327  G4bool relative = pRelativeSearch;
328  std::vector<G4ITNavigator*>::iterator pNavIter
329  = pTransportManager->GetActiveNavigatorsIterator();
330 
331  if( pDirection ) { direction = *pDirection; }
332 
333 #ifdef G4DEBUG_NAVIGATION
334  if( fVerbose > 2 )
335  {
336  G4cout << " Entered G4ITMultiNavigator::LocateGlobalPointAndSetup() "
337  << G4endl;
338  G4cout << " Locating at position: " << position
339  << ", with direction: " << direction << G4endl
340  << " Relative: " << relative
341  << ", ignore direction: " << ignoreDirection << G4endl;
342  G4cout << " Number of active navigators: " << fNoActiveNavigators
343  << G4endl;
344  }
345 #endif
346 
347  for ( G4int num=0; num< fNoActiveNavigators ; ++pNavIter,++num )
348  {
349  if( fWasLimitedByGeometry && fLimitTruth[num] )
350  {
351  (*pNavIter)->SetGeometricallyLimitedStep();
352  }
353 
354  G4VPhysicalVolume *pLocated
355  = (*pNavIter)->LocateGlobalPointAndSetup( position, &direction,
356  relative, ignoreDirection );
357  // Set the state related to the location
358  //
359  fLocatedVolume[num] = pLocated;
360 
361  // Clear state related to the step
362  //
363  fLimitedStep[num] = kDoNot;
364  fCurrentStepSize[num] = 0.0;
365  fLimitTruth[ num ] = false; // Always clear on locating (see Navigator)
366 
367 #ifdef G4DEBUG_NAVIGATION
368  if( fVerbose > 2 )
369  {
370  G4cout << " Located in world: " << num << ", at: " << position << G4endl
371  << " Used geomLimStp: " << fLimitTruth[num]
372  << ", found in volume: " << pLocated << G4endl;
373  G4cout << " Name = '" ;
374  if( pLocated )
375  {
376  G4cout << pLocated->GetName() << "'";
377  G4cout << " - CopyNo= " << pLocated->GetCopyNo();
378  }
379  else
380  {
381  G4cout << "Null' Id: Not-Set ";
382  }
383  G4cout << G4endl;
384  }
385 #endif
386  }
387 
388  fWasLimitedByGeometry = false; // Clear on locating
389  G4VPhysicalVolume* volMassLocated= fLocatedVolume[0];
390 
391  return volMassLocated;
392 }
#define fWasLimitedByGeometry
#define fLocatedVolume
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
#define fLimitTruth
const G4String & GetName() const
bool G4bool
Definition: G4Types.hh:79
#define fLimitedStep
std::vector< G4ITNavigator * >::iterator GetActiveNavigatorsIterator()
virtual G4int GetCopyNo() const =0
#define fCurrentStepSize
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ITMultiNavigator::LocateGlobalPointWithinVolume ( const G4ThreeVector position)

Definition at line 397 of file G4ITMultiNavigator.cc.

398 {
399  // Relocate the point in each geometry
400 
401  std::vector<G4ITNavigator*>::iterator pNavIter
402  = pTransportManager->GetActiveNavigatorsIterator();
403 
404 #ifdef G4DEBUG_NAVIGATION
405  if( fVerbose > 2 )
406  {
407  G4cout << " Entered G4ITMultiNavigator::ReLocate() " << G4endl
408  << " Re-locating at position: " << position << G4endl;
409  }
410 #endif
411 
412  for ( G4int num=0; num< fNoActiveNavigators ; ++pNavIter,++num )
413  {
414  // ... none limited the step
415 
416  (*pNavIter)->LocateGlobalPointWithinVolume( position );
417 
418  // Clear state related to the step
419  //
420  fLimitedStep[num] = kDoNot;
421  fCurrentStepSize[num] = 0.0;
422 
423  fLimitTruth[ num ] = false; // Always clear on locating (see Navigator)
424  }
425  fWasLimitedByGeometry = false; // Clear on locating
427 }
#define fWasLimitedByGeometry
int G4int
Definition: G4Types.hh:78
#define fLastLocatedPosition
#define position
Definition: xmlparse.cc:622
G4GLOB_DLL std::ostream G4cout
#define fLimitTruth
#define fLimitedStep
std::vector< G4ITNavigator * >::iterator GetActiveNavigatorsIterator()
#define fCurrentStepSize
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

G4double G4ITMultiNavigator::ObtainFinalStep ( G4int  navigatorId,
G4double pNewSafety,
G4double minStepLast,
ELimited limitedStep 
)

Definition at line 195 of file G4ITMultiNavigator.cc.

199 {
200  if( navigatorId > fNoActiveNavigators )
201  {
202  std::ostringstream message;
203  message << "Bad Navigator Id!" << G4endl
204  << " Navigator Id = " << navigatorId
205  << " No Active = " << fNoActiveNavigators << ".";
206  G4Exception("G4ITMultiNavigator::ObtainFinalStep()", "GeomNav0002",
207  FatalException, message);
208  }
209 
210  // Prepare the information to return
211  //
212  pNewSafety = fNewSafety[ navigatorId ];
213  limitedStep = fLimitedStep[ navigatorId ];
214  minStep= fMinStep;
215 
216 #ifdef G4DEBUG_NAVIGATION
217  if( fVerbose > 1 )
218  {
219  G4cout << " G4ITMultiNavigator::ComputeStep returns "
220  << fCurrentStepSize[ navigatorId ]
221  << " for Navigator " << navigatorId
222  << " Limited step = " << limitedStep
223  << " Safety(mm) = " << pNewSafety / mm << G4endl;
224  }
225 #endif
226 
227  return fCurrentStepSize[ navigatorId ];
228 }
static constexpr double mm
Definition: G4SIunits.hh:115
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define fMinStep
#define fLimitedStep
#define fCurrentStepSize
#define G4endl
Definition: G4ios.hh:61
#define fNewSafety

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ITMultiNavigator::PrepareNavigators ( )

Definition at line 253 of file G4ITMultiNavigator.cc.

254 {
255  // Key purposes:
256  // - Check and cache set of active navigators
257  // - Reset state for new track
258 
259 #ifdef G4DEBUG_NAVIGATION
260  if( fVerbose > 1 )
261  {
262  G4cout << " Entered G4ITMultiNavigator::PrepareNavigators() " << G4endl;
263  }
264 #endif
265 
266  // Message the transportation-manager to find active navigators
267 
268  std::vector<G4ITNavigator*>::iterator pNavigatorIter;
269  fNoActiveNavigators= pTransportManager-> GetNoActiveNavigators();
270 
271  if( fNoActiveNavigators > fMaxNav )
272  {
273  std::ostringstream message;
274  message << "Too many active Navigators / worlds !" << G4endl
275  << " Active Navigators (worlds): "
276  << fNoActiveNavigators << G4endl
277  << " which is more than the number allowed: "
278  << fMaxNav << " !";
279  G4Exception("G4ITMultiNavigator::PrepareNavigators()", "GeomNav0002",
280  FatalException, message);
281  }
282 
283  pNavigatorIter= pTransportManager-> GetActiveNavigatorsIterator();
284  for( G4int num=0; num< fNoActiveNavigators; ++pNavigatorIter,++num )
285  {
286  fpNavigator[num] = *pNavigatorIter;
287  fLimitTruth[num] = false;
288  fLimitedStep[num] = kDoNot;
289  fCurrentStepSize[num] = 0.0;
290  fLocatedVolume[num] = 0;
291  }
292  fWasLimitedByGeometry = false;
293 
294  // Check the world volume of the mass navigator
295  // in case a call to SetWorldVolume() changed it
296 
297  G4VPhysicalVolume* massWorld = GetWorldVolume();
298 
299  if( (massWorld != fLastMassWorld) && (massWorld!=0) )
300  {
301  // Pass along change to Mass Navigator
302  fpNavigator[0] -> SetWorldVolume( massWorld );
303 
304 #ifdef G4DEBUG_NAVIGATION
305  if( fVerbose > 0 )
306  {
307  G4cout << " G4ITMultiNavigator::PrepareNavigators() changed world volume "
308  << " for mass geometry to " << massWorld->GetName() << G4endl;
309  }
310 #endif
311 
312  fLastMassWorld = massWorld;
313  }
314 }
#define fWasLimitedByGeometry
#define fLocatedVolume
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
#define fLimitTruth
const G4String & GetName() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define fLimitedStep
#define fCurrentStepSize
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ITMultiNavigator::PrepareNewTrack ( const G4ThreeVector  position,
const G4ThreeVector  direction 
)

Definition at line 232 of file G4ITMultiNavigator.cc.

234 {
235 #ifdef G4DEBUG_NAVIGATION
236  if( fVerbose > 1 )
237  {
238  G4cout << " Entered G4ITMultiNavigator::PrepareNewTrack() " << G4endl;
239  }
240 #endif
241 
243 
244  LocateGlobalPointAndSetup( position, &direction, false, false );
245  //
246  // The first location for each Navigator must be non-relative
247  // or else call ResetStackAndState() for each Navigator
248  // Use direction to get correct side of boundary (ignore dir= false)
249 }
G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=0, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void G4ITMultiNavigator::PrintLimited ( )
protected

Definition at line 544 of file G4ITMultiNavigator.cc.

545 {
546  // Report results -- for checking
547 
548  static const G4String StrDoNot("DoNot"), StrUnique("Unique"),
549  StrUndefined("Undefined"),
550  StrSharedTransport("SharedTransport"),
551  StrSharedOther("SharedOther");
552  G4cout << "### G4ITMultiNavigator::PrintLimited() reports: " << G4endl;
553  G4cout << " Minimum step (true): " << fTrueMinStep
554  << ", reported min: " << fMinStep << G4endl;
555 
556 #ifdef G4DEBUG_NAVIGATION
557  if(fVerbose>=2)
558  {
559  G4cout << std::setw(5) << " NavId" << " "
560  << std::setw(12) << " step-size " << " "
561  << std::setw(12) << " raw-size " << " "
562  << std::setw(12) << " pre-safety " << " "
563  << std::setw(15) << " Limited / flag" << " "
564  << std::setw(15) << " World " << " "
565  << G4endl;
566  }
567 #endif
568 
569  for ( G4int num= 0; num < fNoActiveNavigators; num++ )
570  {
571  G4double rawStep = fCurrentStepSize[num];
572  G4double stepLen = fCurrentStepSize[num];
573  if( stepLen > fTrueMinStep )
574  {
575  stepLen = fTrueMinStep; // did not limit (went as far as asked)
576  }
577  G4int oldPrec= G4cout.precision(9);
578 
579  G4cout << std::setw(5) << num << " "
580  << std::setw(12) << stepLen << " "
581  << std::setw(12) << rawStep << " "
582  << std::setw(12) << fNewSafety[num] << " "
583  << std::setw(5) << (fLimitTruth[num] ? "YES" : " NO") << " ";
584  G4String limitedStr;
585  switch ( fLimitedStep[num] )
586  {
587  case kDoNot : limitedStr= StrDoNot; break;
588  case kUnique : limitedStr = StrUnique; break;
589  case kSharedTransport: limitedStr= StrSharedTransport; break;
590  case kSharedOther : limitedStr = StrSharedOther; break;
591  default : limitedStr = StrUndefined; break;
592  }
593  G4cout << " " << std::setw(15) << limitedStr << " ";
594  G4cout.precision(oldPrec);
595 
596  G4ITNavigator *pNav= fpNavigator[ num ];
597  G4String WorldName( "Not-Set" );
598  if (pNav)
599  {
600  G4VPhysicalVolume *pWorld= pNav->GetWorldVolume();
601  if( pWorld )
602  {
603  WorldName = pWorld->GetName();
604  }
605  }
606  G4cout << " " << WorldName ;
607  G4cout << G4endl;
608  }
609 }
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
#define fLimitTruth
const G4String & GetName() const
#define fTrueMinStep
#define fMinStep
#define fLimitedStep
#define fCurrentStepSize
#define G4endl
Definition: G4ios.hh:61
#define fNewSafety
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Definition at line 657 of file G4ITMultiNavigator.cc.

660 {
661  // Reset geometry for all -- and use the touchable for the mass history
662 
663  G4VPhysicalVolume* massVolume=0;
664  G4ITNavigator* pMassNavigator= fpNavigator[0];
665 
666  if( pMassNavigator )
667  {
668  massVolume= pMassNavigator->ResetHierarchyAndLocate( point, direction,
669  MassHistory);
670  }
671  else
672  {
673  G4Exception("G4ITMultiNavigator::ResetHierarchyAndLocate()",
674  "GeomNav0002", FatalException,
675  "Cannot reset hierarchy before navigators are initialised.");
676  }
677 
678  std::vector<G4ITNavigator*>::iterator pNavIter=
679  pTransportManager->GetActiveNavigatorsIterator();
680 
681  for ( G4int num=0; num< fNoActiveNavigators ; ++pNavIter,++num )
682  {
683  G4bool relativeSearch, ignoreDirection;
684 
685  (*pNavIter)-> LocateGlobalPointAndSetup( point,
686  &direction,
687  relativeSearch=false,
688  ignoreDirection=false);
689  }
690  return massVolume;
691 }
G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=0, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
int G4int
Definition: G4Types.hh:78
bool G4bool
Definition: G4Types.hh:79
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::vector< G4ITNavigator * >::iterator GetActiveNavigatorsIterator()

Here is the call graph for this function:

void G4ITMultiNavigator::ResetState ( )
protected

Definition at line 614 of file G4ITMultiNavigator.cc.

615 {
616  fWasLimitedByGeometry= false;
617 
618  G4Exception("G4ITMultiNavigator::ResetState()", "GeomNav0001",
620  "Cannot reset state for navigators of G4ITMultiNavigator.");
621 
622  std::vector<G4ITNavigator*>::iterator pNavigatorIter;
623  pNavigatorIter= pTransportManager-> GetActiveNavigatorsIterator();
624  for( G4int num=0; num< fNoActiveNavigators; ++pNavigatorIter,++num )
625  {
626  // (*pNavigatorIter)->ResetState(); // KEEP THIS comment !!!
627  }
628 }
#define fWasLimitedByGeometry
int G4int
Definition: G4Types.hh:78
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

void G4ITMultiNavigator::SetupHierarchy ( )
protected

Definition at line 632 of file G4ITMultiNavigator.cc.

633 {
634  G4Exception( "G4ITMultiNavigator::SetupHierarchy()",
635  "GeomNav0001", FatalException,
636  "Cannot setup hierarchy for navigators of G4ITMultiNavigator.");
637 }
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

void G4ITMultiNavigator::WhichLimited ( )
protected

Definition at line 486 of file G4ITMultiNavigator.cc.

487 {
488  // Flag which processes limited the step
489 
490  G4int last=-1;
491  const G4int IdTransport= 0; // Id of Mass Navigator !!
492  G4int noLimited=0;
493  ELimited shared= kSharedOther;
494 
495 #ifdef G4DEBUG_NAVIGATION
496  if( fVerbose > 2 )
497  {
498  G4cout << " Entered G4ITMultiNavigator::WhichLimited() " << G4endl;
499  }
500 #endif
501 
502  // Assume that [IdTransport] is Mass / Transport
503  //
504  G4bool transportLimited = (fCurrentStepSize[IdTransport] == fMinStep)
505  && ( fMinStep!= kInfinity);
506  if( transportLimited )
507  {
508  shared= kSharedTransport;
509  }
510 
511  for ( G4int num= 0; num < fNoActiveNavigators; num++ )
512  {
513  G4bool limitedStep;
514 
515  G4double step= fCurrentStepSize[num];
516 
517  limitedStep = ( step == fMinStep ) && ( step != kInfinity);
518 
519  fLimitTruth[ num ] = limitedStep;
520  if( limitedStep )
521  {
522  noLimited++;
523  fLimitedStep[num] = shared;
524  last= num;
525  }
526  else
527  {
528  fLimitedStep[num] = kDoNot;
529  }
530  }
531  if( (last > -1) && (noLimited == 1 ) )
532  {
533  fLimitedStep[ last ] = kUnique;
534  }
535 
536  fNoLimitingStep= noLimited;
537 
538  return;
539 }
static const G4double kInfinity
Definition: geomdefs.hh:42
ELimited
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
#define fLimitTruth
bool G4bool
Definition: G4Types.hh:79
#define fMinStep
#define fLimitedStep
#define fCurrentStepSize
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define fNoLimitingStep

Here is the caller graph for this function:

Friends And Related Function Documentation

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

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