G4Navigator.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id: G4Navigator.cc 2012-11-01 japost Exp $
00028 // GEANT4 tag $ Name:  $
00029 // 
00030 // class G4Navigator Implementation
00031 //
00032 // Original author: Paul Kent, July 95/96
00033 // Responsible 1996-present: John Apostolakis
00034 // Co-maintainer:            Gabriele Cosmo
00035 // Additional revisions by: Pedro Arce, Vladimir Grichine
00036 // --------------------------------------------------------------------
00037 
00038 #include <iomanip>
00039 
00040 #include "G4Navigator.hh"
00041 #include "G4ios.hh"
00042 #include "G4SystemOfUnits.hh"
00043 #include "G4GeometryTolerance.hh"
00044 #include "G4VPhysicalVolume.hh"
00045 
00046 #include "G4VoxelSafety.hh"
00047 
00048 // ********************************************************************
00049 // Constructor
00050 // ********************************************************************
00051 //
00052 G4Navigator::G4Navigator()
00053   : fWasLimitedByGeometry(false), fVerbose(0),
00054     fTopPhysical(0), fCheck(false), fPushed(false), fWarnPush(true)
00055 {
00056   fActive= false; 
00057   fLastTriedStepComputation= false;
00058 
00059   ResetStackAndState();
00060     // Initialises also all 
00061     // - exit / entry flags
00062     // - flags & variables for exit normals
00063     // - zero step counters
00064     // - blocked volume 
00065 
00066   fActionThreshold_NoZeroSteps  = 10; 
00067   fAbandonThreshold_NoZeroSteps = 25; 
00068 
00069   kCarTolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
00070   fregularNav.SetNormalNavigation( &fnormalNav );
00071 
00072   fStepEndPoint = G4ThreeVector( kInfinity, kInfinity, kInfinity ); 
00073   fLastStepEndPointLocal = G4ThreeVector( kInfinity, kInfinity, kInfinity ); 
00074 }
00075 
00076 // ********************************************************************
00077 // Destructor
00078 // ********************************************************************
00079 //
00080 G4Navigator::~G4Navigator()
00081 {;}
00082 
00083 // ********************************************************************
00084 // ResetHierarchyAndLocate
00085 // ********************************************************************
00086 //
00087 G4VPhysicalVolume*
00088 G4Navigator::ResetHierarchyAndLocate(const G4ThreeVector &p,
00089                                      const G4ThreeVector &direction,
00090                                      const G4TouchableHistory &h)
00091 {
00092   ResetState();
00093   fHistory = *h.GetHistory();
00094   SetupHierarchy();
00095   fLastTriedStepComputation= false;  // Redundant, but best
00096   return LocateGlobalPointAndSetup(p, &direction, true, false);
00097 }
00098 
00099 // ********************************************************************
00100 // LocateGlobalPointAndSetup
00101 //
00102 // Locate the point in the hierarchy return 0 if outside
00103 // The direction is required 
00104 //    - if on an edge shared by more than two surfaces 
00105 //      (to resolve likely looping in tracking)
00106 //    - at initial location of a particle
00107 //      (to resolve potential ambiguity at boundary)
00108 // 
00109 // Flags on exit: (comments to be completed)
00110 // fEntering         - True if entering `daughter' volume (or replica)
00111 //                     whether daughter of last mother directly 
00112 //                     or daughter of that volume's ancestor.
00113 // ********************************************************************
00114 //
00115 G4VPhysicalVolume* 
00116 G4Navigator::LocateGlobalPointAndSetup( const G4ThreeVector& globalPoint,
00117                                         const G4ThreeVector* pGlobalDirection,
00118                                         const G4bool relativeSearch,
00119                                         const G4bool ignoreDirection )
00120 {
00121   G4bool notKnownContained=true, noResult;
00122   G4VPhysicalVolume *targetPhysical;
00123   G4LogicalVolume *targetLogical;
00124   G4VSolid *targetSolid=0;
00125   G4ThreeVector localPoint, globalDirection;
00126   EInside insideCode;
00127   
00128   G4bool considerDirection = (!ignoreDirection) || fLocatedOnEdge;
00129 
00130   fLastTriedStepComputation=   false;   
00131   fChangedGrandMotherRefFrame= false;  // For local exit normal
00132    
00133   if( considerDirection && pGlobalDirection != 0 )
00134   {
00135     globalDirection=*pGlobalDirection;
00136   }
00137 
00138 
00139 #ifdef G4VERBOSE
00140   if( fVerbose > 2 )
00141   {
00142     G4int oldcoutPrec = G4cout.precision(8);
00143     G4cout << "*** G4Navigator::LocateGlobalPointAndSetup: ***" << G4endl; 
00144     G4cout << "    Called with arguments: " << G4endl
00145            << "    Globalpoint = " << globalPoint << G4endl
00146            << "    RelativeSearch = " << relativeSearch  << G4endl;
00147     if( fVerbose == 4 )
00148     {
00149       G4cout << "    ----- Upon entering:" << G4endl;
00150       PrintState();
00151     }
00152     G4cout.precision(oldcoutPrec);
00153   }
00154 #endif
00155 
00156   if ( !relativeSearch )
00157   {
00158     ResetStackAndState();
00159   }
00160   else
00161   {
00162     if ( fWasLimitedByGeometry )
00163     {
00164       fWasLimitedByGeometry = false;
00165       fEnteredDaughter = fEntering;   // Remember
00166       fExitedMother = fExiting;       // Remember
00167       if ( fExiting )
00168       {
00169         if ( fHistory.GetDepth() )
00170         {
00171           fBlockedPhysicalVolume = fHistory.GetTopVolume();
00172           fBlockedReplicaNo = fHistory.GetTopReplicaNo();
00173           fHistory.BackLevel();
00174         }
00175         else
00176         {
00177           fLastLocatedPointLocal = localPoint;
00178           fLocatedOutsideWorld = true;
00179           return 0;           // Have exited world volume
00180         }
00181         // A fix for the case where a volume is "entered" at an edge
00182         // and a coincident surface exists outside it.
00183         //  - This stops it from exiting further volumes and cycling
00184         //  - However ReplicaNavigator treats this case itself
00185         //
00186         if ( fLocatedOnEdge && (VolumeType(fBlockedPhysicalVolume)!=kReplica ))
00187         { 
00188           fExiting= false;
00189         }
00190       }
00191       else
00192         if ( fEntering )
00193         {
00194           switch (VolumeType(fBlockedPhysicalVolume))
00195           {
00196             case kNormal:
00197               fHistory.NewLevel(fBlockedPhysicalVolume, kNormal,
00198                                 fBlockedPhysicalVolume->GetCopyNo());
00199               break;
00200             case kReplica:
00201               freplicaNav.ComputeTransformation(fBlockedReplicaNo,
00202                                                 fBlockedPhysicalVolume);
00203               fHistory.NewLevel(fBlockedPhysicalVolume, kReplica,
00204                                 fBlockedReplicaNo);
00205               fBlockedPhysicalVolume->SetCopyNo(fBlockedReplicaNo);
00206               break;
00207             case kParameterised:
00208               if( fBlockedPhysicalVolume->GetRegularStructureId() == 0 )
00209               {
00210                 G4VSolid *pSolid;
00211                 G4VPVParameterisation *pParam;
00212                 G4TouchableHistory parentTouchable( fHistory );
00213                 pParam = fBlockedPhysicalVolume->GetParameterisation();
00214                 pSolid = pParam->ComputeSolid(fBlockedReplicaNo,
00215                                               fBlockedPhysicalVolume);
00216                 pSolid->ComputeDimensions(pParam, fBlockedReplicaNo,
00217                                           fBlockedPhysicalVolume);
00218                 pParam->ComputeTransformation(fBlockedReplicaNo,
00219                                               fBlockedPhysicalVolume);
00220                 fHistory.NewLevel(fBlockedPhysicalVolume, kParameterised,
00221                                   fBlockedReplicaNo);
00222                 fBlockedPhysicalVolume->SetCopyNo(fBlockedReplicaNo);
00223                 //
00224                 // Set the correct solid and material in Logical Volume
00225                 //
00226                 G4LogicalVolume *pLogical;
00227                 pLogical = fBlockedPhysicalVolume->GetLogicalVolume();
00228                 pLogical->SetSolid( pSolid );
00229                 pLogical->UpdateMaterial(pParam ->
00230                   ComputeMaterial(fBlockedReplicaNo,
00231                                   fBlockedPhysicalVolume, 
00232                                   &parentTouchable));
00233               }
00234               break;
00235           }
00236           fEntering = false;
00237           fBlockedPhysicalVolume = 0;
00238           localPoint = fHistory.GetTopTransform().TransformPoint(globalPoint);
00239           notKnownContained = false;
00240         }
00241     }
00242     else
00243     {
00244       fBlockedPhysicalVolume = 0;
00245       fEntering = false;
00246       fEnteredDaughter = false;  // Full Step was not taken, did not enter
00247       fExiting = false;
00248       fExitedMother = false;     // Full Step was not taken, did not exit
00249     }
00250   }
00251   //
00252   // Search from top of history up through geometry until
00253   // containing volume found:
00254   // If on 
00255   // o OUTSIDE - Back up level, not/no longer exiting volumes
00256   // o SURFACE and EXITING - Back up level, setting new blocking no.s
00257   // else
00258   // o containing volume found
00259   //
00260   G4int noLevelsExited=0 ;
00261 
00262   while (notKnownContained)
00263   {
00264     if ( fHistory.GetTopVolumeType()!=kReplica )
00265     {
00266       targetSolid = fHistory.GetTopVolume()->GetLogicalVolume()->GetSolid();
00267       localPoint = fHistory.GetTopTransform().TransformPoint(globalPoint);
00268       insideCode = targetSolid->Inside(localPoint);
00269 #ifdef G4VERBOSE
00270       if(( fVerbose == 1 ) && ( fCheck ))
00271       {
00272          G4String solidResponse = "-kInside-";
00273          if (insideCode == kOutside)
00274            solidResponse = "-kOutside-";
00275          else if (insideCode == kSurface)
00276            solidResponse = "-kSurface-";
00277          G4cout << "*** G4Navigator::LocateGlobalPointAndSetup(): ***" << G4endl
00278                 << "    Invoked Inside() for solid: " << targetSolid->GetName()
00279                 << ". Solid replied: " << solidResponse << G4endl
00280                 << "    For local point p: " << localPoint << G4endl;
00281       }
00282 #endif
00283     }
00284     else
00285     {
00286       insideCode = freplicaNav.BackLocate(fHistory, globalPoint, localPoint,
00287                                           fExiting, notKnownContained);
00288       // !CARE! if notKnownContained returns false then the point is within
00289       // the containing placement volume of the replica(s). If insidecode
00290       // will result in the history being backed up one level, then the
00291       // local point returned is the point in the system of this new level
00292     }
00293 
00294 
00295     if ( insideCode==kOutside )
00296     {
00297       noLevelsExited++; 
00298       if ( fHistory.GetDepth() )
00299       {
00300         fBlockedPhysicalVolume = fHistory.GetTopVolume();
00301         fBlockedReplicaNo = fHistory.GetTopReplicaNo();
00302         fHistory.BackLevel();
00303         fExiting = false;
00304 
00305         if( noLevelsExited > 1 )
00306         {
00307           // The first transformation was done by the sub-navigator
00308           //
00309           const G4RotationMatrix* mRot = fBlockedPhysicalVolume->GetRotation();
00310           if( mRot )
00311           { 
00312             fGrandMotherExitNormal *= (*mRot).inverse();
00313             fChangedGrandMotherRefFrame= true;
00314           }
00315         }
00316       }
00317       else
00318       {
00319         fLastLocatedPointLocal = localPoint;
00320         fLocatedOutsideWorld = true;
00321           // No extra transformation for ExitNormal - is in frame of Top Volume
00322         return 0;         // Have exited world volume
00323       }
00324     }
00325     else
00326       if ( insideCode==kSurface )
00327       {
00328         G4bool isExiting = fExiting;
00329         if( (!fExiting)&&considerDirection )
00330         {
00331           // Figure out whether we are exiting this level's volume
00332           // by using the direction
00333           //
00334           G4bool directionExiting = false;
00335           G4ThreeVector localDirection =
00336               fHistory.GetTopTransform().TransformAxis(globalDirection);
00337 
00338           // Make sure localPoint in correct reference frame
00339           //     ( Was it already correct ? How ? )
00340           //
00341           localPoint= fHistory.GetTopTransform().TransformPoint(globalPoint);
00342           if ( fHistory.GetTopVolumeType()!=kReplica )
00343           {
00344             G4ThreeVector normal = targetSolid->SurfaceNormal(localPoint);
00345             directionExiting = normal.dot(localDirection) > 0.0;
00346             isExiting = isExiting || directionExiting;
00347           }
00348         }
00349         if( isExiting )
00350         {
00351           noLevelsExited++; 
00352           if ( fHistory.GetDepth() )
00353           {
00354             fBlockedPhysicalVolume = fHistory.GetTopVolume();
00355             fBlockedReplicaNo = fHistory.GetTopReplicaNo();
00356             fHistory.BackLevel();
00357             //
00358             // Still on surface but exited volume not necessarily convex
00359             //
00360             fValidExitNormal = false;
00361 
00362             if( noLevelsExited > 1 )
00363             {
00364               // The first transformation was done by the sub-navigator
00365               //
00366               const G4RotationMatrix* mRot=fBlockedPhysicalVolume->GetRotation();
00367               if( mRot )
00368               { 
00369                 fGrandMotherExitNormal *= (*mRot).inverse();
00370                 fChangedGrandMotherRefFrame= true;
00371               }
00372             }
00373           } 
00374           else
00375           {
00376             fLastLocatedPointLocal = localPoint;
00377             fLocatedOutsideWorld = true;
00378               // No extra transformation for ExitNormal, is in frame of Top Vol
00379             return 0;          // Have exited world volume
00380           }
00381         }
00382         else
00383         {
00384           notKnownContained=false;
00385         }
00386       }
00387       else
00388       {
00389         notKnownContained=false;
00390       }
00391   }  // END while (notKnownContained)
00392   //
00393   // Search downwards until deepest containing volume found,
00394   // blocking fBlockedPhysicalVolume/BlockedReplicaNum
00395   //
00396   // 3 Cases:
00397   //
00398   // o Parameterised daughters
00399   //   =>Must be one G4PVParameterised daughter & voxels
00400   // o Positioned daughters & voxels
00401   // o Positioned daughters & no voxels
00402 
00403   noResult = true;  // noResult should be renamed to 
00404                     // something like enteredLevel, as that is its meaning.
00405   do
00406   {
00407     // Determine `type' of current mother volume
00408     //
00409     targetPhysical = fHistory.GetTopVolume();
00410     if (!targetPhysical) { break; }
00411     targetLogical = targetPhysical->GetLogicalVolume();
00412     switch( CharacteriseDaughters(targetLogical) )
00413     {
00414       case kNormal:
00415         if ( targetLogical->GetVoxelHeader() )  // use optimised navigation
00416         {
00417           noResult = fvoxelNav.LevelLocate(fHistory,
00418                                            fBlockedPhysicalVolume,
00419                                            fBlockedReplicaNo,
00420                                            globalPoint,
00421                                            pGlobalDirection,
00422                                            considerDirection,
00423                                            localPoint);
00424         }
00425         else                       // do not use optimised navigation
00426         {
00427           noResult = fnormalNav.LevelLocate(fHistory,
00428                                             fBlockedPhysicalVolume,
00429                                             fBlockedReplicaNo,
00430                                             globalPoint,
00431                                             pGlobalDirection,
00432                                             considerDirection,
00433                                             localPoint);
00434         }
00435         break;
00436       case kReplica:
00437         noResult = freplicaNav.LevelLocate(fHistory,
00438                                            fBlockedPhysicalVolume,
00439                                            fBlockedReplicaNo,
00440                                            globalPoint,
00441                                            pGlobalDirection,
00442                                            considerDirection,
00443                                            localPoint);
00444         break;
00445       case kParameterised:
00446         if( GetDaughtersRegularStructureId(targetLogical) != 1 )
00447         {
00448           noResult = fparamNav.LevelLocate(fHistory,
00449                                            fBlockedPhysicalVolume,
00450                                            fBlockedReplicaNo,
00451                                            globalPoint,
00452                                            pGlobalDirection,
00453                                            considerDirection,
00454                                            localPoint);
00455         }
00456         else  // Regular structure
00457         {
00458           noResult = fregularNav.LevelLocate(fHistory,
00459                                              fBlockedPhysicalVolume,
00460                                              fBlockedReplicaNo,
00461                                              globalPoint,
00462                                              pGlobalDirection,
00463                                              considerDirection,
00464                                              localPoint);
00465         }
00466         break;
00467     }
00468 
00469     // LevelLocate returns true if it finds a daughter volume 
00470     // in which globalPoint is inside (or on the surface).
00471 
00472     if ( noResult )
00473     {
00474       // Entering a daughter after ascending
00475       //
00476       // The blocked volume is no longer valid - it was for another level
00477       //
00478       fBlockedPhysicalVolume = 0;
00479       fBlockedReplicaNo = -1;
00480 
00481       // fEntering should be false -- else blockedVolume is assumed good.
00482       // fEnteredDaughter is used for ExitNormal
00483       //
00484       fEntering = false;
00485       fEnteredDaughter = true;
00486 
00487       if( fExitedMother )
00488       {
00489         G4VPhysicalVolume* enteredPhysical = fHistory.GetTopVolume();
00490         const G4RotationMatrix* mRot = enteredPhysical->GetRotation();
00491         if( mRot )
00492         { 
00493           fGrandMotherExitNormal *= (*mRot).inverse();
00494         }
00495       }
00496 
00497 #ifdef G4DEBUG_NAVIGATION
00498       if( fVerbose > 2 )
00499       { 
00500          G4VPhysicalVolume* enteredPhysical = fHistory.GetTopVolume();
00501          G4cout << "*** G4Navigator::LocateGlobalPointAndSetup() ***" << G4endl;
00502          G4cout << "    Entering volume: " << enteredPhysical->GetName()
00503                 << G4endl;
00504       }
00505 #endif
00506     }
00507   } while (noResult);
00508 
00509   fLastLocatedPointLocal = localPoint;
00510 
00511 #ifdef G4VERBOSE
00512   if( fVerbose >= 4 )
00513   {
00514     G4int oldcoutPrec = G4cout.precision(8);
00515     G4String curPhysVol_Name("None");
00516     if (targetPhysical)  { curPhysVol_Name = targetPhysical->GetName(); }
00517     G4cout << "    Return value = new volume = " << curPhysVol_Name << G4endl;
00518     G4cout << "    ----- Upon exiting:" << G4endl;
00519     PrintState();
00520     if( fVerbose == 5 )
00521     {
00522       G4cout << "Upon exiting LocateGlobalPointAndSetup():" << G4endl;
00523       G4cout << "    History = " << G4endl << fHistory << G4endl << G4endl;
00524     }
00525     G4cout.precision(oldcoutPrec);
00526   }
00527 #endif
00528 
00529   fLocatedOutsideWorld= false;
00530 
00531   return targetPhysical;
00532 }
00533 
00534 // ********************************************************************
00535 // LocateGlobalPointWithinVolume
00536 //
00537 // -> the state information of this Navigator and its subNavigators
00538 //    is updated in order to start the next step at pGlobalpoint
00539 // -> no check is performed whether pGlobalpoint is inside the 
00540 //    original volume (this must be the case).
00541 //
00542 // Note: a direction could be added to the arguments, to aid in future
00543 //       optional checking (via the old code below, flagged by OLD_LOCATE). 
00544 //       [ This would be done only in verbose mode ]
00545 // ********************************************************************
00546 //
00547 void
00548 G4Navigator::LocateGlobalPointWithinVolume(const G4ThreeVector& pGlobalpoint)
00549 {  
00550    fLastLocatedPointLocal = ComputeLocalPoint(pGlobalpoint);
00551    fLastTriedStepComputation= false;
00552    fChangedGrandMotherRefFrame= false;  //  Frame for Exit Normal
00553 
00554 #ifdef G4DEBUG_NAVIGATION
00555    if( fVerbose > 2 )
00556    { 
00557      G4cout << "Entering LocateGlobalWithinVolume(): History = " << G4endl;
00558      G4cout << fHistory << G4endl;
00559    }
00560 #endif
00561 
00562    // For the case of Voxel (or Parameterised) volume the respective 
00563    // Navigator must be messaged to update its voxel information etc
00564 
00565    // Update the state of the Sub Navigators 
00566    // - in particular any voxel information they store/cache
00567    //
00568    G4VPhysicalVolume*  motherPhysical = fHistory.GetTopVolume();
00569    G4LogicalVolume*    motherLogical  = motherPhysical->GetLogicalVolume();
00570    G4SmartVoxelHeader* pVoxelHeader   = motherLogical->GetVoxelHeader();
00571 
00572    if ( fHistory.GetTopVolumeType()!=kReplica )
00573    {
00574      switch( CharacteriseDaughters(motherLogical) )
00575      {
00576        case kNormal:
00577          if ( pVoxelHeader )
00578          {
00579            fvoxelNav.VoxelLocate( pVoxelHeader, fLastLocatedPointLocal );
00580          }
00581          break;
00582        case kParameterised:
00583          if( GetDaughtersRegularStructureId(motherLogical) != 1 )
00584          {
00585            // Resets state & returns voxel node
00586            //
00587            fparamNav.ParamVoxelLocate( pVoxelHeader, fLastLocatedPointLocal );
00588          }
00589          break;
00590        case kReplica:
00591          G4Exception("G4Navigator::LocateGlobalPointWithinVolume()",
00592                      "GeomNav0001", FatalException,
00593                      "Not applicable for replicated volumes.");
00594          break;
00595      }
00596    }
00597 
00598    // Reset the state variables 
00599    //   - which would have been affected
00600    //     by the 'equivalent' call to LocateGlobalPointAndSetup
00601    //   - who's values have been invalidated by the 'move'.
00602    //
00603    fBlockedPhysicalVolume = 0; 
00604    fBlockedReplicaNo = -1;
00605    fEntering = false;
00606    fEnteredDaughter = false;  // Boundary not encountered, did not enter
00607    fExiting = false;
00608    fExitedMother = false;     // Boundary not encountered, did not exit
00609 }
00610 
00611 // ********************************************************************
00612 // SetSavedState
00613 //
00614 // Save the state, in case this is a parasitic call
00615 // Save fValidExitNormal, fExitNormal, fExiting, fEntering, 
00616 //      fBlockedPhysicalVolume, fBlockedReplicaNo, fLastStepWasZero; 
00617 // ********************************************************************
00618 //
00619 void G4Navigator::SetSavedState()
00620 {
00621   fSaveState.sExitNormal = fExitNormal;
00622   fSaveState.sValidExitNormal = fValidExitNormal;
00623   fSaveState.sExiting = fExiting;
00624   fSaveState.sEntering = fEntering;
00625 
00626   fSaveState.spBlockedPhysicalVolume = fBlockedPhysicalVolume;
00627   fSaveState.sBlockedReplicaNo = fBlockedReplicaNo, 
00628 
00629   fSaveState.sLastStepWasZero = fLastStepWasZero;
00630   
00631   fSaveState.sLocatedOutsideWorld = fLocatedOutsideWorld;
00632   fSaveState.sLastLocatedPointLocal= fLastLocatedPointLocal;
00633   fSaveState.sEnteredDaughter= fEnteredDaughter;
00634   fSaveState.sExitedMother= fExitedMother;
00635 
00636   // Even the safety sphere - if you want to change it do it explicitly!
00637   //
00638   fSaveState.sPreviousSftOrigin= fPreviousSftOrigin;
00639   fSaveState.sPreviousSafety= fPreviousSafety;
00640 }
00641 
00642 // ********************************************************************
00643 // RestoreSavedState
00644 //
00645 // Restore the state (in Compute Step), in case this is a parasitic call
00646 // ********************************************************************
00647 //
00648 void G4Navigator::RestoreSavedState()
00649 {
00650   fExitNormal = fSaveState.sExitNormal;
00651   fValidExitNormal = fSaveState.sValidExitNormal;
00652   fExiting = fSaveState.sExiting;
00653   fEntering = fSaveState.sEntering;
00654 
00655   fBlockedPhysicalVolume = fSaveState.spBlockedPhysicalVolume;
00656   fBlockedReplicaNo = fSaveState.sBlockedReplicaNo, 
00657 
00658   fLastStepWasZero = fSaveState.sLastStepWasZero;
00659   
00660   fLocatedOutsideWorld = fSaveState.sLocatedOutsideWorld;
00661   fLastLocatedPointLocal= fSaveState.sLastLocatedPointLocal;
00662   fEnteredDaughter= fSaveState.sEnteredDaughter;
00663   fExitedMother= fSaveState.sExitedMother;
00664   fSaveState.sPreviousSftOrigin= fPreviousSftOrigin;
00665   fSaveState.sPreviousSafety= fPreviousSafety;
00666 }
00667 
00668 // ********************************************************************
00669 // ComputeStep
00670 //
00671 // Computes the next geometric Step: intersections with current
00672 // mother and `daughter' volumes.
00673 //
00674 // NOTE:
00675 //
00676 // Flags on entry:
00677 // --------------
00678 // fValidExitNormal  - Normal of exited volume is valid (convex, not a 
00679 //                     coincident boundary)
00680 // fExitNormal       - Surface normal of exited volume
00681 // fExiting          - True if have exited solid
00682 //
00683 // fBlockedPhysicalVolume - Ptr to exited volume (or 0)
00684 // fBlockedReplicaNo - Replication no of exited volume
00685 // fLastStepWasZero  - True if last Step size was zero.
00686 //
00687 // Flags on exit:
00688 // -------------
00689 // fValidExitNormal  - True if surface normal of exited volume is valid
00690 // fExitNormal       - Surface normal of exited volume rotated to mothers
00691 //                    reference system
00692 // fExiting          - True if exiting mother
00693 // fEntering         - True if entering `daughter' volume (or replica)
00694 // fBlockedPhysicalVolume - Ptr to candidate (entered) volume
00695 // fBlockedReplicaNo - Replication no of candidate (entered) volume
00696 // fLastStepWasZero  - True if this Step size was zero.
00697 // ********************************************************************
00698 //
00699 G4double G4Navigator::ComputeStep( const G4ThreeVector &pGlobalpoint,
00700                                    const G4ThreeVector &pDirection,
00701                                    const G4double pCurrentProposedStepLength,
00702                                          G4double &pNewSafety)
00703 {
00704   G4ThreeVector localDirection = ComputeLocalAxis(pDirection);
00705   G4double Step = kInfinity;
00706   G4VPhysicalVolume  *motherPhysical = fHistory.GetTopVolume();
00707   G4LogicalVolume *motherLogical = motherPhysical->GetLogicalVolume();
00708 
00709   // All state relating to exiting normals must be reset
00710   //
00711   fExitNormalGlobalFrame= G4ThreeVector( 0., 0., 0.);
00712     // Reset value - to erase its memory
00713   fChangedGrandMotherRefFrame= false;
00714     // Reset - used for local exit normal
00715   fGrandMotherExitNormal= G4ThreeVector( 0., 0., 0.); 
00716   fCalculatedExitNormal  = false;
00717     // Reset for new step
00718 
00719   static G4int sNavCScalls=0;
00720   sNavCScalls++;
00721 
00722   fLastTriedStepComputation= true; 
00723 
00724 #ifdef G4VERBOSE
00725   if( fVerbose > 0 )
00726   {
00727     G4cout << "*** G4Navigator::ComputeStep: ***" << G4endl; 
00728     G4cout << "    Volume = " << motherPhysical->GetName() 
00729            << " - Proposed step length = " << pCurrentProposedStepLength
00730            << G4endl; 
00731 #ifdef G4DEBUG_NAVIGATION
00732     if( fVerbose >= 2 )
00733     {
00734       G4cout << "  Called with the arguments: " << G4endl
00735              << "  Globalpoint = " << std::setw(25) << pGlobalpoint << G4endl
00736              << "  Direction   = " << std::setw(25) << pDirection << G4endl;
00737       if( fVerbose >= 4 )
00738       {
00739         G4cout << "  ---- Upon entering : State" << G4endl;
00740         PrintState();
00741       }
00742     }
00743 #endif
00744   }
00745 #endif
00746 
00747   G4ThreeVector newLocalPoint = ComputeLocalPoint(pGlobalpoint);
00748   if( newLocalPoint != fLastLocatedPointLocal )
00749   {
00750     // Check whether the relocation is within safety
00751     //
00752     G4ThreeVector oldLocalPoint = fLastLocatedPointLocal;
00753     G4double moveLenSq = (newLocalPoint-oldLocalPoint).mag2();
00754 
00755     if ( moveLenSq >= kCarTolerance*kCarTolerance )
00756     {
00757 #ifdef G4VERBOSE
00758       ComputeStepLog(pGlobalpoint, moveLenSq);
00759 #endif
00760       // Relocate the point within the same volume
00761       //
00762       LocateGlobalPointWithinVolume( pGlobalpoint );
00763       fLastTriedStepComputation= true;     // Ensure that this is set again !!
00764     }
00765   }
00766   if ( fHistory.GetTopVolumeType()!=kReplica )
00767   {
00768     switch( CharacteriseDaughters(motherLogical) )
00769     {
00770       case kNormal:
00771         if ( motherLogical->GetVoxelHeader() )
00772         {
00773           Step = fvoxelNav.ComputeStep(fLastLocatedPointLocal,
00774                                        localDirection,
00775                                        pCurrentProposedStepLength,
00776                                        pNewSafety,
00777                                        fHistory,
00778                                        fValidExitNormal,
00779                                        fExitNormal,
00780                                        fExiting,
00781                                        fEntering,
00782                                        &fBlockedPhysicalVolume,
00783                                        fBlockedReplicaNo);
00784       
00785         }
00786         else
00787         {
00788           if( motherPhysical->GetRegularStructureId() == 0 )
00789           {
00790             Step = fnormalNav.ComputeStep(fLastLocatedPointLocal,
00791                                           localDirection,
00792                                           pCurrentProposedStepLength,
00793                                           pNewSafety,
00794                                           fHistory,
00795                                           fValidExitNormal,
00796                                           fExitNormal,
00797                                           fExiting,
00798                                           fEntering,
00799                                           &fBlockedPhysicalVolume,
00800                                           fBlockedReplicaNo);
00801           }
00802           else  // Regular (non-voxelised) structure
00803           {
00804             LocateGlobalPointAndSetup( pGlobalpoint, &pDirection, true, true );
00805             fLastTriedStepComputation= true; // Ensure that this is set again !!
00806             //
00807             // if physical process limits the step, the voxel will not be the
00808             // one given by ComputeStepSkippingEqualMaterials() and the local
00809             // point will be wrongly calculated.
00810 
00811             // There is a problem: when msc limits the step and the point is
00812             // assigned wrongly to phantom in previous step (while it is out
00813             // of the container volume). Then LocateGlobalPointAndSetup() has
00814             // reset the history topvolume to world.
00815             //
00816             if(fHistory.GetTopVolume()->GetRegularStructureId() == 0 )
00817             { 
00818               G4Exception("G4Navigator::ComputeStep()",
00819                           "GeomNav1001", JustWarning,
00820                 "Point is relocated in voxels, while it should be outside!");
00821               Step = fnormalNav.ComputeStep(fLastLocatedPointLocal,
00822                                             localDirection,
00823                                             pCurrentProposedStepLength,
00824                                             pNewSafety,
00825                                             fHistory,
00826                                             fValidExitNormal,
00827                                             fExitNormal,
00828                                             fExiting,
00829                                             fEntering,
00830                                             &fBlockedPhysicalVolume,
00831                                             fBlockedReplicaNo);
00832             }
00833             else
00834             {
00835               Step = fregularNav.
00836                    ComputeStepSkippingEqualMaterials(fLastLocatedPointLocal,
00837                                                      localDirection,
00838                                                      pCurrentProposedStepLength,
00839                                                      pNewSafety,
00840                                                      fHistory,
00841                                                      fValidExitNormal,
00842                                                      fExitNormal,
00843                                                      fExiting,
00844                                                      fEntering,
00845                                                      &fBlockedPhysicalVolume,
00846                                                      fBlockedReplicaNo,
00847                                                      motherPhysical);
00848             }
00849           }
00850         }
00851         break;
00852       case kParameterised:
00853         if( GetDaughtersRegularStructureId(motherLogical) != 1 )
00854         {
00855           Step = fparamNav.ComputeStep(fLastLocatedPointLocal,
00856                                        localDirection,
00857                                        pCurrentProposedStepLength,
00858                                        pNewSafety,
00859                                        fHistory,
00860                                        fValidExitNormal,
00861                                        fExitNormal,
00862                                        fExiting,
00863                                        fEntering,
00864                                        &fBlockedPhysicalVolume,
00865                                        fBlockedReplicaNo);
00866         }
00867         else  // Regular structure
00868         {
00869           Step = fregularNav.ComputeStep(fLastLocatedPointLocal,
00870                                          localDirection,
00871                                          pCurrentProposedStepLength,
00872                                          pNewSafety,
00873                                          fHistory,
00874                                          fValidExitNormal,
00875                                          fExitNormal,
00876                                          fExiting,
00877                                          fEntering,
00878                                          &fBlockedPhysicalVolume,
00879                                          fBlockedReplicaNo);
00880         }
00881         break;
00882       case kReplica:
00883         G4Exception("G4Navigator::ComputeStep()", "GeomNav0001",
00884                     FatalException, "Not applicable for replicated volumes.");
00885         break;
00886     }
00887   }
00888   else
00889   {
00890     // In the case of a replica, it must handle the exiting
00891     // edge/corner problem by itself
00892     //
00893     G4bool exitingReplica = fExitedMother;
00894     G4bool calculatedExitNormal;
00895     Step = freplicaNav.ComputeStep(pGlobalpoint,
00896                                    pDirection,
00897                                    fLastLocatedPointLocal,
00898                                    localDirection,
00899                                    pCurrentProposedStepLength,
00900                                    pNewSafety,
00901                                    fHistory,
00902                                    fValidExitNormal,
00903                                    calculatedExitNormal,
00904                                    fExitNormal,
00905                                    exitingReplica,
00906                                    fEntering,
00907                                    &fBlockedPhysicalVolume,
00908                                    fBlockedReplicaNo);
00909     fExiting= exitingReplica;                          
00910     fCalculatedExitNormal= calculatedExitNormal;
00911   }
00912 
00913   // Remember last safety origin & value.
00914   //
00915   fPreviousSftOrigin = pGlobalpoint;
00916   fPreviousSafety = pNewSafety; 
00917 
00918   // Count zero steps - one can occur due to changing momentum at a boundary
00919   //                  - one, two (or a few) can occur at common edges between
00920   //                    volumes
00921   //                  - more than two is likely a problem in the geometry
00922   //                    description or the Navigation 
00923 
00924   // Rule of thumb: likely at an Edge if two consecutive steps are zero,
00925   //                because at least two candidate volumes must have been
00926   //                checked
00927   //
00928   fLocatedOnEdge   = fLastStepWasZero && (Step==0.0);
00929   fLastStepWasZero = (Step==0.0);
00930   if (fPushed)  { fPushed = fLastStepWasZero; }
00931 
00932   // Handle large number of consecutive zero steps
00933   //
00934   if ( fLastStepWasZero )
00935   {
00936     fNumberZeroSteps++;
00937 #ifdef G4DEBUG_NAVIGATION
00938     if( fNumberZeroSteps > 1 )
00939     {
00940        G4cout << "G4Navigator::ComputeStep(): another zero step, # "
00941               << fNumberZeroSteps
00942               << " at " << pGlobalpoint
00943               << " in volume " << motherPhysical->GetName()
00944               << " nav-comp-step calls # " << sNavCScalls
00945               << G4endl;
00946     }
00947 #endif
00948     if( fNumberZeroSteps > fActionThreshold_NoZeroSteps-1 )
00949     {
00950        // Act to recover this stuck track. Pushing it along direction
00951        //
00952        Step += 100*kCarTolerance;
00953 #ifdef G4VERBOSE
00954        if ((!fPushed) && (fWarnPush))
00955        {
00956          std::ostringstream message;
00957          message << "Track stuck or not moving." << G4endl
00958                  << "          Track stuck, not moving for " 
00959                  << fNumberZeroSteps << " steps" << G4endl
00960                  << "          in volume -" << motherPhysical->GetName()
00961                  << "- at point " << pGlobalpoint << G4endl
00962                  << "          direction: " << pDirection << "." << G4endl
00963                  << "          Potential geometry or navigation problem !"
00964                  << G4endl
00965                  << "          Trying pushing it of " << Step << " mm ...";
00966          G4Exception("G4Navigator::ComputeStep()", "GeomNav1002",
00967                      JustWarning, message, "Potential overlap in geometry!");
00968        }
00969 #endif
00970        fPushed = true;
00971     }
00972     if( fNumberZeroSteps > fAbandonThreshold_NoZeroSteps-1 )
00973     {
00974       // Must kill this stuck track
00975       //
00976       std::ostringstream message;
00977       message << "Stuck Track: potential geometry or navigation problem."
00978               << G4endl
00979               << "        Track stuck, not moving for " 
00980               << fNumberZeroSteps << " steps" << G4endl
00981               << "        in volume -" << motherPhysical->GetName()
00982               << "- at point " << pGlobalpoint << G4endl
00983               << "        direction: " << pDirection << ".";
00984       motherPhysical->CheckOverlaps(5000, false);
00985       G4Exception("G4Navigator::ComputeStep()", "GeomNav0003",
00986                   EventMustBeAborted, message);
00987     }
00988   }
00989   else
00990   {
00991     if (!fPushed)  fNumberZeroSteps = 0;
00992   }
00993 
00994   fEnteredDaughter = fEntering;   // I expect to enter a volume in this Step
00995   fExitedMother = fExiting;
00996 
00997   fStepEndPoint = pGlobalpoint + Step * pDirection; 
00998   fLastStepEndPointLocal = fLastLocatedPointLocal + Step * localDirection; 
00999 
01000   if( fExiting )
01001   {
01002 #ifdef G4DEBUG_NAVIGATION
01003     if( fVerbose > 2 )
01004     { 
01005       G4cout << " At G4Nav CompStep End - if(exiting) - fExiting= " << fExiting 
01006              << " fValidExitNormal = " << fValidExitNormal  << G4endl;
01007       G4cout << " fExitNormal= " << fExitNormal << G4endl;
01008     }
01009 #endif
01010 
01011     if(fValidExitNormal || fCalculatedExitNormal)
01012     {
01013       if (  fHistory.GetTopVolumeType()!=kReplica )
01014       {
01015         // Convention: fExitNormal is in the 'grand-mother' coordinate system
01016         //
01017         fGrandMotherExitNormal= fExitNormal;
01018         fCalculatedExitNormal= true;
01019       }
01020       else
01021       {
01022         fGrandMotherExitNormal = fExitNormal;
01023       }
01024     }
01025     else
01026     {  
01027       // We must calculate the normal anyway (in order to have it if requested)
01028       //
01029       G4ThreeVector finalLocalPoint =
01030             fLastLocatedPointLocal + localDirection*Step;
01031 
01032       if (  fHistory.GetTopVolumeType()!=kReplica )
01033       {
01034         // Find normal in the 'mother' coordinate system
01035         //
01036         G4ThreeVector exitNormalMotherFrame=
01037            motherLogical->GetSolid()->SurfaceNormal(finalLocalPoint);
01038         
01039         // Transform it to the 'grand-mother' coordinate system
01040         //
01041         const G4RotationMatrix* mRot = motherPhysical->GetRotation();
01042         if( mRot )
01043         {
01044           fChangedGrandMotherRefFrame= true;           
01045           fGrandMotherExitNormal = (*mRot).inverse() * exitNormalMotherFrame;
01046         }
01047         else
01048         {
01049           fGrandMotherExitNormal = exitNormalMotherFrame;
01050         }
01051 
01052         // Do not set fValidExitNormal -- this signifies
01053         // that the solid is convex!
01054         //
01055         fCalculatedExitNormal= true;
01056       }
01057       else
01058       {
01059         fCalculatedExitNormal = false;
01060         //
01061         // Nothing can be done at this stage currently - to solve this
01062         // Replica Navigation must have calculated the normal for this case
01063         // already.
01064         // Cases: mother is not convex, and exit is at previous replica level
01065 
01066 #ifdef G4DEBUG_NAVIGATION
01067         G4ExceptionDescription desc;
01068 
01069         desc << "Problem in ComputeStep:  Replica Navigation did not provide"
01070              << " valid exit Normal. " << G4endl;
01071         desc << " Do not know how calculate it in this case." << G4endl;
01072         desc << "  Location    = " << finalLocalPoint << G4endl;
01073         desc << "  Volume name = " << motherPhysical->GetName()
01074              << "  copy/replica No = " << motherPhysical->GetCopyNo() << G4endl;
01075         G4Exception("G4Navigator::ComputeStep()", "GeomNav0003",
01076                     JustWarning, desc, "Normal not available for exiting.");
01077 #endif
01078       }
01079     }
01080 
01081     // Now transform it to the global reference frame !!
01082     //
01083     if( fValidExitNormal || fCalculatedExitNormal )
01084     {
01085       G4int depth= fHistory.GetDepth();
01086       if( depth > 0 )
01087       {
01088         G4AffineTransform GrandMotherToGlobalTransf =
01089           fHistory.GetTransform(depth-1).Inverse();
01090         fExitNormalGlobalFrame =
01091           GrandMotherToGlobalTransf.TransformAxis( fGrandMotherExitNormal );
01092       }
01093       else
01094       {
01095         fExitNormalGlobalFrame= fGrandMotherExitNormal;
01096       }
01097     }
01098     else
01099     {
01100       fExitNormalGlobalFrame= G4ThreeVector( 0., 0., 0.);
01101     }
01102   }
01103   fStepEndPoint= pGlobalpoint+Step*pDirection; 
01104 
01105   if( (Step == pCurrentProposedStepLength) && (!fExiting) && (!fEntering) )
01106   {
01107     // This if Step is not really limited by the geometry.
01108     // The Navigator is obliged to return "infinity"
01109     //
01110     Step = kInfinity;
01111   }
01112 
01113 #ifdef G4VERBOSE
01114   if( fVerbose > 1 )
01115   {
01116     if( fVerbose >= 4 )
01117     {
01118       G4cout << "    ----- Upon exiting :" << G4endl;
01119       PrintState();
01120     }
01121     G4cout << "  Returned step= " << Step;
01122     if( fVerbose > 5 )   G4cout << G4endl;
01123     if( Step == kInfinity )
01124     {
01125        G4cout << " Requested step= " << pCurrentProposedStepLength ;
01126        if( fVerbose > 5) G4cout << G4endl;
01127     }
01128     G4cout << "  Safety = " << pNewSafety << G4endl;
01129   }
01130 #endif
01131 
01132   return Step;
01133 }
01134 
01135 // ********************************************************************
01136 // CheckNextStep
01137 //
01138 // Compute the step without altering the navigator state
01139 // ********************************************************************
01140 //
01141 G4double G4Navigator::CheckNextStep( const G4ThreeVector& pGlobalpoint,
01142                                      const G4ThreeVector& pDirection,
01143                                      const G4double pCurrentProposedStepLength,
01144                                            G4double& pNewSafety)
01145 {
01146   G4double step;
01147 
01148   // Save the state, for this parasitic call
01149   //
01150   SetSavedState();
01151 
01152   step = ComputeStep ( pGlobalpoint, 
01153                        pDirection,
01154                        pCurrentProposedStepLength, 
01155                        pNewSafety ); 
01156 
01157   // If a parasitic call, then attempt to restore the key parts of the state
01158   //
01159   RestoreSavedState(); 
01160 
01161   return step; 
01162 }
01163 
01164 // ********************************************************************
01165 // ResetState
01166 //
01167 // Resets stack and minimum of navigator state `machine'
01168 // ********************************************************************
01169 //
01170 void G4Navigator::ResetState()
01171 {
01172   fWasLimitedByGeometry  = false;
01173   fEntering              = false;
01174   fExiting               = false;
01175   fLocatedOnEdge         = false;
01176   fLastStepWasZero       = false;
01177   fEnteredDaughter       = false;
01178   fExitedMother          = false;
01179   fPushed                = false;
01180 
01181   fValidExitNormal       = false;
01182   fChangedGrandMotherRefFrame= false;
01183   fCalculatedExitNormal  = false;
01184 
01185   fExitNormal            = G4ThreeVector(0,0,0);
01186   fGrandMotherExitNormal = G4ThreeVector(0,0,0);
01187   fExitNormalGlobalFrame = G4ThreeVector(0,0,0);
01188 
01189   fPreviousSftOrigin     = G4ThreeVector(0,0,0);
01190   fPreviousSafety        = 0.0; 
01191 
01192   fNumberZeroSteps       = 0;
01193     
01194   fBlockedPhysicalVolume = 0;
01195   fBlockedReplicaNo      = -1;
01196 
01197   fLastLocatedPointLocal = G4ThreeVector( kInfinity, -kInfinity, 0.0 ); 
01198   fLocatedOutsideWorld   = false;
01199 }
01200 
01201 // ********************************************************************
01202 // SetupHierarchy
01203 //
01204 // Renavigates & resets hierarchy described by current history
01205 // o Reset volumes
01206 // o Recompute transforms and/or solids of replicated/parameterised volumes
01207 // ********************************************************************
01208 //
01209 void G4Navigator::SetupHierarchy()
01210 {
01211   G4int i;
01212   const G4int cdepth = fHistory.GetDepth();
01213   G4VPhysicalVolume *current;
01214   G4VSolid *pSolid;
01215   G4VPVParameterisation *pParam;
01216 
01217   for ( i=1; i<=cdepth; i++ )
01218   {
01219     current = fHistory.GetVolume(i);
01220     switch ( fHistory.GetVolumeType(i) )
01221     {
01222       case kNormal:
01223         break;
01224       case kReplica:
01225         freplicaNav.ComputeTransformation(fHistory.GetReplicaNo(i), current);
01226         break;
01227       case kParameterised:
01228         G4int replicaNo;
01229         pParam = current->GetParameterisation();
01230         replicaNo = fHistory.GetReplicaNo(i);
01231         pSolid = pParam->ComputeSolid(replicaNo, current);
01232 
01233         // Set up dimensions & transform in solid/physical volume
01234         //
01235         pSolid->ComputeDimensions(pParam, replicaNo, current);
01236         pParam->ComputeTransformation(replicaNo, current);
01237 
01238         G4TouchableHistory touchable( fHistory );
01239         touchable.MoveUpHistory();  // move up to the parent level
01240       
01241         // Set up the correct solid and material in Logical Volume
01242         //
01243         G4LogicalVolume *pLogical = current->GetLogicalVolume();
01244         pLogical->SetSolid( pSolid );
01245         pLogical->UpdateMaterial( pParam ->
01246           ComputeMaterial(replicaNo, current, &touchable) );
01247         break;
01248     }
01249   }
01250 }
01251 
01252 // ********************************************************************
01253 // GetLocalExitNormal
01254 //
01255 // Obtains the Normal vector to a surface (in local coordinates)
01256 // pointing out of previous volume and into current volume
01257 // ********************************************************************
01258 //
01259 G4ThreeVector G4Navigator::GetLocalExitNormal( G4bool* valid )
01260 {
01261   G4ThreeVector    ExitNormal(0.,0.,0.);
01262   G4VSolid        *currentSolid=0;
01263   G4LogicalVolume *candidateLogical;
01264   if ( fLastTriedStepComputation ) 
01265   {
01266     // use fLastLocatedPointLocal and next candidate volume
01267     //
01268     G4ThreeVector nextSolidExitNormal(0.,0.,0.);
01269 
01270     if( fEntering && (fBlockedPhysicalVolume!=0) ) 
01271     { 
01272       candidateLogical= fBlockedPhysicalVolume->GetLogicalVolume();
01273       if( candidateLogical ) 
01274       {
01275         // fLastStepEndPointLocal is in the coordinates of the mother
01276         // we need it in the daughter's coordinate system.
01277 
01278         // The following code should also work in case of Replica
01279         {
01280           // First transform fLastLocatedPointLocal to the new daughter
01281           // coordinates
01282           //
01283           G4AffineTransform MotherToDaughterTransform=
01284             GetMotherToDaughterTransform( fBlockedPhysicalVolume, 
01285                                           fBlockedReplicaNo,
01286                                           VolumeType(fBlockedPhysicalVolume) ); 
01287           G4ThreeVector daughterPointOwnLocal= 
01288             MotherToDaughterTransform.TransformPoint( fLastStepEndPointLocal ); 
01289 
01290           // OK if it is a parameterised volume
01291           //
01292           EInside  inSideIt; 
01293           G4bool   onSurface;
01294           G4double safety= -1.0; 
01295           currentSolid= candidateLogical->GetSolid(); 
01296           inSideIt  =  currentSolid->Inside(daughterPointOwnLocal); 
01297           onSurface =  (inSideIt == kSurface); 
01298           if( ! onSurface ) 
01299           {
01300             if( inSideIt == kOutside )
01301             {
01302               safety = (currentSolid->DistanceToIn(daughterPointOwnLocal)); 
01303               onSurface = safety < 100.0 * kCarTolerance; 
01304             }
01305             else if (inSideIt == kInside )
01306             {
01307               safety = (currentSolid->DistanceToOut(daughterPointOwnLocal)); 
01308               onSurface = safety < 100.0 * kCarTolerance; 
01309             }
01310           }
01311 
01312           if( onSurface ) 
01313           {
01314             nextSolidExitNormal =
01315               currentSolid->SurfaceNormal(daughterPointOwnLocal); 
01316  
01317             // Entering the solid ==> opposite
01318             //
01319             ExitNormal = -nextSolidExitNormal;
01320             fCalculatedExitNormal= true;
01321           }
01322           else
01323           {
01324 #ifdef G4VERBOSE
01325             if(( fVerbose == 1 ) && ( fCheck ))
01326             {
01327               std::ostringstream message;
01328               message << "Point not on surface ! " << G4endl
01329                       << "  Point           = "
01330                       << daughterPointOwnLocal << G4endl 
01331                       << "  Physical volume = "
01332                       << fBlockedPhysicalVolume->GetName() << G4endl
01333                       << "  Logical volume  = "
01334                       << candidateLogical->GetName() << G4endl
01335                       << "  Solid           = " << currentSolid->GetName() 
01336                       << "  Type            = "
01337                       << currentSolid->GetEntityType() << G4endl
01338                       << *currentSolid << G4endl;
01339               if( inSideIt == kOutside )
01340               { 
01341                 message << "Point is Outside. " << G4endl
01342                         << "  Safety (from outside) = " << safety << G4endl;
01343               }
01344               else // if( inSideIt == kInside ) 
01345               {
01346                 message << "Point is Inside. " << G4endl
01347                         << "  Safety (from inside) = " << safety << G4endl;              
01348               }
01349               G4Exception("G4Navigator::GetLocalExitNormal()", "GeomNav1001",
01350                           JustWarning, message);
01351             }
01352 #endif
01353           }
01354           *valid = onSurface;   //   was =true;
01355         }
01356       }
01357     }
01358     else if ( fExiting ) 
01359     {
01360       ExitNormal = fGrandMotherExitNormal;
01361       *valid = true;
01362       fCalculatedExitNormal= true;  // Should be true already
01363     }
01364     else  // i.e.  ( fBlockedPhysicalVolume == 0 )
01365     {
01366       *valid = false;
01367       G4Exception("G4Navigator::GetLocalExitNormal()",
01368                   "GeomNav0003", JustWarning, 
01369                   "Incorrect call to GetLocalSurfaceNormal." );
01370     }
01371   }
01372   else //  ( ! fLastTriedStepComputation ) ie. last call was to Locate
01373   {
01374     if ( EnteredDaughterVolume() )
01375     {
01376       G4VSolid* daughterSolid =fHistory.GetTopVolume()->GetLogicalVolume()
01377                                                       ->GetSolid();
01378       ExitNormal= -(daughterSolid->SurfaceNormal(fLastLocatedPointLocal));
01379       if( std::fabs(ExitNormal.mag2()-1.0 ) > CLHEP::perMillion )
01380       {
01381         G4ExceptionDescription desc;
01382         desc << " Parameters of solid: " << *daughterSolid
01383              << " Point for surface = " << fLastLocatedPointLocal << std::endl;
01384         G4Exception("G4Navigator::GetLocalExitNormal()",
01385                     "GeomNav0003", FatalException, desc,
01386                     "Surface Normal returned by Solid is not a Unit Vector." );
01387       }
01388       fCalculatedExitNormal= true;
01389       *valid = true;
01390     }
01391     else
01392     {
01393       if( fExitedMother )
01394       {
01395         ExitNormal = fGrandMotherExitNormal;
01396         *valid = true;
01397         fCalculatedExitNormal= true;
01398       }
01399       else  // We are not at a boundary. ExitNormal remains (0,0,0)
01400       { 
01401         *valid = false;
01402         fCalculatedExitNormal= false; 
01403         G4ExceptionDescription message; 
01404         message << "Function called when *NOT* at a Boundary." << G4endl;
01405         G4Exception("G4Navigator::GetLocalExitNormal()",
01406                     "GeomNav0003", JustWarning, message); 
01407       }
01408     }
01409   }
01410   return ExitNormal;
01411 }
01412 
01413 // ********************************************************************
01414 // GetMotherToDaughterTransform
01415 //
01416 // Obtains the mother to daughter affine transformation
01417 // ********************************************************************
01418 //
01419 G4AffineTransform
01420 G4Navigator::GetMotherToDaughterTransform( G4VPhysicalVolume *pEnteringPhysVol,
01421                                            G4int   enteringReplicaNo,
01422                                            EVolume enteringVolumeType ) 
01423 {
01424   switch (enteringVolumeType)
01425   {
01426     case kNormal:  // Nothing is needed to prepare the transformation
01427       break;       // It is stored already in the physical volume (placement)
01428     case kReplica: // Sets the transform in the Replica - tbc
01429       G4Exception("G4Navigator::GetMotherToDaughterTransform()",
01430                   "GeomNav0001", FatalException,
01431                   "Method NOT Implemented yet for replica volumes.");
01432       break;
01433     case kParameterised:
01434       if( pEnteringPhysVol->GetRegularStructureId() == 0 )
01435       {
01436         G4VPVParameterisation *pParam =
01437           pEnteringPhysVol->GetParameterisation();
01438         G4VSolid* pSolid =
01439           pParam->ComputeSolid(enteringReplicaNo, pEnteringPhysVol);
01440         pSolid->ComputeDimensions(pParam, enteringReplicaNo, pEnteringPhysVol);
01441 
01442         // Sets the transform in the Parameterisation
01443         //
01444         pParam->ComputeTransformation(enteringReplicaNo, pEnteringPhysVol);
01445 
01446         // Set the correct solid and material in Logical Volume
01447         //
01448         G4LogicalVolume* pLogical = pEnteringPhysVol->GetLogicalVolume();
01449         pLogical->SetSolid( pSolid );
01450       }
01451       break;
01452   }
01453   return G4AffineTransform(pEnteringPhysVol->GetRotation(), 
01454                            pEnteringPhysVol->GetTranslation()).Invert(); 
01455 }
01456 
01457 // ********************************************************************
01458 // GetLocalExitNormalAndCheck
01459 //
01460 // Obtains the Normal vector to a surface (in local coordinates)
01461 // pointing out of previous volume and into current volume, and
01462 // checks the current point against expected 'local' value.
01463 // ********************************************************************
01464 //
01465 G4ThreeVector G4Navigator::
01466 GetLocalExitNormalAndCheck( 
01467 #ifdef G4DEBUG_NAVIGATION
01468                            const G4ThreeVector& ExpectedBoundaryPointGlobal,
01469 #else
01470                            const G4ThreeVector&,
01471 #endif
01472                                  G4bool*        pValid)
01473 {
01474 #ifdef G4DEBUG_NAVIGATION
01475   // Check Current point against expected 'local' value
01476   //
01477   if ( fLastTriedStepComputation ) 
01478   {
01479     G4ThreeVector ExpectedBoundaryPointLocal;
01480 
01481     const G4AffineTransform& GlobalToLocal= GetGlobalToLocalTransform(); 
01482     ExpectedBoundaryPointLocal =
01483       GlobalToLocal.TransformPoint( ExpectedBoundaryPointGlobal ); 
01484 
01485     // Add here:  Comparison against expected position,
01486     //            i.e. the endpoint of ComputeStep
01487   }
01488 #endif
01489   
01490   return GetLocalExitNormal( pValid); 
01491 }
01492 
01493 // ********************************************************************
01494 // GetGlobalExitNormal
01495 //
01496 // Obtains the Normal vector to a surface (in global coordinates)
01497 // pointing out of previous volume and into current volume
01498 // ********************************************************************
01499 //
01500 G4ThreeVector 
01501 G4Navigator::GetGlobalExitNormal(const G4ThreeVector& IntersectPointGlobal,
01502                                        G4bool*        pNormalCalculated)
01503 {
01504   G4bool         validNormal;
01505   G4ThreeVector  localNormal, globalNormal;
01506 
01507   if( fLastTriedStepComputation && fExiting )  
01508   {
01509     // This was computed in ComputeStep -- and only on arrival at boundary
01510     //
01511     globalNormal = fExitNormalGlobalFrame; 
01512     *pNormalCalculated = true; // ComputeStep always computes it if Exiting
01513                                // (fExiting==true)
01514   }
01515   else
01516   {
01517     localNormal = GetLocalExitNormalAndCheck(IntersectPointGlobal,&validNormal);
01518     *pNormalCalculated = fCalculatedExitNormal;
01519 
01520 #ifdef G4DEBUG_NAVIGATION
01521     if( (!validNormal) && !fCalculatedExitNormal)
01522     {
01523       G4ExceptionDescription edN;
01524       edN << "  Calculated = " << fCalculatedExitNormal << G4endl;
01525       edN << "   Entering= "  << fEntering << G4endl;
01526       G4int oldVerbose= this->GetVerboseLevel();
01527       this->SetVerboseLevel(4);
01528       edN << "   State of Navigator: " << G4endl;
01529       edN << *this << G4endl;
01530       this->SetVerboseLevel( oldVerbose );
01531        
01532       G4Exception("G4Navigator::GetGlobalExitNormal()",
01533                   "GeomNav0003", JustWarning, edN,
01534                   "LocalExitNormalAndCheck() did not calculate Normal.");
01535      }
01536 #endif
01537      
01538      G4double localMag2= localNormal.mag2();
01539      if( validNormal && (std::fabs(localMag2-1.0)) > CLHEP::perMillion )
01540      {
01541        G4ExceptionDescription edN;
01542 
01543        edN << "G4Navigator::GetGlobalExitNormal: "
01544            << "  Using Local Normal - from call to GetLocalExitNormalAndCheck. "
01545            << G4endl
01546            << "  Local  Exit Normal = " << localNormal  << " || = "
01547            << std::sqrt(localMag2) << G4endl
01548            << "  Global Exit Normal = " << globalNormal << " || = "
01549            << globalNormal.mag() << G4endl;
01550        edN << "  Calculated It      = " << fCalculatedExitNormal << G4endl;
01551 
01552        G4Exception("G4Navigator::GetGlobalExitNormal()",
01553                    "GeomNav0003",JustWarning, edN,
01554                    "Value obtained from new local *solid* is incorrect.");
01555        localNormal = localNormal.unit(); // Should we correct it ??
01556      }
01557      G4AffineTransform localToGlobal = GetLocalToGlobalTransform();
01558      globalNormal = localToGlobal.TransformAxis( localNormal );
01559   }
01560 
01561 #ifdef G4DEBUG_NAVIGATION
01562   // Temporary extra checks
01563   if( fLastTriedStepComputation && fExiting)
01564   {
01565     localNormal = GetLocalExitNormalAndCheck( IntersectPointGlobal, &validNormal);
01566     *pNormalCalculated = fCalculatedExitNormal;
01567 
01568     G4AffineTransform localToGlobal = GetLocalToGlobalTransform();
01569     globalNormal = localToGlobal.TransformAxis( localNormal );
01570     
01571     // Check the value computed against fExitNormalGlobalFrame
01572     G4ThreeVector diffNorm = globalNormal - fExitNormalGlobalFrame;
01573     if( diffNorm.mag2() > perMillion*CLHEP::perMillion)
01574     {
01575       G4ExceptionDescription edDfn;
01576       edDfn << "Found difference in normals in case of exiting mother "
01577             << "- when Get is called after ComputingStep " << G4endl;
01578       edDfn << "  Magnitude of diff =      " << diffNorm.mag() << G4endl;
01579       edDfn << "  Normal stored (Global)     = " << fExitNormalGlobalFrame
01580             << G4endl;
01581       edDfn << "  Global Computed from Local = " << globalNormal << G4endl;
01582       G4Exception("G4Navigator::GetGlobalExitNormal()", "GeomNav0003",
01583                   JustWarning, edDfn);
01584     }
01585   }
01586 #endif
01587    
01588   return globalNormal;
01589 }
01590 
01591 // To make the new Voxel Safety the default, uncomment the next line
01592 // #define  G4NEW_SAFETY  1
01593 
01594 // ********************************************************************
01595 // ComputeSafety
01596 //
01597 // It assumes that it will be 
01598 //  i) called at the Point in the same volume as the EndPoint of the
01599 //     ComputeStep.
01600 // ii) after (or at the end of) ComputeStep OR after the relocation.
01601 // ********************************************************************
01602 //
01603 G4double G4Navigator::ComputeSafety( const G4ThreeVector &pGlobalpoint,
01604                                      const G4double pMaxLength,
01605                                      const G4bool keepState)
01606 {
01607   G4double newSafety = 0.0;
01608 
01609 #ifdef G4DEBUG_NAVIGATION
01610   G4int oldcoutPrec = G4cout.precision(8);
01611   if( fVerbose > 0 )
01612   {
01613     G4cout << "*** G4Navigator::ComputeSafety: ***" << G4endl
01614            << "    Called at point: " << pGlobalpoint << G4endl;
01615 
01616     G4VPhysicalVolume  *motherPhysical = fHistory.GetTopVolume();
01617     G4cout << "    Volume = " << motherPhysical->GetName() 
01618            << " - Maximum length = " << pMaxLength << G4endl; 
01619     if( fVerbose >= 4 )
01620     {
01621        G4cout << "    ----- Upon entering Compute Safety:" << G4endl;
01622        PrintState();
01623     }
01624   }
01625 #endif
01626 
01627   if (keepState)  { SetSavedState(); }
01628 
01629   G4double distEndpointSq = (pGlobalpoint-fStepEndPoint).mag2(); 
01630   G4bool   stayedOnEndpoint  = distEndpointSq < kCarTolerance*kCarTolerance; 
01631   G4bool   endpointOnSurface = fEnteredDaughter || fExitedMother;
01632 
01633   if( !(endpointOnSurface && stayedOnEndpoint) )
01634   {
01635     // Pseudo-relocate to this point (updates voxel information only)
01636     //
01637     LocateGlobalPointWithinVolume( pGlobalpoint ); 
01638       // --->> DANGER: Side effects on sub-navigator voxel information <<---
01639       //       Could be replaced again by 'granular' calls to sub-navigator
01640       //       locates (similar side-effects, but faster.  
01641       //       Solutions:
01642       //        1) Re-locate (to where?)
01643       //        2) Insure that the methods using (G4ComputeStep?)
01644       //           does a relocation (if information is disturbed only ?)
01645 
01646 #ifdef G4DEBUG_NAVIGATION
01647     if( fVerbose >= 2 )
01648     {
01649       G4cout << "  G4Navigator::ComputeSafety() relocates-in-volume to point: "
01650              << pGlobalpoint << G4endl;
01651     }
01652 #endif 
01653     G4VPhysicalVolume *motherPhysical = fHistory.GetTopVolume();
01654     G4LogicalVolume *motherLogical = motherPhysical->GetLogicalVolume();
01655     G4SmartVoxelHeader* pVoxelHeader = motherLogical->GetVoxelHeader();
01656     G4ThreeVector localPoint = ComputeLocalPoint(pGlobalpoint);
01657 
01658     if ( fHistory.GetTopVolumeType()!=kReplica )
01659     {
01660       switch(CharacteriseDaughters(motherLogical))
01661       {
01662         case kNormal:
01663           if ( pVoxelHeader )
01664           {
01665 #ifdef G4NEW_SAFETY
01666             static G4VoxelSafety sfVoxelSafety;
01667             G4double safetyTwo = sfVoxelSafety.ComputeSafety(localPoint,
01668                                            *motherPhysical, pMaxLength);
01669 #ifdef G4DEBUG_NAVIGATION
01670             G4double safetyNormal=
01671                fnormalNav.ComputeSafety(localPoint, fHistory, pMaxLength);
01672 
01673             G4double diffSafety= (safetyTwo - safetyNormal);
01674             if( std::fabs(diffSafety) > CLHEP::perMillion
01675                                       * std::fabs(safetyNormal) )
01676             {
01677               G4ExceptionDescription exd;
01678               exd << " ERROR in G4Navigator::ComputeStep " << G4endl;
01679               exd << "  Error>  DIFFERENCE in Safety from G4VoxelSafety "
01680                   << " - compared to value from Normal.  Diff = " << diffSafety << G4endl;
01681               exd << "  New Voxel Safety= " << safetyTwo
01682                   << "  Value from Normal= " << safetyNormal << G4endl;
01683               exd << "  Location:  Local coordinates = " << localPoint
01684                   << "  mother Volume= " << *motherPhysical->GetName()
01685                   << "  Copy# = " << motherPhysical->GetCopyNo() << G4endl;
01686               exd << "          : Global coordinates = " << pGlobalpoint
01687                   << G4endl;
01688               G4Exception("G4Navigator::ComputeSafety()",
01689                           "GeomNav0003", JustWarning, exd);
01690             }
01691             G4double safetyOldVoxel;
01692             safetyOldVoxel =
01693               fvoxelNav.ComputeSafety(localPoint,fHistory,pMaxLength);
01694 
01695             G4double diffBetterToApprox= safetyTwo - safetyOldVoxel;
01696             if( diffBetterToApprox < 0.0 )
01697             {
01698               G4ExceptionDescription exd;
01699               exd << "SMALLER Safety from call G4VoxelSafety compared to value "
01700                   << " from Normal.  Diff = " << diffSafety << G4endl;
01701               exd << " New Voxel Safety= " << safetyTwo
01702                   << "  Old value from Normal= " << newSafety << G4endl;
01703               exd << "  Location:  Local coordinates = " << localPoint
01704                   << "  mother Volume= " << *motherPhysical->GetName()
01705                   << "  Copy# = " << motherPhysical->GetCopyNo() << G4endl;
01706               exd << "          : Global coordinates = " << pGlobalpoint
01707                   << G4endl;
01708               G4Exception("G4Navigator::ComputeSafety()",
01709                           "GeomNav0003", JustWarning, exd);
01710             }
01711 #endif
01712 
01713             // newSafety= safetyOldVoxel;
01714             newSafety= safetyTwo;   // Faster and best available
01715 #else
01716             G4double safetyOldVoxel;
01717             safetyOldVoxel =
01718               fvoxelNav.ComputeSafety(localPoint,fHistory,pMaxLength);
01719             newSafety= safetyOldVoxel;
01720 #endif
01721           }
01722           else
01723           {
01724             newSafety=fnormalNav.ComputeSafety(localPoint,fHistory,pMaxLength);
01725           }
01726           break;
01727         case kParameterised:
01728           if( GetDaughtersRegularStructureId(motherLogical) != 1 )
01729           {
01730             newSafety=fparamNav.ComputeSafety(localPoint,fHistory,pMaxLength);
01731           }
01732           else  // Regular structure
01733           {
01734             newSafety=fregularNav.ComputeSafety(localPoint,fHistory,pMaxLength);
01735           }
01736           break;
01737         case kReplica:
01738           G4Exception("G4Navigator::ComputeSafety()", "GeomNav0001",
01739                       FatalException, "Not applicable for replicated volumes.");
01740           break;
01741       }
01742     }
01743     else
01744     {
01745       newSafety = freplicaNav.ComputeSafety(pGlobalpoint, localPoint,
01746                                             fHistory, pMaxLength);
01747     }
01748   }
01749   else // if( endpointOnSurface && stayedOnEndpoint )
01750   {
01751 #ifdef G4DEBUG_NAVIGATION
01752     if( fVerbose >= 2 )
01753     {
01754       G4cout << "    G4Navigator::ComputeSafety() finds that point - " 
01755              << pGlobalpoint << " - is on surface " << G4endl; 
01756       if( fEnteredDaughter ) { G4cout << "   entered new daughter volume"; }
01757       if( fExitedMother )    { G4cout << "   and exited previous volume."; }
01758       G4cout << G4endl;
01759       G4cout << " EndPoint was = " << fStepEndPoint << G4endl;
01760     } 
01761 #endif
01762     newSafety = 0.0; 
01763   }
01764 
01765   // Remember last safety origin & value
01766   //
01767   fPreviousSftOrigin = pGlobalpoint;
01768   fPreviousSafety = newSafety; 
01769 
01770   if (keepState)  { RestoreSavedState(); }
01771 
01772 #ifdef G4DEBUG_NAVIGATION
01773   if( fVerbose > 1 )
01774   {
01775     G4cout << "   ---- Exiting ComputeSafety  " << G4endl;
01776     if( fVerbose > 2 )  { PrintState(); }
01777     G4cout << "    Returned value of Safety = " << newSafety << G4endl;
01778   }
01779   G4cout.precision(oldcoutPrec);
01780 #endif
01781 
01782   return newSafety;
01783 }
01784 
01785 // ********************************************************************
01786 // CreateTouchableHistoryHandle
01787 // ********************************************************************
01788 //
01789 G4TouchableHistoryHandle G4Navigator::CreateTouchableHistoryHandle() const
01790 {
01791   return G4TouchableHistoryHandle( CreateTouchableHistory() );
01792 }
01793 
01794 // ********************************************************************
01795 // PrintState
01796 // ********************************************************************
01797 //
01798 void  G4Navigator::PrintState() const
01799 {
01800   G4int oldcoutPrec = G4cout.precision(4);
01801   if( fVerbose == 4 )
01802   {
01803     G4cout << "The current state of G4Navigator is: " << G4endl;
01804     G4cout << "  ValidExitNormal= " << fValidExitNormal << G4endl
01805            << "  ExitNormal     = " << fExitNormal      << G4endl
01806            << "  Exiting        = " << fExiting         << G4endl
01807            << "  Entering       = " << fEntering        << G4endl
01808            << "  BlockedPhysicalVolume= " ;
01809     if (fBlockedPhysicalVolume==0)
01810       G4cout << "None";
01811     else
01812       G4cout << fBlockedPhysicalVolume->GetName();
01813     G4cout << G4endl
01814            << "  BlockedReplicaNo     = " <<  fBlockedReplicaNo       << G4endl
01815            << "  LastStepWasZero      = " <<   fLastStepWasZero       << G4endl
01816            << G4endl;   
01817   }
01818   if( ( 1 < fVerbose) && (fVerbose < 4) )
01819   {
01820     G4cout << G4endl; // Make sure to line up
01821     G4cout << std::setw(30) << " ExitNormal "  << " "
01822            << std::setw( 5) << " Valid "       << " "     
01823            << std::setw( 9) << " Exiting "     << " "      
01824            << std::setw( 9) << " Entering"     << " " 
01825            << std::setw(15) << " Blocked:Volume "  << " "   
01826            << std::setw( 9) << " ReplicaNo"        << " "  
01827            << std::setw( 8) << " LastStepZero  "   << " "   
01828            << G4endl;   
01829     G4cout << "( " << std::setw(7) << fExitNormal.x() 
01830            << ", " << std::setw(7) << fExitNormal.y()
01831            << ", " << std::setw(7) << fExitNormal.z() << " ) "
01832            << std::setw( 5)  << fValidExitNormal  << " "   
01833            << std::setw( 9)  << fExiting          << " "
01834            << std::setw( 9)  << fEntering         << " ";
01835     if ( fBlockedPhysicalVolume==0 )
01836       G4cout << std::setw(15) << "None";
01837     else
01838       G4cout << std::setw(15)<< fBlockedPhysicalVolume->GetName();
01839       G4cout << std::setw( 9)  << fBlockedReplicaNo  << " "
01840              << std::setw( 8)  << fLastStepWasZero   << " "
01841              << G4endl;   
01842   }
01843   if( fVerbose > 2 ) 
01844   {
01845     G4cout.precision(8);
01846     G4cout << " Current Localpoint = " << fLastLocatedPointLocal << G4endl;
01847     G4cout << " PreviousSftOrigin  = " << fPreviousSftOrigin << G4endl;
01848     G4cout << " PreviousSafety     = " << fPreviousSafety << G4endl; 
01849   }
01850   G4cout.precision(oldcoutPrec);
01851 }
01852 
01853 // ********************************************************************
01854 // ComputeStepLog
01855 // ********************************************************************
01856 //
01857 void G4Navigator::ComputeStepLog(const G4ThreeVector& pGlobalpoint,
01858                                        G4double moveLenSq) const
01859 {
01860   //  The following checks only make sense if the move is larger
01861   //  than the tolerance.
01862 
01863   static const G4double fAccuracyForWarning   = kCarTolerance,
01864                         fAccuracyForException = 1000*kCarTolerance;
01865 
01866   G4ThreeVector OriginalGlobalpoint = fHistory.GetTopTransform().Inverse().
01867                                       TransformPoint(fLastLocatedPointLocal); 
01868 
01869   G4double shiftOriginSafSq = (fPreviousSftOrigin-pGlobalpoint).mag2();
01870 
01871   // Check that the starting point of this step is 
01872   // within the isotropic safety sphere of the last point
01873   // to a accuracy/precision  given by fAccuracyForWarning.
01874   //   If so give warning.
01875   //   If it fails by more than fAccuracyForException exit with error.
01876   //
01877   if( shiftOriginSafSq >= sqr(fPreviousSafety) )
01878   {
01879     G4double shiftOrigin = std::sqrt(shiftOriginSafSq);
01880     G4double diffShiftSaf = shiftOrigin - fPreviousSafety;
01881 
01882     if( diffShiftSaf > fAccuracyForWarning )
01883     {
01884       G4int oldcoutPrec= G4cout.precision(8);
01885       G4int oldcerrPrec= G4cerr.precision(10);
01886       std::ostringstream message, suggestion;
01887       message << "Accuracy error or slightly inaccurate position shift."
01888               << G4endl
01889               << "     The Step's starting point has moved " 
01890               << std::sqrt(moveLenSq)/mm << " mm " << G4endl
01891               << "     since the last call to a Locate method." << G4endl
01892               << "     This has resulted in moving " 
01893               << shiftOrigin/mm << " mm " 
01894               << " from the last point at which the safety " 
01895               << "     was calculated " << G4endl
01896               << "     which is more than the computed safety= " 
01897               << fPreviousSafety/mm << " mm  at that point." << G4endl
01898               << "     This difference is " 
01899               << diffShiftSaf/mm << " mm." << G4endl
01900               << "     The tolerated accuracy is "
01901               << fAccuracyForException/mm << " mm.";
01902 
01903       suggestion << " ";
01904       static G4int warnNow = 0;
01905       if( ((++warnNow % 100) == 1) )
01906       {
01907         message << G4endl
01908                << "  This problem can be due to either " << G4endl
01909                << "    - a process that has proposed a displacement"
01910                << " larger than the current safety , or" << G4endl
01911                << "    - inaccuracy in the computation of the safety";
01912         suggestion << "We suggest that you " << G4endl
01913                    << "   - find i) what particle is being tracked, and "
01914                    << " ii) through what part of your geometry " << G4endl
01915                    << "      for example by re-running this event with "
01916                    << G4endl
01917                    << "         /tracking/verbose 1 "  << G4endl
01918                    << "    - check which processes you declare for"
01919                    << " this particle (and look at non-standard ones)"
01920                    << G4endl
01921                    << "   - in case, create a detailed logfile"
01922                    << " of this event using:" << G4endl
01923                    << "         /tracking/verbose 6 ";
01924       }
01925       G4Exception("G4Navigator::ComputeStep()",
01926                   "GeomNav1002", JustWarning,
01927                   message, G4String(suggestion.str()));
01928       G4cout.precision(oldcoutPrec);
01929       G4cerr.precision(oldcerrPrec);
01930     }
01931 #ifdef G4DEBUG_NAVIGATION
01932     else
01933     {
01934       G4cerr << "WARNING - G4Navigator::ComputeStep()" << G4endl
01935              << "          The Step's starting point has moved "
01936              << std::sqrt(moveLenSq) << "," << G4endl
01937              << "          which has taken it to the limit of"
01938              << " the current safety. " << G4endl;
01939     }
01940 #endif
01941   }
01942   G4double safetyPlus = fPreviousSafety + fAccuracyForException;
01943   if ( shiftOriginSafSq > sqr(safetyPlus) )
01944   {
01945     std::ostringstream message;
01946     message << "May lead to a crash or unreliable results." << G4endl
01947             << "        Position has shifted considerably without"
01948             << " notifying the navigator !" << G4endl
01949             << "        Tolerated safety: " << safetyPlus << G4endl
01950             << "        Computed shift  : " << shiftOriginSafSq;
01951     G4Exception("G4Navigator::ComputeStep()", "GeomNav1002",
01952                 JustWarning, message);
01953   }
01954 }
01955 
01956 // ********************************************************************
01957 // Operator <<
01958 // ********************************************************************
01959 //
01960 std::ostream& operator << (std::ostream &os,const G4Navigator &n)
01961 {
01962   //  Old version did only the following:
01963   // os << "Current History: " << G4endl << n.fHistory;
01964   //  Old behaviour is recovered for fVerbose = 0
01965   
01966   // Adapted from G4Navigator::PrintState() const
01967 
01968   G4int oldcoutPrec = os.precision(4);
01969   if( n.fVerbose >= 4 )
01970   {
01971     os << "The current state of G4Navigator is: " << G4endl;
01972     os << "  ValidExitNormal= " << n.fValidExitNormal << G4endl
01973     << "  ExitNormal     = " << n.fExitNormal      << G4endl
01974     << "  Exiting        = " << n.fExiting         << G4endl
01975     << "  Entering       = " << n.fEntering        << G4endl
01976     << "  BlockedPhysicalVolume= " ;
01977     if (n.fBlockedPhysicalVolume==0)
01978       os << "None";
01979     else
01980       os << n.fBlockedPhysicalVolume->GetName();
01981     os << G4endl
01982     << "  BlockedReplicaNo     = " <<  n.fBlockedReplicaNo       << G4endl
01983     << "  LastStepWasZero      = " <<   n.fLastStepWasZero       << G4endl
01984     << G4endl;
01985   }
01986   if( ( 1 < n.fVerbose) && (n.fVerbose < 4) )
01987   {
01988     os << G4endl; // Make sure to line up
01989     os << std::setw(30) << " ExitNormal "  << " "
01990     << std::setw( 5) << " Valid "       << " "
01991     << std::setw( 9) << " Exiting "     << " "
01992     << std::setw( 9) << " Entering"     << " "
01993     << std::setw(15) << " Blocked:Volume "  << " "
01994     << std::setw( 9) << " ReplicaNo"        << " "
01995     << std::setw( 8) << " LastStepZero  "   << " "
01996     << G4endl;
01997     os << "( " << std::setw(7) << n.fExitNormal.x()
01998     << ", " << std::setw(7) << n.fExitNormal.y()
01999     << ", " << std::setw(7) << n.fExitNormal.z() << " ) "
02000     << std::setw( 5)  << n.fValidExitNormal  << " "
02001     << std::setw( 9)  << n.fExiting          << " "
02002     << std::setw( 9)  << n.fEntering         << " ";
02003     if ( n.fBlockedPhysicalVolume==0 )
02004       { os << std::setw(15) << "None"; }
02005     else
02006       { os << std::setw(15)<< n.fBlockedPhysicalVolume->GetName(); }
02007     os << std::setw( 9)  << n.fBlockedReplicaNo  << " "
02008     << std::setw( 8)  << n.fLastStepWasZero   << " "
02009     << G4endl;
02010   }
02011   if( n.fVerbose > 2 )
02012   {
02013     os.precision(8);
02014     os << " Current Localpoint = " << n.fLastLocatedPointLocal << G4endl;
02015     os << " PreviousSftOrigin  = " << n.fPreviousSftOrigin << G4endl;
02016     os << " PreviousSafety     = " << n.fPreviousSafety << G4endl;
02017   }
02018   if( n.fVerbose > 3 || n.fVerbose == 0 )
02019   {
02020     os << "Current History: " << G4endl << n.fHistory;
02021   }
02022     
02023   os.precision(oldcoutPrec);
02024   return os;
02025 }

Generated on Mon May 27 17:48:57 2013 for Geant4 by  doxygen 1.4.7