The Geometry Navigator

Navigation through the geometry at tracking time is implemented by the class G4Navigator. The navigator is used to locate points in the geometry and compute distances to geometry boundaries. At tracking time, the navigator is intended to be the only point of interaction with tracking.

Internally, the G4Navigator has several private helper/utility classes:

  • G4NavigationHistory - stores the compounded transformations, replication/parameterisation information, and volume pointers at each level of the hierarchy to the current location. The volume types at each level are also stored - whether normal (placement), replicated or parameterised.

  • G4NormalNavigation - provides location & distance computation functions for geometries containing ‘placement’ volumes, with no voxels.

  • G4VoxelNavigation - provides location and distance computation functions for geometries containing ‘placement’ physical volumes with voxels. Internally a stack of voxel information is maintained. Private functions allow for isotropic distance computation to voxel boundaries and for computation of the ‘next voxel’ in a specified direction.

  • G4ParameterisedNavigation - provides location and distance computation functions for geometries containing parameterised volumes with voxels. Voxel information is maintained similarly to G4VoxelNavigation, but computation can also be simpler by adopting voxels to be one level deep only (unrefined, or 1D optimisation)

  • G4ReplicaNavigation - provides location and distance computation functions for replicated volumes.

  • G4RegularNavigation - provides location and distance computation functions for fast navigation in volumes containing a regular parameterisation. If two contiguous voxels have the same material, navigation does not stop at the surface.

In addition, the navigator maintains a set of flags for exiting/entry optimisation. A navigator is not a singleton class; this is mainly to allow a design extension in future (e.g. geometrical event biasing).

Using the navigator to locate points

More than one navigator object can be created inside an application; these navigators can act independently for different purposes. The main navigator which is activated automatically at the startup of a simulation program is the navigator used for the tracking and attached the world volume of the main tracking (or mass) geometry.

The navigator for tracking can be retrieved at any state of the application by messaging the G4TransportationManager:

G4Navigator* tracking_navigator =
  G4TransportationManager::GetInstance()->GetNavigatorForTracking();

This also allows to retrieve at any time a pointer to the world volume assigned for tracking:

G4VPhysicalVolume* tracking_world = tracking_navigator->GetWorldVolume();

The navigator for tracking also retains all the information of the current history of volumes traversed at a precise moment of the tracking during a run. Therefore, if the navigator for tracking is used during tracking for locating a generic point in the tree of volumes, the actual particle gets also -relocated- in the specified position and tracking will be of course affected !

In order to avoid the problem above and provide information about location of a point without affecting the tracking, it is suggested to either use an alternative G4Navigator object (which can then be assigned to the world-volume), or access the information through the step.

If the user instantiates an alternative G4Navigator, ownership is retained by the user’s code, and the navigator object should be deleted by that code.

Using the ‘step’ to retrieve geometrical information

During the tracking run, geometrical information can be retrieved through the touchable handle associated to the current step. For example, to identify the exact copy-number of a specific physical volume in the mass geometry, one should do the following:

// Given the pointer to the step object ...
//
G4Step* aStep = ..;

// ... retrieve the 'pre-step' point
//
G4StepPoint* preStepPoint = aStep->GetPreStepPoint();

// ... retrieve a touchable handle and access to the information
//
G4TouchableHandle theTouchable = preStepPoint->GetTouchableHandle();
G4int copyNo = theTouchable->GetCopyNumber();
G4int motherCopyNo = theTouchable->GetCopyNumber(1);

To determine the exact position in global coordinates in the mass geometry and convert to local coordinates (local to the current volume):

G4ThreeVector worldPosition = preStepPoint->GetPosition();
G4ThreeVector localPosition = theTouchable->GetHistory()->
              GetTopTransform().TransformPoint(worldPosition);

Using an alternative navigator to locate points

In order to know (when in the idle state of the application) in which physical volume a given point is located in the detector geometry, it is necessary to create an alternative navigator object first and assign it to the world volume:

G4Navigator* aNavigator = new G4Navigator();
aNavigator->SetWorldVolume(worldVolumePointer);

