FAQ.4.  Tracks and steps

Q:

How can I access the track information through the step object and what information am I allowed to access ?

A:

A G4Step object consists of two points:

      G4StepPoint* point1 = step->GetPreStepPoint();
      G4StepPoint* point2 = step->GetPostStepPoint();
  

To get their positions in the global coordinate system:

      G4ThreeVector pos1 = point1->GetPosition();
      G4ThreeVector pos2 = point2->GetPosition();
  

Hereafter we call current volume the volume where the step has just gone through. Geometrical informations are available from preStepPoint. G4VTouchable and its derivates keep these geometrical informations. We retrieve a touchable by creating a handle for it:

      G4TouchableHandle touch1 = point1->GetTouchableHandle();
  

To get the current volume:

       G4VPhysicalVolume* volume = touch1->GetVolume();
  

To get its name:

      G4String name = volume->GetName();
  

To get the physical volume copy number:

      G4int copyNumber = touch1->GetCopyNumber();
  

To get logical volume:

      G4LogicalVolume* lVolume = volume->GetLogicalVolume();
  

To get the associated material: the following statements are equivalent:

      G4Material* material = point1  ->GetMaterial();
      G4Material* material = lVolume ->GetMaterial();
  

To get the geometrical region:

      G4Region* region = lVolume->GetRegion();
  

To get its mother volume:

      G4VPhysicalVolume* mother = touch1->GetVolume(depth=1);
      grandMother: depth=2 ...etc...
  

To get the copy number of the mother volume:

      G4int copyNumber = touch1->GetCopyNumber(depth=1);
      grandMother: depth=2 ...etc...
  

To get the process which has limited the current step:

      G4VProcess* aProcess = point2->GetProcessDefinedStep();
  

To check that the particle has just entered in the current volume (i.e. it is at the first step in the volume; the preStepPoint is at the boundary):

      if (point1->GetStepStatus() == fGeomBoundary) 
  

To check that the particle is leaving the current volume (i.e. it is at the last step in the volume; the postStepPoint is at the boundary):

      if (point2->GetStepStatus() == fGeomBoundary)
  

In the above situation, to get touchable of the next volume:

      G4TouchableHandle touch2 = point2->GetTouchableHandle();
  

From touch2, all informations on the next volume can be retrieved as above.

Physics quantities are available from the step (G4Step) or from the track (G4Track).

To get the energy deposition, step length, displacement and time of flight spent by the current step:

      G4double eDeposit      = step->GetTotalEnergyDeposit();
      G4double sLength       = step->GetStepLength();
      G4ThreeVector displace = step->GetDeltaPosition();
      G4double tof           = step->GetDeltaTime();
  

To get momentum, kinetic energy and global time (time since the beginning of the event) of the track after the completion of the current step:

      G4Track* track         = step->GetTrack();
      G4ThreeVector momentum = track->GetMomentum();
      G4double kinEnergy     = track->GetKineticEnergy();
      G4double globalTime    = track->GetGlobalTime();
      ...etc...
  

Remark

To transform a position from the global coordinate system to the local system of the current volume, use the preStepPoint transformation, as described in the geometry section above.

Q:

How can I get and store (or plot) informations at tracking time from a given volume ?

A:

To get the information at tracking time in a given volume A, one can adopt either one or a combination of the following strategies:

  1. If the geometry is simple enough, and wish to score some commonly used physics quantities (e.g. energy deposition, dose, flux, etc.), just activate G4ScoringManager in your main program, and use the scorer-based UI commands to transform volume A into a scorer.

    See Option 6 below, and the example RE03 in examples/extended/runAndEvent.

  2. Through the SteppingAction, check that the particle is inside volume A and do whatever needed. Hints can be found in the previous section of this FAQ document.

    Usually, the hits containers and histograms are attributes of a Track, Event or Run and can be managed through either a TrackingAction, EventAction and/or RunAction and eventually messaging their pointer to the SteppingAction.

    A similar approach is illustrated in examples/basic B2, B4, extended/electromagnetic, optical, and many others...

  3. In DetectorConstruction, by declaring volume A as a SensitiveDetector. At stepping time, the Geant4 kernel will automatically check that a particle is inside volume A and will handle the control to a specific function G4VSensitiveDetector::ProcessHits(). It is just necessary to instanciate a class inherited from G4VSensitiveDetector, say VolumeA_SD, and do whatever needed by implementing the function VolumeA_SD::ProcessHits(), as described in Option 2 above.

  4. In addition to Option 3 above, should create a HitsCollection to store the information. A HitsCollection can be created in VolumeA_SD::Initialize(). A Hit can be created or filled in VolumeA_SD::ProcessHits(). Additional operations on HitsCollection can be performed in VolumeA_SD::EndOfEvent().

    This approach is illustrated in examples/basic B2, B4 and extended/analysis, extended/runAndEvent RE01, etc...

  5. In DetectorConstruction, volume A can be declared as SensitiveDetector, and one or several pre-defined scorers can be attached to volume A. In this case, neither a SteppingAction nor a spcific VolumeA_SD sensitive detector is needed any longer. It is just necessary to create a dedicated scorer, e.g. MyRunScorer, inherited from G4Run, and handle the HitsCollections within MyRunScorer::RecordEvent(). MyRunScorer itself can be instanciated from RunAction::GenerateRun().

    This approach is illustrated in examples/novice N07, extended/runAndEvent RE02.

  6. A set of build-in scorer-based UI commands allows to perform most possible operations described through the previous Option 5 directly from run-time macros.

    See example extended/runAndEvent RE03.