Geant4 11.3.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 ()=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 *)
 
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 G4ClonedSmoothTrajectory
 

Detailed Description

Definition at line 65 of file G4SmoothTrajectory.hh.

Constructor & Destructor Documentation

◆ G4SmoothTrajectory() [1/3]

G4SmoothTrajectory::G4SmoothTrajectory ( )
default

◆ G4SmoothTrajectory() [2/3]

G4SmoothTrajectory::G4SmoothTrajectory ( const G4Track * aTrack)

Definition at line 62 of file G4SmoothTrajectory.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 //
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
80 //
81 positionRecord->push_back(new G4SmoothTrajectoryPoint(aTrack->GetPosition(), nullptr));
82}
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 108 of file G4SmoothTrajectory.cc.

109{
110 if (positionRecord != nullptr) {
111 for (auto& i : *positionRecord) {
112 delete i;
113 }
114 positionRecord->clear();
115 delete positionRecord;
116 }
117}

◆ G4SmoothTrajectory() [3/3]

G4SmoothTrajectory::G4SmoothTrajectory ( G4SmoothTrajectory & right)

Definition at line 84 of file G4SmoothTrajectory.cc.

85{
86 ParticleName = right.ParticleName;
87 PDGCharge = right.PDGCharge;
88 PDGEncoding = right.PDGEncoding;
89 fTrackID = right.fTrackID;
90 fParentID = right.fParentID;
91 initialKineticEnergy = right.initialKineticEnergy;
92 initialMomentum = right.initialMomentum;
93 positionRecord = new G4TrajectoryPointContainer();
94
95 for (auto& i : *right.positionRecord) {
96 auto rightPoint = (G4SmoothTrajectoryPoint*)i;
97 positionRecord->push_back(new G4SmoothTrajectoryPoint(*rightPoint));
98 }
99}

Member Function Documentation

◆ AppendStep()

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

Implements G4VTrajectory.

Definition at line 202 of file G4SmoothTrajectory.cc.