Then, locate the point myPoint (defined in global coordinates), retrieve a touchable handle and do whatever you need with it:

aNavigator->LocateGlobalPointAndSetup(myPoint);
G4TouchableHistoryHandle aTouchable =
     aNavigator->CreateTouchableHistoryHandle();

   // Do whatever you need with it ...
   // ... convert point in local coordinates (local to the current volume)
   //
   G4ThreeVector localPosition = aTouchable->GetHistory()->
                 GetTopTransform().TransformPoint(myPoint);

   // ... convert back to global coordinates system
   G4ThreeVector globalPosition = aTouchable->GetHistory()->
                 GetTopTransform().Inverse().TransformPoint(localPosition);

If outside of the tracking run and given a generic local position (local to a given volume in the geometry tree), it is -not- possible to determine a priori its global position and convert it to the global coordinates system. The reason for this is rather simple, nobody can guarantee that the given (local) point is located in the right -copy- of the physical volume ! In order to retrieve this information, some extra knowledge related to the absolute position of the physical volume is required first, i.e. one should first determine a global point belonging to that volume, eventually making a dedicated scan of the geometry tree through a dedicated G4Navigator object and then apply the method above after having created the touchable for it.

Fast navigation in regular patterned geometries and phantoms

Since release 9.1 of Geant4, a specialised navigation algorithm has been introduced to allow for optimal memory use and extremely efficient navigation in geometries represented by a regular pattern of volumes and particularly three-dimensional grids of boxes. A typical application of this kind is the case of DICOM phantoms for medical physics studies.

The class G4RegularNavigation is used and automatically activated when such geometries are defined. It is required to the user to implement a parameterisation of the kind G4PhantomParameterisation and place the parameterised volume containing it in a container volume, so that all cells in the three-dimensional grid (voxels) completely fill the container volume. This way the location of a point inside a voxel can be done in a fast way, transforming the position to the coordinate system of the container volume and doing a simple calculation of the kind:

copyNo_x = (localPoint.x()+fVoxelHalfX*fNoVoxelX)/(fVoxelHalfX*2.)

where fVoxelHalfX is the half dimension of the voxel along X and fNoVoxelX is the number of voxels in the X dimension. Voxel 0 will be the one closest to the corner (fVoxelHalfX*fNoVoxelX, fVoxelHalfY*fNoVoxelY, fVoxelHalfZ*fNoVoxelZ).

Having the voxels filling completely the container volume allows to avoid the lengthy computation of ComputeStep() and ComputeSafety methods required in the traditional navigation algorithm. In this case, when a track is inside the parent volume, it has always to be inside one of the voxels and it will be only necessary to calculate the distance to the walls of the current voxel.

Skipping borders of voxels with same material

Another speed optimisation can be provided by skipping the frontiers of two voxels which the same material assigned, so that bigger steps can be done. This optimisation may be not very useful when the number of materials is very big (in which case the probability of having contiguous voxels with same material is reduced), or when the physical step is small compared to the voxel dimensions (very often the case of electrons). The optimisation can be switched off in such cases, by invoking the following method with argument skip = 0:

Phantoms with only one material

If you want to describe a phantom of a unique material, you may spare some memory by not filling the set of indices of materials of each voxel. If the method SetMaterialIndices() is not invoked, the index for all voxels will be 0, that is the first (and unique) material in your list.

G4RegularParameterisation::SetSkipEqualMaterials( G4bool skip );

Example

To use the specialised navigation, it is required to first create an object of type G4PhantomParameterisation:

G4PhantomParameterisation* param = new G4PhantomParameterisation();

Then, fill it with the all the necessary data:

// Voxel dimensions in the three dimensions
//
G4double halfX = ...;
G4double halfY = ...;
G4double halfZ = ...;
param->SetVoxelDimensions( halfX, halfY, halfZ );

// Number of voxels in the three dimensions
//
G4int nVoxelX = ...;
G4int nVoxelY = ...;
G4int nVoxelZ = ...;
param->SetNoVoxel( nVoxelX, nVoxelY, nVoxelZ );

