Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4Track Class Reference

#include <G4Track.hh>

Public Member Functions

 G4Track ()
 
 G4Track (G4DynamicParticle *apValueDynamicParticle, G4double aValueTime, const G4ThreeVector &aValuePosition)
 
 G4Track (const G4Track &)
 
 ~G4Track ()
 
voidoperator new (size_t)
 
void operator delete (void *aTrack)
 
G4bool operator== (const G4Track &)
 
void CopyTrackInfo (const G4Track &)
 
G4int GetTrackID () const
 
void SetTrackID (const G4int aValue)
 
G4int GetParentID () const
 
void SetParentID (const G4int aValue)
 
const G4DynamicParticleGetDynamicParticle () const
 
const G4ParticleDefinitionGetParticleDefinition () const
 
G4ParticleDefinitionGetDefinition () const
 
const G4ThreeVectorGetPosition () const
 
void SetPosition (const G4ThreeVector &aValue)
 
G4double GetGlobalTime () const
 
void SetGlobalTime (const G4double aValue)
 
G4double GetLocalTime () const
 
void SetLocalTime (const G4double aValue)
 
G4double GetProperTime () const
 
void SetProperTime (const G4double aValue)
 
G4VPhysicalVolumeGetVolume () const
 
G4VPhysicalVolumeGetNextVolume () const
 
G4MaterialGetMaterial () const
 
G4MaterialGetNextMaterial () const
 
const G4MaterialCutsCoupleGetMaterialCutsCouple () const
 
const G4MaterialCutsCoupleGetNextMaterialCutsCouple () const
 
const G4VTouchableGetTouchable () const
 
const G4TouchableHandleGetTouchableHandle () const
 
void SetTouchableHandle (const G4TouchableHandle &apValue)
 
const G4VTouchableGetNextTouchable () const
 
const G4TouchableHandleGetNextTouchableHandle () const
 
void SetNextTouchableHandle (const G4TouchableHandle &apValue)
 
const G4VTouchableGetOriginTouchable () const
 
const G4TouchableHandleGetOriginTouchableHandle () const
 
void SetOriginTouchableHandle (const G4TouchableHandle &apValue)
 
G4double GetKineticEnergy () const
 
void SetKineticEnergy (const G4double aValue)
 
G4double GetTotalEnergy () const
 
const G4ThreeVectorGetMomentumDirection () const
 
void SetMomentumDirection (const G4ThreeVector &aValue)
 
G4ThreeVector GetMomentum () const
 
G4double GetVelocity () const
 
void SetVelocity (G4double val)
 
G4double CalculateVelocity () const
 
G4double CalculateVelocityForOpticalPhoton () const
 
G4bool UseGivenVelocity () const
 
void UseGivenVelocity (G4bool val)
 
const G4ThreeVectorGetPolarization () const
 
void SetPolarization (const G4ThreeVector &aValue)
 
G4TrackStatus GetTrackStatus () const
 
void SetTrackStatus (const G4TrackStatus aTrackStatus)
 
G4bool IsBelowThreshold () const
 
void SetBelowThresholdFlag (G4bool value=true)
 
G4bool IsGoodForTracking () const
 
void SetGoodForTrackingFlag (G4bool value=true)
 
G4double GetTrackLength () const
 
void AddTrackLength (const G4double aValue)
 
const G4StepGetStep () const
 
void SetStep (const G4Step *aValue)
 
G4int GetCurrentStepNumber () const
 
void IncrementCurrentStepNumber ()
 
G4double GetStepLength () const
 
void SetStepLength (G4double value)
 
const G4ThreeVectorGetVertexPosition () const
 
void SetVertexPosition (const G4ThreeVector &aValue)
 
const G4ThreeVectorGetVertexMomentumDirection () const
 
void SetVertexMomentumDirection (const G4ThreeVector &aValue)
 
G4double GetVertexKineticEnergy () const
 
void SetVertexKineticEnergy (const G4double aValue)
 
