Geant4 11.4.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
std::shared_ptr< std::vector< G4AttValue > > GetAttValues () const

Friends

class G4ClonedTrajectory

Additional Inherited Members

Protected Member Functions inherited from G4VTrajectory
 G4VTrajectory (const G4VTrajectory &right)=default
G4VTrajectoryoperator= (const G4VTrajectory &right)=default
 G4VTrajectory (G4VTrajectory &&)=default
G4VTrajectoryoperator= (G4VTrajectory &&)=default

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.

80{
81 ParticleName = right.ParticleName;
82 PDGCharge = right.PDGCharge;
83 PDGEncoding = right.PDGEncoding;
84 fTrackID = right.fTrackID;
85 fParentID = right.fParentID;
86 initialKineticEnergy = right.initialKineticEnergy;
87 initialMomentum = right.initialMomentum;
88 positionRecord = new G4TrajectoryPointContainer();
89
90 for (auto& i : *right.positionRecord) {
91 auto rightPoint = (G4TrajectoryPoint*)i;
92 positionRecord->push_back(new G4TrajectoryPoint(*rightPoint));
93 }
94}
G4VTrajectory()=default

◆ ~G4Trajectory()

G4Trajectory::~G4Trajectory ( )
override

Definition at line 103 of file G4Trajectory.cc.

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

Member Function Documentation

◆ AppendStep()

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

Implements G4VTrajectory.

Definition at line 195 of file G4Trajectory.cc.

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

◆ CloneForMaster()

G4VTrajectory * G4Trajectory::CloneForMaster ( ) const
overridevirtual

Reimplemented from G4VTrajectory.

Definition at line 96 of file G4Trajectory.cc.

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

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Definition at line 166 of file G4Trajectory.cc.

167{
168 auto values = new std::vector<G4AttValue>;
169
170 values->push_back(G4AttValue("ID", G4UIcommand::ConvertToString(fTrackID), ""));
171
172 values->push_back(G4AttValue("PID", G4UIcommand::ConvertToString(fParentID), ""));
173
174 values->push_back(G4AttValue("PN", ParticleName, ""));
175
176 values->push_back(G4AttValue("Ch", G4UIcommand::ConvertToString(PDGCharge), ""));
177
178 values->push_back(G4AttValue("PDG", G4UIcommand::ConvertToString(PDGEncoding), ""));
179
180 values->push_back(G4AttValue("IKE", G4BestUnit(initialKineticEnergy, "Energy"), ""));
181
182 values->push_back(G4AttValue("IMom", G4BestUnit(initialMomentum, "Energy"), ""));
183
184 values->push_back(G4AttValue("IMag", G4BestUnit(initialMomentum.mag(), "Energy"), ""));
185
186 values->push_back(G4AttValue("NTP", G4UIcommand::ConvertToString(GetPointEntries()), ""));
187
188#ifdef G4ATTDEBUG
189 G4cout << G4AttCheck(values, GetAttDefs());
190#endif
191
192 return values;
193}
#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 122 of file G4Trajectory.cc.

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

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Definition at line 130 of file G4Trajectory.cc.

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

201{
202 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
203}
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 205 of file G4Trajectory.cc.

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

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

◆ 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: