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

#include <G4Trajectory.hh>

Inheritance diagram for G4Trajectory:
Collaboration diagram for G4Trajectory:

Public Member Functions

 G4Trajectory ()
 
 G4Trajectory (const G4Track *aTrack)
 
 G4Trajectory (G4Trajectory &)
 
virtual ~G4Trajectory ()
 
voidoperator new (size_t)
 
void operator delete (void *)
 
int operator== (const G4Trajectory &right) const
 
G4int GetTrackID () const
 
G4int GetParentID () const
 
G4String GetParticleName () const
 
G4double GetCharge () const
 
G4int GetPDGEncoding () const
 
G4double GetInitialKineticEnergy () const
 
G4ThreeVector GetInitialMomentum () const
 
virtual void ShowTrajectory (std::ostream &os=G4cout) const
 
virtual void DrawTrajectory () const
 
virtual void AppendStep (const G4Step *aStep)
 
virtual int GetPointEntries () const
 
virtual G4VTrajectoryPointGetPoint (G4int i) const
 
virtual void MergeTrajectory (G4VTrajectory *secondTrajectory)
 
G4ParticleDefinitionGetParticleDefinition ()
 
virtual const std::map
< G4String, G4AttDef > * 
GetAttDefs () const
 
virtual std::vector< G4AttValue > * CreateAttValues () const
 
- Public Member Functions inherited from G4VTrajectory
 G4VTrajectory ()
 
virtual ~G4VTrajectory ()
 
G4bool operator== (const G4VTrajectory &right) const
 

Detailed Description

Definition at line 72 of file G4Trajectory.hh.

Constructor & Destructor Documentation

G4Trajectory::G4Trajectory ( )

Definition at line 57 of file G4Trajectory.cc.

58 : positionRecord(0), fTrackID(0), fParentID(0),
59  PDGEncoding( 0 ), PDGCharge(0.0), ParticleName(""),
60  initialKineticEnergy( 0. ), initialMomentum( G4ThreeVector() )
61 {;}
CLHEP::Hep3Vector G4ThreeVector
G4Trajectory::G4Trajectory ( const G4Track aTrack)

Definition at line 63 of file G4Trajectory.cc.

64 {
65  G4ParticleDefinition * fpParticleDefinition = aTrack->GetDefinition();
66  ParticleName = fpParticleDefinition->GetParticleName();
67  PDGCharge = fpParticleDefinition->GetPDGCharge();
68  PDGEncoding = fpParticleDefinition->GetPDGEncoding();
69  fTrackID = aTrack->GetTrackID();
70  fParentID = aTrack->GetParentID();
71  initialKineticEnergy = aTrack->GetKineticEnergy();
72  initialMomentum = aTrack->GetMomentum();
73  positionRecord = new TrajectoryPointContainer();
74  // Following is for the first trajectory point
75  positionRecord->push_back(new G4TrajectoryPoint(aTrack->GetPosition()));
76 }
G4ParticleDefinition * GetDefinition() const
G4int GetParentID() const
std::vector< G4VTrajectoryPoint * > TrajectoryPointContainer
const G4ThreeVector & GetPosition() const
const G4String & GetParticleName() const
G4double GetKineticEnergy() const
G4int GetTrackID() const
G4ThreeVector GetMomentum() const
G4double GetPDGCharge() const

Here is the call graph for this function:

G4Trajectory::G4Trajectory ( G4Trajectory right)

Definition at line 78 of file G4Trajectory.cc.

78  :G4VTrajectory()
79 {
80  ParticleName = right.ParticleName;
81  PDGCharge = right.PDGCharge;
82  PDGEncoding = right.PDGEncoding;
83  fTrackID = right.fTrackID;
84  fParentID = right.fParentID;
85  initialKineticEnergy = right.initialKineticEnergy;
86  initialMomentum = right.initialMomentum;
87  positionRecord = new TrajectoryPointContainer();
88 
89  for(size_t i=0;i<right.positionRecord->size();i++)
90  {
91  G4TrajectoryPoint* rightPoint = (G4TrajectoryPoint*)((*(right.positionRecord))[i]);
92  positionRecord->push_back(new G4TrajectoryPoint(*rightPoint));
93  }
94 }
std::vector< G4VTrajectoryPoint * > TrajectoryPointContainer
G4Trajectory::~G4Trajectory ( )
virtual

