Geant4 10.7.0
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 ()
 
 G4SmoothTrajectory (const G4Track *aTrack)
 
 G4SmoothTrajectory (G4SmoothTrajectory &)
 
virtual ~G4SmoothTrajectory ()
 
G4SmoothTrajectoryoperator= (const G4SmoothTrajectory &)=delete
 
G4int operator== (const G4SmoothTrajectory &r) const
 
void * operator new (size_t)
 
void operator delete (void *)
 
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 63 of file G4SmoothTrajectory.hh.

Constructor & Destructor Documentation

◆ G4SmoothTrajectory() [1/3]

G4SmoothTrajectory::G4SmoothTrajectory ( )

Definition at line 55 of file G4SmoothTrajectory.cc.

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

◆ G4SmoothTrajectory() [2/3]

G4SmoothTrajectory::G4SmoothTrajectory ( const G4Track aTrack)

Definition at line 60 of file G4SmoothTrajectory.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 //
74 positionRecord
75 ->push_back(new G4SmoothTrajectoryPoint(aTrack->GetPosition()));
76
77 // The first point has no auxiliary points, so set the auxiliary
78 // points vector to NULL
79 //
80 positionRecord
81 ->push_back(new G4SmoothTrajectoryPoint(aTrack->GetPosition(), nullptr));
82}
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

◆ G4SmoothTrajectory() [3/3]

G4SmoothTrajectory::G4SmoothTrajectory ( G4SmoothTrajectory right)

Definition at line 84 of file G4SmoothTrajectory.cc.

86{
87 ParticleName = right.ParticleName;
88 PDGCharge = right.PDGCharge;
89 PDGEncoding = right.PDGEncoding;
90 fTrackID = right.fTrackID;
91 fParentID = right.fParentID;
92 initialKineticEnergy = right.initialKineticEnergy;
93 initialMomentum = right.initialMomentum;
94 positionRecord = new G4TrajectoryPointContainer();
95
96 for(std::size_t i=0; i<right.positionRecord->size(); ++i)
97 {
98 G4SmoothTrajectoryPoint* rightPoint
99 = (G4SmoothTrajectoryPoint*)((*(right.positionRecord))[i]);
100 positionRecord->push_back(new G4SmoothTrajectoryPoint(*rightPoint));
101 }
102}

◆ ~G4SmoothTrajectory()

G4SmoothTrajectory::~G4SmoothTrajectory ( )
virtual

Definition at line 104 of file G4SmoothTrajectory.cc.

105{
106 if (positionRecord)
107 {
108 for(std::size_t i=0; i<positionRecord->size(); ++i)
109 {
110 delete (*positionRecord)[i];
111 }
112 positionRecord->clear();
113 delete positionRecord;
114 }
115}

Member Function Documentation

◆ AppendStep()

void G4SmoothTrajectory::AppendStep ( const G4Step aStep)
virtual

Implements G4VTrajectory.

Definition at line 213 of file G4SmoothTrajectory.cc.

