Geant4 9.6.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 ()
 
void * operator new (size_t)
 
void operator delete (void *)
 
int operator== (const G4SmoothTrajectory &right) const
 
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 (G4int i_mode=0) const
 
virtual void AppendStep (const G4Step *aStep)
 
virtual int 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 int GetPointEntries () const =0
 
virtual G4VTrajectoryPointGetPoint (G4int i) const =0
 
virtual void ShowTrajectory (std::ostream &os=G4cout) const
 
virtual void DrawTrajectory (G4int i_mode=0) 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 64 of file G4SmoothTrajectory.hh.

Constructor & Destructor Documentation

◆ G4SmoothTrajectory() [1/3]

G4SmoothTrajectory::G4SmoothTrajectory ( )

Definition at line 58 of file G4SmoothTrajectory.cc.

59: positionRecord(0), fTrackID(0), fParentID(0),
60 PDGEncoding( 0 ), PDGCharge(0.0), ParticleName(""),
61 initialKineticEnergy( 0. ), initialMomentum( G4ThreeVector() )
62{;}
CLHEP::Hep3Vector G4ThreeVector

◆ G4SmoothTrajectory() [2/3]

G4SmoothTrajectory::G4SmoothTrajectory ( const G4Track aTrack)

Definition at line 64 of file G4SmoothTrajectory.cc.

65{
66 G4ParticleDefinition * fpParticleDefinition = aTrack->GetDefinition();
67 ParticleName = fpParticleDefinition->GetParticleName();
68 PDGCharge = fpParticleDefinition->GetPDGCharge();
69 PDGEncoding = fpParticleDefinition->GetPDGEncoding();
70 fTrackID = aTrack->GetTrackID();
71 fParentID = aTrack->GetParentID();
72 initialKineticEnergy = aTrack->GetKineticEnergy();
73 initialMomentum = aTrack->GetMomentum();
74 positionRecord = new TrajectoryPointContainer();
75 // Following is for the first trajectory point
76 positionRecord->push_back(new G4SmoothTrajectoryPoint(aTrack->GetPosition()));
77
78 // The first point has no auxiliary points, so set the auxiliary
79 // points vector to NULL (jacek 31/10/2002)
80 positionRecord->push_back(new G4SmoothTrajectoryPoint(aTrack->GetPosition(), 0));
81}
std::vector< G4VTrajectoryPoint * > TrajectoryPointContainer
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 83 of file G4SmoothTrajectory.cc.

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

◆ ~G4SmoothTrajectory()

G4SmoothTrajectory::~G4SmoothTrajectory ( )
virtual

Definition at line 101 of file G4SmoothTrajectory.cc.

102{
103 if (positionRecord) {
104 // positionRecord->clearAndDestroy();
105 size_t i;
106 for(i=0;i<positionRecord->size();i++){
107 delete (*positionRecord)[i];
108 }
109 positionRecord->clear();
110 delete positionRecord;
111 }
112}

Member Function Documentation

◆ AppendStep()

void G4SmoothTrajectory::AppendStep ( const G4Step aStep)
virtual

Implements G4VTrajectory.

Definition at line 218 of file G4SmoothTrajectory.cc.