const G4LogicalVolumeGetLogicalVolumeAtVertex () const
 
void SetLogicalVolumeAtVertex (const G4LogicalVolume *)
 
const G4VProcessGetCreatorProcess () const
 
void SetCreatorProcess (const G4VProcess *aValue)
 
void SetCreatorModelIndex (G4int idx)
 
const G4StringGetCreatorModelName () const
 
G4int GetCreatorModelID () const
 
G4double GetWeight () const
 
void SetWeight (G4double aValue)
 
G4VUserTrackInformationGetUserInformation () const
 
void SetUserInformation (G4VUserTrackInformation *aValue) const
 
void SetAuxiliaryTrackInformation (G4int idx, G4VAuxiliaryTrackInformation *info) const
 
G4VAuxiliaryTrackInformationGetAuxiliaryTrackInformation (G4int idx) const
 
std::map< G4int,
G4VAuxiliaryTrackInformation * > * 
GetAuxiliaryTrackInformationMap () const
 
void RemoveAuxiliaryTrackInformation (G4int idx)
 
void RemoveAuxiliaryTrackInformation (G4String &name)
 

Static Public Member Functions

static void SetVelocityTableProperties (G4double t_max, G4double t_min, G4int nbin)
 
static G4double GetMaxTOfVelocityTable ()
 
static G4double GetMinTOfVelocityTable ()
 
static G4int GetNbinOfVelocityTable ()
 

Detailed Description

Definition at line 76 of file G4Track.hh.

Constructor & Destructor Documentation

G4Track::G4Track ( )

Definition at line 98 of file G4Track.cc.

100  : fCurrentStepNumber(0),
101  fGlobalTime(0), fLocalTime(0.),
102  fTrackLength(0.),
103  fParentID(0), fTrackID(0),
104  fVelocity(c_light),
105  fpDynamicParticle(0),
106  fTrackStatus(fAlive),
107  fBelowThreshold(false), fGoodForTracking(false),
108  fStepLength(0.0), fWeight(1.0),
109  fpStep(0),
110  fVtxKineticEnergy(0.0),
111  fpLVAtVertex(0), fpCreatorProcess(0),
112  fCreatorModelIndex(-1),
113  fpUserInformation(0),
114  prev_mat(0), groupvel(0),
115  prev_velocity(0.0), prev_momentum(0.0),
116  is_OpticalPhoton(false),
117  useGivenVelocity(false),
118  fpAuxiliaryTrackInformationMap(0)
119 {
120 }
float c_light
Definition: hepunit.py:257
G4Track::G4Track ( G4DynamicParticle apValueDynamicParticle,
G4double  aValueTime,
const G4ThreeVector aValuePosition 
)

Definition at line 57 of file G4Track.cc.

61  : fCurrentStepNumber(0), fPosition(aValuePosition),
62  fGlobalTime(aValueTime), fLocalTime(0.),
63  fTrackLength(0.),
64  fParentID(0), fTrackID(0),
65  fVelocity(c_light),
66  fpDynamicParticle(apValueDynamicParticle),
67  fTrackStatus(fAlive),
68  fBelowThreshold(false), fGoodForTracking(false),
69  fStepLength(0.0), fWeight(1.0),
70  fpStep(0),
71  fVtxKineticEnergy(0.0),
72  fpLVAtVertex(0), fpCreatorProcess(0),
73  fCreatorModelIndex(-1),
74  fpUserInformation(0),
75  prev_mat(0), groupvel(0),
76  prev_velocity(0.0), prev_momentum(0.0),
77  is_OpticalPhoton(false),
78  useGivenVelocity(false),
79  fpAuxiliaryTrackInformationMap(0)
80 {
81  static G4ThreadLocal G4bool isFirstTime = true;
82  static G4ThreadLocal G4ParticleDefinition* fOpticalPhoton =0;
83  if ( isFirstTime ) {
84  isFirstTime = false;
85  // set fOpticalPhoton
86  fOpticalPhoton = G4ParticleTable::GetParticleTable()->FindParticle("opticalphoton");
87  }
88  // check if the particle type is Optical Photon
89  is_OpticalPhoton = (fpDynamicParticle->GetDefinition() == fOpticalPhoton);
90 
91  if (velTable ==0 ) velTable = G4VelocityTable::GetVelocityTable();
92 
93  fVelocity = CalculateVelocity();
94 
95 }
static G4VelocityTable * GetVelocityTable()
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
G4ParticleDefinition * GetDefinition() const
#define G4ThreadLocal
Definition: tls.hh:89
bool G4bool
Definition: G4Types.hh:79
G4double CalculateVelocity() const
Definition: G4Track.cc:222
static G4ParticleTable * GetParticleTable()
float c_light
Definition: hepunit.py:257

