Geant4 11.1.1
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 *)
 
G4int operator== (const G4Trajectory &r) 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 G4int 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 G4int GetPointEntries () const =0
 
virtual G4VTrajectoryPointGetPoint (G4int i) const =0
 
virtual void ShowTrajectory (std::ostream &os=G4cout) const
 
virtual void DrawTrajectory () 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 61 of file G4Trajectory.hh.

Constructor & Destructor Documentation

◆ G4Trajectory() [1/3]

G4Trajectory::G4Trajectory ( )

Definition at line 55 of file G4Trajectory.cc.

56 : initialMomentum( G4ThreeVector() )
57{
58}
CLHEP::Hep3Vector G4ThreeVector

◆ G4Trajectory() [2/3]

G4Trajectory::G4Trajectory ( const G4Track aTrack)

Definition at line 60 of file G4Trajectory.cc.

61{
62 G4ParticleDefinition* fpParticleDefinition = aTrack->GetDefinition();
63 ParticleName = fpParticleDefinition->GetParticleName();
64 PDGCharge = fpParticleDefinition->GetPDGCharge();
65 PDGEncoding = fpParticleDefinition->GetPDGEncoding();
66 fTrackID = aTrack->GetTrackID();
67 fParentID = aTrack->GetParentID();
68 initialKineticEnergy = aTrack->GetKineticEnergy();
69 initialMomentum = aTrack->GetMomentum();
70 positionRecord = new G4TrajectoryPointContainer();
71
72 // Following is for the first trajectory point
73 positionRecord->push_back(new G4TrajectoryPoint(aTrack->GetPosition()));
74}
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 76 of file G4Trajectory.cc.

78{
79 ParticleName = right.ParticleName;
80 PDGCharge = right.PDGCharge;
81 PDGEncoding = right.PDGEncoding;
82 fTrackID = right.fTrackID;
83 fParentID = right.fParentID;
84 initialKineticEnergy = right.initialKineticEnergy;
85 initialMomentum = right.initialMomentum;
86 positionRecord = new G4TrajectoryPointContainer();
87
88 for(std::size_t i=0; i<right.positionRecord->size(); ++i)
89 {
90 G4TrajectoryPoint* rightPoint
91 = (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 != nullptr)
99 {
100 for(std::size_t i=0; i<positionRecord->size(); ++i)
101 {
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 205 of file G4Trajectory.cc.

206{
207 positionRecord->push_back( new G4TrajectoryPoint(aStep->GetPostStepPoint()->
208 GetPosition()) );
209}
G4StepPoint * GetPostStepPoint() const

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 168 of file G4Trajectory.cc.

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

Referenced by G4RichTrajectory::CreateAttValues().

◆ DrawTrajectory()

void G4Trajectory::DrawTrajectory ( ) const
virtual

Reimplemented from G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 117 of file G4Trajectory.cc.

118{
119 // Invoke the default implementation in G4VTrajectory...
121
122 // ... or override with your own code here.
123}
virtual void DrawTrajectory() const

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 125 of file G4Trajectory.cc.

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

90 { return PDGCharge; }

◆ GetInitialKineticEnergy()

G4double G4Trajectory::GetInitialKineticEnergy ( ) const
inline

Definition at line 93 of file G4Trajectory.hh.

94 { return initialKineticEnergy; }

◆ GetInitialMomentum()

G4ThreeVector G4Trajectory::GetInitialMomentum ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 95 of file G4Trajectory.hh.

96 { return initialMomentum; }

◆ GetParentID()

G4int G4Trajectory::GetParentID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 85 of file G4Trajectory.hh.

86 { return fParentID; }

◆ GetParticleDefinition()

G4ParticleDefinition * G4Trajectory::GetParticleDefinition ( )

Definition at line 211 of file G4Trajectory.cc.

212{
213 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
214}
static G4ParticleTable * GetParticleTable()

◆ GetParticleName()

G4String G4Trajectory::GetParticleName ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 87 of file G4Trajectory.hh.

88 { return ParticleName; }

◆ GetPDGEncoding()

G4int G4Trajectory::GetPDGEncoding ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 91 of file G4Trajectory.hh.

92 { return PDGEncoding; }

◆ GetPoint()

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

Implements G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 105 of file G4Trajectory.hh.

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

◆ GetPointEntries()

virtual G4int G4Trajectory::GetPointEntries ( ) const
inlinevirtual

Implements G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 103 of file G4Trajectory.hh.

104 { return G4int(positionRecord->size()); }
int G4int
Definition: G4Types.hh:85

Referenced by CreateAttValues(), and MergeTrajectory().

◆ GetTrackID()

G4int G4Trajectory::GetTrackID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 83 of file G4Trajectory.hh.

84 { return fTrackID; }

◆ MergeTrajectory()

void G4Trajectory::MergeTrajectory ( G4VTrajectory secondTrajectory)
virtual

Implements G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 216 of file G4Trajectory.cc.

217{
218 if(secondTrajectory == nullptr) return;
219
220 G4Trajectory* seco = (G4Trajectory*)secondTrajectory;
221 G4int ent = seco->GetPointEntries();
222 for(G4int i=1; i<ent; ++i) // initial pt of 2nd trajectory shouldn't be merged
223 {
224 positionRecord->push_back((*(seco->positionRecord))[i]);
225 }
226 delete (*seco->positionRecord)[0];
227 seco->positionRecord->clear();
228}

◆ operator delete()

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

Definition at line 137 of file G4Trajectory.hh.

138{
139 aTrajectoryAllocator()->FreeSingle((G4Trajectory*)aTrajectory);
140}
G4TRACKING_DLL G4Allocator< G4Trajectory > *& aTrajectoryAllocator()
Definition: G4Trajectory.cc:49

◆ operator new()

void * G4Trajectory::operator new ( size_t  )
inline

Definition at line 128 of file G4Trajectory.hh.

129{
130 if (aTrajectoryAllocator() == nullptr)
131 {
133 }
134 return (void*)aTrajectoryAllocator()->MallocSingle();
135}

◆ operator==()

G4int G4Trajectory::operator== ( const G4Trajectory r) const
inline

Definition at line 142 of file G4Trajectory.hh.

143{
144 return (this==&r);
145}

◆ ShowTrajectory()

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

Reimplemented from G4VTrajectory.

Reimplemented in G4RichTrajectory.

Definition at line 109 of file G4Trajectory.cc.

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

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