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

Friends

class G4ClonedSmoothTrajectory

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 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 109 of file G4SmoothTrajectory.cc.

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

◆ 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 (auto& i : *right.positionRecord) {
97 auto rightPoint = (G4SmoothTrajectoryPoint*)i;
98 positionRecord->push_back(new G4SmoothTrajectoryPoint(*rightPoint));
99 }
100}
G4VTrajectory()=default

Member Function Documentation

◆ AppendStep()

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

Implements G4VTrajectory.

Definition at line 203 of file G4SmoothTrajectory.cc.

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

◆ CloneForMaster()

G4VTrajectory * G4SmoothTrajectory::CloneForMaster ( ) const
overridevirtual

Reimplemented from G4VTrajectory.

Definition at line 102 of file G4SmoothTrajectory.cc.

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

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Definition at line 174 of file G4SmoothTrajectory.cc.

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

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

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Definition at line 138 of file G4SmoothTrajectory.cc.

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

210{
211 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
212}
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 214 of file G4SmoothTrajectory.cc.

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

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

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