Here is the call graph for this function:

G4Track::G4Track ( const G4Track right)

Definition at line 123 of file G4Track.cc.

125  : fCurrentStepNumber(0),
126  fGlobalTime(0), fLocalTime(0.),
127  fTrackLength(0.),
128  fParentID(0), fTrackID(0),
129  fVelocity(c_light),
130  fpDynamicParticle(0),
131  fTrackStatus(fAlive),
132  fBelowThreshold(false), fGoodForTracking(false),
133  fStepLength(0.0), fWeight(1.0),
134  fpStep(0),
135  fVtxKineticEnergy(0.0),
136  fpLVAtVertex(0), fpCreatorProcess(0),
137  fCreatorModelIndex(-1),
138  fpUserInformation(0),
139  prev_mat(0), groupvel(0),
140  prev_velocity(0.0), prev_momentum(0.0),
141  is_OpticalPhoton(false),
142  useGivenVelocity(false),
143  fpAuxiliaryTrackInformationMap(0)
144 {
145  *this = right;
146 }
float c_light
Definition: hepunit.py:257
G4Track::~G4Track ( )

Definition at line 149 of file G4Track.cc.

151 {
152  delete fpDynamicParticle;
153  delete fpUserInformation;
154  ClearAuxiliaryTrackInformation();
155 }

Member Function Documentation

void G4Track::AddTrackLength ( const G4double  aValue)

Here is the caller graph for this function:

G4double G4Track::CalculateVelocity ( ) const

Definition at line 222 of file G4Track.cc.

224 {
225  if (useGivenVelocity) return fVelocity;
226 
227  G4double velocity = c_light ;
228 
229  G4double mass = fpDynamicParticle->GetMass();
230 
231  // special case for photons
232  if ( is_OpticalPhoton ) return CalculateVelocityForOpticalPhoton();
233 
234  // particles other than optical photon
235  if (mass<DBL_MIN) {
236  // Zero Mass
237  velocity = c_light;
238  } else {
239  G4double T = (fpDynamicParticle->GetKineticEnergy())/mass;
240  if (T > GetMaxTOfVelocityTable()) {
241  velocity = c_light;
242  } else if (T<DBL_MIN) {
243  velocity =0.;
244  } else if (T<GetMinTOfVelocityTable()) {
245  velocity = c_light*std::sqrt(T*(T+2.))/(T+1.0);
246  } else {
247  velocity = velTable->Value(T);
248  }
249 
250  }
251  return velocity ;
252 }
G4double GetKineticEnergy() const
G4double CalculateVelocityForOpticalPhoton() const
Definition: G4Track.cc:255
G4double GetMass() const
static G4double GetMinTOfVelocityTable()
Definition: G4Track.cc:314
#define DBL_MIN
Definition: templates.hh:75
static G4double GetMaxTOfVelocityTable()
Definition: G4Track.cc:309
double G4double
Definition: G4Types.hh:76
G4double Value(G4double theEnergy)
float c_light
Definition: hepunit.py:257

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4Track::CalculateVelocityForOpticalPhoton ( ) const

Definition at line 255 of file G4Track.cc.