203{
204 positionRecord->push_back(new G4SmoothTrajectoryPoint(
206}
const G4ThreeVector & GetPosition() const
std::vector< G4ThreeVector > * GetPointerToVectorOfAuxiliaryPoints() const
G4StepPoint * GetPostStepPoint() const

◆ CloneForMaster()

G4VTrajectory * G4SmoothTrajectory::CloneForMaster ( ) const
overridevirtual

Reimplemented from G4VTrajectory.

Definition at line 101 of file G4SmoothTrajectory.cc.

102{
103 G4AutoLock lock(&CloneSmoothTrajectoryMutex);
104 auto* cloned = new G4ClonedSmoothTrajectory(*this);
105 return cloned;
106}
G4TemplateAutoLock< G4Mutex > G4AutoLock
friend class G4ClonedSmoothTrajectory

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Definition at line 173 of file G4SmoothTrajectory.cc.

174{
175 auto values = new std::vector<G4AttValue>;
176
177 values->push_back(G4AttValue("ID", G4UIcommand::ConvertToString(fTrackID), ""));
178
179 values->push_back(G4AttValue("PID", G4UIcommand::ConvertToString(fParentID), ""));
180
181 values->push_back(G4AttValue("PN", ParticleName, ""));
182
183 values->push_back(G4AttValue("Ch", G4UIcommand::ConvertToString(PDGCharge), ""));
184
185 values->push_back(G4AttValue("PDG", G4UIcommand::ConvertToString(PDGEncoding), ""));
186
187 values->push_back(G4AttValue("IKE", G4BestUnit(initialKineticEnergy, "Energy"), ""));
188
189 values->push_back(G4AttValue("IMom", G4BestUnit(initialMomentum, "Energy"), ""));
190
191 values->push_back(G4AttValue("IMag", G4BestUnit(initialMomentum.mag(), "Energy"), ""));
192
193 values->push_back(G4AttValue("NTP", G4UIcommand::ConvertToString(GetPointEntries()), ""));
194
195#ifdef G4ATTDEBUG
196 G4cout << G4AttCheck(values, GetAttDefs());
197#endif
198
199 return values;
200}
#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 G4SmoothTrajectory::DrawTrajectory ( ) const
overridevirtual

Reimplemented from G4VTrajectory.

Definition at line 128 of file G4SmoothTrajectory.cc.

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

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Definition at line 137 of file G4SmoothTrajectory.cc.

138{
139 G4bool isNew;
140 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("G4SmoothTrajectory", isNew);
141 if (isNew) {
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", "Physics", "G4BestUnit", "G4double");
159
160 G4String IMom("IMom");
161 (*store)[IMom] = G4AttDef(IMom, "Initial momentum", "Physics", "G4BestUnit", "G4ThreeVector");
162
163 G4String IMag("IMag");
164 (*store)[IMag] =
165 G4AttDef(IMag, "Initial momentum magnitude", "Physics", "G4BestUnit", "G4double");
166
167 G4String NTP("NTP");
168 (*store)[NTP] = G4AttDef(NTP, "No. of points", "Physics", "", "G4int");
169 }
170 return store;
171}
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 94 of file G4SmoothTrajectory.hh.

94{ return PDGCharge; }

◆ GetInitialKineticEnergy()

G4double G4SmoothTrajectory::GetInitialKineticEnergy ( ) const
inline

Definition at line 96 of file G4SmoothTrajectory.hh.

96{ return initialKineticEnergy; }

◆ GetInitialMomentum()

G4ThreeVector G4SmoothTrajectory::GetInitialMomentum ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 97 of file G4SmoothTrajectory.hh.

97{ return initialMomentum; }

◆ GetParentID()

G4int G4SmoothTrajectory::GetParentID ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 92 of file G4SmoothTrajectory.hh.

92{ return fParentID; }

◆ GetParticleDefinition()

G4ParticleDefinition * G4SmoothTrajectory::GetParticleDefinition ( )

Definition at line 208 of file G4SmoothTrajectory.cc.

209{
210 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
211}
static G4ParticleTable * GetParticleTable()

◆ GetParticleName()

G4String G4SmoothTrajectory::GetParticleName ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 93 of file G4SmoothTrajectory.hh.

93{ return ParticleName; }

◆ GetPDGEncoding()

G4int G4SmoothTrajectory::GetPDGEncoding ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 95 of file G4SmoothTrajectory.hh.

95{ return PDGEncoding; }

◆ GetPoint()

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

Implements G4VTrajectory.

Definition at line 105 of file G4SmoothTrajectory.hh.

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

◆ GetPointEntries()

G4int G4SmoothTrajectory::GetPointEntries ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 104 of file G4SmoothTrajectory.hh.

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

Referenced by CreateAttValues().

◆ GetTrackID()

G4int G4SmoothTrajectory::GetTrackID ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 91 of file G4SmoothTrajectory.hh.

91{ return fTrackID; }

◆ MergeTrajectory()

void G4SmoothTrajectory::MergeTrajectory ( G4VTrajectory * secondTrajectory)
overridevirtual

Implements G4VTrajectory.

Definition at line 213 of file G4SmoothTrajectory.cc.

214{
215 if (secondTrajectory == nullptr) return;
216
217 auto seco = (G4SmoothTrajectory*)secondTrajectory;
218 G4int ent = seco->GetPointEntries();
219
220 // initial point of the second trajectory should not be merged
221 //
222 for (G4int i = 1; i < ent; ++i) {
223 positionRecord->push_back((*(seco->positionRecord))[i]);
224 }
225 delete (*seco->positionRecord)[0];
226 seco->positionRecord->clear();
227}
G4SmoothTrajectory()=default

◆ operator delete()

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

Definition at line 136 of file G4SmoothTrajectory.hh.

137{
138 aSmoothTrajectoryAllocator()->FreeSingle((G4SmoothTrajectory*)aTrajectory);
139}
G4TRACKING_DLL G4Allocator< G4SmoothTrajectory > *& aSmoothTrajectoryAllocator()

◆ operator new()

void * G4SmoothTrajectory::operator new ( size_t )
inline

Definition at line 128 of file G4SmoothTrajectory.hh.

129{
130 if (aSmoothTrajectoryAllocator() == nullptr) {
131 aSmoothTrajectoryAllocator() = new G4Allocator<G4SmoothTrajectory>;
132 }
133 return (void*)aSmoothTrajectoryAllocator()->MallocSingle();
134}

◆ operator=()

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

◆ operator==()

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

Definition at line 141 of file G4SmoothTrajectory.hh.

142{
143 return (this == &r);
144}

◆ ShowTrajectory()

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

Reimplemented from G4VTrajectory.

Definition at line 119 of file G4SmoothTrajectory.cc.

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

Friends And Related Symbol Documentation

◆ G4ClonedSmoothTrajectory

friend class G4ClonedSmoothTrajectory
friend

Definition at line 69 of file G4SmoothTrajectory.hh.

Referenced by CloneForMaster(), and G4ClonedSmoothTrajectory.


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