214{
215 positionRecord->push_back(
218}
const G4ThreeVector & GetPosition() const
std::vector< G4ThreeVector > * GetPointerToVectorOfAuxiliaryPoints() const
G4StepPoint * GetPostStepPoint() const

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Definition at line 176 of file G4SmoothTrajectory.cc.

177{
178 std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
179
180 values->push_back
181 (G4AttValue("ID",G4UIcommand::ConvertToString(fTrackID),""));
182
183 values->push_back
184 (G4AttValue("PID",G4UIcommand::ConvertToString(fParentID),""));
185
186 values->push_back(G4AttValue("PN",ParticleName,""));
187
188 values->push_back
189 (G4AttValue("Ch",G4UIcommand::ConvertToString(PDGCharge),""));
190
191 values->push_back
192 (G4AttValue("PDG",G4UIcommand::ConvertToString(PDGEncoding),""));
193
194 values->push_back
195 (G4AttValue("IKE",G4BestUnit(initialKineticEnergy,"Energy"),""));
196
197 values->push_back
198 (G4AttValue("IMom",G4BestUnit(initialMomentum,"Energy"),""));
199
200 values->push_back
201 (G4AttValue("IMag",G4BestUnit(initialMomentum.mag(),"Energy"),""));
202
203 values->push_back
205
206#ifdef G4ATTDEBUG
207 G4cout << G4AttCheck(values,GetAttDefs());
208#endif
209
210 return values;
211}
#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:430

◆ DrawTrajectory()

void G4SmoothTrajectory::DrawTrajectory ( ) const
virtual

Reimplemented from G4VTrajectory.

Definition at line 126 of file G4SmoothTrajectory.cc.

127{
128 // Invoke the default implementation in G4VTrajectory...
129 //
131
132 // ... or override with your own code here.
133}
virtual void DrawTrajectory() const

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Definition at line 135 of file G4SmoothTrajectory.cc.

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

Implements G4VTrajectory.

Definition at line 91 of file G4SmoothTrajectory.hh.

92 { return PDGCharge; }

◆ GetInitialKineticEnergy()

G4double G4SmoothTrajectory::GetInitialKineticEnergy ( ) const
inline

Definition at line 95 of file G4SmoothTrajectory.hh.

96 { return initialKineticEnergy; }

◆ GetInitialMomentum()

G4ThreeVector G4SmoothTrajectory::GetInitialMomentum ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 97 of file G4SmoothTrajectory.hh.

98 { return initialMomentum; }

◆ GetParentID()

G4int G4SmoothTrajectory::GetParentID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 87 of file G4SmoothTrajectory.hh.

88 { return fParentID; }

◆ GetParticleDefinition()

G4ParticleDefinition * G4SmoothTrajectory::GetParticleDefinition ( )

Definition at line 220 of file G4SmoothTrajectory.cc.

221{
222 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
223}
static G4ParticleTable * GetParticleTable()

◆ GetParticleName()

G4String G4SmoothTrajectory::GetParticleName ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 89 of file G4SmoothTrajectory.hh.

90 { return ParticleName; }

◆ GetPDGEncoding()

G4int G4SmoothTrajectory::GetPDGEncoding ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 93 of file G4SmoothTrajectory.hh.

94 { return PDGEncoding; }

◆ GetPoint()

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

Implements G4VTrajectory.

Definition at line 107 of file G4SmoothTrajectory.hh.

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

◆ GetPointEntries()

virtual G4int G4SmoothTrajectory::GetPointEntries ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 105 of file G4SmoothTrajectory.hh.

106 { return positionRecord->size(); }

Referenced by CreateAttValues(), and MergeTrajectory().

◆ GetTrackID()

G4int G4SmoothTrajectory::GetTrackID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 85 of file G4SmoothTrajectory.hh.

86 { return fTrackID; }

◆ MergeTrajectory()

void G4SmoothTrajectory::MergeTrajectory ( G4VTrajectory secondTrajectory)
virtual

Implements G4VTrajectory.

Definition at line 225 of file G4SmoothTrajectory.cc.

226{
227 if(secondTrajectory == nullptr) return;
228
229 G4SmoothTrajectory* seco = (G4SmoothTrajectory*)secondTrajectory;
230 G4int ent = seco->GetPointEntries();
231
232 // initial point of the second trajectory should not be merged
233 //
234 for(G4int i=1; i<ent; ++i)
235 {
236 positionRecord->push_back((*(seco->positionRecord))[i]);
237 }
238 delete (*seco->positionRecord)[0];
239 seco->positionRecord->clear();
240}
int G4int
Definition: G4Types.hh:85

◆ operator delete()

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

Definition at line 142 of file G4SmoothTrajectory.hh.

143{
144 aSmoothTrajectoryAllocator()->FreeSingle((G4SmoothTrajectory*)aTrajectory);
145}
G4TRACKING_DLL G4Allocator< G4SmoothTrajectory > *& aSmoothTrajectoryAllocator()

◆ operator new()

void * G4SmoothTrajectory::operator new ( size_t  )
inline

Definition at line 133 of file G4SmoothTrajectory.hh.

134{
135 if (aSmoothTrajectoryAllocator() == nullptr)
136 {
138 }
139 return (void*)aSmoothTrajectoryAllocator()->MallocSingle();
140}

◆ operator=()

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

◆ operator==()

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

Definition at line 147 of file G4SmoothTrajectory.hh.

148{
149 return (this==&r);
150}

◆ ShowTrajectory()

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

Reimplemented from G4VTrajectory.

Definition at line 117 of file G4SmoothTrajectory.cc.

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

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