Definition at line 96 of file G4Trajectory.cc.

97 {
98  if (positionRecord) {
99  size_t i;
100  for(i=0;i<positionRecord->size();i++){
101  delete (*positionRecord)[i];
102  }
103  positionRecord->clear();
104  delete positionRecord;
105  }
106 }

Member Function Documentation

void G4Trajectory::AppendStep ( const G4Step aStep)
virtual

Implements G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 202 of file G4Trajectory.cc.

203 {
204  positionRecord->push_back( new G4TrajectoryPoint(aStep->GetPostStepPoint()->
205  GetPosition() ));
206 }
G4StepPoint * GetPostStepPoint() const

Here is the call graph for this function:

std::vector< G4AttValue > * G4Trajectory::CreateAttValues ( ) const
virtual

Reimplemented from G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 165 of file G4Trajectory.cc.

166 {
167  std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
168 
169  values->push_back
170  (G4AttValue("ID",G4UIcommand::ConvertToString(fTrackID),""));
171 
172  values->push_back
173  (G4AttValue("PID",G4UIcommand::ConvertToString(fParentID),""));
174 
175  values->push_back(G4AttValue("PN",ParticleName,""));
176 
177  values->push_back
178  (G4AttValue("Ch",G4UIcommand::ConvertToString(PDGCharge),""));
179 
180  values->push_back
181  (G4AttValue("PDG",G4UIcommand::ConvertToString(PDGEncoding),""));
182 
183  values->push_back
184  (G4AttValue("IKE",G4BestUnit(initialKineticEnergy,"Energy"),""));
185 
186  values->push_back
187  (G4AttValue("IMom",G4BestUnit(initialMomentum,"Energy"),""));
188 
189  values->push_back
190  (G4AttValue("IMag",G4BestUnit(initialMomentum.mag(),"Energy"),""));
191 
192  values->push_back
194 
195 #ifdef G4ATTDEBUG
196  G4cout << G4AttCheck(values,GetAttDefs());
197 #endif
198 
199  return values;
200 }
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:372
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual int GetPointEntries() const
G4GLOB_DLL std::ostream G4cout
double mag() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Trajectory::DrawTrajectory ( ) const
virtual

Reimplemented from G4VTrajectory.

Reimplemented in G4RichTrajectory, and LXeTrajectory.

Definition at line 115 of file G4Trajectory.cc.

116 {
117  // Invoke the default implementation in G4VTrajectory...
119  // ... or override with your own code here.
120 }
virtual void DrawTrajectory() const

Here is the call graph for this function:

const std::map< G4String, G4AttDef > * G4Trajectory::GetAttDefs ( ) const
virtual

Reimplemented from G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 122 of file G4Trajectory.cc.

123 {
124  G4bool isNew;
125  std::map<G4String,G4AttDef>* store
126  = G4AttDefStore::GetInstance("G4Trajectory",isNew);
127  if (isNew) {
128 
129  G4String ID("ID");
130  (*store)[ID] = G4AttDef(ID,"Track ID","Physics","","G4int");
131 
132  G4String PID("PID");
133  (*store)[PID] = G4AttDef(PID,"Parent ID","Physics","","G4int");
134 
135  G4String PN("PN");
136  (*store)[PN] = G4AttDef(PN,"Particle Name","Physics","","G4String");
137 
138  G4String Ch("Ch");
139  (*store)[Ch] = G4AttDef(Ch,"Charge","Physics","e+","G4double");
140 
141  G4String PDG("PDG");
142  (*store)[PDG] = G4AttDef(PDG,"PDG Encoding","Physics","","G4int");
143 
144  G4String IKE("IKE");
145  (*store)[IKE] =
146  G4AttDef(IKE, "Initial kinetic energy",
147  "Physics","G4BestUnit","G4double");
148 
149  G4String IMom("IMom");
150  (*store)[IMom] = G4AttDef(IMom, "Initial momentum",
151  "Physics","G4BestUnit","G4ThreeVector");
152 
153  G4String IMag("IMag");
154  (*store)[IMag] =
155  G4AttDef(IMag, "Initial momentum magnitude",
156  "Physics","G4BestUnit","G4double");
157 
158  G4String NTP("NTP");
159  (*store)[NTP] = G4AttDef(NTP,"No. of points","Physics","","G4int");
160 
161  }
162  return store;
163 }
bool G4bool
Definition: G4Types.hh:79
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4Trajectory::GetCharge ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 101 of file G4Trajectory.hh.

