How can I get the information of the original primary particle?


G4Track is a transient class and there is no given way to store the generation of the particles. There is a user class named G4VUserTrackInformation which can be attached to G4Track. That is, by creating a concrete class object derived from this abstract base class and attaching it to G4Track, you can keep primary information with any secondary particle.

A sample implementation of G4VUserTrackInformation is the following.

Header file (T01TrackInformation.hh)

#ifndef T01TrackInformation_h
#define T01TrackInformation_h 1

#include "globals.hh"
#include "G4ThreeVector.hh"
#include "G4ParticleDefinition.hh"
#include "G4Track.hh"
#include "G4Allocator.hh"
#include "G4VUserTrackInformation.hh"

class T01TrackInformation : public G4VUserTrackInformation 
{
  public:
    T01TrackInformation();
    T01TrackInformation(const G4Track* aTrack);
    T01TrackInformation(const T01TrackInformation* aTrackInfo);
    virtual ~T01TrackInformation();
   
    inline void *operator new(size_t);
    inline void operator delete(void *aTrackInfo);
    inline int operator ==(const T01TrackInformation& right) const
    {return (this==&right);}

    void Print() const;

  private:
    G4int                 originalTrackID;
    G4ParticleDefinition* particleDefinition;
    G4ThreeVector         originalPosition;
    G4ThreeVector         originalMomentum;
    G4double              originalEnergy;
    G4double              originalTime;

  public:
    inline G4int GetOriginalTrackID() const {return originalTrackID;}
    inline G4ParticleDefinition* GetOriginalParticle() const {return particleDefinition;}
    inline G4ThreeVector GetOriginalPosition() const {return originalPosition;}
    inline G4ThreeVector GetOriginalMomentum() const {return originalMomentum;}
    inline G4double GetOriginalEnergy() const {return originalEnergy;}
    inline G4double GetOriginalTime() const {return originalTime;}
};

extern G4Allocator<T01TrackInformation> aTrackInformationAllocator;

inline void* T01TrackInformation::operator new(size_t)
{ void* aTrackInfo;
  aTrackInfo = (void*)aTrackInformationAllocator.MallocSingle();
  return aTrackInfo;
}

inline void T01TrackInformation::operator delete(void *aTrackInfo)
{ aTrackInformationAllocator.FreeSingle((T01TrackInformation*)aTrackInfo);}

#endif

Source file (T01TrackInformation.cc)

#include "T01TrackInformation.hh"
#include "G4ios.hh"

G4Allocator<T01TrackInformation> aTrackInformationAllocator;

T01TrackInformation::T01TrackInformation()
{
    originalTrackID = 0;
    particleDefinition = 0;
    originalPosition = G4ThreeVector(0.,0.,0.);
    originalMomentum = G4ThreeVector(0.,0.,0.);
    originalEnergy = 0.;
    originalTime = 0.;
}

T01TrackInformation::T01TrackInformation(const G4Track* aTrack)
{
    originalTrackID = aTrack->GetTrackID();
    particleDefinition = aTrack->GetDefinition();
    originalPosition = aTrack->GetPosition();
    originalMomentum = aTrack->GetMomentum();
    originalEnergy = aTrack->GetTotalEnergy();
    originalTime = aTrack->GetGlobalTime();
}

T01TrackInformation::T01TrackInformation(const T01TrackInformation* aTrackInfo)
{
    originalTrackID = aTrackInfo->originalTrackID;
    particleDefinition = aTrackInfo->particleDefinition;
    originalPosition = aTrackInfo->originalPosition;
    originalMomentum = aTrackInfo->originalMomentum;
    originalEnergy = aTrackInfo->originalEnergy;
    originalTime = aTrackInfo->originalTime;
}

T01TrackInformation::~T01TrackInformation(){;}

void T01TrackInformation::Print() const
{
    G4cout 
     << "Original track ID " << originalTrackID 
     << " at " << originalPosition << G4endl;
}

Once you create above, the class you have to implement for creating this class objects is G4UserTrackingAction. Here is the sample implementation.

Source file (T01TrackingAction.cc)

#include "T01TrackingAction.hh"
#include "G4TrackingManager.hh"
#include "G4Track.hh"
#include "G4TrackVector.hh"
#include "T01TrackInformation.hh"

T01TrackingAction::T01TrackingAction()
{;}

T01TrackingAction::~T01TrackingAction()
{;}

void T01TrackingAction::PreUserTrackingAction(const G4Track* aTrack)
{
  if(aTrack->GetParentID()==0 && aTrack->GetUserInformation()==0)
  {
    T01TrackInformation* anInfo = new T01TrackInformation(aTrack);
    G4Track* theTrack = (G4Track*)aTrack;
    theTrack->SetUserInformation(anInfo);
  }
}

void T01TrackingAction::PostUserTrackingAction(const G4Track* aTrack)
{
  G4TrackVector* secondaries = fpTrackingManager->GimmeSecondaries();
  if(secondaries)
  {
    T01TrackInformation* info = (T01TrackInformation*)(aTrack->GetUserInformation());
    size_t nSeco = secondaries->size();
    if(nSeco>0)
    {
      for(size_t i=0;i<nSeco;i++)
      { 
        T01TrackInformation* infoNew = new T01TrackInformation(info);
        (*secondaries)[i]->SetUserInformation(infoNew);
      }
    }
  }
}

Now you can get the information of the primary particle associating with seconday G4Track object. Thus, for example in your sensitive detector implementation or in your G4UserSteppingAction, you can retrieve it as following.

Source file for the case of a sensitive detector
#include "T01TrackInformation.hh"

...

G4bool MySensitiveDetector::ProcessHits(G4Step*aStep,G4TouchableHistory*ROhist)
{
  T01TrackInformation* info = (T01TrackInformation*)(aStep->GetTrack()->GetUserInformation());
  G4cout << " Original Track ID " << info->GetOriginalTrackID() << G4endl;

  ...
}


Makoto Asai