219{
220 // (jacek 30/10/2002)
221 positionRecord->push_back(
224}
const G4ThreeVector & GetPosition() const
std::vector< G4ThreeVector > * GetPointerToVectorOfAuxiliaryPoints() const
Definition: G4Step.hh:240
G4StepPoint * GetPostStepPoint() const

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Definition at line 181 of file G4SmoothTrajectory.cc.

182{
183 std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
184
185 values->push_back
186 (G4AttValue("ID",G4UIcommand::ConvertToString(fTrackID),""));
187
188 values->push_back
189 (G4AttValue("PID",G4UIcommand::ConvertToString(fParentID),""));
190
191 values->push_back(G4AttValue("PN",ParticleName,""));
192
193 values->push_back
194 (G4AttValue("Ch",G4UIcommand::ConvertToString(PDGCharge),""));
195
196 values->push_back
197 (G4AttValue("PDG",G4UIcommand::ConvertToString(PDGEncoding),""));
198
199 values->push_back
200 (G4AttValue("IKE",G4BestUnit(initialKineticEnergy,"Energy"),""));
201
202 values->push_back
203 (G4AttValue("IMom",G4BestUnit(initialMomentum,"Energy"),""));
204
205 values->push_back
206 (G4AttValue("IMag",G4BestUnit(initialMomentum.mag(),"Energy"),""));
207
208 values->push_back
210
211#ifdef G4ATTDEBUG
212 G4cout << G4AttCheck(values,GetAttDefs());
213#endif
214
215 return values;
216}
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4DLLIMPORT std::ostream G4cout
double mag() const
virtual int GetPointEntries() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:349

◆ DrawTrajectory()

void G4SmoothTrajectory::DrawTrajectory ( G4int  i_mode = 0) const
virtual

Reimplemented from G4VTrajectory.

Definition at line 130 of file G4SmoothTrajectory.cc.

131{
132 // Invoke the default implementation in G4VTrajectory...
134 // ... or override with your own code here.
135}
virtual void DrawTrajectory(G4int i_mode=0) const

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Definition at line 137 of file G4SmoothTrajectory.cc.

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

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

◆ GetCharge()

G4double G4SmoothTrajectory::GetCharge ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 96 of file G4SmoothTrajectory.hh.

97 { return PDGCharge; }

◆ GetInitialKineticEnergy()

G4double G4SmoothTrajectory::GetInitialKineticEnergy ( ) const
inline

Definition at line 100 of file G4SmoothTrajectory.hh.

101 { return initialKineticEnergy; }

◆ GetInitialMomentum()

G4ThreeVector G4SmoothTrajectory::GetInitialMomentum ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 102 of file G4SmoothTrajectory.hh.

103 { return initialMomentum; }

◆ GetParentID()

G4int G4SmoothTrajectory::GetParentID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 92 of file G4SmoothTrajectory.hh.

93 { return fParentID; }

◆ GetParticleDefinition()

G4ParticleDefinition * G4SmoothTrajectory::GetParticleDefinition ( )

Definition at line 226 of file G4SmoothTrajectory.cc.

227{
228 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
229}
static G4ParticleTable * GetParticleTable()

◆ GetParticleName()

G4String G4SmoothTrajectory::GetParticleName ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 94 of file G4SmoothTrajectory.hh.

95 { return ParticleName; }

◆ GetPDGEncoding()

G4int G4SmoothTrajectory::GetPDGEncoding ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 98 of file G4SmoothTrajectory.hh.

99 { return PDGEncoding; }

◆ GetPoint()

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

Implements G4VTrajectory.

Definition at line 111 of file G4SmoothTrajectory.hh.

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

◆ GetPointEntries()

virtual int G4SmoothTrajectory::GetPointEntries ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 110 of file G4SmoothTrajectory.hh.

110{ return positionRecord->size(); }

Referenced by CreateAttValues(), and MergeTrajectory().

◆ GetTrackID()

G4int G4SmoothTrajectory::GetTrackID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 90 of file G4SmoothTrajectory.hh.

91 { return fTrackID; }

◆ MergeTrajectory()

void G4SmoothTrajectory::MergeTrajectory ( G4VTrajectory secondTrajectory)
virtual

Implements G4VTrajectory.

Definition at line 231 of file G4SmoothTrajectory.cc.

232{
233 if(!secondTrajectory) return;
234
235 G4SmoothTrajectory* seco = (G4SmoothTrajectory*)secondTrajectory;
236 G4int ent = seco->GetPointEntries();
237 for(G4int i=1;i<ent;i++) // initial point of the second trajectory should not be merged
238 {
239 positionRecord->push_back((*(seco->positionRecord))[i]);
240 // positionRecord->push_back(seco->positionRecord->removeAt(1));
241 }
242 delete (*seco->positionRecord)[0];
243 seco->positionRecord->clear();
244}
int G4int
Definition: G4Types.hh:66

◆ operator delete()

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

Definition at line 149 of file G4SmoothTrajectory.hh.

150{
151 aSmoothTrajectoryAllocator.FreeSingle((G4SmoothTrajectory*)aTrajectory);
152}
G4DLLIMPORT G4Allocator< G4SmoothTrajectory > aSmoothTrajectoryAllocator

◆ operator new()

void * G4SmoothTrajectory::operator new ( size_t  )
inline

Definition at line 142 of file G4SmoothTrajectory.hh.

143{
144 void* aTrajectory;
145 aTrajectory = (void*)aSmoothTrajectoryAllocator.MallocSingle();
146 return aTrajectory;
147}

◆ operator==()

int G4SmoothTrajectory::operator== ( const G4SmoothTrajectory right) const
inline

Definition at line 86 of file G4SmoothTrajectory.hh.

87 {return (this==&right);}

◆ ShowTrajectory()

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

Reimplemented from G4VTrajectory.

Definition at line 114 of file G4SmoothTrajectory.cc.

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

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