// Vector of materials of the voxels
//
std::vector < G4Material* > theMaterials;
theMaterials.push_back( new G4Material( ...
theMaterials.push_back( new G4Material( ...
param->SetMaterials( theMaterials );

// List of material indices
// For each voxel it is a number that correspond to the index of its
// material in the vector of materials defined above;
//
size_t* mateIDs = new size_t[nVoxelX*nVoxelY*nVoxelZ];
mateIDs[0] = n0;
mateIDs[1] = n1;
...
param->SetMaterialIndices( mateIDs );

Then, define the volume that contains all the voxels:

G4Box* cont_solid = new G4Box("PhantomContainer",nVoxelX*halfX.,nVoxelY*halfY.,nVoxelZ*halfZ);
G4LogicalVolume* cont_logic =
  new G4LogicalVolume( cont_solid,
           matePatient,        // material is not relevant here...
           "PhantomContainer",
           0, 0, 0 );
G4VPhysicalVolume * cont_phys =
  new G4PVPlacement(rotm,                  // rotation
            pos,                   // translation
            cont_logic,            // logical volume
            "PhantomContainer",    // name
            world_logic,           // mother volume
            false,                 // No op. bool.
            1);                    // Copy number

The physical volume should be assigned as the container volume of the parameterisation:

param->BuildContainerSolid(cont_phys);

// Assure that the voxels are completely filling the container volume
//
param->CheckVoxelsFillContainer( cont_solid->GetXHalfLength(),
                                    cont_solid->GetyHalfLength(),
                                    cont_solid->GetzHalfLength() );

// The parameterised volume which uses this parameterisation is placed
// in the container logical volume
//
G4PVParameterised * patient_phys =
  new G4PVParameterised("Patient",               // name
                        patient_logic,           // logical volume
                        cont_logic,              // mother volume
            kXAxis,                  // optimisation hint
                        nVoxelX*nVoxelY*nVoxelZ, // number of voxels
                        param);                  // parameterisation

// Indicate that this physical volume is having a regular structure
//
patient_phys->SetRegularStructureId(1);

An example showing the application of the optimised navigation algorithm for phantoms geometries is available in examples/extended/medical/DICOM. It implements a real application for reading DICOM images and convert them to Geant4 geometries with defined materials and densities, allowing for different implementation solutions to be chosen (non-optimised, classical 3D optimisation, nested parameterisations and use of G4PhantomParameterisation).

Run-time commands

When running in verbose mode (i.e. the default, G4VERBOSE set while installing the Geant4 kernel libraries), the navigator provides a few commands to control its behavior. It is possible to select different verbosity levels (up to 5), with the command:

geometry/navigator/verbose [verbose_level]

or to force the navigator to run in check mode:

geometry/navigator/check_mode [true/false]

The latter will force more strict and less tolerant checks in step/safety computation to verify the correctness of the solids’ response in the geometry.

By combining check_mode with verbosity level-1, additional verbosity checks on the response from the solids can be activated.

Setting Geometry Tolerance to be relative

The tolerance value defining the accuracy of tracking on the surfaces is by default set to a reasonably small value of 10E-9 mm. Such accuracy may be however redundant for use on simulation of detectors of big size or macroscopic dimensions. Since release 9.0, it is possible to specify the surface tolerance to be relative to the extent of the world volume defined for containing the geometry setup.

The class G4GeometryManager can be used to activate the computation of the surface tolerance to be relative to the geometry setup which has been defined. It can be done this way:

G4GeometryManager::GetInstance()->SetWorldMaximumExtent(WorldExtent);

where, WorldExtent is the actual maximum extent of the world volume used for placing the whole geometry setup.

Such call to G4GeometryManager must be done before defining any geometrical component of the setup (solid shape or volume), and can be done only once!

The class G4GeometryTolerance is to be used for retrieving the actual values defined for tolerances, surface (Cartesian), angular or radial respectively:

G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
G4GeometryTolerance::GetInstance()->GetAngularTolerance();
G4GeometryTolerance::GetInstance()->GetRadialTolerance();