257 {
258 
259  G4double velocity = c_light ;
260 
261 
262  G4Material* mat=0;
263  G4bool update_groupvel = false;
264  if ( fpStep !=0 ){
265  mat= this->GetMaterial(); // Fix for repeated volumes
266  }else{
267  if (fpTouchable!=0){
268  mat=fpTouchable->GetVolume()->GetLogicalVolume()->GetMaterial();
269  }
270  }
271  // check if previous step is in the same volume
272  // and get new GROUPVELOCITY table if necessary
273  if ((mat != 0) && ((mat != prev_mat)||(groupvel==0))) {
274  groupvel = 0;
275  if(mat->GetMaterialPropertiesTable() != 0)
276  groupvel = mat->GetMaterialPropertiesTable()->GetProperty("GROUPVEL");
277  update_groupvel = true;
278  }
279  prev_mat = mat;
280 
281  if (groupvel != 0 ) {
282  // light velocity = c/(rindex+d(rindex)/d(log(E_phot)))
283  // values stored in GROUPVEL material properties vector
284  velocity = prev_velocity;
285 
286  // check if momentum is same as in the previous step
287  // and calculate group velocity if necessary
288  G4double current_momentum = fpDynamicParticle->GetTotalMomentum();
289  if( update_groupvel || (current_momentum != prev_momentum) ) {
290  velocity =
291  groupvel->Value(current_momentum);
292  prev_velocity = velocity;
293  prev_momentum = current_momentum;
294  }
295  }
296 
297  return velocity ;
298 }
G4Material * GetMaterial() const
G4double GetTotalMomentum() const
bool G4bool
Definition: G4Types.hh:79
G4double Value(G4double theEnergy, size_t &lastidx) const
G4Material * GetMaterial() const
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
double G4double
Definition: G4Types.hh:76
G4MaterialPropertyVector * GetProperty(const char *key)
float c_light
Definition: hepunit.py:257

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Track::CopyTrackInfo ( const G4Track right)

Definition at line 215 of file G4Track.cc.

217 {
218  *this = right;
219 }
G4VAuxiliaryTrackInformation * G4Track::GetAuxiliaryTrackInformation ( G4int  idx) const

Definition at line 340 of file G4Track.cc.

342 {
343  if(!fpAuxiliaryTrackInformationMap) return 0;
344  std::map<G4int,G4VAuxiliaryTrackInformation*>::iterator itr
345  = fpAuxiliaryTrackInformationMap->find(idx);
346  if(itr==fpAuxiliaryTrackInformationMap->end()) return 0;
347  else return (*itr).second;
348 }
std::map<G4int,G4VAuxiliaryTrackInformation*>* G4Track::GetAuxiliaryTrackInformationMap ( ) const
inline

Definition at line 338 of file G4Track.hh.

339  { return fpAuxiliaryTrackInformationMap; }
G4int G4Track::GetCreatorModelID ( ) const
inline

Here is the caller graph for this function:

const G4String& G4Track::GetCreatorModelName ( ) const
inline
const G4VProcess* G4Track::GetCreatorProcess ( ) const

Here is the caller graph for this function:

G4int G4Track::GetCurrentStepNumber ( ) const
G4ParticleDefinition* G4Track::GetDefinition ( ) const
const G4DynamicParticle* G4Track::GetDynamicParticle ( ) const
G4double G4Track::GetGlobalTime ( ) const
G4double G4Track::GetKineticEnergy ( ) const
G4double G4Track::GetLocalTime ( ) const

Here is the caller graph for this function:

const G4LogicalVolume* G4Track::GetLogicalVolumeAtVertex ( ) const

Here is the caller graph for this function:

G4Material* G4Track::GetMaterial ( ) const
const G4MaterialCutsCouple* G4Track::GetMaterialCutsCouple ( ) const
G4double G4Track::GetMaxTOfVelocityTable ( )
static

Definition at line 309 of file G4Track.cc.

static G4double GetMaxTOfVelocityTable()

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4Track::GetMinTOfVelocityTable ( )
static

Definition at line 314 of file G4Track.cc.

