299 lines
9.6 KiB
Plaintext
299 lines
9.6 KiB
Plaintext
//
|
|
// ********************************************************************
|
|
// * DISCLAIMER *
|
|
// * *
|
|
// * The following disclaimer summarizes all the specific disclaimers *
|
|
// * of contributors to this software. The specific disclaimers,which *
|
|
// * govern, are listed with their locations in: *
|
|
// * http://cern.ch/geant4/license *
|
|
// * *
|
|
// * Neither the authors of this software system, nor their employing *
|
|
// * institutes,nor the agencies providing financial support for this *
|
|
// * work make any representation or warranty, express or implied, *
|
|
// * regarding this software system or assume any liability for its *
|
|
// * use. *
|
|
// * *
|
|
// * This code implementation is the intellectual property of the *
|
|
// * GEANT4 collaboration. *
|
|
// * By copying, distributing or modifying the Program (or any work *
|
|
// * based on the Program) you indicate your acceptance of this *
|
|
// * statement, and all its terms. *
|
|
// ********************************************************************
|
|
//
|
|
//
|
|
// $Id: G4Track.icc,v 1.10 2002/12/16 11:59:12 gcosmo Exp $
|
|
// GEANT4 tag $Name: geant4-05-02-patch-01 $
|
|
//
|
|
//-----------------------------------------------------------------
|
|
// Definitions of inline functions
|
|
//-----------------------------------------------------------------
|
|
// change GetMaterial 16 Feb. 2000 H.Kurashige
|
|
|
|
// Operators
|
|
extern G4Allocator<G4Track> aTrackAllocator;
|
|
inline void* G4Track::operator new(size_t)
|
|
{ void *aTrack;
|
|
aTrack = (void *) aTrackAllocator.MallocSingle();
|
|
return aTrack;
|
|
}
|
|
// Override "new" for "G4Allocator".
|
|
|
|
inline void G4Track::operator delete(void *aTrack)
|
|
{ aTrackAllocator.FreeSingle((G4Track *) aTrack);}
|
|
// Override "delete" for "G4Allocator".
|
|
|
|
inline G4bool G4Track::operator==( const G4Track& s)
|
|
{ return (this==&s); }
|
|
// Define "==" operator because "G4TrackVector" uses
|
|
// "RWPtrOrderdVector" which requires this.
|
|
|
|
// Get/Set functions
|
|
// dynamic particle
|
|
inline const G4DynamicParticle* G4Track::GetDynamicParticle() const
|
|
{ return fpDynamicParticle; }
|
|
|
|
// particle definition
|
|
inline G4ParticleDefinition* G4Track::GetDefinition() const
|
|
{ return fpDynamicParticle->GetDefinition(); }
|
|
|
|
// parent track ID
|
|
inline G4int G4Track::GetParentID() const
|
|
{ return fParentID; }
|
|
|
|
inline void G4Track::SetParentID(const G4int aValue)
|
|
{ fParentID = aValue; }
|
|
|
|
// current track ID
|
|
inline G4int G4Track::GetTrackID() const
|
|
{ return fTrackID; }
|
|
|
|
inline void G4Track::SetTrackID(const G4int aValue)
|
|
{ fTrackID = aValue; }
|
|
|
|
// position
|
|
inline const G4ThreeVector& G4Track::GetPosition() const
|
|
{ return fPosition; }
|
|
|
|
inline void G4Track::SetPosition(const G4ThreeVector& aValue)
|
|
{ fPosition = aValue; }
|
|
|
|
// global time
|
|
inline G4double G4Track::GetGlobalTime() const
|
|
{ return fGlobalTime; }
|
|
|
|
inline void G4Track::SetGlobalTime(const G4double aValue)
|
|
{ fGlobalTime = aValue; }
|
|
// Time since the event in which the track belongs is created.
|
|
|
|
// local time
|
|
inline G4double G4Track::GetLocalTime() const
|
|
{ return fLocalTime; }
|
|
|
|
inline void G4Track::SetLocalTime(const G4double aValue)
|
|
{ fLocalTime = aValue; }
|
|
// Time since the current track is created.
|
|
|
|
// proper time
|
|
inline G4double G4Track::GetProperTime() const
|
|
{ return fpDynamicParticle->GetProperTime(); }
|
|
|
|
inline void G4Track::SetProperTime(const G4double aValue)
|
|
{ fpDynamicParticle->SetProperTime(aValue); }
|
|
// Proper time of the current track
|
|
|
|
// volume
|
|
inline G4VPhysicalVolume* G4Track::GetVolume() const
|
|
{ return fpTouchable->GetVolume(); }
|
|
|
|
inline G4VPhysicalVolume* G4Track::GetNextVolume() const
|
|
{ return fpNextTouchable->GetVolume(); }
|
|
|
|
// material
|
|
inline
|
|
const G4MaterialCutsCouple* G4Track::GetMaterialCutsCouple() const
|
|
{ return fpStep->GetPreStepPoint()->GetMaterialCutsCouple(); }
|
|
|
|
inline
|
|
const G4MaterialCutsCouple* G4Track::GetNextMaterialCutsCouple() const
|
|
{ return fpStep->GetPostStepPoint()->GetMaterialCutsCouple(); }
|
|
|
|
// material
|
|
inline G4Material* G4Track::GetMaterial() const
|
|
{ return fpStep->GetPreStepPoint()->GetMaterial(); }
|
|
|
|
inline G4Material* G4Track::GetNextMaterial() const
|
|
{ return fpStep->GetPostStepPoint()->GetMaterial(); }
|
|
|
|
|
|
// touchable
|
|
inline const G4VTouchable* G4Track::GetTouchable() const
|
|
{
|
|
return fpTouchable();
|
|
}
|
|
|
|
inline const G4TouchableHandle& G4Track::GetTouchableHandle() const
|
|
{
|
|
return fpTouchable;
|
|
}
|
|
|
|
inline void G4Track::SetTouchableHandle( const G4TouchableHandle& apValue)
|
|
{
|
|
fpTouchable = apValue;
|
|
}
|
|
|
|
inline const G4VTouchable* G4Track::GetNextTouchable() const
|
|
{
|
|
return fpNextTouchable();
|
|
}
|
|
|
|
inline const G4TouchableHandle& G4Track::GetNextTouchableHandle() const
|
|
{
|
|
return fpNextTouchable;
|
|
}
|
|
|
|
inline void G4Track::SetNextTouchableHandle( const G4TouchableHandle& apValue)
|
|
{
|
|
fpNextTouchable = apValue;
|
|
}
|
|
|
|
|
|
// kinetic energy
|
|
inline G4double G4Track::GetKineticEnergy() const
|
|
{ return fpDynamicParticle->GetKineticEnergy(); }
|
|
|
|
inline void G4Track::SetKineticEnergy(const G4double aValue)
|
|
{ fpDynamicParticle->SetKineticEnergy(aValue); }
|
|
|
|
// total energy
|
|
inline G4double G4Track::GetTotalEnergy() const
|
|
{ return fpDynamicParticle->GetTotalEnergy(); }
|
|
|
|
// momentum
|
|
inline G4ThreeVector G4Track::GetMomentum() const
|
|
{ return fpDynamicParticle->GetMomentum(); }
|
|
|
|
// momentum (direction)
|
|
inline const G4ThreeVector& G4Track::GetMomentumDirection() const
|
|
{ return fpDynamicParticle->GetMomentumDirection(); }
|
|
|
|
inline void G4Track::SetMomentumDirection(const G4ThreeVector& aValue)
|
|
{ fpDynamicParticle->SetMomentumDirection(aValue) ;}
|
|
|
|
// polarization
|
|
inline const G4ThreeVector& G4Track::GetPolarization() const
|
|
{ return fpDynamicParticle->GetPolarization(); }
|
|
|
|
inline void G4Track::SetPolarization(const G4ThreeVector& aValue)
|
|
{ fpDynamicParticle->SetPolarization(aValue.x(),
|
|
aValue.y(),
|
|
aValue.z()); }
|
|
|
|
// track status
|
|
inline G4TrackStatus G4Track::GetTrackStatus() const
|
|
{ return fTrackStatus; }
|
|
|
|
inline void G4Track::SetTrackStatus(const G4TrackStatus aTrackStatus)
|
|
{ fTrackStatus = aTrackStatus; }
|
|
|
|
// track length
|
|
inline G4double G4Track::GetTrackLength() const
|
|
{ return fTrackLength; }
|
|
|
|
inline void G4Track::AddTrackLength(const G4double aValue)
|
|
{ fTrackLength += aValue; }
|
|
// Accumulated track length
|
|
|
|
// step number
|
|
inline G4int G4Track::GetCurrentStepNumber() const
|
|
{ return fCurrentStepNumber; }
|
|
|
|
inline void G4Track::IncrementCurrentStepNumber()
|
|
{ fCurrentStepNumber++; }
|
|
|
|
// step length
|
|
inline G4double G4Track::GetStepLength() const
|
|
{ return fStepLength; }
|
|
|
|
inline void G4Track::SetStepLength(G4double value)
|
|
{ fStepLength = value; }
|
|
|
|
// vertex (where this track was created) information
|
|
inline const G4ThreeVector& G4Track::GetVertexPosition() const
|
|
{ return fVtxPosition; }
|
|
|
|
inline void G4Track::SetVertexPosition(const G4ThreeVector& aValue)
|
|
{ fVtxPosition = aValue; }
|
|
|
|
inline const G4ThreeVector& G4Track::GetVertexMomentumDirection() const
|
|
{ return fVtxMomentumDirection; }
|
|
inline void G4Track::SetVertexMomentumDirection(const G4ThreeVector& aValue)
|
|
{ fVtxMomentumDirection = aValue ;}
|
|
|
|
inline G4double G4Track::GetVertexKineticEnergy() const
|
|
{ return fVtxKineticEnergy; }
|
|
|
|
inline void G4Track::SetVertexKineticEnergy(const G4double aValue)
|
|
{ fVtxKineticEnergy = aValue; }
|
|
|
|
inline G4LogicalVolume* G4Track::GetLogicalVolumeAtVertex() const
|
|
{ return fpLVAtVertex; }
|
|
|
|
inline void G4Track::SetLogicalVolumeAtVertex(G4LogicalVolume* aValue)
|
|
{ fpLVAtVertex = aValue; }
|
|
|
|
inline const G4VProcess* G4Track::GetCreatorProcess() const
|
|
{ return fpCreatorProcess; }
|
|
// If the pointer is 0, this means the track is created
|
|
// by the event generator, i.e. the primary track.If it is not
|
|
// 0, it points to the process which created this track.
|
|
|
|
inline void G4Track::SetCreatorProcess(G4VProcess* aValue)
|
|
{ fpCreatorProcess = aValue; }
|
|
|
|
// flag for "Below Threshold"
|
|
inline G4bool G4Track::IsBelowThreshold() const
|
|
{ return fBelowThreshold; }
|
|
|
|
inline void G4Track::SetBelowThresholdFlag(G4bool value)
|
|
{ fBelowThreshold = value; }
|
|
|
|
// flag for " Good for Tracking"
|
|
inline G4bool G4Track::IsGoodForTracking() const
|
|
{ return fGoodForTracking; }
|
|
|
|
inline void G4Track::SetGoodForTrackingFlag(G4bool value)
|
|
{ fGoodForTracking = value; }
|
|
|
|
// track weight
|
|
inline void G4Track::SetWeight(G4double aValue)
|
|
{ fWeight = aValue; }
|
|
|
|
inline G4double G4Track::GetWeight() const
|
|
{ return fWeight; }
|
|
|
|
// user information
|
|
inline G4VUserTrackInformation* G4Track::GetUserInformation() const
|
|
{ return fpUserInformation; }
|
|
inline void G4Track::SetUserInformation(G4VUserTrackInformation* aValue)
|
|
{ fpUserInformation = aValue; }
|
|
|
|
|
|
|
|
//-------------------------------------------------------------
|
|
// To implement bi-directional association between G4Step and
|
|
// and G4Track, a combined usage of 'forward declaration' and
|
|
// 'include' is necessary.
|
|
//-------------------------------------------------------------
|
|
#include "G4Step.hh"
|
|
|
|
inline const G4Step* G4Track::GetStep() const
|
|
{ return fpStep; }
|
|
|
|
inline void G4Track::SetStep(const G4Step* aValue)
|
|
{ fpStep = (aValue); }
|
|
|
|
|
|
|
|
|
|
|