Geant4 11.3.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 ()=default
 
 G4Trajectory (const G4Track *aTrack)
 
 G4Trajectory (G4Trajectory &)
 
 ~G4Trajectory () override
 
void * operator new (size_t)
 
void operator delete (void *)
 
G4bool operator== (const G4Trajectory &r) const
 
G4VTrajectoryCloneForMaster () const override
 
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
 

Friends

class G4ClonedTrajectory
 

Detailed Description

Definition at line 64 of file G4Trajectory.hh.

Constructor & Destructor Documentation

◆ G4Trajectory() [1/3]

G4Trajectory::G4Trajectory ( )
default

◆ G4Trajectory() [2/3]

G4Trajectory::G4Trajectory ( const G4Track * aTrack)

Definition at line 62 of file G4Trajectory.cc.

63{
64 G4ParticleDefinition* fpParticleDefinition = aTrack->GetDefinition();
65 ParticleName = fpParticleDefinition->GetParticleName();
66 PDGCharge = fpParticleDefinition->GetPDGCharge();
67 PDGEncoding = fpParticleDefinition->GetPDGEncoding();
68 fTrackID = aTrack->GetTrackID();
69 fParentID = aTrack->GetParentID();
70 initialKineticEnergy = aTrack->GetKineticEnergy();
71 initialMomentum = aTrack->GetMomentum();
72 positionRecord = new G4TrajectoryPointContainer();
73
74 // Following is for the first trajectory point
75 positionRecord->push_back(new G4TrajectoryPoint(aTrack->GetPosition()));
76}
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 G4TrajectoryPointContainer();
88
89 for (auto& i : *right.positionRecord) {
90 auto rightPoint = (G4TrajectoryPoint*)i;
91 positionRecord->push_back(new G4TrajectoryPoint(*rightPoint));
92 }
93}

◆ ~G4Trajectory()

G4Trajectory::~G4Trajectory ( )
override

Definition at line 102 of file G4Trajectory.cc.

103{
104 if (positionRecord != nullptr) {
105 for (auto& i : *positionRecord) {
106 delete i;
107 }
108 positionRecord->clear();
109 delete positionRecord;
110 }
111}

Member Function Documentation

◆ AppendStep()

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

Implements G4VTrajectory.

Definition at line 194 of file G4Trajectory.cc.

195{
196 positionRecord->push_back(new G4TrajectoryPoint(aStep->GetPostStepPoint()->GetPosition()));
197}
const G4ThreeVector & GetPosition() const
G4StepPoint * GetPostStepPoint() const

◆ CloneForMaster()

G4VTrajectory * G4Trajectory::CloneForMaster ( ) const
overridevirtual

Reimplemented from G4VTrajectory.

Definition at line 95 of file G4Trajectory.cc.

96{
97 G4AutoLock lock(&CloneTrajectoryMutex);
98 auto* cloned = new G4ClonedTrajectory(*this);
99 return cloned;
100}
G4TemplateAutoLock< G4Mutex > G4AutoLock
friend class G4ClonedTrajectory

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Definition at line 165 of file G4Trajectory.cc.

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

◆ DrawTrajectory()

void G4Trajectory::DrawTrajectory ( ) const
overridevirtual

Reimplemented from G4VTrajectory.

Definition at line 121 of file G4Trajectory.cc.

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

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Definition at line 129 of file G4Trajectory.cc.

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

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

◆ GetCharge()

G4double G4Trajectory::GetCharge ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 92 of file G4Trajectory.hh.

92{ return PDGCharge; }

◆ GetInitialKineticEnergy()

G4double G4Trajectory::GetInitialKineticEnergy ( ) const
inline

Definition at line 94 of file G4Trajectory.hh.

94{ return initialKineticEnergy; }

◆ GetInitialMomentum()

G4ThreeVector G4Trajectory::GetInitialMomentum ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 95 of file G4Trajectory.hh.

95{ return initialMomentum; }

◆ GetParentID()

G4int G4Trajectory::GetParentID ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 90 of file G4Trajectory.hh.

90{ return fParentID; }

◆ GetParticleDefinition()

G4ParticleDefinition * G4Trajectory::GetParticleDefinition ( )

Definition at line 199 of file G4Trajectory.cc.

200{
201 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
202}
static G4ParticleTable * GetParticleTable()

◆ GetParticleName()

G4String G4Trajectory::GetParticleName ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 91 of file G4Trajectory.hh.

91{ return ParticleName; }

◆ GetPDGEncoding()

G4int G4Trajectory::GetPDGEncoding ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 93 of file G4Trajectory.hh.

93{ return PDGEncoding; }

◆ GetPoint()

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

Implements G4VTrajectory.

Definition at line 103 of file G4Trajectory.hh.

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

◆ GetPointEntries()

G4int G4Trajectory::GetPointEntries ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 102 of file G4Trajectory.hh.

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

Referenced by CreateAttValues().

◆ GetTrackID()

G4int G4Trajectory::GetTrackID ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 89 of file G4Trajectory.hh.

89{ return fTrackID; }

◆ MergeTrajectory()

void G4Trajectory::MergeTrajectory ( G4VTrajectory * secondTrajectory)
overridevirtual

Implements G4VTrajectory.

Definition at line 204 of file G4Trajectory.cc.

205{
206 if (secondTrajectory == nullptr) return;
207
208 auto seco = (G4Trajectory*)secondTrajectory;
209 G4int ent = seco->GetPointEntries();
210 for (G4int i = 1; i < ent; ++i) // initial pt of 2nd trajectory shouldn't be merged
211 {
212 positionRecord->push_back((*(seco->positionRecord))[i]);
213 }
214 delete (*seco->positionRecord)[0];
215 seco->positionRecord->clear();
216}
G4Trajectory()=default

◆ operator delete()

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

Definition at line 132 of file G4Trajectory.hh.

133{
134 aTrajectoryAllocator()->FreeSingle((G4Trajectory*)aTrajectory);
135}
G4TRACKING_DLL G4Allocator< G4Trajectory > *& aTrajectoryAllocator()

◆ operator new()

void * G4Trajectory::operator new ( size_t )
inline

Definition at line 124 of file G4Trajectory.hh.

125{
126 if (aTrajectoryAllocator() == nullptr) {
127 aTrajectoryAllocator() = new G4Allocator<G4Trajectory>;
128 }
129 return (void*)aTrajectoryAllocator()->MallocSingle();
130}

◆ operator==()

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

Definition at line 137 of file G4Trajectory.hh.

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

◆ ShowTrajectory()

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

Reimplemented from G4VTrajectory.

Definition at line 113 of file G4Trajectory.cc.

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

Friends And Related Symbol Documentation

◆ G4ClonedTrajectory

friend class G4ClonedTrajectory
friend

Definition at line 68 of file G4Trajectory.hh.

Referenced by CloneForMaster(), and G4ClonedTrajectory.


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