Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Trajectory Class Reference

#include <G4Trajectory.hh>

+ Inheritance diagram for G4Trajectory:

Public Member Functions

 G4Trajectory ()
 
 G4Trajectory (const G4Track *aTrack)
 
 G4Trajectory (G4Trajectory &)
 
virtual ~G4Trajectory ()
 
void * operator 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 (G4int i_mode=0) 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
 
virtual G4int GetTrackID () const =0
 
virtual G4int GetParentID () const =0
 
virtual G4String GetParticleName () const =0
 
virtual G4double GetCharge () const =0
 
virtual G4int GetPDGEncoding () const =0
 
virtual G4ThreeVector GetInitialMomentum () const =0
 
virtual int GetPointEntries () const =0
 
virtual G4VTrajectoryPointGetPoint (G4int i) const =0
 
virtual void ShowTrajectory (std::ostream &os=G4cout) const
 
virtual void DrawTrajectory (G4int i_mode=0) const
 
virtual const std::map< G4String, G4AttDef > * GetAttDefs () const
 
virtual std::vector< G4AttValue > * CreateAttValues () const
 
virtual void AppendStep (const G4Step *aStep)=0
 
virtual void MergeTrajectory (G4VTrajectory *secondTrajectory)=0
 

Detailed Description

Definition at line 69 of file G4Trajectory.hh.

Constructor & Destructor Documentation

◆ G4Trajectory() [1/3]

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() [2/3]

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}
std::vector< G4VTrajectoryPoint * > TrajectoryPointContainer
Definition: G4Trajectory.hh:67
G4double GetPDGCharge() const
const G4String & GetParticleName() const
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
G4ThreeVector GetMomentum() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4int GetParentID() const

◆ G4Trajectory() [3/3]

G4Trajectory::G4Trajectory ( G4Trajectory right)

Definition at line 78 of file G4Trajectory.cc.

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}

◆ ~G4Trajectory()

G4Trajectory::~G4Trajectory ( )
virtual

Definition at line 96 of file G4Trajectory.cc.

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

Member Function Documentation

◆ AppendStep()

void G4Trajectory::AppendStep ( const G4Step aStep)
virtual

Implements G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 212 of file G4Trajectory.cc.

213{
214 positionRecord->push_back( new G4TrajectoryPoint(aStep->GetPostStepPoint()->
215 GetPosition() ));
216}
G4StepPoint * GetPostStepPoint() const

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 175 of file G4Trajectory.cc.

176{
177 std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
178
179 values->push_back
180 (G4AttValue("ID",G4UIcommand::ConvertToString(fTrackID),""));
181
182 values->push_back
183 (G4AttValue("PID",G4UIcommand::ConvertToString(fParentID),""));
184
185 values->push_back(G4AttValue("PN",ParticleName,""));
186
187 values->push_back
188 (G4AttValue("Ch",G4UIcommand::ConvertToString(PDGCharge),""));
189
190 values->push_back
191 (G4AttValue("PDG",G4UIcommand::ConvertToString(PDGEncoding),""));
192
193 values->push_back
194 (G4AttValue("IKE",G4BestUnit(initialKineticEnergy,"Energy"),""));
195
196 values->push_back
197 (G4AttValue("IMom",G4BestUnit(initialMomentum,"Energy"),""));
198
199 values->push_back
200 (G4AttValue("IMag",G4BestUnit(initialMomentum.mag(),"Energy"),""));
201
202 values->push_back
204
205#ifdef G4ATTDEBUG
206 G4cout << G4AttCheck(values,GetAttDefs());
207#endif
208
209 return values;
210}
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4DLLIMPORT std::ostream G4cout
double mag() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual int GetPointEntries() const
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:349

Referenced by G4RichTrajectory::CreateAttValues().

◆ DrawTrajectory()

void G4Trajectory::DrawTrajectory ( G4int  i_mode = 0) const
virtual

Reimplemented from G4VTrajectory.

Definition at line 125 of file G4Trajectory.cc.