static G4double GetMinTOfVelocityTable()

Here is the call graph for this function:

Here is the caller graph for this function:

G4ThreeVector G4Track::GetMomentum ( ) const

Here is the caller graph for this function:

const G4ThreeVector& G4Track::GetMomentumDirection ( ) const
G4int G4Track::GetNbinOfVelocityTable ( )
static

Definition at line 319 of file G4Track.cc.

static G4int GetNbinOfVelocityTable()

Here is the call graph for this function:

G4Material* G4Track::GetNextMaterial ( ) const
const G4MaterialCutsCouple* G4Track::GetNextMaterialCutsCouple ( ) const
const G4VTouchable* G4Track::GetNextTouchable ( ) const
const G4TouchableHandle& G4Track::GetNextTouchableHandle ( ) const

Here is the caller graph for this function:

G4VPhysicalVolume* G4Track::GetNextVolume ( ) const

Here is the caller graph for this function:

const G4VTouchable* G4Track::GetOriginTouchable ( ) const
const G4TouchableHandle& G4Track::GetOriginTouchableHandle ( ) const
G4int G4Track::GetParentID ( ) const
const G4ParticleDefinition* G4Track::GetParticleDefinition ( ) const

Here is the caller graph for this function:

const G4ThreeVector& G4Track::GetPolarization ( ) const

Here is the caller graph for this function:

const G4ThreeVector& G4Track::GetPosition ( ) const
G4double G4Track::GetProperTime ( ) const

Here is the caller graph for this function:

const G4Step* G4Track::GetStep ( ) const

Here is the caller graph for this function:

G4double G4Track::GetStepLength ( ) const

Here is the caller graph for this function:

G4double G4Track::GetTotalEnergy ( ) const

Here is the caller graph for this function:

const G4VTouchable* G4Track::GetTouchable ( ) const

Here is the caller graph for this function:

const G4TouchableHandle& G4Track::GetTouchableHandle ( ) const

Here is the caller graph for this function:

G4int G4Track::GetTrackID ( ) const
G4double G4Track::GetTrackLength ( ) const

Here is the caller graph for this function:

G4TrackStatus G4Track::GetTrackStatus ( ) const
G4VUserTrackInformation* G4Track::GetUserInformation ( ) const

Here is the caller graph for this function:

G4double G4Track::GetVelocity ( ) const

Here is the caller graph for this function:

G4double G4Track::GetVertexKineticEnergy ( ) const

Here is the caller graph for this function:

const G4ThreeVector& G4Track::GetVertexMomentumDirection ( ) const

Here is the caller graph for this function:

const G4ThreeVector& G4Track::GetVertexPosition ( ) const

Here is the caller graph for this function:

G4VPhysicalVolume* G4Track::GetVolume ( ) const
G4double G4Track::GetWeight ( ) const

Here is the caller graph for this function:

void G4Track::IncrementCurrentStepNumber ( )

Here is the caller graph for this function:

G4bool G4Track::IsBelowThreshold ( ) const
G4bool G4Track::IsGoodForTracking ( ) const

Here is the caller graph for this function:

void G4Track::operator delete ( void aTrack)
inline
void* G4Track::operator new ( size_t  )
inline
G4bool G4Track::operator== ( const G4Track )
void G4Track::RemoveAuxiliaryTrackInformation ( G4int  idx)

Definition at line 351 of file G4Track.cc.

353 {
354  if(fpAuxiliaryTrackInformationMap && idx>=0 && idx<G4PhysicsModelCatalog::Entries())
355  fpAuxiliaryTrackInformationMap->erase(idx);
356 }

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Track::RemoveAuxiliaryTrackInformation ( G4String name)

Definition at line 359 of file G4Track.cc.

361 {
362  if(fpAuxiliaryTrackInformationMap)
363  {
366  }
367 }
void RemoveAuxiliaryTrackInformation(G4int idx)
Definition: G4Track.cc:351
int G4int
Definition: G4Types.hh:78
static G4int GetIndex(const G4String &)

