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

Detailed Description

Definition at line 62 of file G4Trajectory.hh.

Constructor & Destructor Documentation

◆ G4Trajectory() [1/3]

G4Trajectory::G4Trajectory ( )
default

◆ G4Trajectory() [2/3]

G4Trajectory::G4Trajectory ( const G4Track * aTrack)

Definition at line 56 of file G4Trajectory.cc.

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

73{
74 ParticleName = right.ParticleName;
75 PDGCharge = right.PDGCharge;
76 PDGEncoding = right.PDGEncoding;
77 fTrackID = right.fTrackID;
78 fParentID = right.fParentID;
79 initialKineticEnergy = right.initialKineticEnergy;
80 initialMomentum = right.initialMomentum;
81 positionRecord = new G4TrajectoryPointContainer();
82
83 for (auto& i : *right.positionRecord) {
84 auto rightPoint = (G4TrajectoryPoint*)i;
85 positionRecord->push_back(new G4TrajectoryPoint(*rightPoint));
86 }
87}

◆ ~G4Trajectory()

G4Trajectory::~G4Trajectory ( )
override

Definition at line 89 of file G4Trajectory.cc.

90{
91 if (positionRecord != nullptr) {
92 for (auto& i : *positionRecord) {
93 delete i;
94 }
95 positionRecord->clear();
96 delete positionRecord;
97 }
98}

Member Function Documentation

◆ AppendStep()

void G4Trajectory::AppendStep ( const G4Step * aStep)
overridevirtual

Implements G4VTrajectory.

Definition at line 181 of file G4Trajectory.cc.

182{
183 positionRecord->push_back(new G4TrajectoryPoint(aStep->GetPostStepPoint()->GetPosition()));
184}
const G4ThreeVector & GetPosition() const
G4StepPoint * GetPostStepPoint() const

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Definition at line 152 of file G4Trajectory.cc.

153{
154 auto values = new std::vector<G4AttValue>;
155
156 values->push_back(G4AttValue("ID", G4UIcommand::ConvertToString(fTrackID), ""));
157
158 values->push_back(G4AttValue("PID", G4UIcommand::ConvertToString(fParentID), ""));
159
160 values->push_back(G4AttValue("PN", ParticleName, ""));
161
162 values->push_back(G4AttValue("Ch", G4UIcommand::ConvertToString(PDGCharge), ""));
163
164 values->push_back(G4AttValue("PDG", G4UIcommand::ConvertToString(PDGEncoding), ""));
165
166 values->push_back(G4AttValue("IKE", G4BestUnit(initialKineticEnergy, "Energy"), ""));
167
168 values->push_back(G4AttValue("IMom", G4BestUnit(initialMomentum, "Energy"), ""));
169
170 values->push_back(G4AttValue("IMag", G4BestUnit(initialMomentum.mag(), "Energy"), ""));
171
172 values->push_back(G4AttValue("NTP", G4UIcommand::ConvertToString(GetPointEntries()), ""));
173
174#ifdef G4ATTDEBUG
175 G4cout << G4AttCheck(values, GetAttDefs());
176#endif
177
178 return values;
179}
#define G4BestUnit(a, b)
G4GLOB_DLL std::ostream G4cout
double mag() const
G4int GetPointEntries() const override
const std::map< G4String, G4AttDef > * GetAttDefs() const override
static G4String ConvertToString(G4bool boolVal)

Referenced by G4RichTrajectory::CreateAttValues().

◆ DrawTrajectory()

void G4Trajectory::DrawTrajectory ( ) const
overridevirtual

Reimplemented from G4VTrajectory.

Definition at line 108 of file G4Trajectory.cc.

109{
110 // Invoke the default implementation in G4VTrajectory...
112
113 // ... or override with your own code here.
114}
virtual void DrawTrajectory() const

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Definition at line 116 of file G4Trajectory.cc.