126{
127 // Invoke the default implementation in G4VTrajectory...
129 // ... or override with your own code here.
130}
virtual void DrawTrajectory(G4int i_mode=0) const

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 132 of file G4Trajectory.cc.

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

Referenced by CreateAttValues(), G4RichTrajectory::GetAttDefs(), G4VisCommandList::SetNewValue(), and G4VisCommandSceneAddTrajectories::SetNewValue().

◆ GetCharge()

G4double G4Trajectory::GetCharge ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 98 of file G4Trajectory.hh.

99 { return PDGCharge; }

◆ GetInitialKineticEnergy()

G4double G4Trajectory::GetInitialKineticEnergy ( ) const
inline

Definition at line 102 of file G4Trajectory.hh.

103 { return initialKineticEnergy; }

◆ GetInitialMomentum()

G4ThreeVector G4Trajectory::GetInitialMomentum ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 104 of file G4Trajectory.hh.

105 { return initialMomentum; }

◆ GetParentID()

G4int G4Trajectory::GetParentID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 94 of file G4Trajectory.hh.

95 { return fParentID; }

◆ GetParticleDefinition()

G4ParticleDefinition * G4Trajectory::GetParticleDefinition ( )

Definition at line 218 of file G4Trajectory.cc.

219{
220 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
221}
static G4ParticleTable * GetParticleTable()

◆ GetParticleName()

G4String G4Trajectory::GetParticleName ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 96 of file G4Trajectory.hh.

97 { return ParticleName; }

◆ GetPDGEncoding()

G4int G4Trajectory::GetPDGEncoding ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 100 of file G4Trajectory.hh.

101 { return PDGEncoding; }

◆ GetPoint()

virtual G4VTrajectoryPoint * G4Trajectory::GetPoint ( G4int  i) const
inlinevirtual

Implements G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 113 of file G4Trajectory.hh.

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

◆ GetPointEntries()

virtual int G4Trajectory::GetPointEntries ( ) const
inlinevirtual

Implements G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 112 of file G4Trajectory.hh.

112{ return positionRecord->size(); }

Referenced by CreateAttValues(), and MergeTrajectory().

◆ GetTrackID()

G4int G4Trajectory::GetTrackID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 92 of file G4Trajectory.hh.

93 { return fTrackID; }

◆ MergeTrajectory()

void G4Trajectory::MergeTrajectory ( G4VTrajectory secondTrajectory)
virtual

Implements G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 223 of file G4Trajectory.cc.

224{
225 if(!secondTrajectory) return;
226
227 G4Trajectory* seco = (G4Trajectory*)secondTrajectory;
228 G4int ent = seco->GetPointEntries();
229 for(G4int i=1;i<ent;i++) // initial point of the second trajectory should not be merged
230 {
231 positionRecord->push_back((*(seco->positionRecord))[i]);
232 // positionRecord->push_back(seco->positionRecord->removeAt(1));
233 }
234 delete (*seco->positionRecord)[0];
235 seco->positionRecord->clear();
236}
int G4int
Definition: G4Types.hh:66

◆ operator delete()

void G4Trajectory::operator delete ( void *  aTrajectory)
inline

Definition at line 150 of file G4Trajectory.hh.

151{
152 aTrajectoryAllocator.FreeSingle((G4Trajectory*)aTrajectory);
153}
G4DLLIMPORT G4Allocator< G4Trajectory > aTrajectoryAllocator
Definition: G4Trajectory.cc:55

◆ operator new()

void * G4Trajectory::operator new ( size_t  )
inline

Definition at line 143 of file G4Trajectory.hh.

144{
145 void* aTrajectory;
146 aTrajectory = (void*)aTrajectoryAllocator.MallocSingle();
147 return aTrajectory;
148}

◆ operator==()

int G4Trajectory::operator== ( const G4Trajectory right) const
inline

Definition at line 88 of file G4Trajectory.hh.

89 {return (this==&right);}

◆ ShowTrajectory()

void G4Trajectory::ShowTrajectory ( std::ostream &  os = G4cout) const
virtual

Reimplemented from G4VTrajectory.

Definition at line 109 of file G4Trajectory.cc.

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

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