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

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