102  { return PDGCharge; }
G4double G4Trajectory::GetInitialKineticEnergy ( ) const
inline

Definition at line 105 of file G4Trajectory.hh.

106  { return initialKineticEnergy; }
G4ThreeVector G4Trajectory::GetInitialMomentum ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 107 of file G4Trajectory.hh.

108  { return initialMomentum; }
G4int G4Trajectory::GetParentID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 97 of file G4Trajectory.hh.

98  { return fParentID; }
G4ParticleDefinition * G4Trajectory::GetParticleDefinition ( void  )

Definition at line 208 of file G4Trajectory.cc.

209 {
210  return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
211 }
static G4ParticleTable * GetParticleTable()

Here is the call graph for this function:

G4String G4Trajectory::GetParticleName ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 99 of file G4Trajectory.hh.

100  { return ParticleName; }

Here is the caller graph for this function:

G4int G4Trajectory::GetPDGEncoding ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 103 of file G4Trajectory.hh.

104  { return PDGEncoding; }
virtual G4VTrajectoryPoint* G4Trajectory::GetPoint ( G4int  i) const
inlinevirtual

Implements G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 115 of file G4Trajectory.hh.

116  { return (*positionRecord)[i]; }

Here is the caller graph for this function:

virtual int G4Trajectory::GetPointEntries ( ) const
inlinevirtual

Implements G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 114 of file G4Trajectory.hh.

114 { return positionRecord->size(); }

Here is the caller graph for this function:

G4int G4Trajectory::GetTrackID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 95 of file G4Trajectory.hh.

96  { return fTrackID; }
void G4Trajectory::MergeTrajectory ( G4VTrajectory secondTrajectory)
virtual

Implements G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 213 of file G4Trajectory.cc.

214 {
215  if(!secondTrajectory) return;
216 
217  G4Trajectory* seco = (G4Trajectory*)secondTrajectory;
218  G4int ent = seco->GetPointEntries();
219  for(G4int i=1;i<ent;i++) // initial point of the second trajectory should not be merged
220  {
221  positionRecord->push_back((*(seco->positionRecord))[i]);
222  }
223  delete (*seco->positionRecord)[0];
224  seco->positionRecord->clear();
225 }
int G4int
Definition: G4Types.hh:78
virtual int GetPointEntries() const

Here is the call graph for this function:

void G4Trajectory::operator delete ( void aTrajectory)
inline

Definition at line 149 of file G4Trajectory.hh.

150 {
151  aTrajectoryAllocator->FreeSingle((G4Trajectory*)aTrajectory);
152 }
G4TRACKING_DLL G4ThreadLocal G4Allocator< G4Trajectory > * aTrajectoryAllocator
Definition: G4Trajectory.cc:55
void * G4Trajectory::operator new ( size_t  )
inline

Definition at line 142 of file G4Trajectory.hh.

143 {
146  return (void*)aTrajectoryAllocator->MallocSingle();
147 }
G4TRACKING_DLL G4ThreadLocal G4Allocator< G4Trajectory > * aTrajectoryAllocator
Definition: G4Trajectory.cc:55
int G4Trajectory::operator== ( const G4Trajectory right) const
inline

Definition at line 91 of file G4Trajectory.hh.

92  {return (this==&right);}
void G4Trajectory::ShowTrajectory ( std::ostream &  os = G4cout) const
virtual

Reimplemented from G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 108 of file G4Trajectory.cc.

109 {
110  // Invoke the default implementation in G4VTrajectory...
112  // ... or override with your own code here.
113 }
virtual void ShowTrajectory(std::ostream &os=G4cout) const

Here is the call graph for this function:


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