117{
118 G4bool isNew;
119 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("G4Trajectory", isNew);
120 if (isNew) {
121 G4String ID("ID");
122 (*store)[ID] = G4AttDef(ID, "Track ID", "Physics", "", "G4int");
123
124 G4String PID("PID");
125 (*store)[PID] = G4AttDef(PID, "Parent ID", "Physics", "", "G4int");
126
127 G4String PN("PN");
128 (*store)[PN] = G4AttDef(PN, "Particle Name", "Physics", "", "G4String");
129
130 G4String Ch("Ch");
131 (*store)[Ch] = G4AttDef(Ch, "Charge", "Physics", "e+", "G4double");
132
133 G4String PDG("PDG");
134 (*store)[PDG] = G4AttDef(PDG, "PDG Encoding", "Physics", "", "G4int");
135
136 G4String IKE("IKE");
137 (*store)[IKE] = G4AttDef(IKE, "Initial kinetic energy", "Physics", "G4BestUnit", "G4double");
138
139 G4String IMom("IMom");
140 (*store)[IMom] = G4AttDef(IMom, "Initial momentum", "Physics", "G4BestUnit", "G4ThreeVector");
141
142 G4String IMag("IMag");
143 (*store)[IMag] =
144 G4AttDef(IMag, "Initial momentum magnitude", "Physics", "G4BestUnit", "G4double");
145
146 G4String NTP("NTP");
147 (*store)[NTP] = G4AttDef(NTP, "No. of points", "Physics", "", "G4int");
148 }
149 return store;
150}
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
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 85 of file G4Trajectory.hh.

85{ return PDGCharge; }

◆ GetInitialKineticEnergy()

G4double G4Trajectory::GetInitialKineticEnergy ( ) const
inline

Definition at line 87 of file G4Trajectory.hh.

87{ return initialKineticEnergy; }

◆ GetInitialMomentum()

G4ThreeVector G4Trajectory::GetInitialMomentum ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 88 of file G4Trajectory.hh.

88{ return initialMomentum; }

◆ GetParentID()

G4int G4Trajectory::GetParentID ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 83 of file G4Trajectory.hh.

83{ return fParentID; }

◆ GetParticleDefinition()

G4ParticleDefinition * G4Trajectory::GetParticleDefinition ( )

Definition at line 186 of file G4Trajectory.cc.

187{
188 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
189}
static G4ParticleTable * GetParticleTable()

◆ GetParticleName()

G4String G4Trajectory::GetParticleName ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 84 of file G4Trajectory.hh.

84{ return ParticleName; }

◆ GetPDGEncoding()

G4int G4Trajectory::GetPDGEncoding ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 86 of file G4Trajectory.hh.

86{ return PDGEncoding; }

◆ GetPoint()

G4VTrajectoryPoint * G4Trajectory::GetPoint ( G4int i) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 96 of file G4Trajectory.hh.

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

◆ GetPointEntries()

G4int G4Trajectory::GetPointEntries ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 95 of file G4Trajectory.hh.

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

Referenced by CreateAttValues().

◆ GetTrackID()

G4int G4Trajectory::GetTrackID ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 82 of file G4Trajectory.hh.

82{ return fTrackID; }

◆ MergeTrajectory()

void G4Trajectory::MergeTrajectory ( G4VTrajectory * secondTrajectory)
overridevirtual

Implements G4VTrajectory.

Definition at line 191 of file G4Trajectory.cc.

192{
193 if (secondTrajectory == nullptr) return;
194
195 auto seco = (G4Trajectory*)secondTrajectory;
196 G4int ent = seco->GetPointEntries();
197 for (G4int i = 1; i < ent; ++i) // initial pt of 2nd trajectory shouldn't be merged
198 {
199 positionRecord->push_back((*(seco->positionRecord))[i]);
200 }
201 delete (*seco->positionRecord)[0];
202 seco->positionRecord->clear();
203}

◆ operator delete()

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

Definition at line 125 of file G4Trajectory.hh.

126{
127 aTrajectoryAllocator()->FreeSingle((G4Trajectory*)aTrajectory);
128}
G4TRACKING_DLL G4Allocator< G4Trajectory > *& aTrajectoryAllocator()

◆ operator new()

void * G4Trajectory::operator new ( size_t )
inline

Definition at line 117 of file G4Trajectory.hh.

118{
119 if (aTrajectoryAllocator() == nullptr) {
121 }
122 return (void*)aTrajectoryAllocator()->MallocSingle();
123}

◆ operator==()

G4bool G4Trajectory::operator== ( const G4Trajectory & r) const
inline

Definition at line 130 of file G4Trajectory.hh.

130{ return (this == &r); }

◆ ShowTrajectory()

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

Reimplemented from G4VTrajectory.

Definition at line 100 of file G4Trajectory.cc.

101{
102 // Invoke the default implementation in G4VTrajectory...
104
105 // ... or override with your own code here.
106}
virtual void ShowTrajectory(std::ostream &os=G4cout) const

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