Here is the call graph for this function:

void G4Track::SetAuxiliaryTrackInformation ( G4int  idx,
G4VAuxiliaryTrackInformation info 
) const

Definition at line 324 of file G4Track.cc.

326 {
327  if(!fpAuxiliaryTrackInformationMap)
328  { fpAuxiliaryTrackInformationMap = new std::map<G4int,G4VAuxiliaryTrackInformation*>; }
329  if(idx<0 || idx>=G4PhysicsModelCatalog::Entries())
330  {
332  ED << "Process/model index <" << idx << "> is invalid.";
333  G4Exception("G4VAuxiliaryTrackInformation::G4VAuxiliaryTrackInformation()","TRACK0982",
334  FatalException, ED);
335  }
336  (*fpAuxiliaryTrackInformationMap)[idx] = info;
337 }
const XML_Char XML_Encoding * info
Definition: expat.h:530
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Track::SetBelowThresholdFlag ( G4bool  value = true)
void G4Track::SetCreatorModelIndex ( G4int  idx)
inline

Here is the caller graph for this function:

void G4Track::SetCreatorProcess ( const G4VProcess aValue)

Here is the caller graph for this function:

void G4Track::SetGlobalTime ( const G4double  aValue)
void G4Track::SetGoodForTrackingFlag ( G4bool  value = true)

Here is the caller graph for this function:

void G4Track::SetKineticEnergy ( const G4double  aValue)

Here is the caller graph for this function:

void G4Track::SetLocalTime ( const G4double  aValue)
void G4Track::SetLogicalVolumeAtVertex ( const G4LogicalVolume )

Here is the caller graph for this function:

void G4Track::SetMomentumDirection ( const G4ThreeVector aValue)

Here is the caller graph for this function:

void G4Track::SetNextTouchableHandle ( const G4TouchableHandle apValue)

Here is the caller graph for this function:

void G4Track::SetOriginTouchableHandle ( const G4TouchableHandle apValue)

Here is the caller graph for this function:

void G4Track::SetParentID ( const G4int  aValue)

Here is the caller graph for this function:

void G4Track::SetPolarization ( const G4ThreeVector aValue)
void G4Track::SetPosition ( const G4ThreeVector aValue)
void G4Track::SetProperTime ( const G4double  aValue)
void G4Track::SetStep ( const G4Step aValue)

Here is the caller graph for this function:

void G4Track::SetStepLength ( G4double  value)

Here is the caller graph for this function:

void G4Track::SetTouchableHandle ( const G4TouchableHandle apValue)

Here is the caller graph for this function:

void G4Track::SetTrackID ( const G4int  aValue)

Here is the caller graph for this function:

void G4Track::SetTrackStatus ( const G4TrackStatus  aTrackStatus)

Here is the caller graph for this function:

void G4Track::SetUserInformation ( G4VUserTrackInformation aValue) const

Here is the caller graph for this function:

void G4Track::SetVelocity ( G4double  val)

Here is the caller graph for this function:

void G4Track::SetVelocityTableProperties ( G4double  t_max,
G4double  t_min,
G4int  nbin 
)
static

Definition at line 301 of file G4Track.cc.

303 {
306 }
static G4VelocityTable * GetVelocityTable()
static void SetVelocityTableProperties(G4double t_max, G4double t_min, G4int nbin)

Here is the call graph for this function:

void G4Track::SetVertexKineticEnergy ( const G4double  aValue)

Here is the caller graph for this function:

void G4Track::SetVertexMomentumDirection ( const G4ThreeVector aValue)

Here is the caller graph for this function:

void G4Track::SetVertexPosition ( const G4ThreeVector aValue)

Here is the caller graph for this function:

void G4Track::SetWeight ( G4double  aValue)

Here is the caller graph for this function:

G4bool G4Track::UseGivenVelocity ( ) const

Here is the caller graph for this function:

void G4Track::UseGivenVelocity ( G4bool  val)

The documentation for this class was generated from the following files: