Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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;
57  fLastTriedStepComputation= false;
58 
60  // Initialises also all
61  // - exit / entry flags
62  // - flags & variables for exit normals
63  // - zero step counters
64  // - blocked volume
65 
66  fActionThreshold_NoZeroSteps = 10;
67  fAbandonThreshold_NoZeroSteps = 25;
68 
70  fregularNav.SetNormalNavigation( &fnormalNav );
71 
72  fStepEndPoint = G4ThreeVector( kInfinity, kInfinity, kInfinity );
73  fLastStepEndPointLocal = G4ThreeVector( kInfinity, kInfinity, kInfinity );
74 }
75 
76 // ********************************************************************
77 // Destructor
78 // ********************************************************************
79 //
81 {;}
82 
83 // ********************************************************************
84 // ResetHierarchyAndLocate
85 // ********************************************************************
86 //
89  const G4ThreeVector &direction,
90  const G4TouchableHistory &h)
91 {
92  ResetState();
93  fHistory = *h.GetHistory();
95  fLastTriedStepComputation= false; // Redundant, but best
96  return LocateGlobalPointAndSetup(p, &direction, true, false);
97 }
98 
99 // ********************************************************************
100 // LocateGlobalPointAndSetup
101 //
102 // Locate the point in the hierarchy return 0 if outside
103 // The direction is required
104 // - if on an edge shared by more than two surfaces
105 // (to resolve likely looping in tracking)
106 // - at initial location of a particle
107 // (to resolve potential ambiguity at boundary)
108 //
109 // Flags on exit: (comments to be completed)
110 // fEntering - True if entering `daughter' volume (or replica)
111 // whether daughter of last mother directly
112 // or daughter of that volume's ancestor.
113 // ********************************************************************
114 //
117  const G4ThreeVector* pGlobalDirection,
118  const G4bool relativeSearch,
119  const G4bool ignoreDirection )
120 {
121  G4bool notKnownContained=true, noResult;
122  G4VPhysicalVolume *targetPhysical;
123  G4LogicalVolume *targetLogical;
124  G4VSolid *targetSolid=0;
125  G4ThreeVector localPoint, globalDirection;
126  EInside insideCode;
127 
128  G4bool considerDirection = (!ignoreDirection) || fLocatedOnEdge;
129 
130  fLastTriedStepComputation= false;
131  fChangedGrandMotherRefFrame= false; // For local exit normal
132 
133  if( considerDirection && pGlobalDirection != 0 )
134  {
135  globalDirection=*pGlobalDirection;
136  }
137 
138 
139 #ifdef G4VERBOSE
140  if( fVerbose > 2 )
141  {
142  G4int oldcoutPrec = G4cout.precision(8);
143  G4cout << "*** G4Navigator::LocateGlobalPointAndSetup: ***" << G4endl;
144  G4cout << " Called with arguments: " << G4endl
145  << " Globalpoint = " << globalPoint << G4endl
146  << " RelativeSearch = " << relativeSearch << G4endl;
147  if( fVerbose == 4 )
148  {
149  G4cout << " ----- Upon entering:" << G4endl;
150  PrintState();
151  }
152  G4cout.precision(oldcoutPrec);
153  }
154 #endif
155 
156  if ( !relativeSearch )
157  {
159  }
160  else
161  {
162  if ( fWasLimitedByGeometry )
163  {
164  fWasLimitedByGeometry = false;
165  fEnteredDaughter = fEntering; // Remember
166  fExitedMother = fExiting; // Remember
167  if ( fExiting )
168  {
169  if ( fHistory.GetDepth() )
170  {
171  fBlockedPhysicalVolume = fHistory.GetTopVolume();
172  fBlockedReplicaNo = fHistory.GetTopReplicaNo();
174  }
175  else
176  {
177  fLastLocatedPointLocal = localPoint;
178  fLocatedOutsideWorld = true;
179  return 0; // Have exited world volume
180  }
181  // A fix for the case where a volume is "entered" at an edge
182  // and a coincident surface exists outside it.
183  // - This stops it from exiting further volumes and cycling
184  // - However ReplicaNavigator treats this case itself
185  //
186  if ( fLocatedOnEdge && (VolumeType(fBlockedPhysicalVolume)!=kReplica ))
187  {
188  fExiting= false;
189  }
190  }
191  else
192  if ( fEntering )
193  {
194  switch (VolumeType(fBlockedPhysicalVolume))
195  {
196  case kNormal:
197  fHistory.NewLevel(fBlockedPhysicalVolume, kNormal,
198  fBlockedPhysicalVolume->GetCopyNo());
199  break;
200  case kReplica:
201  freplicaNav.ComputeTransformation(fBlockedReplicaNo,
202  fBlockedPhysicalVolume);
203  fHistory.NewLevel(fBlockedPhysicalVolume, kReplica,
204  fBlockedReplicaNo);
205  fBlockedPhysicalVolume->SetCopyNo(fBlockedReplicaNo);
206  break;
207  case kParameterised:
208  if( fBlockedPhysicalVolume->GetRegularStructureId() == 0 )
209  {
210  G4VSolid *pSolid;
211  G4VPVParameterisation *pParam;
212  G4TouchableHistory parentTouchable( fHistory );
213  pParam = fBlockedPhysicalVolume->GetParameterisation();
214  pSolid = pParam->ComputeSolid(fBlockedReplicaNo,
215  fBlockedPhysicalVolume);
216  pSolid->ComputeDimensions(pParam, fBlockedReplicaNo,
217  fBlockedPhysicalVolume);
218  pParam->ComputeTransformation(fBlockedReplicaNo,
219  fBlockedPhysicalVolume);
220  fHistory.NewLevel(fBlockedPhysicalVolume, kParameterised,
221  fBlockedReplicaNo);
222  fBlockedPhysicalVolume->SetCopyNo(fBlockedReplicaNo);
223  //
224  // Set the correct solid and material in Logical Volume
225  //
226  G4LogicalVolume *pLogical;
227  pLogical = fBlockedPhysicalVolume->GetLogicalVolume();
228  pLogical->SetSolid( pSolid );
229  pLogical->UpdateMaterial(pParam ->
230  ComputeMaterial(fBlockedReplicaNo,
231  fBlockedPhysicalVolume,
232  &parentTouchable));
233  }
234  break;
235  }
236  fEntering = false;
237  fBlockedPhysicalVolume = 0;
238  localPoint = fHistory.GetTopTransform().TransformPoint(globalPoint);
239  notKnownContained = false;
240  }
241  }
242  else
243  {
244  fBlockedPhysicalVolume = 0;
245  fEntering = false;
246  fEnteredDaughter = false; // Full Step was not taken, did not enter
247  fExiting = false;
248  fExitedMother = false; // Full Step was not taken, did not exit
249  }
250  }
251  //
252  // Search from top of history up through geometry until
253  // containing volume found:
254  // If on
255  // o OUTSIDE - Back up level, not/no longer exiting volumes
256  // o SURFACE and EXITING - Back up level, setting new blocking no.s
257  // else
258  // o containing volume found
259  //
260  G4int noLevelsExited=0 ;
261 
262  while (notKnownContained)
263  {
265  {
266  targetSolid = fHistory.GetTopVolume()->GetLogicalVolume()->GetSolid();
267  localPoint = fHistory.GetTopTransform().TransformPoint(globalPoint);
268  insideCode = targetSolid->Inside(localPoint);
269 #ifdef G4VERBOSE
270  if(( fVerbose == 1 ) && ( fCheck ))
271  {
272  G4String solidResponse = "-kInside-";
273  if (insideCode == kOutside)
274  solidResponse = "-kOutside-";
275  else if (insideCode == kSurface)
276  solidResponse = "-kSurface-";
277  G4cout << "*** G4Navigator::LocateGlobalPointAndSetup(): ***" << G4endl
278  << " Invoked Inside() for solid: " << targetSolid->GetName()
279  << ". Solid replied: " << solidResponse << G4endl
280  << " For local point p: " << localPoint << G4endl;
281  }
282 #endif
283  }
284  else
285  {
286  insideCode = freplicaNav.BackLocate(fHistory, globalPoint, localPoint,
287  fExiting, notKnownContained);
288  // !CARE! if notKnownContained returns false then the point is within
289  // the containing placement volume of the replica(s). If insidecode
290  // will result in the history being backed up one level, then the
291  // local point returned is the point in the system of this new level
292  }
293 
294 
295  if ( insideCode==kOutside )
296  {
297  noLevelsExited++;
298  if ( fHistory.GetDepth() )
299  {
300  fBlockedPhysicalVolume = fHistory.GetTopVolume();
301  fBlockedReplicaNo = fHistory.GetTopReplicaNo();
303  fExiting = false;
304 
305  if( noLevelsExited > 1 )
306  {
307  // The first transformation was done by the sub-navigator
308  //
309  const G4RotationMatrix* mRot = fBlockedPhysicalVolume->GetRotation();
310  if( mRot )
311  {
312  fGrandMotherExitNormal *= (*mRot).inverse();
313  fChangedGrandMotherRefFrame= true;
314  }
315  }
316  }
317  else
318  {
319  fLastLocatedPointLocal = localPoint;
320  fLocatedOutsideWorld = true;
321  // No extra transformation for ExitNormal - is in frame of Top Volume
322  return 0; // Have exited world volume
323  }
324  }
325  else
326  if ( insideCode==kSurface )
327  {
328  G4bool isExiting = fExiting;
329  if( (!fExiting)&&considerDirection )
330  {
331  // Figure out whether we are exiting this level's volume
332  // by using the direction
333  //
334  G4bool directionExiting = false;
335  G4ThreeVector localDirection =
336  fHistory.GetTopTransform().TransformAxis(globalDirection);
337 
338  // Make sure localPoint in correct reference frame
339  // ( Was it already correct ? How ? )
340  //
341  localPoint= fHistory.GetTopTransform().TransformPoint(globalPoint);
343  {
344  G4ThreeVector normal = targetSolid->SurfaceNormal(localPoint);
345  directionExiting = normal.dot(localDirection) > 0.0;
346  isExiting = isExiting || directionExiting;
347  }
348  }
349  if( isExiting )
350  {
351  noLevelsExited++;
352  if ( fHistory.GetDepth() )
353  {
354  fBlockedPhysicalVolume = fHistory.GetTopVolume();
355  fBlockedReplicaNo = fHistory.GetTopReplicaNo();
357  //
358  // Still on surface but exited volume not necessarily convex
359  //
360  fValidExitNormal = false;
361 
362  if( noLevelsExited > 1 )
363  {
364  // The first transformation was done by the sub-navigator
365  //
366  const G4RotationMatrix* mRot=fBlockedPhysicalVolume->GetRotation();
367  if( mRot )
368  {
369  fGrandMotherExitNormal *= (*mRot).inverse();
370  fChangedGrandMotherRefFrame= true;
371  }
372  }
373  }
374  else
375  {
376  fLastLocatedPointLocal = localPoint;
377  fLocatedOutsideWorld = true;
378  // No extra transformation for ExitNormal, is in frame of Top Vol
379  return 0; // Have exited world volume
380  }
381  }
382  else
383  {
384  notKnownContained=false;
385  }
386  }
387  else
388  {
389  notKnownContained=false;
390  }
391  } // END while (notKnownContained)
392  //
393  // Search downwards until deepest containing volume found,
394  // blocking fBlockedPhysicalVolume/BlockedReplicaNum
395  //
396  // 3 Cases:
397  //
398  // o Parameterised daughters
399  // =>Must be one G4PVParameterised daughter & voxels
400  // o Positioned daughters & voxels
401  // o Positioned daughters & no voxels
402 
403  noResult = true; // noResult should be renamed to
404  // something like enteredLevel, as that is its meaning.
405  do
406  {
407  // Determine `type' of current mother volume
408  //
409  targetPhysical = fHistory.GetTopVolume();
410  if (!targetPhysical) { break; }
411  targetLogical = targetPhysical->GetLogicalVolume();
412  switch( CharacteriseDaughters(targetLogical) )
413  {
414  case kNormal:
415  if ( targetLogical->GetVoxelHeader() ) // use optimised navigation
416  {
417  noResult = fvoxelNav.LevelLocate(fHistory,
418  fBlockedPhysicalVolume,
419  fBlockedReplicaNo,
420  globalPoint,
421  pGlobalDirection,
422  considerDirection,
423  localPoint);
424  }
425  else // do not use optimised navigation
426  {
427  noResult = fnormalNav.LevelLocate(fHistory,
428  fBlockedPhysicalVolume,
429  fBlockedReplicaNo,
430  globalPoint,
431  pGlobalDirection,
432  considerDirection,
433  localPoint);
434  }
435  break;
436  case kReplica:
437  noResult = freplicaNav.LevelLocate(fHistory,
438  fBlockedPhysicalVolume,
439  fBlockedReplicaNo,
440  globalPoint,
441  pGlobalDirection,
442  considerDirection,
443  localPoint);
444  break;
445  case kParameterised:
446  if( GetDaughtersRegularStructureId(targetLogical) != 1 )
447  {
448  noResult = fparamNav.LevelLocate(fHistory,
449  fBlockedPhysicalVolume,
450  fBlockedReplicaNo,
451  globalPoint,
452  pGlobalDirection,
453  considerDirection,
454  localPoint);
455  }
456  else // Regular structure
457  {
458  noResult = fregularNav.LevelLocate(fHistory,
459  fBlockedPhysicalVolume,
460  fBlockedReplicaNo,
461  globalPoint,
462  pGlobalDirection,
463  considerDirection,
464  localPoint);
465  }
466  break;
467  }
468 
469  // LevelLocate returns true if it finds a daughter volume
470  // in which globalPoint is inside (or on the surface).
471 
472  if ( noResult )
473  {
474  // Entering a daughter after ascending
475  //
476  // The blocked volume is no longer valid - it was for another level
477  //
478  fBlockedPhysicalVolume = 0;
479  fBlockedReplicaNo = -1;
480 
481  // fEntering should be false -- else blockedVolume is assumed good.
482  // fEnteredDaughter is used for ExitNormal
483  //
484  fEntering = false;
485  fEnteredDaughter = true;
486 
487  if( fExitedMother )
488  {
489  G4VPhysicalVolume* enteredPhysical = fHistory.GetTopVolume();
490  const G4RotationMatrix* mRot = enteredPhysical->GetRotation();
491  if( mRot )
492  {
493  fGrandMotherExitNormal *= (*mRot).inverse();
494  }
495  }
496 
497 #ifdef G4DEBUG_NAVIGATION
498  if( fVerbose > 2 )
499  {
500  G4VPhysicalVolume* enteredPhysical = fHistory.GetTopVolume();
501  G4cout << "*** G4Navigator::LocateGlobalPointAndSetup() ***" << G4endl;
502  G4cout << " Entering volume: " << enteredPhysical->GetName()
503  << G4endl;
504  }
505 #endif
506  }
507  } while (noResult);
508 
509  fLastLocatedPointLocal = localPoint;
510 
511 #ifdef G4VERBOSE
512  if( fVerbose >= 4 )
513  {
514  G4int oldcoutPrec = G4cout.precision(8);
515  G4String curPhysVol_Name("None");
516  if (targetPhysical) { curPhysVol_Name = targetPhysical->GetName(); }
517  G4cout << " Return value = new volume = " << curPhysVol_Name << G4endl;
518  G4cout << " ----- Upon exiting:" << G4endl;
519  PrintState();
520  if( fVerbose == 5 )
521  {
522  G4cout << "Upon exiting LocateGlobalPointAndSetup():" << G4endl;
523  G4cout << " History = " << G4endl << fHistory << G4endl << G4endl;
524  }
525  G4cout.precision(oldcoutPrec);
526  }
527 #endif
528 
529  fLocatedOutsideWorld= false;
530 
531  return targetPhysical;
532 }
533 
534 // ********************************************************************
535 // LocateGlobalPointWithinVolume
536 //
537 // -> the state information of this Navigator and its subNavigators
538 // is updated in order to start the next step at pGlobalpoint
539 // -> no check is performed whether pGlobalpoint is inside the
540 // original volume (this must be the case).
541 //
542 // Note: a direction could be added to the arguments, to aid in future
543 // optional checking (via the old code below, flagged by OLD_LOCATE).
544 // [ This would be done only in verbose mode ]
545 // ********************************************************************
546 //
547 void
549 {
550  fLastLocatedPointLocal = ComputeLocalPoint(pGlobalpoint);
551  fLastTriedStepComputation= false;
552  fChangedGrandMotherRefFrame= false; // Frame for Exit Normal
553 
554 #ifdef G4DEBUG_NAVIGATION
555  if( fVerbose > 2 )
556  {
557  G4cout << "Entering LocateGlobalWithinVolume(): History = " << G4endl;
558  G4cout << fHistory << G4endl;
559  }
560 #endif
561 
562  // For the case of Voxel (or Parameterised) volume the respective
563  // Navigator must be messaged to update its voxel information etc
564 
565  // Update the state of the Sub Navigators
566  // - in particular any voxel information they store/cache
567  //
568  G4VPhysicalVolume* motherPhysical = fHistory.GetTopVolume();
569  G4LogicalVolume* motherLogical = motherPhysical->GetLogicalVolume();
570  G4SmartVoxelHeader* pVoxelHeader = motherLogical->GetVoxelHeader();
571 
573  {
574  switch( CharacteriseDaughters(motherLogical) )
575  {
576  case kNormal:
577  if ( pVoxelHeader )
578  {
579  fvoxelNav.VoxelLocate( pVoxelHeader, fLastLocatedPointLocal );
580  }
581  break;
582  case kParameterised:
583  if( GetDaughtersRegularStructureId(motherLogical) != 1 )
584  {
585  // Resets state & returns voxel node
586  //
587  fparamNav.ParamVoxelLocate( pVoxelHeader, fLastLocatedPointLocal );
588  }
589  break;
590  case kReplica:
591  G4Exception("G4Navigator::LocateGlobalPointWithinVolume()",
592  "GeomNav0001", FatalException,
593  "Not applicable for replicated volumes.");
594  break;
595  }
596  }
597 
598  // Reset the state variables
599  // - which would have been affected
600  // by the 'equivalent' call to LocateGlobalPointAndSetup
601  // - who's values have been invalidated by the 'move'.
602  //
603  fBlockedPhysicalVolume = 0;
604  fBlockedReplicaNo = -1;
605  fEntering = false;
606  fEnteredDaughter = false; // Boundary not encountered, did not enter
607  fExiting = false;
608  fExitedMother = false; // Boundary not encountered, did not exit
609 }
610 
611 // ********************************************************************
612 // SetSavedState
613 //
614 // Save the state, in case this is a parasitic call
615 // Save fValidExitNormal, fExitNormal, fExiting, fEntering,
616 // fBlockedPhysicalVolume, fBlockedReplicaNo, fLastStepWasZero;
617 // ********************************************************************
618 //
620 {
621  fSaveState.sExitNormal = fExitNormal;
622  fSaveState.sValidExitNormal = fValidExitNormal;
623  fSaveState.sExiting = fExiting;
624  fSaveState.sEntering = fEntering;
625 
626  fSaveState.spBlockedPhysicalVolume = fBlockedPhysicalVolume;
627  fSaveState.sBlockedReplicaNo = fBlockedReplicaNo,
628 
629  fSaveState.sLastStepWasZero = fLastStepWasZero;
630 
631  fSaveState.sLocatedOutsideWorld = fLocatedOutsideWorld;
632  fSaveState.sLastLocatedPointLocal= fLastLocatedPointLocal;
633  fSaveState.sEnteredDaughter= fEnteredDaughter;
634  fSaveState.sExitedMother= fExitedMother;
635 
636  // Even the safety sphere - if you want to change it do it explicitly!
637  //
638  fSaveState.sPreviousSftOrigin= fPreviousSftOrigin;
639  fSaveState.sPreviousSafety= fPreviousSafety;
640 }
641 
642 // ********************************************************************
643 // RestoreSavedState
644 //
645 // Restore the state (in Compute Step), in case this is a parasitic call
646 // ********************************************************************
647 //
649 {
650  fExitNormal = fSaveState.sExitNormal;
651  fValidExitNormal = fSaveState.sValidExitNormal;
652  fExiting = fSaveState.sExiting;
653  fEntering = fSaveState.sEntering;
654 
655  fBlockedPhysicalVolume = fSaveState.spBlockedPhysicalVolume;
656  fBlockedReplicaNo = fSaveState.sBlockedReplicaNo,
657 
658  fLastStepWasZero = fSaveState.sLastStepWasZero;
659 
660  fLocatedOutsideWorld = fSaveState.sLocatedOutsideWorld;
661  fLastLocatedPointLocal= fSaveState.sLastLocatedPointLocal;
662  fEnteredDaughter= fSaveState.sEnteredDaughter;
663  fExitedMother= fSaveState.sExitedMother;
664  fSaveState.sPreviousSftOrigin= fPreviousSftOrigin;
665  fSaveState.sPreviousSafety= fPreviousSafety;
666 }
667 
668 // ********************************************************************
669 // ComputeStep
670 //
671 // Computes the next geometric Step: intersections with current
672 // mother and `daughter' volumes.
673 //
674 // NOTE:
675 //
676 // Flags on entry:
677 // --------------
678 // fValidExitNormal - Normal of exited volume is valid (convex, not a
679 // coincident boundary)
680 // fExitNormal - Surface normal of exited volume
681 // fExiting - True if have exited solid
682 //
683 // fBlockedPhysicalVolume - Ptr to exited volume (or 0)
684 // fBlockedReplicaNo - Replication no of exited volume
685 // fLastStepWasZero - True if last Step size was zero.
686 //
687 // Flags on exit:
688 // -------------
689 // fValidExitNormal - True if surface normal of exited volume is valid
690 // fExitNormal - Surface normal of exited volume rotated to mothers
691 // reference system
692 // fExiting - True if exiting mother
693 // fEntering - True if entering `daughter' volume (or replica)
694 // fBlockedPhysicalVolume - Ptr to candidate (entered) volume
695 // fBlockedReplicaNo - Replication no of candidate (entered) volume
696 // fLastStepWasZero - True if this Step size was zero.
697 // ********************************************************************
698 //
700  const G4ThreeVector &pDirection,
701  const G4double pCurrentProposedStepLength,
702  G4double &pNewSafety)
703 {
704  G4ThreeVector localDirection = ComputeLocalAxis(pDirection);
705  G4double Step = kInfinity;
706  G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
707  G4LogicalVolume *motherLogical = motherPhysical->GetLogicalVolume();
708 
709  // All state relating to exiting normals must be reset
710  //
711  fExitNormalGlobalFrame= G4ThreeVector( 0., 0., 0.);
712  // Reset value - to erase its memory
713  fChangedGrandMotherRefFrame= false;
714  // Reset - used for local exit normal
715  fGrandMotherExitNormal= G4ThreeVector( 0., 0., 0.);
716  fCalculatedExitNormal = false;
717  // Reset for new step
718 
719  static G4int sNavCScalls=0;
720  sNavCScalls++;
721 
722  fLastTriedStepComputation= true;
723 
724 #ifdef G4VERBOSE
725  if( fVerbose > 0 )
726  {
727  G4cout << "*** G4Navigator::ComputeStep: ***" << G4endl;
728  G4cout << " Volume = " << motherPhysical->GetName()
729  << " - Proposed step length = " << pCurrentProposedStepLength
730  << G4endl;
731 #ifdef G4DEBUG_NAVIGATION
732  if( fVerbose >= 2 )
733  {
734  G4cout << " Called with the arguments: " << G4endl
735  << " Globalpoint = " << std::setw(25) << pGlobalpoint << G4endl
736  << " Direction = " << std::setw(25) << pDirection << G4endl;
737  if( fVerbose >= 4 )
738  {
739  G4cout << " ---- Upon entering : State" << G4endl;
740  PrintState();
741  }
742  }
743 #endif
744  }
745 #endif
746 
747  G4ThreeVector newLocalPoint = ComputeLocalPoint(pGlobalpoint);
748  if( newLocalPoint != fLastLocatedPointLocal )
749  {
750  // Check whether the relocation is within safety
751  //
752  G4ThreeVector oldLocalPoint = fLastLocatedPointLocal;
753  G4double moveLenSq = (newLocalPoint-oldLocalPoint).mag2();
754 
755  if ( moveLenSq >= kCarTolerance*kCarTolerance )
756  {
757 #ifdef G4VERBOSE
758  ComputeStepLog(pGlobalpoint, moveLenSq);
759 #endif
760  // Relocate the point within the same volume
761  //
762  LocateGlobalPointWithinVolume( pGlobalpoint );
763  fLastTriedStepComputation= true; // Ensure that this is set again !!
764  }
765  }
767  {
768  switch( CharacteriseDaughters(motherLogical) )
769  {
770  case kNormal:
771  if ( motherLogical->GetVoxelHeader() )
772  {
773  Step = fvoxelNav.ComputeStep(fLastLocatedPointLocal,
774  localDirection,
775  pCurrentProposedStepLength,
776  pNewSafety,
777  fHistory,
778  fValidExitNormal,
779  fExitNormal,
780  fExiting,
781  fEntering,
782  &fBlockedPhysicalVolume,
783  fBlockedReplicaNo);
784 
785  }
786  else
787  {
788  if( motherPhysical->GetRegularStructureId() == 0 )
789  {
790  Step = fnormalNav.ComputeStep(fLastLocatedPointLocal,
791  localDirection,
792  pCurrentProposedStepLength,
793  pNewSafety,
794  fHistory,
795  fValidExitNormal,
796  fExitNormal,
797  fExiting,
798  fEntering,
799  &fBlockedPhysicalVolume,
800  fBlockedReplicaNo);
801  }
802  else // Regular (non-voxelised) structure
803  {
804  LocateGlobalPointAndSetup( pGlobalpoint, &pDirection, true, true );
805  fLastTriedStepComputation= true; // Ensure that this is set again !!
806  //
807  // if physical process limits the step, the voxel will not be the
808  // one given by ComputeStepSkippingEqualMaterials() and the local
809  // point will be wrongly calculated.
810 
811  // There is a problem: when msc limits the step and the point is
812  // assigned wrongly to phantom in previous step (while it is out
813  // of the container volume). Then LocateGlobalPointAndSetup() has
814  // reset the history topvolume to world.
815  //
817  {
818  G4Exception("G4Navigator::ComputeStep()",
819  "GeomNav1001", JustWarning,
820  "Point is relocated in voxels, while it should be outside!");
821  Step = fnormalNav.ComputeStep(fLastLocatedPointLocal,
822  localDirection,
823  pCurrentProposedStepLength,
824  pNewSafety,
825  fHistory,
826  fValidExitNormal,
827  fExitNormal,
828  fExiting,
829  fEntering,
830  &fBlockedPhysicalVolume,
831  fBlockedReplicaNo);
832  }
833  else
834  {
835  Step = fregularNav.
836  ComputeStepSkippingEqualMaterials(fLastLocatedPointLocal,
837  localDirection,
838  pCurrentProposedStepLength,
839  pNewSafety,
840  fHistory,
841  fValidExitNormal,
842  fExitNormal,
843  fExiting,
844  fEntering,
845  &fBlockedPhysicalVolume,
846  fBlockedReplicaNo,
847  motherPhysical);
848  }
849  }
850  }
851  break;
852  case kParameterised:
853  if( GetDaughtersRegularStructureId(motherLogical) != 1 )
854  {
855  Step = fparamNav.ComputeStep(fLastLocatedPointLocal,
856  localDirection,
857  pCurrentProposedStepLength,
858  pNewSafety,
859  fHistory,
860  fValidExitNormal,
861  fExitNormal,
862  fExiting,
863  fEntering,
864  &fBlockedPhysicalVolume,
865  fBlockedReplicaNo);
866  }
867  else // Regular structure
868  {
869  Step = fregularNav.ComputeStep(fLastLocatedPointLocal,
870  localDirection,
871  pCurrentProposedStepLength,
872  pNewSafety,
873  fHistory,
874  fValidExitNormal,
875  fExitNormal,
876  fExiting,
877  fEntering,
878  &fBlockedPhysicalVolume,
879  fBlockedReplicaNo);
880  }
881  break;
882  case kReplica:
883  G4Exception("G4Navigator::ComputeStep()", "GeomNav0001",
884  FatalException, "Not applicable for replicated volumes.");
885  break;
886  }
887  }
888  else
889  {
890  // In the case of a replica, it must handle the exiting
891  // edge/corner problem by itself
892  //
893  G4bool exitingReplica = fExitedMother;
894  G4bool calculatedExitNormal;
895  Step = freplicaNav.ComputeStep(pGlobalpoint,
896  pDirection,
897  fLastLocatedPointLocal,
898  localDirection,
899  pCurrentProposedStepLength,
900  pNewSafety,
901  fHistory,
902  fValidExitNormal,
903  calculatedExitNormal,
904  fExitNormal,
905  exitingReplica,
906  fEntering,
907  &fBlockedPhysicalVolume,
908  fBlockedReplicaNo);
909  fExiting= exitingReplica;
910  fCalculatedExitNormal= calculatedExitNormal;
911  }
912 
913  // Remember last safety origin & value.
914  //
915  fPreviousSftOrigin = pGlobalpoint;
916  fPreviousSafety = pNewSafety;
917 
918  // Count zero steps - one can occur due to changing momentum at a boundary
919  // - one, two (or a few) can occur at common edges between
920  // volumes
921  // - more than two is likely a problem in the geometry
922  // description or the Navigation
923 
924  // Rule of thumb: likely at an Edge if two consecutive steps are zero,
925  // because at least two candidate volumes must have been
926  // checked
927  //
928  fLocatedOnEdge = fLastStepWasZero && (Step==0.0);
929  fLastStepWasZero = (Step==0.0);
930  if (fPushed) { fPushed = fLastStepWasZero; }
931 
932  // Handle large number of consecutive zero steps
933  //
934  if ( fLastStepWasZero )
935  {
936  fNumberZeroSteps++;
937 #ifdef G4DEBUG_NAVIGATION
938  if( fNumberZeroSteps > 1 )
939  {
940  G4cout << "G4Navigator::ComputeStep(): another zero step, # "
941  << fNumberZeroSteps
942  << " at " << pGlobalpoint
943  << " in volume " << motherPhysical->GetName()
944  << " nav-comp-step calls # " << sNavCScalls
945  << G4endl;
946  }
947 #endif
948  if( fNumberZeroSteps > fActionThreshold_NoZeroSteps-1 )
949  {
950  // Act to recover this stuck track. Pushing it along direction
951  //
952  Step += 100*kCarTolerance;
953 #ifdef G4VERBOSE
954  if ((!fPushed) && (fWarnPush))
955  {
956  std::ostringstream message;
957  message << "Track stuck or not moving." << G4endl
958  << " Track stuck, not moving for "
959  << fNumberZeroSteps << " steps" << G4endl
960  << " in volume -" << motherPhysical->GetName()
961  << "- at point " << pGlobalpoint << G4endl
962  << " direction: " << pDirection << "." << G4endl
963  << " Potential geometry or navigation problem !"
964  << G4endl
965  << " Trying pushing it of " << Step << " mm ...";
966  G4Exception("G4Navigator::ComputeStep()", "GeomNav1002",
967  JustWarning, message, "Potential overlap in geometry!");
968  }
969 #endif
970  fPushed = true;
971  }
972  if( fNumberZeroSteps > fAbandonThreshold_NoZeroSteps-1 )
973  {
974  // Must kill this stuck track
975  //
976  std::ostringstream message;
977  message << "Stuck Track: potential geometry or navigation problem."
978  << G4endl
979  << " Track stuck, not moving for "
980  << fNumberZeroSteps << " steps" << G4endl
981  << " in volume -" << motherPhysical->GetName()
982  << "- at point " << pGlobalpoint << G4endl
983  << " direction: " << pDirection << ".";
984  motherPhysical->CheckOverlaps(5000, false);
985  G4Exception("G4Navigator::ComputeStep()", "GeomNav0003",
986  EventMustBeAborted, message);
987  }
988  }
989  else
990  {
991  if (!fPushed) fNumberZeroSteps = 0;
992  }
993 
994  fEnteredDaughter = fEntering; // I expect to enter a volume in this Step
995  fExitedMother = fExiting;
996 
997  fStepEndPoint = pGlobalpoint + Step * pDirection;
998  fLastStepEndPointLocal = fLastLocatedPointLocal + Step * localDirection;
999 
1000  if( fExiting )
1001  {
1002 #ifdef G4DEBUG_NAVIGATION
1003  if( fVerbose > 2 )
1004  {
1005  G4cout << " At G4Nav CompStep End - if(exiting) - fExiting= " << fExiting
1006  << " fValidExitNormal = " << fValidExitNormal << G4endl;
1007  G4cout << " fExitNormal= " << fExitNormal << G4endl;
1008  }
1009 #endif
1010 
1011  if(fValidExitNormal || fCalculatedExitNormal)
1012  {
1014  {
1015  // Convention: fExitNormal is in the 'grand-mother' coordinate system
1016  //
1017  fGrandMotherExitNormal= fExitNormal;
1018  fCalculatedExitNormal= true;
1019  }
1020  else
1021  {
1022  fGrandMotherExitNormal = fExitNormal;
1023  }
1024  }
1025  else
1026  {
1027  // We must calculate the normal anyway (in order to have it if requested)
1028  //
1029  G4ThreeVector finalLocalPoint =
1030  fLastLocatedPointLocal + localDirection*Step;
1031 
1033  {
1034  // Find normal in the 'mother' coordinate system
1035  //
1036  G4ThreeVector exitNormalMotherFrame=
1037  motherLogical->GetSolid()->SurfaceNormal(finalLocalPoint);
1038 
1039  // Transform it to the 'grand-mother' coordinate system
1040  //
1041  const G4RotationMatrix* mRot = motherPhysical->GetRotation();
1042  if( mRot )
1043  {
1044  fChangedGrandMotherRefFrame= true;
1045  fGrandMotherExitNormal = (*mRot).inverse() * exitNormalMotherFrame;
1046  }
1047  else
1048  {
1049  fGrandMotherExitNormal = exitNormalMotherFrame;
1050  }
1051 
1052  // Do not set fValidExitNormal -- this signifies
1053  // that the solid is convex!
1054  //
1055  fCalculatedExitNormal= true;
1056  }
1057  else
1058  {
1059  fCalculatedExitNormal = false;
1060  //
1061  // Nothing can be done at this stage currently - to solve this
1062  // Replica Navigation must have calculated the normal for this case
1063  // already.
1064  // Cases: mother is not convex, and exit is at previous replica level
1065 
1066 #ifdef G4DEBUG_NAVIGATION
1068 
1069  desc << "Problem in ComputeStep: Replica Navigation did not provide"
1070  << " valid exit Normal. " << G4endl;
1071  desc << " Do not know how calculate it in this case." << G4endl;
1072  desc << " Location = " << finalLocalPoint << G4endl;
1073  desc << " Volume name = " << motherPhysical->GetName()
1074  << " copy/replica No = " << motherPhysical->GetCopyNo() << G4endl;
1075  G4Exception("G4Navigator::ComputeStep()", "GeomNav0003",
1076  JustWarning, desc, "Normal not available for exiting.");
1077 #endif
1078  }
1079  }
1080 
1081  // Now transform it to the global reference frame !!
1082  //
1083  if( fValidExitNormal || fCalculatedExitNormal )
1084  {
1085  G4int depth= fHistory.GetDepth();
1086  if( depth > 0 )
1087  {
1088  G4AffineTransform GrandMotherToGlobalTransf =
1089  fHistory.GetTransform(depth-1).Inverse();
1090  fExitNormalGlobalFrame =
1091  GrandMotherToGlobalTransf.TransformAxis( fGrandMotherExitNormal );
1092  }
1093  else
1094  {
1095  fExitNormalGlobalFrame= fGrandMotherExitNormal;
1096  }
1097  }
1098  else
1099  {
1100  fExitNormalGlobalFrame= G4ThreeVector( 0., 0., 0.);
1101  }
1102  }
1103  fStepEndPoint= pGlobalpoint+Step*pDirection;
1104 
1105  if( (Step == pCurrentProposedStepLength) && (!fExiting) && (!fEntering) )
1106  {
1107  // This if Step is not really limited by the geometry.
1108  // The Navigator is obliged to return "infinity"
1109  //
1110  Step = kInfinity;
1111  }
1112 
1113 #ifdef G4VERBOSE
1114  if( fVerbose > 1 )
1115  {
1116  if( fVerbose >= 4 )
1117  {
1118  G4cout << " ----- Upon exiting :" << G4endl;
1119  PrintState();
1120  }
1121  G4cout << " Returned step= " << Step;
1122  if( fVerbose > 5 ) G4cout << G4endl;
1123  if( Step == kInfinity )
1124  {
1125  G4cout << " Requested step= " << pCurrentProposedStepLength ;
1126  if( fVerbose > 5) G4cout << G4endl;
1127  }
1128  G4cout << " Safety = " << pNewSafety << G4endl;
1129  }
1130 #endif
1131 
1132  return Step;
1133 }
1134 
1135 // ********************************************************************
1136 // CheckNextStep
1137 //
1138 // Compute the step without altering the navigator state
1139 // ********************************************************************
1140 //
1142  const G4ThreeVector& pDirection,
1143  const G4double pCurrentProposedStepLength,
1144  G4double& pNewSafety)
1145 {
1146  G4double step;
1147 
1148  // Save the state, for this parasitic call
1149  //
1150  SetSavedState();
1151 
1152  step = ComputeStep ( pGlobalpoint,
1153  pDirection,
1154  pCurrentProposedStepLength,
1155  pNewSafety );
1156 
1157  // If a parasitic call, then attempt to restore the key parts of the state
1158  //
1159  RestoreSavedState();
1160 
1161  return step;
1162 }
1163 
1164 // ********************************************************************
1165 // ResetState
1166 //
1167 // Resets stack and minimum of navigator state `machine'
1168 // ********************************************************************
1169 //
1171 {
1172  fWasLimitedByGeometry = false;
1173  fEntering = false;
1174  fExiting = false;
1175  fLocatedOnEdge = false;
1176  fLastStepWasZero = false;
1177  fEnteredDaughter = false;
1178  fExitedMother = false;
1179  fPushed = false;
1180 
1181  fValidExitNormal = false;
1182  fChangedGrandMotherRefFrame= false;
1183  fCalculatedExitNormal = false;
1184 
1185  fExitNormal = G4ThreeVector(0,0,0);
1186  fGrandMotherExitNormal = G4ThreeVector(0,0,0);
1187  fExitNormalGlobalFrame = G4ThreeVector(0,0,0);
1188 
1189  fPreviousSftOrigin = G4ThreeVector(0,0,0);
1190  fPreviousSafety = 0.0;
1191 
1192  fNumberZeroSteps = 0;
1193 
1194  fBlockedPhysicalVolume = 0;
1195  fBlockedReplicaNo = -1;
1196 
1197  fLastLocatedPointLocal = G4ThreeVector( kInfinity, -kInfinity, 0.0 );
1198  fLocatedOutsideWorld = false;
1199 }
1200 
1201 // ********************************************************************
1202 // SetupHierarchy
1203 //
1204 // Renavigates & resets hierarchy described by current history
1205 // o Reset volumes
1206 // o Recompute transforms and/or solids of replicated/parameterised volumes
1207 // ********************************************************************
1208 //
1210 {
1211  G4int i;
1212  const G4int cdepth = fHistory.GetDepth();
1213  G4VPhysicalVolume *current;
1214  G4VSolid *pSolid;
1215  G4VPVParameterisation *pParam;
1216 
1217  for ( i=1; i<=cdepth; i++ )
1218  {
1219  current = fHistory.GetVolume(i);
1220  switch ( fHistory.GetVolumeType(i) )
1221  {
1222  case kNormal:
1223  break;
1224  case kReplica:
1225  freplicaNav.ComputeTransformation(fHistory.GetReplicaNo(i), current);
1226  break;
1227  case kParameterised:
1228  G4int replicaNo;
1229  pParam = current->GetParameterisation();
1230  replicaNo = fHistory.GetReplicaNo(i);
1231  pSolid = pParam->ComputeSolid(replicaNo, current);
1232 
1233  // Set up dimensions & transform in solid/physical volume
1234  //
1235  pSolid->ComputeDimensions(pParam, replicaNo, current);
1236  pParam->ComputeTransformation(replicaNo, current);
1237 
1238  G4TouchableHistory touchable( fHistory );
1239  touchable.MoveUpHistory(); // move up to the parent level
1240 
1241  // Set up the correct solid and material in Logical Volume
1242  //
1243  G4LogicalVolume *pLogical = current->GetLogicalVolume();
1244  pLogical->SetSolid( pSolid );
1245  pLogical->UpdateMaterial( pParam ->
1246  ComputeMaterial(replicaNo, current, &touchable) );
1247  break;
1248  }
1249  }
1250 }
1251 
1252 // ********************************************************************
1253 // GetLocalExitNormal
1254 //
1255 // Obtains the Normal vector to a surface (in local coordinates)
1256 // pointing out of previous volume and into current volume
1257 // ********************************************************************
1258 //
1260 {
1261  G4ThreeVector ExitNormal(0.,0.,0.);
1262  G4VSolid *currentSolid=0;
1263  G4LogicalVolume *candidateLogical;
1264  if ( fLastTriedStepComputation )
1265  {
1266  // use fLastLocatedPointLocal and next candidate volume
1267  //
1268  G4ThreeVector nextSolidExitNormal(0.,0.,0.);
1269 
1270  if( fEntering && (fBlockedPhysicalVolume!=0) )
1271  {
1272  candidateLogical= fBlockedPhysicalVolume->GetLogicalVolume();
1273  if( candidateLogical )
1274  {
1275  // fLastStepEndPointLocal is in the coordinates of the mother
1276  // we need it in the daughter's coordinate system.
1277 
1278  // The following code should also work in case of Replica
1279  {
1280  // First transform fLastLocatedPointLocal to the new daughter
1281  // coordinates
1282  //
1283  G4AffineTransform MotherToDaughterTransform=
1284  GetMotherToDaughterTransform( fBlockedPhysicalVolume,
1285  fBlockedReplicaNo,
1286  VolumeType(fBlockedPhysicalVolume) );
1287  G4ThreeVector daughterPointOwnLocal=
1288  MotherToDaughterTransform.TransformPoint( fLastStepEndPointLocal );
1289 
1290  // OK if it is a parameterised volume
1291  //
1292  EInside inSideIt;
1293  G4bool onSurface;
1294  G4double safety= -1.0;
1295  currentSolid= candidateLogical->GetSolid();
1296  inSideIt = currentSolid->Inside(daughterPointOwnLocal);
1297  onSurface = (inSideIt == kSurface);
1298  if( ! onSurface )
1299  {
1300  if( inSideIt == kOutside )
1301  {
1302  safety = (currentSolid->DistanceToIn(daughterPointOwnLocal));
1303  onSurface = safety < 100.0 * kCarTolerance;
1304  }
1305  else if (inSideIt == kInside )
1306  {
1307  safety = (currentSolid->DistanceToOut(daughterPointOwnLocal));
1308  onSurface = safety < 100.0 * kCarTolerance;
1309  }
1310  }
1311 
1312  if( onSurface )
1313  {
1314  nextSolidExitNormal =
1315  currentSolid->SurfaceNormal(daughterPointOwnLocal);
1316 
1317  // Entering the solid ==> opposite
1318  //
1319  ExitNormal = -nextSolidExitNormal;
1320  fCalculatedExitNormal= true;
1321  }
1322  else
1323  {
1324 #ifdef G4VERBOSE
1325  if(( fVerbose == 1 ) && ( fCheck ))
1326  {
1327  std::ostringstream message;
1328  message << "Point not on surface ! " << G4endl
1329  << " Point = "
1330  << daughterPointOwnLocal << G4endl
1331  << " Physical volume = "
1332  << fBlockedPhysicalVolume->GetName() << G4endl
1333  << " Logical volume = "
1334  << candidateLogical->GetName() << G4endl
1335  << " Solid = " << currentSolid->GetName()
1336  << " Type = "
1337  << currentSolid->GetEntityType() << G4endl
1338  << *currentSolid << G4endl;
1339  if( inSideIt == kOutside )
1340  {
1341  message << "Point is Outside. " << G4endl
1342  << " Safety (from outside) = " << safety << G4endl;
1343  }
1344  else // if( inSideIt == kInside )
1345  {
1346  message << "Point is Inside. " << G4endl
1347  << " Safety (from inside) = " << safety << G4endl;
1348  }
1349  G4Exception("G4Navigator::GetLocalExitNormal()", "GeomNav1001",
1350  JustWarning, message);
1351  }
1352 #endif
1353  }
1354  *valid = onSurface; // was =true;
1355  }
1356  }
1357  }
1358  else if ( fExiting )
1359  {
1360  ExitNormal = fGrandMotherExitNormal;
1361  *valid = true;
1362  fCalculatedExitNormal= true; // Should be true already
1363  }
1364  else // i.e. ( fBlockedPhysicalVolume == 0 )
1365  {
1366  *valid = false;
1367  G4Exception("G4Navigator::GetLocalExitNormal()",
1368  "GeomNav0003", JustWarning,
1369  "Incorrect call to GetLocalSurfaceNormal." );
1370  }
1371  }
1372  else // ( ! fLastTriedStepComputation ) ie. last call was to Locate
1373  {
1374  if ( EnteredDaughterVolume() )
1375  {
1376  G4VSolid* daughterSolid =fHistory.GetTopVolume()->GetLogicalVolume()
1377  ->GetSolid();
1378  ExitNormal= -(daughterSolid->SurfaceNormal(fLastLocatedPointLocal));
1379  if( std::fabs(ExitNormal.mag2()-1.0 ) > CLHEP::perMillion )
1380  {
1382  desc << " Parameters of solid: " << *daughterSolid
1383  << " Point for surface = " << fLastLocatedPointLocal << std::endl;
1384  G4Exception("G4Navigator::GetLocalExitNormal()",
1385  "GeomNav0003", FatalException, desc,
1386  "Surface Normal returned by Solid is not a Unit Vector." );
1387  }
1388  fCalculatedExitNormal= true;
1389  *valid = true;
1390  }
1391  else
1392  {
1393  if( fExitedMother )
1394  {
1395  ExitNormal = fGrandMotherExitNormal;
1396  *valid = true;
1397  fCalculatedExitNormal= true;
1398  }
1399  else // We are not at a boundary. ExitNormal remains (0,0,0)
1400  {
1401  *valid = false;
1402  fCalculatedExitNormal= false;
1403  G4ExceptionDescription message;
1404  message << "Function called when *NOT* at a Boundary." << G4endl;
1405  G4Exception("G4Navigator::GetLocalExitNormal()",
1406  "GeomNav0003", JustWarning, message);
1407  }
1408  }
1409  }
1410  return ExitNormal;
1411 }
1412 
1413 // ********************************************************************
1414 // GetMotherToDaughterTransform
1415 //
1416 // Obtains the mother to daughter affine transformation
1417 // ********************************************************************
1418 //
1421  G4int enteringReplicaNo,
1422  EVolume enteringVolumeType )
1423 {
1424  switch (enteringVolumeType)
1425  {
1426  case kNormal: // Nothing is needed to prepare the transformation
1427  break; // It is stored already in the physical volume (placement)
1428  case kReplica: // Sets the transform in the Replica - tbc
1429  G4Exception("G4Navigator::GetMotherToDaughterTransform()",
1430  "GeomNav0001", FatalException,
1431  "Method NOT Implemented yet for replica volumes.");
1432  break;
1433  case kParameterised:
1434  if( pEnteringPhysVol->GetRegularStructureId() == 0 )
1435  {
1436  G4VPVParameterisation *pParam =
1437  pEnteringPhysVol->GetParameterisation();
1438  G4VSolid* pSolid =
1439  pParam->ComputeSolid(enteringReplicaNo, pEnteringPhysVol);
1440  pSolid->ComputeDimensions(pParam, enteringReplicaNo, pEnteringPhysVol);
1441 
1442  // Sets the transform in the Parameterisation
1443  //
1444  pParam->ComputeTransformation(enteringReplicaNo, pEnteringPhysVol);
1445 
1446  // Set the correct solid and material in Logical Volume
1447  //
1448  G4LogicalVolume* pLogical = pEnteringPhysVol->GetLogicalVolume();
1449  pLogical->SetSolid( pSolid );
1450  }
1451  break;
1452  }
1453  return G4AffineTransform(pEnteringPhysVol->GetRotation(),
1454  pEnteringPhysVol->GetTranslation()).Invert();
1455 }
1456 
1457 // ********************************************************************
1458 // GetLocalExitNormalAndCheck
1459 //
1460 // Obtains the Normal vector to a surface (in local coordinates)
1461 // pointing out of previous volume and into current volume, and
1462 // checks the current point against expected 'local' value.
1463 // ********************************************************************
1464 //
1467 #ifdef G4DEBUG_NAVIGATION
1468  const G4ThreeVector& ExpectedBoundaryPointGlobal,
1469 #else
1470  const G4ThreeVector&,
1471 #endif
1472  G4bool* pValid)
1473 {
1474 #ifdef G4DEBUG_NAVIGATION
1475  // Check Current point against expected 'local' value
1476  //
1477  if ( fLastTriedStepComputation )
1478  {
1479  G4ThreeVector ExpectedBoundaryPointLocal;
1480 
1481  const G4AffineTransform& GlobalToLocal= GetGlobalToLocalTransform();
1482  ExpectedBoundaryPointLocal =
1483  GlobalToLocal.TransformPoint( ExpectedBoundaryPointGlobal );
1484 
1485  // Add here: Comparison against expected position,
1486  // i.e. the endpoint of ComputeStep
1487  }
1488 #endif
1489 
1490  return GetLocalExitNormal( pValid);
1491 }
1492 
1493 // ********************************************************************
1494 // GetGlobalExitNormal
1495 //
1496 // Obtains the Normal vector to a surface (in global coordinates)
1497 // pointing out of previous volume and into current volume
1498 // ********************************************************************
1499 //
1502  G4bool* pNormalCalculated)
1503 {
1504  G4bool validNormal;
1505  G4ThreeVector localNormal, globalNormal;
1506 
1507  if( fLastTriedStepComputation && fExiting )
1508  {
1509  // This was computed in ComputeStep -- and only on arrival at boundary
1510  //
1511  globalNormal = fExitNormalGlobalFrame;
1512  *pNormalCalculated = true; // ComputeStep always computes it if Exiting
1513  // (fExiting==true)
1514  }
1515  else
1516  {
1517  localNormal = GetLocalExitNormalAndCheck(IntersectPointGlobal,&validNormal);
1518  *pNormalCalculated = fCalculatedExitNormal;
1519 
1520 #ifdef G4DEBUG_NAVIGATION
1521  if( (!validNormal) && !fCalculatedExitNormal)
1522  {
1524  edN << " Calculated = " << fCalculatedExitNormal << G4endl;
1525  edN << " Entering= " << fEntering << G4endl;
1526  G4int oldVerbose= this->GetVerboseLevel();
1527  this->SetVerboseLevel(4);
1528  edN << " State of Navigator: " << G4endl;
1529  edN << *this << G4endl;
1530  this->SetVerboseLevel( oldVerbose );
1531 
1532  G4Exception("G4Navigator::GetGlobalExitNormal()",
1533  "GeomNav0003", JustWarning, edN,
1534  "LocalExitNormalAndCheck() did not calculate Normal.");
1535  }
1536 #endif
1537 
1538  G4double localMag2= localNormal.mag2();
1539  if( validNormal && (std::fabs(localMag2-1.0)) > CLHEP::perMillion )
1540  {
1542 
1543  edN << "G4Navigator::GetGlobalExitNormal: "
1544  << " Using Local Normal - from call to GetLocalExitNormalAndCheck. "
1545  << G4endl
1546  << " Local Exit Normal = " << localNormal << " || = "
1547  << std::sqrt(localMag2) << G4endl
1548  << " Global Exit Normal = " << globalNormal << " || = "
1549  << globalNormal.mag() << G4endl;
1550  edN << " Calculated It = " << fCalculatedExitNormal << G4endl;
1551 
1552  G4Exception("G4Navigator::GetGlobalExitNormal()",
1553  "GeomNav0003",JustWarning, edN,
1554  "Value obtained from new local *solid* is incorrect.");
1555  localNormal = localNormal.unit(); // Should we correct it ??
1556  }
1557  G4AffineTransform localToGlobal = GetLocalToGlobalTransform();
1558  globalNormal = localToGlobal.TransformAxis( localNormal );
1559  }
1560 
1561 #ifdef G4DEBUG_NAVIGATION
1562  // Temporary extra checks
1563  if( fLastTriedStepComputation && fExiting)
1564  {
1565  localNormal = GetLocalExitNormalAndCheck( IntersectPointGlobal, &validNormal);
1566  *pNormalCalculated = fCalculatedExitNormal;
1567 
1568  G4AffineTransform localToGlobal = GetLocalToGlobalTransform();
1569  globalNormal = localToGlobal.TransformAxis( localNormal );
1570 
1571  // Check the value computed against fExitNormalGlobalFrame
1572  G4ThreeVector diffNorm = globalNormal - fExitNormalGlobalFrame;
1573  if( diffNorm.mag2() > perMillion*CLHEP::perMillion)
1574  {
1575  G4ExceptionDescription edDfn;
1576  edDfn << "Found difference in normals in case of exiting mother "
1577  << "- when Get is called after ComputingStep " << G4endl;
1578  edDfn << " Magnitude of diff = " << diffNorm.mag() << G4endl;
1579  edDfn << " Normal stored (Global) = " << fExitNormalGlobalFrame
1580  << G4endl;
1581  edDfn << " Global Computed from Local = " << globalNormal << G4endl;
1582  G4Exception("G4Navigator::GetGlobalExitNormal()", "GeomNav0003",
1583  JustWarning, edDfn);
1584  }
1585  }
1586 #endif
1587 
1588  return globalNormal;
1589 }
1590 
1591 // To make the new Voxel Safety the default, uncomment the next line
1592 // #define G4NEW_SAFETY 1
1593 
1594 // ********************************************************************
1595 // ComputeSafety
1596 //
1597 // It assumes that it will be
1598 // i) called at the Point in the same volume as the EndPoint of the
1599 // ComputeStep.
1600 // ii) after (or at the end of) ComputeStep OR after the relocation.
1601 // ********************************************************************
1602 //
1604  const G4double pMaxLength,
1605  const G4bool keepState)
1606 {
1607  G4double newSafety = 0.0;
1608 
1609 #ifdef G4DEBUG_NAVIGATION
1610  G4int oldcoutPrec = G4cout.precision(8);
1611  if( fVerbose > 0 )
1612  {
1613  G4cout << "*** G4Navigator::ComputeSafety: ***" << G4endl
1614  << " Called at point: " << pGlobalpoint << G4endl;
1615 
1616  G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
1617  G4cout << " Volume = " << motherPhysical->GetName()
1618  << " - Maximum length = " << pMaxLength << G4endl;
1619  if( fVerbose >= 4 )
1620  {
1621  G4cout << " ----- Upon entering Compute Safety:" << G4endl;
1622  PrintState();
1623  }
1624  }
1625 #endif
1626 
1627  if (keepState) { SetSavedState(); }
1628 
1629  G4double distEndpointSq = (pGlobalpoint-fStepEndPoint).mag2();
1630  G4bool stayedOnEndpoint = distEndpointSq < kCarTolerance*kCarTolerance;
1631  G4bool endpointOnSurface = fEnteredDaughter || fExitedMother;
1632 
1633  if( !(endpointOnSurface && stayedOnEndpoint) )
1634  {
1635  // Pseudo-relocate to this point (updates voxel information only)
1636  //
1637  LocateGlobalPointWithinVolume( pGlobalpoint );
1638  // --->> DANGER: Side effects on sub-navigator voxel information <<---
1639  // Could be replaced again by 'granular' calls to sub-navigator
1640  // locates (similar side-effects, but faster.
1641  // Solutions:
1642  // 1) Re-locate (to where?)
1643  // 2) Insure that the methods using (G4ComputeStep?)
1644  // does a relocation (if information is disturbed only ?)
1645 
1646 #ifdef G4DEBUG_NAVIGATION
1647  if( fVerbose >= 2 )
1648  {
1649  G4cout << " G4Navigator::ComputeSafety() relocates-in-volume to point: "
1650  << pGlobalpoint << G4endl;
1651  }
1652 #endif
1653  G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
1654  G4LogicalVolume *motherLogical = motherPhysical->GetLogicalVolume();
1655  G4SmartVoxelHeader* pVoxelHeader = motherLogical->GetVoxelHeader();
1656  G4ThreeVector localPoint = ComputeLocalPoint(pGlobalpoint);
1657 
1659  {
1660  switch(CharacteriseDaughters(motherLogical))
1661  {
1662  case kNormal:
1663  if ( pVoxelHeader )
1664  {
1665 #ifdef G4NEW_SAFETY
1666  static G4VoxelSafety sfVoxelSafety;
1667  G4double safetyTwo = sfVoxelSafety.ComputeSafety(localPoint,
1668  *motherPhysical, pMaxLength);
1669 #ifdef G4DEBUG_NAVIGATION
1670  G4double safetyNormal=
1671  fnormalNav.ComputeSafety(localPoint, fHistory, pMaxLength);
1672 
1673  G4double diffSafety= (safetyTwo - safetyNormal);
1674  if( std::fabs(diffSafety) > CLHEP::perMillion
1675  * std::fabs(safetyNormal) )
1676  {
1678  exd << " ERROR in G4Navigator::ComputeStep " << G4endl;
1679  exd << " Error> DIFFERENCE in Safety from G4VoxelSafety "
1680  << " - compared to value from Normal. Diff = " << diffSafety << G4endl;
1681  exd << " New Voxel Safety= " << safetyTwo
1682  << " Value from Normal= " << safetyNormal << G4endl;
1683  exd << " Location: Local coordinates = " << localPoint
1684  << " mother Volume= " << *motherPhysical->GetName()
1685  << " Copy# = " << motherPhysical->GetCopyNo() << G4endl;
1686  exd << " : Global coordinates = " << pGlobalpoint
1687  << G4endl;
1688  G4Exception("G4Navigator::ComputeSafety()",
1689  "GeomNav0003", JustWarning, exd);
1690  }
1691  G4double safetyOldVoxel;
1692  safetyOldVoxel =
1693  fvoxelNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1694 
1695  G4double diffBetterToApprox= safetyTwo - safetyOldVoxel;
1696  if( diffBetterToApprox < 0.0 )
1697  {
1699  exd << "SMALLER Safety from call G4VoxelSafety compared to value "
1700  << " from Normal. Diff = " << diffSafety << G4endl;
1701  exd << " New Voxel Safety= " << safetyTwo
1702  << " Old value from Normal= " << newSafety << G4endl;
1703  exd << " Location: Local coordinates = " << localPoint
1704  << " mother Volume= " << *motherPhysical->GetName()
1705  << " Copy# = " << motherPhysical->GetCopyNo() << G4endl;
1706  exd << " : Global coordinates = " << pGlobalpoint
1707  << G4endl;
1708  G4Exception("G4Navigator::ComputeSafety()",
1709  "GeomNav0003", JustWarning, exd);
1710  }
1711 #endif
1712 
1713  // newSafety= safetyOldVoxel;
1714  newSafety= safetyTwo; // Faster and best available
1715 #else
1716  G4double safetyOldVoxel;
1717  safetyOldVoxel =
1718  fvoxelNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1719  newSafety= safetyOldVoxel;
1720 #endif
1721  }
1722  else
1723  {
1724  newSafety=fnormalNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1725  }
1726  break;
1727  case kParameterised:
1728  if( GetDaughtersRegularStructureId(motherLogical) != 1 )
1729  {
1730  newSafety=fparamNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1731  }
1732  else // Regular structure
1733  {
1734  newSafety=fregularNav.ComputeSafety(localPoint,fHistory,pMaxLength);
1735  }
1736  break;
1737  case kReplica:
1738  G4Exception("G4Navigator::ComputeSafety()", "GeomNav0001",
1739  FatalException, "Not applicable for replicated volumes.");
1740  break;
1741  }
1742  }
1743  else
1744  {
1745  newSafety = freplicaNav.ComputeSafety(pGlobalpoint, localPoint,
1746  fHistory, pMaxLength);
1747  }
1748  }
1749  else // 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  }
1764 
1765  // Remember last safety origin & value
1766  //
1767  fPreviousSftOrigin = pGlobalpoint;
1768  fPreviousSafety = newSafety;
1769 
1770  if (keepState) { RestoreSavedState(); }
1771 
1772 #ifdef G4DEBUG_NAVIGATION
1773  if( fVerbose > 1 )
1774  {
1775  G4cout << " ---- Exiting ComputeSafety " << G4endl;
1776  if( fVerbose > 2 ) { PrintState(); }
1777  G4cout << " Returned value of Safety = " << newSafety << G4endl;
1778  }
1779  G4cout.precision(oldcoutPrec);
1780 #endif
1781 
1782  return newSafety;
1783 }
1784 
1785 // ********************************************************************
1786 // CreateTouchableHistoryHandle
1787 // ********************************************************************
1788 //
1790 {
1792 }
1793 
1794 // ********************************************************************
1795 // PrintState
1796 // ********************************************************************
1797 //
1799 {
1800  G4int oldcoutPrec = G4cout.precision(4);
1801  if( fVerbose == 4 )
1802  {
1803  G4cout << "The current state of G4Navigator is: " << G4endl;
1804  G4cout << " ValidExitNormal= " << fValidExitNormal << G4endl
1805  << " ExitNormal = " << fExitNormal << G4endl
1806  << " Exiting = " << fExiting << G4endl
1807  << " Entering = " << fEntering << G4endl
1808  << " BlockedPhysicalVolume= " ;
1809  if (fBlockedPhysicalVolume==0)
1810  G4cout << "None";
1811  else
1812  G4cout << fBlockedPhysicalVolume->GetName();
1813  G4cout << G4endl
1814  << " BlockedReplicaNo = " << fBlockedReplicaNo << G4endl
1815  << " LastStepWasZero = " << fLastStepWasZero << G4endl
1816  << G4endl;
1817  }
1818  if( ( 1 < fVerbose) && (fVerbose < 4) )
1819  {
1820  G4cout << G4endl; // Make sure to line up
1821  G4cout << std::setw(30) << " ExitNormal " << " "
1822  << std::setw( 5) << " Valid " << " "
1823  << std::setw( 9) << " Exiting " << " "
1824  << std::setw( 9) << " Entering" << " "
1825  << std::setw(15) << " Blocked:Volume " << " "
1826  << std::setw( 9) << " ReplicaNo" << " "
1827  << std::setw( 8) << " LastStepZero " << " "
1828  << G4endl;
1829  G4cout << "( " << std::setw(7) << fExitNormal.x()
1830  << ", " << std::setw(7) << fExitNormal.y()
1831  << ", " << std::setw(7) << fExitNormal.z() << " ) "
1832  << std::setw( 5) << fValidExitNormal << " "
1833  << std::setw( 9) << fExiting << " "
1834  << std::setw( 9) << fEntering << " ";
1835  if ( fBlockedPhysicalVolume==0 )
1836  G4cout << std::setw(15) << "None";
1837  else
1838  G4cout << std::setw(15)<< fBlockedPhysicalVolume->GetName();
1839  G4cout << std::setw( 9) << fBlockedReplicaNo << " "
1840  << std::setw( 8) << fLastStepWasZero << " "
1841  << G4endl;
1842  }
1843  if( fVerbose > 2 )
1844  {
1845  G4cout.precision(8);
1846  G4cout << " Current Localpoint = " << fLastLocatedPointLocal << G4endl;
1847  G4cout << " PreviousSftOrigin = " << fPreviousSftOrigin << G4endl;
1848  G4cout << " PreviousSafety = " << fPreviousSafety << G4endl;
1849  }
1850  G4cout.precision(oldcoutPrec);
1851 }
1852 
1853 // ********************************************************************
1854 // ComputeStepLog
1855 // ********************************************************************
1856 //
1857 void G4Navigator::ComputeStepLog(const G4ThreeVector& pGlobalpoint,
1858  G4double moveLenSq) const
1859 {
1860  // The following checks only make sense if the move is larger
1861  // than the tolerance.
1862 
1863  static const G4double fAccuracyForWarning = kCarTolerance,
1864  fAccuracyForException = 1000*kCarTolerance;
1865 
1866  G4ThreeVector OriginalGlobalpoint = fHistory.GetTopTransform().Inverse().
1867  TransformPoint(fLastLocatedPointLocal);
1868 
1869  G4double shiftOriginSafSq = (fPreviousSftOrigin-pGlobalpoint).mag2();
1870 
1871  // Check that the starting point of this step is
1872  // within the isotropic safety sphere of the last point
1873  // to a accuracy/precision given by fAccuracyForWarning.
1874  // If so give warning.
1875  // If it fails by more than fAccuracyForException exit with error.
1876  //
1877  if( shiftOriginSafSq >= sqr(fPreviousSafety) )
1878  {
1879  G4double shiftOrigin = std::sqrt(shiftOriginSafSq);
1880  G4double diffShiftSaf = shiftOrigin - fPreviousSafety;
1881 
1882  if( diffShiftSaf > fAccuracyForWarning )
1883  {
1884  G4int oldcoutPrec= G4cout.precision(8);
1885  G4int oldcerrPrec= G4cerr.precision(10);
1886  std::ostringstream message, suggestion;
1887  message << "Accuracy error or slightly inaccurate position shift."
1888  << G4endl
1889  << " The Step's starting point has moved "
1890  << std::sqrt(moveLenSq)/mm << " mm " << G4endl
1891  << " since the last call to a Locate method." << G4endl
1892  << " This has resulted in moving "
1893  << shiftOrigin/mm << " mm "
1894  << " from the last point at which the safety "
1895  << " was calculated " << G4endl
1896  << " which is more than the computed safety= "
1897  << fPreviousSafety/mm << " mm at that point." << G4endl
1898  << " This difference is "
1899  << diffShiftSaf/mm << " mm." << G4endl
1900  << " The tolerated accuracy is "
1901  << fAccuracyForException/mm << " mm.";
1902 
1903  suggestion << " ";
1904  static G4int warnNow = 0;
1905  if( ((++warnNow % 100) == 1) )
1906  {
1907  message << G4endl
1908  << " This problem can be due to either " << G4endl
1909  << " - a process that has proposed a displacement"
1910  << " larger than the current safety , or" << G4endl
1911  << " - inaccuracy in the computation of the safety";
1912  suggestion << "We suggest that you " << G4endl
1913  << " - find i) what particle is being tracked, and "
1914  << " ii) through what part of your geometry " << G4endl
1915  << " for example by re-running this event with "
1916  << G4endl
1917  << " /tracking/verbose 1 " << G4endl
1918  << " - check which processes you declare for"
1919  << " this particle (and look at non-standard ones)"
1920  << G4endl
1921  << " - in case, create a detailed logfile"
1922  << " of this event using:" << G4endl
1923  << " /tracking/verbose 6 ";
1924  }
1925  G4Exception("G4Navigator::ComputeStep()",
1926  "GeomNav1002", JustWarning,
1927  message, G4String(suggestion.str()));
1928  G4cout.precision(oldcoutPrec);
1929  G4cerr.precision(oldcerrPrec);
1930  }
1931 #ifdef G4DEBUG_NAVIGATION
1932  else
1933  {
1934  G4cerr << "WARNING - G4Navigator::ComputeStep()" << G4endl
1935  << " The Step's starting point has moved "
1936  << std::sqrt(moveLenSq) << "," << G4endl
1937  << " which has taken it to the limit of"
1938  << " the current safety. " << G4endl;
1939  }
1940 #endif
1941  }
1942  G4double safetyPlus = fPreviousSafety + fAccuracyForException;
1943  if ( shiftOriginSafSq > sqr(safetyPlus) )
1944  {
1945  std::ostringstream message;
1946  message << "May lead to a crash or unreliable results." << G4endl
1947  << " Position has shifted considerably without"
1948  << " notifying the navigator !" << G4endl
1949  << " Tolerated safety: " << safetyPlus << G4endl
1950  << " Computed shift : " << shiftOriginSafSq;
1951  G4Exception("G4Navigator::ComputeStep()", "GeomNav1002",
1952  JustWarning, message);
1953  }
1954 }
1955 
1956 // ********************************************************************
1957 // Operator <<
1958 // ********************************************************************
1959 //
1960 std::ostream& operator << (std::ostream &os,const G4Navigator &n)
1961 {
1962  // Old version did only the following:
1963  // os << "Current History: " << G4endl << n.fHistory;
1964  // Old behaviour is recovered for fVerbose = 0
1965 
1966  // Adapted from G4Navigator::PrintState() const
1967 
1968  G4int oldcoutPrec = os.precision(4);
1969  if( n.fVerbose >= 4 )
1970  {
1971  os << "The current state of G4Navigator is: " << G4endl;
1972  os << " ValidExitNormal= " << n.fValidExitNormal << G4endl
1973  << " ExitNormal = " << n.fExitNormal << G4endl
1974  << " Exiting = " << n.fExiting << G4endl
1975  << " Entering = " << n.fEntering << G4endl
1976  << " BlockedPhysicalVolume= " ;
1977  if (n.fBlockedPhysicalVolume==0)
1978  os << "None";
1979  else
1980  os << n.fBlockedPhysicalVolume->GetName();
1981  os << G4endl
1982  << " BlockedReplicaNo = " << n.fBlockedReplicaNo << G4endl
1983  << " LastStepWasZero = " << n.fLastStepWasZero << G4endl
1984  << G4endl;
1985  }
1986  if( ( 1 < n.fVerbose) && (n.fVerbose < 4) )
1987  {
1988  os << G4endl; // Make sure to line up
1989  os << std::setw(30) << " ExitNormal " << " "
1990  << std::setw( 5) << " Valid " << " "
1991  << std::setw( 9) << " Exiting " << " "
1992  << std::setw( 9) << " Entering" << " "
1993  << std::setw(15) << " Blocked:Volume " << " "
1994  << std::setw( 9) << " ReplicaNo" << " "
1995  << std::setw( 8) << " LastStepZero " << " "
1996  << G4endl;
1997  os << "( " << std::setw(7) << n.fExitNormal.x()
1998  << ", " << std::setw(7) << n.fExitNormal.y()
1999  << ", " << std::setw(7) << n.fExitNormal.z() << " ) "
2000  << std::setw( 5) << n.fValidExitNormal << " "
2001  << std::setw( 9) << n.fExiting << " "
2002  << std::setw( 9) << n.fEntering << " ";
2003  if ( n.fBlockedPhysicalVolume==0 )
2004  { os << std::setw(15) << "None"; }
2005  else
2006  { os << std::setw(15)<< n.fBlockedPhysicalVolume->GetName(); }
2007  os << std::setw( 9) << n.fBlockedReplicaNo << " "
2008  << std::setw( 8) << n.fLastStepWasZero << " "
2009  << G4endl;
2010  }
2011  if( n.fVerbose > 2 )
2012  {
2013  os.precision(8);
2014  os << " Current Localpoint = " << n.fLastLocatedPointLocal << G4endl;
2015  os << " PreviousSftOrigin = " << n.fPreviousSftOrigin << G4endl;
2016  os << " PreviousSafety = " << n.fPreviousSafety << G4endl;
2017  }
2018  if( n.fVerbose > 3 || n.fVerbose == 0 )
2019  {
2020  os << "Current History: " << G4endl << n.fHistory;
2021  }
2022 
2023  os.precision(oldcoutPrec);
2024  return os;
2025 }