Geant4  10.01.p03
G4Navigator.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4Navigator.cc 2012-11-01 japost Exp $
28 // GEANT4 tag $ Name: $
29 //
30 // class G4Navigator Implementation
31 //
32 // Original author: Paul Kent, July 95/96
33 // Responsible 1996-present: John Apostolakis
34 // Co-maintainer: Gabriele Cosmo
35 // Additional revisions by: Pedro Arce, Vladimir Grichine
36 // --------------------------------------------------------------------
37 
38 #include <iomanip>
39 
40 #include "G4Navigator.hh"
41 #include "G4ios.hh"
42 #include "G4SystemOfUnits.hh"
43 #include "G4GeometryTolerance.hh"
44 #include "G4VPhysicalVolume.hh"
45 
46 #include "G4VoxelSafety.hh"
47 
48 // ********************************************************************
49 // Constructor
50 // ********************************************************************
51 //
53  : fWasLimitedByGeometry(false), fVerbose(0),
54  fTopPhysical(0), fCheck(false), fPushed(false), fWarnPush(true)
55 {
56  fActive= false;
58 
60  // Initialises also all
61  // - exit / entry flags
62  // - flags & variables for exit normals
63  // - zero step counters
64  // - blocked volume
65 
68 
71 
74 
76 }
77 
78 // ********************************************************************
79 // Destructor
80 // ********************************************************************
81 //
83 {
84  delete fpVoxelSafety;
85 }
86 
87 // ********************************************************************
88 // ResetHierarchyAndLocate
89 // ********************************************************************
90 //
93  const G4ThreeVector &direction,
94  const G4TouchableHistory &h)
95 {
96  ResetState();
97  fHistory = *h.GetHistory();
99  fLastTriedStepComputation= false; // Redundant, but best
100  return LocateGlobalPointAndSetup(p, &direction, true, false);
101 }
102 
103 // ********************************************************************
104 // LocateGlobalPointAndSetup
105 //
106 // Locate the point in the hierarchy return 0 if outside
107 // The direction is required
108 // - if on an edge shared by more than two surfaces
109 // (to resolve likely looping in tracking)
110 // - at initial location of a particle
111 // (to resolve potential ambiguity at boundary)
112 //
113 // Flags on exit: (comments to be completed)
114 // fEntering - True if entering `daughter' volume (or replica)
115 // whether daughter of last mother directly
116 // or daughter of that volume's ancestor.
117 // fExiting - True if exited 'mother' volume
118 // (always ? - how about if going back down ? - tbc)
119 // ********************************************************************
120 //
123  const G4ThreeVector* pGlobalDirection,
124  const G4bool relativeSearch,
125  const G4bool ignoreDirection )
126 {
127  G4bool notKnownContained=true, noResult;
128  G4VPhysicalVolume *targetPhysical;
129  G4LogicalVolume *targetLogical;
130  G4VSolid *targetSolid=0;
131  G4ThreeVector localPoint, globalDirection;
132  EInside insideCode;
133 
134  G4bool considerDirection = (!ignoreDirection) || fLocatedOnEdge;
135 
137  fChangedGrandMotherRefFrame= false; // For local exit normal
138 
139  if( considerDirection && pGlobalDirection != 0 )
140  {
141  globalDirection=*pGlobalDirection;
142  }
143 
144 #ifdef G4VERBOSE
145  if( fVerbose > 2 )
146  {
147  G4int oldcoutPrec = G4cout.precision(8);
148  G4cout << "*** G4Navigator::LocateGlobalPointAndSetup: ***" << G4endl;
149  G4cout << " Called with arguments: " << G4endl
150  << " Globalpoint = " << globalPoint << G4endl
151  << " RelativeSearch = " << relativeSearch << G4endl;
152  if( fVerbose >= 4 )
153  {
154  G4cout << " ----- Upon entering:" << G4endl;
155  PrintState();
156  }
157  G4cout.precision(oldcoutPrec);
158  }
159 #endif
160 
161  G4int noLevelsExited=0 ;
162  G4int noLevelsEntered= 0;
163 
164  if ( !relativeSearch )
165  {
167  }
168  else
169  {
170  if ( fWasLimitedByGeometry )
171  {
172  fWasLimitedByGeometry = false;
173  fEnteredDaughter = fEntering; // Remember
174  fExitedMother = fExiting; // Remember
175  if ( fExiting )
176  {
177  noLevelsExited++; // count this first level entered too
178 
179  if ( fHistory.GetDepth() )
180  {
184  }
185  else
186  {
187  fLastLocatedPointLocal = localPoint;
188  fLocatedOutsideWorld = true;
189  fBlockedPhysicalVolume = 0; // to be sure
190  fBlockedReplicaNo = -1;
191  fEntering = false; // No longer
192  fEnteredDaughter = false;
193  fExitedMother = true; // ??
194 
195  return 0; // Have exited world volume
196  }
197  // A fix for the case where a volume is "entered" at an edge
198  // and a coincident surface exists outside it.
199  // - This stops it from exiting further volumes and cycling
200  // - However ReplicaNavigator treats this case itself
201  //
202  // assert( fBlockedPhysicalVolume!=0 );
203 
204  // Expect to be on edge => on surface
205  //
207  {
208  fExiting= false;
209  // Consider effect on Exit Normal !?
210  }
211  }
212  else
213  if ( fEntering )
214  {
215  // assert( fBlockedPhysicalVolume!=0 );
216 
217  noLevelsEntered++; // count the first level entered too
218 
220  {
221  case kNormal:
224  break;
225  case kReplica:
231  break;
232  case kParameterised:
234  {
235  G4VSolid *pSolid;
236  G4VPVParameterisation *pParam;
237  G4TouchableHistory parentTouchable( fHistory );
239  pSolid = pParam->ComputeSolid(fBlockedReplicaNo,
241  pSolid->ComputeDimensions(pParam, fBlockedReplicaNo,
248  //
249  // Set the correct solid and material in Logical Volume
250  //
251  G4LogicalVolume *pLogical;
253  pLogical->SetSolid( pSolid );
254  pLogical->UpdateMaterial(pParam ->
255  ComputeMaterial(fBlockedReplicaNo,
257  &parentTouchable));
258  }
259  break;
260  }
261  fEntering = false;
263  localPoint = fHistory.GetTopTransform().TransformPoint(globalPoint);
264  notKnownContained = false;
265  }
266  }
267  else
268  {
270  fEntering = false;
271  fEnteredDaughter = false; // Full Step was not taken, did not enter
272  fExiting = false;
273  fExitedMother = false; // Full Step was not taken, did not exit
274  }
275  }
276  //
277  // Search from top of history up through geometry until
278  // containing volume found:
279  // If on
280  // o OUTSIDE - Back up level, not/no longer exiting volumes
281  // o SURFACE and EXITING - Back up level, setting new blocking no.s
282  // else
283  // o containing volume found
284  //
285 
286  while (notKnownContained)
287  {
289  {
290  targetSolid = fHistory.GetTopVolume()->GetLogicalVolume()->GetSolid();
291  localPoint = fHistory.GetTopTransform().TransformPoint(globalPoint);
292  insideCode = targetSolid->Inside(localPoint);
293 #ifdef G4VERBOSE
294  if(( fVerbose == 1 ) && ( fCheck ))
295  {
296  G4String solidResponse = "-kInside-";
297  if (insideCode == kOutside)
298  solidResponse = "-kOutside-";
299  else if (insideCode == kSurface)
300  solidResponse = "-kSurface-";
301  G4cout << "*** G4Navigator::LocateGlobalPointAndSetup(): ***" << G4endl
302  << " Invoked Inside() for solid: " << targetSolid->GetName()
303  << ". Solid replied: " << solidResponse << G4endl
304  << " For local point p: " << localPoint << G4endl;
305  }
306 #endif
307  }
308  else
309  {
310  insideCode = freplicaNav.BackLocate(fHistory, globalPoint, localPoint,
311  fExiting, notKnownContained);
312  // !CARE! if notKnownContained returns false then the point is within
313  // the containing placement volume of the replica(s). If insidecode
314  // will result in the history being backed up one level, then the
315  // local point returned is the point in the system of this new level
316  }
317 
318 
319  if ( insideCode==kOutside )
320  {
321  noLevelsExited++;
322  if ( fHistory.GetDepth() )
323  {
327  fExiting = false;
328 
329  if( noLevelsExited > 1 )
330  {
331  // The first transformation was done by the sub-navigator
332  //
334  if( mRot )
335  {
336  fGrandMotherExitNormal *= (*mRot).inverse();
338  }
339  }
340  }
341  else
342  {
343  fLastLocatedPointLocal = localPoint;
344  fLocatedOutsideWorld = true;
345  // No extra transformation for ExitNormal - is in frame of Top Volume
346  return 0; // Have exited world volume
347  }
348  }
349  else
350  if ( insideCode==kSurface )
351  {
352  G4bool isExiting = fExiting;
353  if( (!fExiting)&&considerDirection )
354  {
355  // Figure out whether we are exiting this level's volume
356  // by using the direction
357  //
358  G4bool directionExiting = false;
359  G4ThreeVector localDirection =
360  fHistory.GetTopTransform().TransformAxis(globalDirection);
361 
362  // Make sure localPoint in correct reference frame
363  // ( Was it already correct ? How ? )
364  //
365  localPoint= fHistory.GetTopTransform().TransformPoint(globalPoint);
367  {
368  G4ThreeVector normal = targetSolid->SurfaceNormal(localPoint);
369  directionExiting = normal.dot(localDirection) > 0.0;
370  isExiting = isExiting || directionExiting;
371  }
372  }
373  if( isExiting )
374  {
375  noLevelsExited++;
376  if ( fHistory.GetDepth() )
377  {
381  //
382  // Still on surface but exited volume not necessarily convex
383  //
384  fValidExitNormal = false;
385 
386  if( noLevelsExited > 1 )
387  {
388  // The first transformation was done by the sub-navigator
389  //
390  const G4RotationMatrix* mRot =
392  if( mRot )
393  {
394  fGrandMotherExitNormal *= (*mRot).inverse();
396  }
397  }
398  }
399  else
400  {
401  fLastLocatedPointLocal = localPoint;
402  fLocatedOutsideWorld = true;
403  // No extra transformation for ExitNormal, is in frame of Top Vol
404  return 0; // Have exited world volume
405  }
406  }
407  else
408  {
409  notKnownContained=false;
410  }
411  }
412  else
413  {
414  notKnownContained=false;
415  }
416  } // END while (notKnownContained)
417  //
418  // Search downwards until deepest containing volume found,
419  // blocking fBlockedPhysicalVolume/BlockedReplicaNum
420  //
421  // 3 Cases:
422  //
423  // o Parameterised daughters
424  // =>Must be one G4PVParameterised daughter & voxels
425  // o Positioned daughters & voxels
426  // o Positioned daughters & no voxels
427 
428  noResult = true; // noResult should be renamed to
429  // something like enteredLevel, as that is its meaning.
430  do
431  {
432  // Determine `type' of current mother volume
433  //
434  targetPhysical = fHistory.GetTopVolume();
435  if (!targetPhysical) { break; }
436  targetLogical = targetPhysical->GetLogicalVolume();
437  switch( CharacteriseDaughters(targetLogical) )
438  {
439  case kNormal:
440  if ( targetLogical->GetVoxelHeader() ) // use optimised navigation
441  {
442  noResult = fvoxelNav.LevelLocate(fHistory,
445  globalPoint,
446  pGlobalDirection,
447  considerDirection,
448  localPoint);
449  }
450  else // do not use optimised navigation
451  {
452  noResult = fnormalNav.LevelLocate(fHistory,
455  globalPoint,
456  pGlobalDirection,
457  considerDirection,
458  localPoint);
459  }
460  break;
461  case kReplica:
462  noResult = freplicaNav.LevelLocate(fHistory,
465  globalPoint,
466  pGlobalDirection,
467  considerDirection,
468  localPoint);
469  break;
470  case kParameterised:
471  if( GetDaughtersRegularStructureId(targetLogical) != 1 )
472  {
473  noResult = fparamNav.LevelLocate(fHistory,
476  globalPoint,
477  pGlobalDirection,
478  considerDirection,
479  localPoint);
480  }
481  else // Regular structure
482  {
483  noResult = fregularNav.LevelLocate(fHistory,
486  globalPoint,
487  pGlobalDirection,
488  considerDirection,
489  localPoint);
490  }
491  break;
492  }
493 
494  // LevelLocate returns true if it finds a daughter volume
495  // in which globalPoint is inside (or on the surface).
496 
497  if ( noResult )
498  {
499  noLevelsEntered++;
500 
501  // Entering a daughter after ascending
502  //
503  // The blocked volume is no longer valid - it was for another level
504  //
506  fBlockedReplicaNo = -1;
507 
508  // fEntering should be false -- else blockedVolume is assumed good.
509  // fEnteredDaughter is used for ExitNormal
510  //
511  fEntering = false;
512  fEnteredDaughter = true;
513 
514  if( fExitedMother )
515  {
516  G4VPhysicalVolume* enteredPhysical = fHistory.GetTopVolume();
517  const G4RotationMatrix* mRot = enteredPhysical->GetRotation();
518  if( mRot )
519  {
520  // Go deeper, i.e. move 'down' in the hierarchy
521  // Apply direct rotation, not inverse
522  //
523  fGrandMotherExitNormal *= (*mRot);
525  }
526  }
527 
528 #ifdef G4DEBUG_NAVIGATION
529  if( fVerbose > 2 )
530  {
531  G4VPhysicalVolume* enteredPhysical = fHistory.GetTopVolume();
532  G4cout << "*** G4Navigator::LocateGlobalPointAndSetup() ***" << G4endl;
533  G4cout << " Entering volume: " << enteredPhysical->GetName()
534  << G4endl;
535  }
536 #endif
537  }
538  } while (noResult);
539 
540  fLastLocatedPointLocal = localPoint;
541 
542 #ifdef G4VERBOSE
543  if( fVerbose >= 4 )
544  {
545  G4int oldcoutPrec = G4cout.precision(8);
546  G4String curPhysVol_Name("None");
547  if (targetPhysical) { curPhysVol_Name = targetPhysical->GetName(); }
548  G4cout << " Return value = new volume = " << curPhysVol_Name << G4endl;
549  G4cout << " ----- Upon exiting:" << G4endl;
550  PrintState();
551  if( fVerbose >= 5 )
552  {
553  G4cout << "Upon exiting LocateGlobalPointAndSetup():" << G4endl;
554  G4cout << " History = " << G4endl << fHistory << G4endl << G4endl;
555  }
556  G4cout.precision(oldcoutPrec);
557  }
558 #endif
559 
560  fLocatedOutsideWorld= false;
561 
562  return targetPhysical;
563 }
564 
565 // ********************************************************************
566 // LocateGlobalPointWithinVolume
567 //
568 // -> the state information of this Navigator and its subNavigators
569 // is updated in order to start the next step at pGlobalpoint
570 // -> no check is performed whether pGlobalpoint is inside the
571 // original volume (this must be the case).
572 //
573 // Note: a direction could be added to the arguments, to aid in future
574 // optional checking (via the old code below, flagged by OLD_LOCATE).
575 // [ This would be done only in verbose mode ]
576 // ********************************************************************
577 //
578 void
580 {
581 #ifdef G4DEBUG_NAVIGATION
582  // Check: Either step was not limited by a boundary
583  // or else the full step is no longer being taken
584  assert( !fWasLimitedByGeometry );
585 #endif
586 
589  fChangedGrandMotherRefFrame= false; // Frame for Exit Normal
590 
591 #ifdef G4DEBUG_NAVIGATION
592  if( fVerbose > 2 )
593  {
594  G4cout << "Entering LocateGlobalWithinVolume(): History = " << G4endl;
595  G4cout << fHistory << G4endl;
596  }
597 #endif
598 
599  // For the case of Voxel (or Parameterised) volume the respective
600  // Navigator must be messaged to update its voxel information etc
601 
602  // Update the state of the Sub Navigators
603  // - in particular any voxel information they store/cache
604  //
605  G4VPhysicalVolume* motherPhysical = fHistory.GetTopVolume();
606  G4LogicalVolume* motherLogical = motherPhysical->GetLogicalVolume();
607  G4SmartVoxelHeader* pVoxelHeader = motherLogical->GetVoxelHeader();
608 
610  {
611  switch( CharacteriseDaughters(motherLogical) )
612  {
613  case kNormal:
614  if ( pVoxelHeader )
615  {
617  }
618  break;
619  case kParameterised:
620  if( GetDaughtersRegularStructureId(motherLogical) != 1 )
621  {
622  // Resets state & returns voxel node
623  //
625  }
626  break;
627  case kReplica:
628  G4Exception("G4Navigator::LocateGlobalPointWithinVolume()",
629  "GeomNav0001", FatalException,
630  "Not applicable for replicated volumes.");
631  break;
632  }
633  }
634 
635  // Reset the state variables
636  // - which would have been affected
637  // by the 'equivalent' call to LocateGlobalPointAndSetup
638  // - who's values have been invalidated by the 'move'.
639  //
641  fBlockedReplicaNo = -1;
642  fEntering = false;
643  fEnteredDaughter = false; // Boundary not encountered, did not enter
644  fExiting = false;
645  fExitedMother = false; // Boundary not encountered, did not exit
646 }
647 
648 // ********************************************************************
649 // SetSavedState
650 //
651 // Save the state, in case this is a parasitic call
652 // Save fValidExitNormal, fExitNormal, fExiting, fEntering,
653 // fBlockedPhysicalVolume, fBlockedReplicaNo, fLastStepWasZero;
654 // ********************************************************************
655 //
657 {
658  // Note: the state of dependent objects is not currently saved.
659  // ( This means that the full state is changed by calls between
660  // SetSavedState() and RestoreSavedState();
661 
666 
669 
671 
677 
678  // Even the safety sphere - if you want to change it do it explicitly!
679  //
682 }
683 
684 // ********************************************************************
685 // RestoreSavedState
686 //
687 // Restore the state (in Compute Step), in case this is a parasitic call
688 // ********************************************************************
689 //
691 {
696 
699 
701 
707 
708  // The 'expected' behaviour is to restore these too (fix 2014.05.26)
711 }
712 
713 // ********************************************************************
714 // ComputeStep
715 //
716 // Computes the next geometric Step: intersections with current
717 // mother and `daughter' volumes.
718 //
719 // NOTE:
720 //
721 // Flags on entry:
722 // --------------
723 // fValidExitNormal - Normal of exited volume is valid (convex, not a
724 // coincident boundary)
725 // fExitNormal - Surface normal of exited volume
726 // fExiting - True if have exited solid
727 //
728 // fBlockedPhysicalVolume - Ptr to exited volume (or 0)
729 // fBlockedReplicaNo - Replication no of exited volume
730 // fLastStepWasZero - True if last Step size was almost zero.
731 //
732 // Flags on exit:
733 // -------------
734 // fValidExitNormal - True if surface normal of exited volume is valid
735 // fExitNormal - Surface normal of exited volume rotated to mothers
736 // reference system
737 // fExiting - True if exiting mother
738 // fEntering - True if entering `daughter' volume (or replica)
739 // fBlockedPhysicalVolume - Ptr to candidate (entered) volume
740 // fBlockedReplicaNo - Replication no of candidate (entered) volume
741 // fLastStepWasZero - True if this Step size was almost zero.
742 // ********************************************************************
743 //
745  const G4ThreeVector &pDirection,
746  const G4double pCurrentProposedStepLength,
747  G4double &pNewSafety)
748 {
749  G4ThreeVector localDirection = ComputeLocalAxis(pDirection);
751  G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
752  G4LogicalVolume *motherLogical = motherPhysical->GetLogicalVolume();
753 
754  // All state relating to exiting normals must be reset
755  //
757  // Reset value - to erase its memory
759  // Reset - used for local exit normal
760  fGrandMotherExitNormal= G4ThreeVector( 0., 0., 0.);
761  fCalculatedExitNormal = false;
762  // Reset for new step
763 
764  static const G4double minStep = 0.05*kCarTolerance;
765  static G4ThreadLocal G4int sNavCScalls=0;
766  sNavCScalls++;
767 
769 
770 #ifdef G4VERBOSE
771  if( fVerbose > 0 )
772  {
773  G4cout << "*** G4Navigator::ComputeStep: ***" << G4endl;
774  G4cout << " Volume = " << motherPhysical->GetName()
775  << " - Proposed step length = " << pCurrentProposedStepLength
776  << G4endl;
777 #ifdef G4DEBUG_NAVIGATION
778  if( fVerbose >= 2 )
779  {
780  G4cout << " Called with the arguments: " << G4endl
781  << " Globalpoint = " << std::setw(25) << pGlobalpoint << G4endl
782  << " Direction = " << std::setw(25) << pDirection << G4endl;
783  if( fVerbose >= 4 )
784  {
785  G4cout << " ---- Upon entering : State" << G4endl;
786  PrintState();
787  }
788  }
789 #endif
790  }
791 #endif
792 
793  G4ThreeVector newLocalPoint = ComputeLocalPoint(pGlobalpoint);
794  if( newLocalPoint != fLastLocatedPointLocal )
795  {
796  // Check whether the relocation is within safety
797  //
798  G4ThreeVector oldLocalPoint = fLastLocatedPointLocal;
799  G4double moveLenSq = (newLocalPoint-oldLocalPoint).mag2();
800 
801  if ( moveLenSq >= kCarTolerance*kCarTolerance )
802  {
803 #ifdef G4VERBOSE
804  ComputeStepLog(pGlobalpoint, moveLenSq);
805 #endif
806  // Relocate the point within the same volume
807  //
808  LocateGlobalPointWithinVolume( pGlobalpoint );
809  fLastTriedStepComputation= true; // Ensure that this is set again !!
810  }
811  }
813  {
814  switch( CharacteriseDaughters(motherLogical) )
815  {
816  case kNormal:
817  if ( motherLogical->GetVoxelHeader() )
818  {
820  localDirection,
821  pCurrentProposedStepLength,
822  pNewSafety,
823  fHistory,
825  fExitNormal,
826  fExiting,
827  fEntering,
830 
831  }
832  else
833  {
834  if( motherPhysical->GetRegularStructureId() == 0 )
835  {
837  localDirection,
838  pCurrentProposedStepLength,
839  pNewSafety,
840  fHistory,
842  fExitNormal,
843  fExiting,
844  fEntering,
847  }
848  else // Regular (non-voxelised) structure
849  {
850  LocateGlobalPointAndSetup( pGlobalpoint, &pDirection, true, true );
851  fLastTriedStepComputation= true; // Ensure that this is set again !!
852  //
853  // if physical process limits the step, the voxel will not be the
854  // one given by ComputeStepSkippingEqualMaterials() and the local
855  // point will be wrongly calculated.
856 
857  // There is a problem: when msc limits the step and the point is
858  // assigned wrongly to phantom in previous step (while it is out
859  // of the container volume). Then LocateGlobalPointAndSetup() has
860  // reset the history topvolume to world.
861  //
863  {
864  G4Exception("G4Navigator::ComputeStep()",
865  "GeomNav1001", JustWarning,
866  "Point is relocated in voxels, while it should be outside!");
868  localDirection,
869  pCurrentProposedStepLength,
870  pNewSafety,
871  fHistory,
873  fExitNormal,
874  fExiting,
875  fEntering,
878  }
879  else
880  {
881  Step = fregularNav.
882  ComputeStepSkippingEqualMaterials(fLastLocatedPointLocal,
883  localDirection,
884  pCurrentProposedStepLength,
885  pNewSafety,
886  fHistory,
888  fExitNormal,
889  fExiting,
890  fEntering,
893  motherPhysical);
894  }
895  }
896  }
897  break;
898  case kParameterised:
899  if( GetDaughtersRegularStructureId(motherLogical) != 1 )
900  {
902  localDirection,
903  pCurrentProposedStepLength,
904  pNewSafety,
905  fHistory,
907  fExitNormal,
908  fExiting,
909  fEntering,
912  }
913  else // Regular structure
914  {
916  localDirection,
917  pCurrentProposedStepLength,
918  pNewSafety,
919  fHistory,
921  fExitNormal,
922  fExiting,
923  fEntering,
926  }
927  break;
928  case kReplica:
929  G4Exception("G4Navigator::ComputeStep()", "GeomNav0001",
930  FatalException, "Not applicable for replicated volumes.");
931  break;
932  }
933  }
934  else
935  {
936  // In the case of a replica, it must handle the exiting
937  // edge/corner problem by itself
938  //
939  G4bool exitingReplica = fExitedMother;
940  G4bool calculatedExitNormal;
941  Step = freplicaNav.ComputeStep(pGlobalpoint,
942  pDirection,
944  localDirection,
945  pCurrentProposedStepLength,
946  pNewSafety,
947  fHistory,
949  calculatedExitNormal,
950  fExitNormal,
951  exitingReplica,
952  fEntering,
955  fExiting= exitingReplica;
956  fCalculatedExitNormal= calculatedExitNormal;
957  }
958 
959  // Remember last safety origin & value.
960  //
961  fPreviousSftOrigin = pGlobalpoint;
962  fPreviousSafety = pNewSafety;
963 
964  // Count zero steps - one can occur due to changing momentum at a boundary
965  // - one, two (or a few) can occur at common edges between
966  // volumes
967  // - more than two is likely a problem in the geometry
968  // description or the Navigation
969 
970  // Rule of thumb: likely at an Edge if two consecutive steps are zero,
971  // because at least two candidate volumes must have been
972  // checked
973  //
974  fLocatedOnEdge = fLastStepWasZero && (Step==0.0);
975  fLastStepWasZero = (Step<minStep);
976  if (fPushed) { fPushed = fLastStepWasZero; }
977 
978  // Handle large number of consecutive zero steps
979  //
980  if ( fLastStepWasZero )
981  {
983 #ifdef G4DEBUG_NAVIGATION
984  if( fNumberZeroSteps > 1 )
985  {
986  G4cout << "G4Navigator::ComputeStep(): another 'zero' step, # "
988  << ", at " << pGlobalpoint
989  << ", in volume " << motherPhysical->GetName()
990  << ", nav-comp-step calls # " << sNavCScalls
991  << ", Step= " << Step
992  << G4endl;
993  }
994 #endif
996  {
997  // Act to recover this stuck track. Pushing it along direction
998  //
999  Step += 100*kCarTolerance;
1000 #ifdef G4VERBOSE
1001  if ((!fPushed) && (fWarnPush))
1002  {
1003  std::ostringstream message;
1004  message << "Track stuck or not moving." << G4endl
1005  << " Track stuck, not moving for "
1006  << fNumberZeroSteps << " steps" << G4endl
1007  << " in volume -" << motherPhysical->GetName()
1008  << "- at point " << pGlobalpoint << G4endl
1009  << " direction: " << pDirection << "." << G4endl
1010  << " Potential geometry or navigation problem !"
1011  << G4endl
1012  << " Trying pushing it of " << Step << " mm ...";
1013  G4Exception("G4Navigator::ComputeStep()", "GeomNav1002",
1014  JustWarning, message, "Potential overlap in geometry!");
1015  }
1016 #endif
1017  fPushed = true;
1018  }
1020  {
1021  // Must kill this stuck track
1022  //
1023  std::ostringstream message;
1024  message << "Stuck Track: potential geometry or navigation problem."
1025  << G4endl
1026  << " Track stuck, not moving for "
1027  << fNumberZeroSteps << " steps" << G4endl
1028  << " in volume -" << motherPhysical->GetName()
1029  << "- at point " << pGlobalpoint << G4endl
1030  << " direction: " << pDirection << ".";
1031 #ifdef G4VERBOSE
1032  motherPhysical->CheckOverlaps(5000, false);
1033 #endif
1034  G4Exception("G4Navigator::ComputeStep()", "GeomNav0003",
1035  EventMustBeAborted, message);
1036  }
1037  }
1038  else
1039  {
1040  if (!fPushed) fNumberZeroSteps = 0;
1041  }
1042 
1043  fEnteredDaughter = fEntering; // I expect to enter a volume in this Step
1045 
1046  fStepEndPoint = pGlobalpoint
1047  + std::min(Step,pCurrentProposedStepLength) * pDirection;
1048  fLastStepEndPointLocal = fLastLocatedPointLocal + Step * localDirection;
1049 
1050  if( fExiting )
1051  {
1052 #ifdef G4DEBUG_NAVIGATION
1053  if( fVerbose > 2 )
1054  {
1055  G4cout << " At G4Nav CompStep End - if(exiting) - fExiting= " << fExiting
1056  << " fValidExitNormal = " << fValidExitNormal << G4endl;
1057  G4cout << " fExitNormal= " << fExitNormal << G4endl;
1058  }
1059 #endif
1060 
1062  {
1064  {
1065  // Convention: fExitNormal is in the 'grand-mother' coordinate system
1066  //
1068  fCalculatedExitNormal= true;
1069  }
1070  else
1071  {
1073  }
1074  }
1075  else
1076  {
1077  // We must calculate the normal anyway (in order to have it if requested)
1078  //
1079  G4ThreeVector finalLocalPoint =
1080  fLastLocatedPointLocal + localDirection*Step;
1081 
1083  {
1084  // Find normal in the 'mother' coordinate system
1085  //
1086  G4ThreeVector exitNormalMotherFrame=
1087  motherLogical->GetSolid()->SurfaceNormal(finalLocalPoint);
1088 
1089  // Transform it to the 'grand-mother' coordinate system
1090  //
1091  const G4RotationMatrix* mRot = motherPhysical->GetRotation();
1092  if( mRot )
1093  {
1095  fGrandMotherExitNormal = (*mRot).inverse() * exitNormalMotherFrame;
1096  }
1097  else
1098  {
1099  fGrandMotherExitNormal = exitNormalMotherFrame;
1100  }
1101 
1102  // Do not set fValidExitNormal -- this signifies
1103  // that the solid is convex!
1104  //
1105  fCalculatedExitNormal= true;
1106  }
1107  else
1108  {
1109  fCalculatedExitNormal = false;
1110  //
1111  // Nothing can be done at this stage currently - to solve this
1112  // Replica Navigation must have calculated the normal for this case
1113  // already.
1114  // Cases: mother is not convex, and exit is at previous replica level
1115 
1116 #ifdef G4DEBUG_NAVIGATION
1118 
1119  desc << "Problem in ComputeStep: Replica Navigation did not provide"
1120  << " valid exit Normal. " << G4endl;
1121  desc << " Do not know how calculate it in this case." << G4endl;
1122  desc << " Location = " << finalLocalPoint << G4endl;
1123  desc << " Volume name = " << motherPhysical->GetName()
1124  << " copy/replica No = " << motherPhysical->GetCopyNo() << G4endl;
1125  G4Exception("G4Navigator::ComputeStep()", "GeomNav0003",
1126  JustWarning, desc, "Normal not available for exiting.");
1127 #endif
1128  }
1129  }
1130 
1131  // Now transform it to the global reference frame !!
1132  //
1134  {
1135  G4int depth= fHistory.GetDepth();
1136  if( depth > 0 )
1137  {
1138  G4AffineTransform GrandMotherToGlobalTransf =
1139  fHistory.GetTransform(depth-1).Inverse();
1141  GrandMotherToGlobalTransf.TransformAxis( fGrandMotherExitNormal );
1142  }
1143  else
1144  {
1146  }
1147  }
1148  else
1149  {
1150  fExitNormalGlobalFrame= G4ThreeVector( 0., 0., 0.);
1151  }
1152  }
1153 
1154  if( (Step == pCurrentProposedStepLength) && (!fExiting) && (!fEntering) )
1155  {
1156  // This if Step is not really limited by the geometry.
1157  // The Navigator is obliged to return "infinity"
1158  //
1159  Step = kInfinity;
1160  }
1161 
1162 #ifdef G4VERBOSE
1163  if( fVerbose > 1 )
1164  {
1165  if( fVerbose >= 4 )
1166  {
1167  G4cout << " ----- Upon exiting :" << G4endl;
1168  PrintState();
1169  }
1170  G4cout << " Returned step= " << Step;
1171  if( fVerbose > 5 ) G4cout << G4endl;
1172  if( Step == kInfinity )
1173  {
1174  G4cout << " Requested step= " << pCurrentProposedStepLength ;
1175  if( fVerbose > 5) G4cout << G4endl;
1176  }
1177  G4cout << " Safety = " << pNewSafety << G4endl;
1178  }
1179 #endif
1180 
1181  return Step;
1182 }
1183 
1184 // ********************************************************************
1185 // CheckNextStep
1186 //
1187 // Compute the step without altering the navigator state
1188 // ********************************************************************
1189 //
1191  const G4ThreeVector& pDirection,
1192  const G4double pCurrentProposedStepLength,
1193  G4double& pNewSafety)
1194 {
1195  G4double step;
1196 
1197  // Save the state, for this parasitic call
1198  //
1199  SetSavedState();
1200 
1201  step = ComputeStep ( pGlobalpoint,
1202  pDirection,
1203  pCurrentProposedStepLength,
1204  pNewSafety );
1205 
1206  // It is a parasitic call, so attempt to restore the key parts of the state
1207  //
1208  RestoreSavedState();
1209  // NOTE: the state of the current subnavigator is NOT restored.
1210  // ***> TODO: restore subnavigator state
1211  // if( last_located) Need Position of last location
1212  // if( last_computed step) Need Endposition of last step
1213 
1214  return step;
1215 }
1216 
1217 // ********************************************************************
1218 // ResetState
1219 //
1220 // Resets stack and minimum of navigator state `machine'
1221 // ********************************************************************
1222 //
1224 {
1225  fWasLimitedByGeometry = false;
1226  fEntering = false;
1227  fExiting = false;
1228  fLocatedOnEdge = false;
1229  fLastStepWasZero = false;
1230  fEnteredDaughter = false;
1231  fExitedMother = false;
1232  fPushed = false;
1233 
1234  fValidExitNormal = false;
1236  fCalculatedExitNormal = false;
1237 
1238  fExitNormal = G4ThreeVector(0,0,0);
1241 
1243  fPreviousSafety = 0.0;
1244 
1245  fNumberZeroSteps = 0;
1246 
1248  fBlockedReplicaNo = -1;
1249 
1251  fLocatedOutsideWorld = false;
1252 }
1253 
1254 // ********************************************************************
1255 // SetupHierarchy
1256 //
1257 // Renavigates & resets hierarchy described by current history
1258 // o Reset volumes
1259 // o Recompute transforms and/or solids of replicated/parameterised volumes
1260 // ********************************************************************
1261 //
1263 {
1264  G4int i;
1265  const G4int cdepth = fHistory.GetDepth();
1266  G4VPhysicalVolume *current;
1267  G4VSolid *pSolid;
1268  G4VPVParameterisation *pParam;
1269 
1270  for ( i=1; i<=cdepth; i++ )
1271  {
1272  current = fHistory.GetVolume(i);
1273  switch ( fHistory.GetVolumeType(i) )
1274  {
1275  case kNormal:
1276  break;
1277  case kReplica:
1279  break;
1280  case kParameterised:
1281  G4int replicaNo;
1282  pParam = current->GetParameterisation();
1283  replicaNo = fHistory.GetReplicaNo(i);
1284  pSolid = pParam->ComputeSolid(replicaNo, current);
1285 
1286  // Set up dimensions & transform in solid/physical volume
1287  //
1288  pSolid->ComputeDimensions(pParam, replicaNo, current);
1289  pParam->ComputeTransformation(replicaNo, current);
1290 
1291  G4TouchableHistory *pTouchable= 0;
1292  if( pParam->IsNested() )
1293  {
1294  pTouchable= new G4TouchableHistory( fHistory );
1295  pTouchable->MoveUpHistory(); // Move up to the parent level
1296  // Adequate only if Nested at the Branch level (last)
1297  // To extend to other cases:
1298  // pTouchable->MoveUpHistory(cdepth-i-1);
1299  // Move to the parent level of *Current* level
1300  // Could replace this line and constructor with a revised
1301  // c-tor for History(levels to drop)
1302  }
1303  // Set up the correct solid and material in Logical Volume
1304  //
1305  G4LogicalVolume *pLogical = current->GetLogicalVolume();
1306  pLogical->SetSolid( pSolid );
1307  pLogical->UpdateMaterial( pParam ->
1308  ComputeMaterial(replicaNo, current, pTouchable) );
1309  delete pTouchable;
1310  break;
1311  }
1312  }
1313 }
1314 
1315 // ********************************************************************
1316 // GetLocalExitNormal
1317 //
1318 // Obtains the Normal vector to a surface (in local coordinates)
1319 // pointing out of previous volume and into current volume
1320 // ********************************************************************
1321 //
1323 {
1324  G4ThreeVector ExitNormal(0.,0.,0.);
1325  G4VSolid *currentSolid=0;
1326  G4LogicalVolume *candidateLogical;
1328  {
1329  // use fLastLocatedPointLocal and next candidate volume
1330  //
1331  G4ThreeVector nextSolidExitNormal(0.,0.,0.);
1332 
1333  if( fEntering && (fBlockedPhysicalVolume!=0) )
1334  {
1335  candidateLogical= fBlockedPhysicalVolume->GetLogicalVolume();
1336  if( candidateLogical )
1337  {
1338  // fLastStepEndPointLocal is in the coordinates of the mother
1339  // we need it in the daughter's coordinate system.
1340 
1341  // The following code should also work in case of Replica
1342  {
1343  // First transform fLastLocatedPointLocal to the new daughter
1344  // coordinates
1345  //
1346  G4AffineTransform MotherToDaughterTransform=
1350  G4ThreeVector daughterPointOwnLocal=
1351  MotherToDaughterTransform.TransformPoint( fLastStepEndPointLocal );
1352 
1353  // OK if it is a parameterised volume
1354  //
1355  EInside inSideIt;
1356  G4bool onSurface;
1357  G4double safety= -1.0;
1358  currentSolid= candidateLogical->GetSolid();
1359  inSideIt = currentSolid->Inside(daughterPointOwnLocal);
1360  onSurface = (inSideIt == kSurface);
1361  if( ! onSurface )
1362  {
1363  if( inSideIt == kOutside )
1364  {
1365  safety = (currentSolid->DistanceToIn(daughterPointOwnLocal));
1366  onSurface = safety < 100.0 * kCarTolerance;
1367  }
1368  else if (inSideIt == kInside )
1369  {
1370  safety = (currentSolid->DistanceToOut(daughterPointOwnLocal));
1371  onSurface = safety < 100.0 * kCarTolerance;
1372  }
1373  }
1374 
1375  if( onSurface )
1376  {
1377  nextSolidExitNormal =
1378  currentSolid->SurfaceNormal(daughterPointOwnLocal);
1379 
1380  // Entering the solid ==> opposite
1381  //
1382  ExitNormal = -nextSolidExitNormal;
1383  fCalculatedExitNormal= true;
1384  }
1385  else
1386  {
1387 #ifdef G4VERBOSE
1388  if(( fVerbose == 1 ) && ( fCheck ))
1389  {
1390  std::ostringstream message;
1391  message << "Point not on surface ! " << G4endl
1392  << " Point = "
1393  << daughterPointOwnLocal << G4endl
1394  << " Physical volume = "
1396  << " Logical volume = "
1397  << candidateLogical->GetName() << G4endl
1398  << " Solid = " << currentSolid->GetName()
1399  << " Type = "
1400  << currentSolid->GetEntityType() << G4endl
1401  << *currentSolid << G4endl;
1402  if( inSideIt == kOutside )
1403  {
1404  message << "Point is Outside. " << G4endl
1405  << " Safety (from outside) = " << safety << G4endl;
1406  }
1407  else // if( inSideIt == kInside )
1408  {
1409  message << "Point is Inside. " << G4endl
1410  << " Safety (from inside) = " << safety << G4endl;
1411  }
1412  G4Exception("G4Navigator::GetLocalExitNormal()", "GeomNav1001",
1413  JustWarning, message);
1414  }
1415 #endif
1416  }
1417  *valid = onSurface; // was =true;
1418  }
1419  }
1420  }
1421  else if ( fExiting )
1422  {
1423  ExitNormal = fGrandMotherExitNormal;
1424  *valid = true;
1425  fCalculatedExitNormal= true; // Should be true already
1426  }
1427  else // i.e. ( fBlockedPhysicalVolume == 0 )
1428  {
1429  *valid = false;
1430  G4Exception("G4Navigator::GetLocalExitNormal()",
1431  "GeomNav0003", JustWarning,
1432  "Incorrect call to GetLocalSurfaceNormal." );
1433  }
1434  }
1435  else // ( ! fLastTriedStepComputation ) ie. last call was to Locate
1436  {
1437  if ( EnteredDaughterVolume() )
1438  {
1439  G4VSolid* daughterSolid =fHistory.GetTopVolume()->GetLogicalVolume()
1440  ->GetSolid();
1441  ExitNormal= -(daughterSolid->SurfaceNormal(fLastLocatedPointLocal));
1442  if( std::fabs(ExitNormal.mag2()-1.0 ) > CLHEP::perMillion )
1443  {
1445  desc << " Parameters of solid: " << *daughterSolid
1446  << " Point for surface = " << fLastLocatedPointLocal << std::endl;
1447  G4Exception("G4Navigator::GetLocalExitNormal()",
1448  "GeomNav0003", FatalException, desc,
1449  "Surface Normal returned by Solid is not a Unit Vector." );
1450  }
1451  fCalculatedExitNormal= true;
1452  *valid = true;
1453  }
1454  else
1455  {
1456  if( fExitedMother )
1457  {
1458  ExitNormal = fGrandMotherExitNormal;
1459  *valid = true;
1460  fCalculatedExitNormal= true;
1461  }
1462  else // We are not at a boundary. ExitNormal remains (0,0,0)
1463  {
1464  *valid = false;
1465  fCalculatedExitNormal= false;
1466  G4ExceptionDescription message;
1467  message << "Function called when *NOT* at a Boundary." << G4endl;
1468  message << "Exit Normal not calculated." << G4endl;
1469  G4Exception("G4Navigator::GetLocalExitNormal()",
1470  "GeomNav0003", JustWarning, message);
1471  }
1472  }
1473  }
1474  return ExitNormal;
1475 }
1476 
1477 // ********************************************************************
1478 // GetMotherToDaughterTransform
1479 //
1480 // Obtains the mother to daughter affine transformation
1481 // ********************************************************************
1482 //
1485  G4int enteringReplicaNo,
1486  EVolume enteringVolumeType )
1487 {
1488  switch (enteringVolumeType)
1489  {
1490  case kNormal: // Nothing is needed to prepare the transformation
1491  break; // It is stored already in the physical volume (placement)
1492  case kReplica: // Sets the transform in the Replica - tbc
1493  G4Exception("G4Navigator::GetMotherToDaughterTransform()",
1494  "GeomNav0001", FatalException,
1495  "Method NOT Implemented yet for replica volumes.");
1496  break;
1497  case kParameterised:
1498  if( pEnteringPhysVol->GetRegularStructureId() == 0 )
1499  {
1500  G4VPVParameterisation *pParam =
1501  pEnteringPhysVol->GetParameterisation();
1502  G4VSolid* pSolid =
1503  pParam->ComputeSolid(enteringReplicaNo, pEnteringPhysVol);
1504  pSolid->ComputeDimensions(pParam, enteringReplicaNo, pEnteringPhysVol);
1505 
1506  // Sets the transform in the Parameterisation
1507  //
1508  pParam->ComputeTransformation(enteringReplicaNo, pEnteringPhysVol);
1509 
1510  // Set the correct solid and material in Logical Volume
1511  //
1512  G4LogicalVolume* pLogical = pEnteringPhysVol->GetLogicalVolume();
1513  pLogical->SetSolid( pSolid );
1514  }
1515  break;
1516  }
1517  return G4AffineTransform(pEnteringPhysVol->GetRotation(),
1518  pEnteringPhysVol->GetTranslation()).Invert();
1519 }
1520 
1521 
1522 // ********************************************************************
1523 // GetLocalExitNormalAndCheck
1524 //
1525 // Obtains the Normal vector to a surface (in local coordinates)
1526 // pointing out of previous volume and into current volume, and
1527 // checks the current point against expected 'local' value.
1528 // ********************************************************************
1529 //
1532 #ifdef G4DEBUG_NAVIGATION
1533  const G4ThreeVector& ExpectedBoundaryPointGlobal,
1534 #else
1535  const G4ThreeVector&,
1536 #endif
1537  G4bool* pValid)
1538 {
1539 #ifdef G4DEBUG_NAVIGATION
1540  // Check Current point against expected 'local' value
1541  //
1543  {
1544  G4ThreeVector ExpectedBoundaryPointLocal;
1545 
1546  const G4AffineTransform& GlobalToLocal= GetGlobalToLocalTransform();
1547  ExpectedBoundaryPointLocal =
1548  GlobalToLocal.TransformPoint( ExpectedBoundaryPointGlobal );
1549 
1550  // Add here: Comparison against expected position,
1551  // i.e. the endpoint of ComputeStep
1552  }
1553 #endif
1554 
1555  return GetLocalExitNormal( pValid);
1556 }
1557 
1558 
1559 // ********************************************************************
1560 // GetGlobalExitNormal
1561 //
1562 // Obtains the Normal vector to a surface (in global coordinates)
1563 // pointing out of previous volume and into current volume
1564 // ********************************************************************
1565 //
1568  G4bool* pNormalCalculated)
1569 {
1570  G4bool validNormal;
1571  G4ThreeVector localNormal, globalNormal;
1572  G4bool usingStored;
1573 
1574  usingStored=
1576  ( ( fLastTriedStepComputation && fExiting) // Just calculated it
1577  || // No locate in between
1579  && (IntersectPointGlobal-fStepEndPoint).mag2()
1580  < (10.0*kCarTolerance*kCarTolerance)
1581  ) // Calculated it 'just' before & then called locate
1582  // but it did not move position
1583  );
1584 
1585  if( usingStored )
1586  {
1587  // This was computed in last call to ComputeStep
1588  // and only if it arrived at boundary
1589  //
1590  globalNormal = fExitNormalGlobalFrame;
1591  G4double normMag2 = globalNormal.mag2();
1592  if( std::fabs ( normMag2 - 1.0 ) < perMillion ) // Value is good
1593  {
1594  *pNormalCalculated = true; // ComputeStep always computes it if Exiting
1595  // (fExiting==true)
1596  }
1597  else
1598  {
1599  G4ExceptionDescription message;
1600  message.precision(10);
1601  message << " WARNING> Expected normal-global-frame to be valid, "
1602  << " i.e. a unit vector!" << G4endl
1603  << " - but |normal| = " << std::sqrt(normMag2)
1604  << " - and |normal|^2 = " << normMag2 << G4endl
1605  << " which differs from 1.0 by " << normMag2 - 1.0 << G4endl
1606  << " n = " << fExitNormalGlobalFrame << G4endl;
1607  message << "============================================================"
1608  << G4endl;
1609  G4int oldVerbose = fVerbose;
1610  fVerbose=4;
1611  message << " State of Navigator: " << G4endl;
1612  message << *this << G4endl;
1613  fVerbose = oldVerbose;
1614  message << "============================================================"
1615  << G4endl;
1616 
1617  G4Exception("G4Navigator::GetGlobalExitNormal()",
1618  "GeomNav0003",JustWarning, message,
1619  "Value obtained from stored global-normal is not a unit vector.");
1620 
1621  // (Re)Compute it now -- as either it was not computed, or it is wrong.
1622  //
1623  localNormal = GetLocalExitNormalAndCheck(IntersectPointGlobal,
1624  &validNormal);
1625  *pNormalCalculated = fCalculatedExitNormal;
1626 
1627  G4AffineTransform localToGlobal = GetLocalToGlobalTransform();
1628  globalNormal = localToGlobal.TransformAxis( localNormal );
1629  }
1630  }
1631  else
1632  {
1633  localNormal = GetLocalExitNormalAndCheck(IntersectPointGlobal,&validNormal);
1634  *pNormalCalculated = fCalculatedExitNormal;
1635 
1636 #ifdef G4DEBUG_NAVIGATION
1637  usingStored= false;
1638 
1639  if( (!validNormal) && !fCalculatedExitNormal)
1640  {
1642  edN << " Calculated = " << fCalculatedExitNormal << G4endl;
1643  edN << " Entering= " << fEntering << G4endl;
1644  G4int oldVerbose= this->GetVerboseLevel();
1645  this->SetVerboseLevel(4);
1646  edN << " State of Navigator: " << G4endl;
1647  edN << *this << G4endl;
1648  this->SetVerboseLevel( oldVerbose );
1649 
1650  G4Exception("G4Navigator::GetGlobalExitNormal()",
1651  "GeomNav0003", JustWarning, edN,
1652  "LocalExitNormalAndCheck() did not calculate Normal.");
1653  }
1654 #endif
1655 
1656  G4double localMag2= localNormal.mag2();
1657  if( validNormal && (std::fabs(localMag2-1.0)) > CLHEP::perMillion )
1658  {
1660 
1661  edN << "G4Navigator::GetGlobalExitNormal: "
1662  << " Using Local Normal - from call to GetLocalExitNormalAndCheck. "
1663  << G4endl
1664  << " Local Exit Normal : " << " || = " << std::sqrt(localMag2)
1665  << " vec = " << localNormal << G4endl
1666  << " Global Exit Normal : " << " || = " << globalNormal.mag()
1667  << " vec = " << globalNormal << G4endl;
1668  edN << " Calculated It = " << fCalculatedExitNormal << G4endl;
1669 
1670  G4Exception("G4Navigator::GetGlobalExitNormal()",
1671  "GeomNav0003",JustWarning, edN,
1672  "Value obtained from new local *solid* is incorrect.");
1673  localNormal = localNormal.unit(); // Should we correct it ??
1674  }
1675  G4AffineTransform localToGlobal = GetLocalToGlobalTransform();
1676  globalNormal = localToGlobal.TransformAxis( localNormal );
1677  }
1678 
1679 #ifdef G4DEBUG_NAVIGATION
1680  if( usingStored )
1681  {
1682  G4ThreeVector globalNormAgn;
1683 
1684  localNormal= GetLocalExitNormalAndCheck(IntersectPointGlobal, &validNormal);
1685 
1686  G4AffineTransform localToGlobal = GetLocalToGlobalTransform();
1687  globalNormAgn = localToGlobal.TransformAxis( localNormal );
1688 
1689  // Check the value computed against fExitNormalGlobalFrame
1690  G4ThreeVector diffNorm = globalNormAgn - fExitNormalGlobalFrame;
1691  if( diffNorm.mag2() > perMillion*CLHEP::perMillion)
1692  {
1693  G4ExceptionDescription edDfn;
1694  edDfn << "Found difference in normals in case of exiting mother "
1695  << "- when Get is called after ComputingStep " << G4endl;
1696  edDfn << " Magnitude of diff = " << diffNorm.mag() << G4endl;
1697  edDfn << " Normal stored (Global) = " << fExitNormalGlobalFrame
1698  << G4endl;
1699  edDfn << " Global Computed from Local = " << globalNormAgn << G4endl;
1700  G4Exception("G4Navigator::GetGlobalExitNormal()", "GeomNav0003",
1701  JustWarning, edDfn);
1702  }
1703  }
1704 #endif
1705 
1706  return globalNormal;
1707 }
1708 
1709 // To make the new Voxel Safety the default, uncomment the next line
1710 #define G4NEW_SAFETY 1
1711 
1712 // ********************************************************************
1713 // ComputeSafety
1714 //
1715 // It assumes that it will be
1716 // i) called at the Point in the same volume as the EndPoint of the
1717 // ComputeStep.
1718 // ii) after (or at the end of) ComputeStep OR after the relocation.
1719 // ********************************************************************
1720 //
1722  const G4double pMaxLength,
1723  const G4bool keepState)
1724 {
1725  G4double newSafety = 0.0;
1726 
1727 #ifdef G4DEBUG_NAVIGATION
1728  G4int oldcoutPrec = G4cout.precision(8);
1729  if( fVerbose > 0 )
1730  {
1731  G4cout << "*** G4Navigator::ComputeSafety: ***" << G4endl
1732  << " Called at point: " << pGlobalpoint << G4endl;
1733 
1734  G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
1735  G4cout << " Volume = " << motherPhysical->GetName()
1736  << " - Maximum length = " << pMaxLength << G4endl;
1737  if( fVerbose >= 4 )
1738  {
1739  G4cout << " ----- Upon entering Compute Safety:" << G4endl;
1740  PrintState();
1741  }
1742  }
1743 #endif
1744 
1745  G4double distEndpointSq = (pGlobalpoint-fStepEndPoint).mag2();
1746  G4bool stayedOnEndpoint = distEndpointSq < kCarTolerance*kCarTolerance;
1747  G4bool endpointOnSurface = fEnteredDaughter || fExitedMother;
1748 
1749  if( endpointOnSurface && stayedOnEndpoint )
1750  {
1751 #ifdef G4DEBUG_NAVIGATION
1752  if( fVerbose >= 2 )
1753  {
1754  G4cout << " G4Navigator::ComputeSafety() finds that point - "
1755  << pGlobalpoint << " - is on surface " << G4endl;
1756  if( fEnteredDaughter ) { G4cout << " entered new daughter volume"; }
1757  if( fExitedMother ) { G4cout << " and exited previous volume."; }
1758  G4cout << G4endl;
1759  G4cout << " EndPoint was = " << fStepEndPoint << G4endl;
1760  }
1761 #endif
1762  newSafety = 0.0;
1763  // return newSafety;
1764  }
1765  else // if( !(endpointOnSurface && stayedOnEndpoint) )
1766  {
1767  if (keepState) { SetSavedState(); }
1768 
1769  // Pseudo-relocate to this point (updates voxel information only)
1770  //
1771  LocateGlobalPointWithinVolume( pGlobalpoint );
1772  // --->> DANGER: Side effects on sub-navigator voxel information <<---
1773  // Could be replaced again by 'granular' calls to sub-navigator
1774  // locates (similar side-effects, but faster.
1775  // Solutions:
1776  // 1) Re-locate (to where?)
1777  // 2) Insure that the methods using (G4ComputeStep?)
1778  // does a relocation (if information is disturbed only ?)
1779 
1780 #ifdef G4DEBUG_NAVIGATION
1781  if( fVerbose >= 2 )
1782  {
1783  G4cout << " G4Navigator::ComputeSafety() relocates-in-volume to point: "
1784  << pGlobalpoint << G4endl;
1785  }
1786 #endif
1787  G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
1788  G4LogicalVolume *motherLogical = motherPhysical->GetLogicalVolume();
1789  G4SmartVoxelHeader* pVoxelHeader = motherLogical->GetVoxelHeader();
1790  G4ThreeVector localPoint = ComputeLocalPoint(pGlobalpoint);
1791 
1793  {
1794  switch(CharacteriseDaughters(motherLogical))
1795  {
1796  case kNormal:
1797  if ( pVoxelHeader )
1798  {
1799 #ifdef G4NEW_SAFETY
1800  G4double safetyTwo = fpVoxelSafety->ComputeSafety(localPoint,
1801  *motherPhysical, pMaxLength);
1802  newSafety= safetyTwo; // Faster and best available
1803 #else
1804  G4double safetyOldVoxel;
1805  safetyOldVoxel =
1806  fvoxelNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1807  newSafety= safetyOldVoxel;
1808 #endif
1809  }
1810  else
1811  {
1812  newSafety=fnormalNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1813  }
1814  break;
1815  case kParameterised:
1816  if( GetDaughtersRegularStructureId(motherLogical) != 1 )
1817  {
1818  newSafety=fparamNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1819  }
1820  else // Regular structure
1821  {
1822  newSafety=fregularNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1823  }
1824  break;
1825  case kReplica:
1826  G4Exception("G4Navigator::ComputeSafety()", "GeomNav0001",
1827  FatalException, "Not applicable for replicated volumes.");
1828  break;
1829  }
1830  }
1831  else
1832  {
1833  newSafety = freplicaNav.ComputeSafety(pGlobalpoint, localPoint,
1834  fHistory, pMaxLength);
1835  }
1836 
1837  if (keepState)
1838  {
1840  // This now overwrites the values of the Safety 'sphere' (correction)
1841  }
1842 
1843  // Remember last safety origin & value
1844  //
1845  // We overwrite the Safety 'sphere' - keeping old behaviour
1846  fPreviousSftOrigin = pGlobalpoint;
1847  fPreviousSafety = newSafety;
1848  }
1849 
1850 #ifdef G4DEBUG_NAVIGATION
1851  if( fVerbose > 1 )
1852  {
1853  G4cout << " ---- Exiting ComputeSafety " << G4endl;
1854  if( fVerbose > 2 ) { PrintState(); }
1855  G4cout << " Returned value of Safety = " << newSafety << G4endl;
1856  }
1857  G4cout.precision(oldcoutPrec);
1858 #endif
1859 
1860  return newSafety;
1861 }
1862 
1863 
1864 // ********************************************************************
1865 // RecheckDistanceToCurrentBoundary
1866 //
1867 // Trial method for checking potential displacement for MS
1868 // Check position aDisplacedGlobalpoint, to see whether it is in the
1869 // current volume (mother).
1870 // If in mother, check distance to boundary along aNewDirection.
1871 // If in entering daughter, check distance back to boundary.
1872 // NOTE:
1873 // Can be called only after ComputeStep() is called - before ReLocation
1874 // Deals only with current volume (and potentially entered)
1875 // Caveats
1876 // First VERSION: Does not consider daughter volumes if it remained in mother
1877 // neither for checking whether it has exited current (mother) volume,
1878 // nor for determining distance to exit this (mother) volume.
1879 // ********************************************************************
1880 //
1882  const G4ThreeVector &aDisplacedGlobalPoint,
1883  const G4ThreeVector &aNewDirection,
1884  const G4double ProposedMove,
1885  G4double *prDistance,
1886  G4double *prNewSafety) const
1887 {
1888  G4ThreeVector localPosition = ComputeLocalPoint(aDisplacedGlobalPoint);
1889  G4ThreeVector localDirection = ComputeLocalAxis(aNewDirection);
1890  // G4double Step = kInfinity;
1891 
1892  G4bool validExitNormal;
1893  G4ThreeVector exitNormal;
1894  // Check against mother solid
1895  G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
1896  G4LogicalVolume *motherLogical = motherPhysical->GetLogicalVolume();
1897 
1898 #ifdef CHECK_ORDER_OF_METHODS
1900  {
1901  G4Exception("G4Navigator::RecheckDistanceToCurrentBoundary()",
1902  "GeomNav0001", FatalException,
1903  "Method must be called after ComputeStep(), before call to LocateMethod.");
1904  }
1905 #endif
1906 
1907  EInside locatedDaug; // = kUndefined;
1908  G4double daughterStep= DBL_MAX;
1909  G4double daughterSafety= DBL_MAX;
1910 
1911  if( fEnteredDaughter )
1912  {
1913  if( motherLogical->CharacteriseDaughters() ==kReplica ) { return false; }
1914 
1915  // Track arrived at boundary of a daughter volume at
1916  // the last call of ComputeStep().
1917  // In case the proposed displaced point is inside this daughter,
1918  // it must backtrack at least to the entry point.
1919  // NOTE: No check is made against other daughter volumes. It is
1920  // assumed that the proposed displacement is small enough that
1921  // this is not needed.
1922 
1923  // Must check boundary of current daughter
1924  //
1925  G4VPhysicalVolume *candPhysical= fBlockedPhysicalVolume;
1926  G4LogicalVolume *candLogical= candPhysical->GetLogicalVolume();
1927  G4VSolid *candSolid= candLogical->GetSolid();
1928 
1929  G4AffineTransform nextLevelTrf(candPhysical->GetRotation(),
1930  candPhysical->GetTranslation());
1931 
1932  G4ThreeVector dgPosition= nextLevelTrf.TransformPoint(localPosition);
1933  G4ThreeVector dgDirection= nextLevelTrf.TransformAxis(localDirection);
1934  locatedDaug = candSolid->Inside(dgPosition);
1935 
1936  if( locatedDaug == kInside )
1937  {
1938  // Reverse direction - and find first exit. ( Is it valid?)
1939  // Must backtrack
1940  G4double distanceBackOut =
1941  candSolid->DistanceToOut(dgPosition,
1942  - dgDirection, // Reverse direction
1943  true, &validExitNormal, &exitNormal);
1944  daughterStep= - distanceBackOut;
1945  // No check is made whether the particle could have arrived at
1946  // at this location without encountering another volume or
1947  // a different psurface of the current volume
1948  if( prNewSafety )
1949  {
1950  daughterSafety= candSolid->DistanceToOut(dgPosition);
1951  }
1952  }
1953  else
1954  {
1955  if( locatedDaug == kOutside )
1956  {
1957  // See whether it still intersects it
1958  //
1959  daughterStep= candSolid->DistanceToIn(dgPosition,
1960  dgDirection);
1961  if( prNewSafety )
1962  {
1963  daughterSafety= candSolid->DistanceToIn(dgPosition);
1964  }
1965  }
1966  else
1967  {
1968  // The point remains on the surface of candidate solid
1969  //
1970  daughterStep= 0.0;
1971  daughterSafety= 0.0;
1972  }
1973  }
1974 
1975  // If trial point is in daughter (or on its surface) we have the
1976  // answer, the rest is not relevant
1977  //
1978  if( locatedDaug != kOutside )
1979  {
1980  *prDistance= daughterStep;
1981  if( prNewSafety ) { *prNewSafety= daughterSafety; }
1982  return true;
1983  }
1984  // If ever extended, so that some type of mother cut daughter,
1985  // this would change
1986  }
1987 
1988  G4VSolid *motherSolid= motherLogical->GetSolid();
1989 
1990  G4double motherStep= DBL_MAX, motherSafety= DBL_MAX;
1991 
1992  // Check distance to boundary of mother
1993  //
1995  {
1996  EInside locatedMoth = motherSolid->Inside(localPosition);
1997 
1998  if( locatedMoth == kInside )
1999  {
2000  motherSafety= motherSolid->DistanceToOut(localPosition);
2001  if( ProposedMove >= motherSafety )
2002  {
2003  motherStep= motherSolid->DistanceToOut(localPosition,
2004  localDirection,
2005  true, &validExitNormal, &exitNormal);
2006  }
2007  else
2008  {
2009  motherStep= ProposedMove;
2010  }
2011  }
2012  else if( locatedMoth == kOutside)
2013  {
2014  motherSafety= motherSolid->DistanceToIn(localPosition);
2015  if( ProposedMove >= motherSafety )
2016  {
2017  motherStep= - motherSolid->DistanceToIn(localPosition,
2018  -localDirection);
2019  }
2020  }
2021  else
2022  {
2023  motherSafety= 0.0;
2024  *prDistance= 0.0; // On surface - no move // motherStep;
2025  if( prNewSafety ) { *prNewSafety= motherSafety; }
2026  return false;
2027  }
2028  }
2029  else
2030  {
2031  return false;
2032  }
2033 
2034  *prDistance= std::min( motherStep, daughterStep );
2035  if( prNewSafety )
2036  {
2037  *prNewSafety= std::min( motherSafety, daughterSafety );
2038  }
2039  return true;
2040 }
2041 
2042 
2043 // ********************************************************************
2044 // CreateTouchableHistoryHandle
2045 // ********************************************************************
2046 //
2048 {
2050 }
2051 
2052 // ********************************************************************
2053 // PrintState
2054 // ********************************************************************
2055 //
2057 {
2058  G4int oldcoutPrec = G4cout.precision(4);
2059  if( fVerbose >= 4 )
2060  {
2061  G4cout << "The current state of G4Navigator is: " << G4endl;
2062  G4cout << " ValidExitNormal= " << fValidExitNormal // << G4endl
2063  << " ExitNormal = " << fExitNormal // << G4endl
2064  << " Exiting = " << fExiting // << G4endl
2065  << " Entering = " << fEntering // << G4endl
2066  << " BlockedPhysicalVolume= " ;
2067  if (fBlockedPhysicalVolume==0)
2068  {
2069  G4cout << "None";
2070  }
2071  else
2072  {
2074  }
2075  G4cout << G4endl
2076  << " BlockedReplicaNo = " << fBlockedReplicaNo // << G4endl
2077  << " LastStepWasZero = " << fLastStepWasZero // << G4endl
2078  << G4endl;
2079  }
2080  if( ( 1 < fVerbose) && (fVerbose < 4) )
2081  {
2082  G4cout << G4endl; // Make sure to line up
2083  G4cout << std::setw(30) << " ExitNormal " << " "
2084  << std::setw( 5) << " Valid " << " "
2085  << std::setw( 9) << " Exiting " << " "
2086  << std::setw( 9) << " Entering" << " "
2087  << std::setw(15) << " Blocked:Volume " << " "
2088  << std::setw( 9) << " ReplicaNo" << " "
2089  << std::setw( 8) << " LastStepZero " << " "
2090  << G4endl;
2091  G4cout << "( " << std::setw(7) << fExitNormal.x()
2092  << ", " << std::setw(7) << fExitNormal.y()
2093  << ", " << std::setw(7) << fExitNormal.z() << " ) "
2094  << std::setw( 5) << fValidExitNormal << " "
2095  << std::setw( 9) << fExiting << " "
2096  << std::setw( 9) << fEntering << " ";
2097  if ( fBlockedPhysicalVolume==0 )
2098  G4cout << std::setw(15) << "None";
2099  else
2100  G4cout << std::setw(15)<< fBlockedPhysicalVolume->GetName();
2101  G4cout << std::setw( 9) << fBlockedReplicaNo << " "
2102  << std::setw( 8) << fLastStepWasZero << " "
2103  << G4endl;
2104  }
2105  if( fVerbose > 2 )
2106  {
2107  G4cout.precision(8);
2108  G4cout << " Current Localpoint = " << fLastLocatedPointLocal << G4endl;
2109  G4cout << " PreviousSftOrigin = " << fPreviousSftOrigin << G4endl;
2110  G4cout << " PreviousSafety = " << fPreviousSafety << G4endl;
2111  }
2112  G4cout.precision(oldcoutPrec);
2113 }
2114 
2115 // ********************************************************************
2116 // ComputeStepLog
2117 // ********************************************************************
2118 //
2120  G4double moveLenSq) const
2121 {
2122  // The following checks only make sense if the move is larger
2123  // than the tolerance.
2124 
2125  static const G4double fAccuracyForWarning = kCarTolerance,
2126  fAccuracyForException = 1000*kCarTolerance;
2127 
2128  G4ThreeVector OriginalGlobalpoint = fHistory.GetTopTransform().Inverse().
2129  TransformPoint(fLastLocatedPointLocal);
2130 
2131  G4double shiftOriginSafSq = (fPreviousSftOrigin-pGlobalpoint).mag2();
2132 
2133  // Check that the starting point of this step is
2134  // within the isotropic safety sphere of the last point
2135  // to a accuracy/precision given by fAccuracyForWarning.
2136  // If so give warning.
2137  // If it fails by more than fAccuracyForException exit with error.
2138  //
2139  if( shiftOriginSafSq >= sqr(fPreviousSafety) )
2140  {
2141  G4double shiftOrigin = std::sqrt(shiftOriginSafSq);
2142  G4double diffShiftSaf = shiftOrigin - fPreviousSafety;
2143 
2144  if( diffShiftSaf > fAccuracyForWarning )
2145  {
2146  G4int oldcoutPrec= G4cout.precision(8);
2147  G4int oldcerrPrec= G4cerr.precision(10);
2148  std::ostringstream message, suggestion;
2149  message << "Accuracy error or slightly inaccurate position shift."
2150  << G4endl
2151  << " The Step's starting point has moved "
2152  << std::sqrt(moveLenSq)/mm << " mm " << G4endl
2153  << " since the last call to a Locate method." << G4endl
2154  << " This has resulted in moving "
2155  << shiftOrigin/mm << " mm "
2156  << " from the last point at which the safety "
2157  << " was calculated " << G4endl
2158  << " which is more than the computed safety= "
2159  << fPreviousSafety/mm << " mm at that point." << G4endl
2160  << " This difference is "
2161  << diffShiftSaf/mm << " mm." << G4endl
2162  << " The tolerated accuracy is "
2163  << fAccuracyForException/mm << " mm.";
2164 
2165  suggestion << " ";
2166  static G4ThreadLocal G4int warnNow = 0;
2167  if( ((++warnNow % 100) == 1) )
2168  {
2169  message << G4endl
2170  << " This problem can be due to either " << G4endl
2171  << " - a process that has proposed a displacement"
2172  << " larger than the current safety , or" << G4endl
2173  << " - inaccuracy in the computation of the safety";
2174  suggestion << "We suggest that you " << G4endl
2175  << " - find i) what particle is being tracked, and "
2176  << " ii) through what part of your geometry " << G4endl
2177  << " for example by re-running this event with "
2178  << G4endl
2179  << " /tracking/verbose 1 " << G4endl
2180  << " - check which processes you declare for"
2181  << " this particle (and look at non-standard ones)"
2182  << G4endl
2183  << " - in case, create a detailed logfile"
2184  << " of this event using:" << G4endl
2185  << " /tracking/verbose 6 ";
2186  }
2187  G4Exception("G4Navigator::ComputeStep()",
2188  "GeomNav1002", JustWarning,
2189  message, G4String(suggestion.str()));
2190  G4cout.precision(oldcoutPrec);
2191  G4cerr.precision(oldcerrPrec);
2192  }
2193 #ifdef G4DEBUG_NAVIGATION
2194  else
2195  {
2196  G4cerr << "WARNING - G4Navigator::ComputeStep()" << G4endl
2197  << " The Step's starting point has moved "
2198  << std::sqrt(moveLenSq) << "," << G4endl
2199  << " which has taken it to the limit of"
2200  << " the current safety. " << G4endl;
2201  }
2202 #endif
2203  }
2204  G4double safetyPlus = fPreviousSafety + fAccuracyForException;
2205  if ( shiftOriginSafSq > sqr(safetyPlus) )
2206  {
2207  std::ostringstream message;
2208  message << "May lead to a crash or unreliable results." << G4endl
2209  << " Position has shifted considerably without"
2210  << " notifying the navigator !" << G4endl
2211  << " Tolerated safety: " << safetyPlus << G4endl
2212  << " Computed shift : " << shiftOriginSafSq;
2213  G4Exception("G4Navigator::ComputeStep()", "GeomNav1002",
2214  JustWarning, message);
2215  }
2216 }
2217 
2218 // ********************************************************************
2219 // Operator <<
2220 // ********************************************************************
2221 //
2222 std::ostream& operator << (std::ostream &os,const G4Navigator &n)
2223 {
2224  // Old version did only the following:
2225  // os << "Current History: " << G4endl << n.fHistory;
2226  // Old behaviour is recovered for fVerbose = 0
2227 
2228  // Adapted from G4Navigator::PrintState() const
2229 
2230  G4int oldcoutPrec = os.precision(4);
2231  if( n.fVerbose >= 4 )
2232  {
2233  os << "The current state of G4Navigator is: " << G4endl;
2234  os << " ValidExitNormal= " << n.fValidExitNormal << G4endl
2235  << " ExitNormal = " << n.fExitNormal << G4endl
2236  << " Exiting = " << n.fExiting << G4endl
2237  << " Entering = " << n.fEntering << G4endl
2238  << " BlockedPhysicalVolume= " ;
2239  if (n.fBlockedPhysicalVolume==0)
2240  os << "None";
2241  else
2242  os << n.fBlockedPhysicalVolume->GetName();
2243  os << G4endl
2244  << " BlockedReplicaNo = " << n.fBlockedReplicaNo << G4endl
2245  << " LastStepWasZero = " << n.fLastStepWasZero << G4endl
2246  << G4endl;
2247  }
2248  if( ( 1 < n.fVerbose) && (n.fVerbose < 4) )
2249  {
2250  os << G4endl; // Make sure to line up
2251  os << std::setw(30) << " ExitNormal " << " "
2252  << std::setw( 5) << " Valid " << " "
2253  << std::setw( 9) << " Exiting " << " "
2254  << std::setw( 9) << " Entering" << " "
2255  << std::setw(15) << " Blocked:Volume " << " "
2256  << std::setw( 9) << " ReplicaNo" << " "
2257  << std::setw( 8) << " LastStepZero " << " "
2258  << G4endl;
2259  os << "( " << std::setw(7) << n.fExitNormal.x()
2260  << ", " << std::setw(7) << n.fExitNormal.y()
2261  << ", " << std::setw(7) << n.fExitNormal.z() << " ) "
2262  << std::setw( 5) << n.fValidExitNormal << " "
2263  << std::setw( 9) << n.fExiting << " "
2264  << std::setw( 9) << n.fEntering << " ";
2265  if ( n.fBlockedPhysicalVolume==0 )
2266  { os << std::setw(15) << "None"; }
2267  else
2268  { os << std::setw(15)<< n.fBlockedPhysicalVolume->GetName(); }
2269  os << std::setw( 9) << n.fBlockedReplicaNo << " "
2270  << std::setw( 8) << n.fLastStepWasZero << " "
2271  << G4endl;
2272  }
2273  if( n.fVerbose > 2 )
2274  {
2275  os.precision(8);
2276  os << " Current Localpoint = " << n.fLastLocatedPointLocal << G4endl;
2277  os << " PreviousSftOrigin = " << n.fPreviousSftOrigin << G4endl;
2278  os << " PreviousSafety = " << n.fPreviousSafety << G4endl;
2279  }
2280  if( n.fVerbose > 3 || n.fVerbose == 0 )
2281  {
2282  os << "Current History: " << G4endl << n.fHistory;
2283  }
2284 
2285  os.precision(oldcoutPrec);
2286  return os;
2287 }
G4String GetName() const
G4ReplicaNavigation freplicaNav
Definition: G4Navigator.hh:503
G4bool fWarnPush
Definition: G4Navigator.hh:495
G4ParameterisedNavigation fparamNav
Definition: G4Navigator.hh:502
G4SmartVoxelHeader * GetVoxelHeader() const
G4bool fValidExitNormal
Definition: G4Navigator.hh:423
G4VPhysicalVolume * GetTopVolume() const
G4VoxelNavigation fvoxelNav
Definition: G4Navigator.hh:501
G4VPhysicalVolume * fBlockedPhysicalVolume
Definition: G4Navigator.hh:415
virtual void ResetState()
G4bool fLastTriedStepComputation
Definition: G4Navigator.hh:400
virtual G4bool LevelLocate(G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint)
G4VoxelSafety * fpVoxelSafety
Definition: G4Navigator.hh:505
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)
G4bool fCheck
Definition: G4Navigator.hh:493
const G4ThreeVector & GetTranslation() const
static const G4double kInfinity
Definition: geomdefs.hh:42
virtual G4bool IsNested() const
#define fWasLimitedByGeometry
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void UpdateMaterial(G4Material *pMaterial)
G4String GetName() const
CLHEP::Hep3Vector G4ThreeVector
CLHEP::HepRotation G4RotationMatrix
G4ThreeVector fLastStepEndPointLocal
Definition: G4Navigator.hh:388
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:137
G4AffineTransform Inverse() const
G4int GetVerboseLevel() const
G4ThreeVector fGrandMotherExitNormal
Definition: G4Navigator.hh:427
G4ThreeVector ComputeLocalAxis(const G4ThreeVector &pVec) const
G4bool fPushed
Definition: G4Navigator.hh:495
virtual G4TouchableHistoryHandle CreateTouchableHistoryHandle() const
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
void ResetStackAndState()
G4int GetDepth() const
G4int fNumberZeroSteps
Definition: G4Navigator.hh:447
G4double GetSurfaceTolerance() const
G4VPhysicalVolume * spBlockedPhysicalVolume
Definition: G4Navigator.hh:471
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
Definition: G4Navigator.cc:744
G4int fVerbose
Definition: G4Navigator.hh:392
G4double ComputeSafety(const G4ThreeVector &localPoint, const G4VPhysicalVolume &currentPhysical, G4double maxLength=DBL_MAX)
void PrintState() const
void SetSolid(G4VSolid *pSolid)
G4int fAbandonThreshold_NoZeroSteps
Definition: G4Navigator.hh:451
G4SmartVoxelNode * ParamVoxelLocate(G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint)
G4bool fEnteredDaughter
Definition: G4Navigator.hh:372
EVolume GetVolumeType(G4int n) const
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)
G4bool fLocatedOnEdge
Definition: G4Navigator.hh:445
const G4AffineTransform GetLocalToGlobalTransform() const
virtual G4GeometryType GetEntityType() const =0
#define G4ThreadLocal
Definition: tls.hh:89
G4ReferenceCountedHandle< G4TouchableHistory > G4TouchableHistoryHandle
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
int G4int
Definition: G4Types.hh:78
G4bool fLastStepWasZero
Definition: G4Navigator.hh:442
virtual G4VPhysicalVolume * ResetHierarchyAndLocate(const G4ThreeVector &point, const G4ThreeVector &direction, const G4TouchableHistory &h)
Definition: G4Navigator.cc:92
G4ThreeVector fExitNormalGlobalFrame
Definition: G4Navigator.hh:431
virtual G4ThreeVector GetLocalExitNormalAndCheck(const G4ThreeVector &point, G4bool *valid)
G4int GetTopReplicaNo() const
void SetVerboseLevel(G4int level)
#define fPushed
G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4NavigationHistory &history, const G4double pMaxLength=DBL_MAX)
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
G4bool fWasLimitedByGeometry
Definition: G4Navigator.hh:382
G4NormalNavigation fnormalNav
Definition: G4Navigator.hh:500
G4bool fEntering
Definition: G4Navigator.hh:406
void SetNormalNavigation(G4NormalNavigation *fnormnav)
G4double ComputeSafety(const G4ThreeVector &localPoint, const G4NavigationHistory &history, const G4double pProposedMaxLength=DBL_MAX)
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
virtual EInside Inside(const G4ThreeVector &p) const =0
virtual G4ThreeVector GetLocalExitNormal(G4bool *valid)
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)
G4ThreeVector fPreviousSftOrigin
Definition: G4Navigator.hh:454
G4double CheckNextStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
bool G4bool
Definition: G4Types.hh:79
EVolume GetTopVolumeType() const
virtual G4int GetRegularStructureId() const =0
G4NavigationHistory fHistory
Definition: G4Navigator.hh:368
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4NavigationHistory &history, const G4double pMaxLength=DBL_MAX)
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
virtual G4VPVParameterisation * GetParameterisation() const =0
virtual void SetCopyNo(G4int CopyNo)=0
G4int GetReplicaNo(G4int n) const
G4bool LevelLocate(G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint)
G4bool fExitedMother
Definition: G4Navigator.hh:378
EVolume CharacteriseDaughters() const
virtual G4bool RecheckDistanceToCurrentBoundary(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double CurrentProposedStepLength, G4double *prDistance, G4double *prNewSafety=0) const
G4bool EnteredDaughterVolume() const
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
G4int fBlockedReplicaNo
Definition: G4Navigator.hh:416
#define G4DEBUG_NAVIGATION
const G4int n
EVolume CharacteriseDaughters(const G4LogicalVolume *pLog) const
G4ThreeVector fLastLocatedPointLocal
Definition: G4Navigator.hh:418
G4double fPreviousSafety
Definition: G4Navigator.hh:455
G4bool fLocatedOutsideWorld
Definition: G4Navigator.hh:420
void SetSavedState()
Definition: G4Navigator.cc:656
G4bool fActive
Definition: G4Navigator.hh:397
G4ThreeVector ComputeLocalPoint(const G4ThreeVector &rGlobPoint) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double kCarTolerance
Definition: G4Navigator.hh:361
static const double perMillion
Definition: G4SIunits.hh:298
void NewLevel(G4VPhysicalVolume *pNewMother, EVolume vType=kNormal, G4int nReplica=-1)
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4TouchableHistory * CreateTouchableHistory() const
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)
G4AffineTransform GetMotherToDaughterTransform(G4VPhysicalVolume *dVolume, G4int dReplicaNo, EVolume dVolumeType)
G4LogicalVolume * GetLogicalVolume() const
const G4AffineTransform & GetTransform(G4int n) const
virtual G4ThreeVector GetGlobalExitNormal(const G4ThreeVector &point, G4bool *valid)
G4bool fExiting
Definition: G4Navigator.hh:406
EInside
Definition: geomdefs.hh:58
G4int GetDaughtersRegularStructureId(const G4LogicalVolume *pLog) const
virtual void ComputeTransformation(const G4int, G4VPhysicalVolume *) const =0
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
virtual G4int GetCopyNo() const =0
const G4AffineTransform & GetGlobalToLocalTransform() const
virtual G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=0, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
Definition: G4Navigator.cc:122
struct G4Navigator::G4SaveNavigatorState fSaveState
const G4RotationMatrix * GetRotation() const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
Definition: Step.hh:41
G4SmartVoxelNode * VoxelLocate(G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint)
G4double ComputeSafety(const G4ThreeVector &localPoint, const G4NavigationHistory &history, const G4double pProposedMaxLength=DBL_MAX)
#define G4endl
Definition: G4ios.hh:61
const G4AffineTransform & GetTopTransform() const
std::ostream & operator<<(std::ostream &os, const G4Navigator &n)
G4bool fCalculatedExitNormal
Definition: G4Navigator.hh:433
G4ThreeVector fStepEndPoint
Definition: G4Navigator.hh:385
const G4NavigationHistory * GetHistory() const
G4ThreeVector fExitNormal
Definition: G4Navigator.hh:424
G4int fActionThreshold_NoZeroSteps
Definition: G4Navigator.hh:449
T sqr(const T &x)
Definition: templates.hh:145
G4bool fChangedGrandMotherRefFrame
Definition: G4Navigator.hh:429
EVolume
Definition: geomdefs.hh:68
G4int MoveUpHistory(G4int num_levels=1)
G4double ComputeSafety(const G4ThreeVector &globalPoint, const G4ThreeVector &localPoint, G4NavigationHistory &history, const G4double pProposedMaxLength=DBL_MAX)
G4bool LevelLocate(G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint)
double G4double
Definition: G4Types.hh:76
G4VPhysicalVolume * GetVolume(G4int n) const
G4bool LevelLocate(G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint)
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
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)
G4RegularNavigation fregularNav
Definition: G4Navigator.hh:504
EInside BackLocate(G4NavigationHistory &history, const G4ThreeVector &globalPoint, G4ThreeVector &localPoint, const G4bool &exiting, G4bool &notKnownInside) const
#define DBL_MAX
Definition: templates.hh:83
virtual void SetupHierarchy()
static const double mm
Definition: G4SIunits.hh:102
virtual G4bool CheckOverlaps(G4int res=1000, G4double tol=0., G4bool verbose=true, G4int errMax=1)
EVolume VolumeType(const G4VPhysicalVolume *pVol) const
static G4GeometryTolerance * GetInstance()
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
Definition: G4Navigator.cc:579
G4VSolid * GetSolid() const
void ComputeStepLog(const G4ThreeVector &pGlobalpoint, G4double moveLenSq) const
void RestoreSavedState()
Definition: G4Navigator.cc:690
virtual ~G4Navigator()
Definition: G4Navigator.cc:82
G4GLOB_DLL std::ostream G4cerr