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

#include <G4SmoothTrajectory.hh>

+ Inheritance diagram for G4SmoothTrajectory:

Public Member Functions

 G4SmoothTrajectory ()=default
 
 G4SmoothTrajectory (const G4Track *aTrack)
 
 ~G4SmoothTrajectory () override
 
 G4SmoothTrajectory (G4SmoothTrajectory &)
 
G4SmoothTrajectoryoperator= (const G4SmoothTrajectory &)=delete
 
G4bool operator== (const G4SmoothTrajectory &r) const
 
void * operator new (size_t)
 
void operator delete (void *)
 
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 64 of file G4SmoothTrajectory.hh.

Constructor & Destructor Documentation

◆ G4SmoothTrajectory() [1/3]

G4SmoothTrajectory::G4SmoothTrajectory ( )
default

◆ G4SmoothTrajectory() [2/3]

G4SmoothTrajectory::G4SmoothTrajectory ( const G4Track * aTrack)

Definition at line 56 of file G4SmoothTrajectory.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 //
70 positionRecord->push_back(new G4SmoothTrajectoryPoint(aTrack->GetPosition()));
71
72 // The first point has no auxiliary points, so set the auxiliary
73 // points vector to NULL
74 //
75 positionRecord->push_back(new G4SmoothTrajectoryPoint(aTrack->GetPosition(), nullptr));
76}
const G4String & GetParticleName() const
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
G4ThreeVector GetMomentum() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4int GetParentID() const

◆ ~G4SmoothTrajectory()

G4SmoothTrajectory::~G4SmoothTrajectory ( )
override

Definition at line 95 of file G4SmoothTrajectory.cc.

96{
97 if (positionRecord != nullptr) {
98 for (auto& i : *positionRecord) {
99 delete i;
100 }
101 positionRecord->clear();
102 delete positionRecord;
103 }
104}

◆ G4SmoothTrajectory() [3/3]

G4SmoothTrajectory::G4SmoothTrajectory ( G4SmoothTrajectory & right)

Definition at line 78 of file G4SmoothTrajectory.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 = (G4SmoothTrajectoryPoint*)i;
91 positionRecord->push_back(new G4SmoothTrajectoryPoint(*rightPoint));
92 }
93}

Member Function Documentation

◆ AppendStep()

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

Implements G4VTrajectory.

Definition at line 189 of file G4SmoothTrajectory.cc.

190{
191 positionRecord->push_back(new G4SmoothTrajectoryPoint(
193}
const G4ThreeVector & GetPosition() const
std::vector< G4ThreeVector > * GetPointerToVectorOfAuxiliaryPoints() const
G4StepPoint * GetPostStepPoint() const

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Definition at line 160 of file G4SmoothTrajectory.cc.

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

◆ DrawTrajectory()

void G4SmoothTrajectory::DrawTrajectory ( ) const
overridevirtual

Reimplemented from G4VTrajectory.

Definition at line 115 of file G4SmoothTrajectory.cc.

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

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Definition at line 124 of file G4SmoothTrajectory.cc.

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

Implements G4VTrajectory.

Definition at line 88 of file G4SmoothTrajectory.hh.

88{ return PDGCharge; }

◆ GetInitialKineticEnergy()

G4double G4SmoothTrajectory::GetInitialKineticEnergy ( ) const
inline

Definition at line 90 of file G4SmoothTrajectory.hh.

90{ return initialKineticEnergy; }

◆ GetInitialMomentum()

G4ThreeVector G4SmoothTrajectory::GetInitialMomentum ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 91 of file G4SmoothTrajectory.hh.

91{ return initialMomentum; }

◆ GetParentID()

G4int G4SmoothTrajectory::GetParentID ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 86 of file G4SmoothTrajectory.hh.

86{ return fParentID; }

◆ GetParticleDefinition()

G4ParticleDefinition * G4SmoothTrajectory::GetParticleDefinition ( )

Definition at line 195 of file G4SmoothTrajectory.cc.

196{
197 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
198}
static G4ParticleTable * GetParticleTable()

◆ GetParticleName()

G4String G4SmoothTrajectory::GetParticleName ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 87 of file G4SmoothTrajectory.hh.

87{ return ParticleName; }

◆ GetPDGEncoding()

G4int G4SmoothTrajectory::GetPDGEncoding ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 89 of file G4SmoothTrajectory.hh.

89{ return PDGEncoding; }

◆ GetPoint()

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

Implements G4VTrajectory.

Definition at line 99 of file G4SmoothTrajectory.hh.

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

◆ GetPointEntries()

G4int G4SmoothTrajectory::GetPointEntries ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 98 of file G4SmoothTrajectory.hh.

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

Referenced by CreateAttValues().

◆ GetTrackID()

G4int G4SmoothTrajectory::GetTrackID ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 85 of file G4SmoothTrajectory.hh.

85{ return fTrackID; }

◆ MergeTrajectory()

void G4SmoothTrajectory::MergeTrajectory ( G4VTrajectory * secondTrajectory)
overridevirtual

Implements G4VTrajectory.

Definition at line 200 of file G4SmoothTrajectory.cc.

201{
202 if (secondTrajectory == nullptr) return;
203
204 auto seco = (G4SmoothTrajectory*)secondTrajectory;
205 G4int ent = seco->GetPointEntries();
206
207 // initial point of the second trajectory should not be merged
208 //
209 for (G4int i = 1; i < ent; ++i) {
210 positionRecord->push_back((*(seco->positionRecord))[i]);
211 }
212 delete (*seco->positionRecord)[0];
213 seco->positionRecord->clear();
214}

◆ operator delete()

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

Definition at line 130 of file G4SmoothTrajectory.hh.

131{
132 aSmoothTrajectoryAllocator()->FreeSingle((G4SmoothTrajectory*)aTrajectory);
133}
G4TRACKING_DLL G4Allocator< G4SmoothTrajectory > *& aSmoothTrajectoryAllocator()

◆ operator new()

void * G4SmoothTrajectory::operator new ( size_t )
inline

Definition at line 122 of file G4SmoothTrajectory.hh.

123{
124 if (aSmoothTrajectoryAllocator() == nullptr) {
126 }
127 return (void*)aSmoothTrajectoryAllocator()->MallocSingle();
128}

◆ operator=()

G4SmoothTrajectory & G4SmoothTrajectory::operator= ( const G4SmoothTrajectory & )
delete

◆ operator==()

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

Definition at line 135 of file G4SmoothTrajectory.hh.

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

◆ ShowTrajectory()

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

Reimplemented from G4VTrajectory.

Definition at line 106 of file G4SmoothTrajectory.cc.

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

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