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

#include <G4ClonedTrajectory.hh>

+ Inheritance diagram for G4ClonedTrajectory:

Public Member Functions

 G4ClonedTrajectory ()=default
 
 G4ClonedTrajectory (const G4Trajectory &)
 
 ~G4ClonedTrajectory () override
 
void * operator new (size_t)
 
void operator delete (void *)
 
G4bool operator== (const G4ClonedTrajectory &r) const
 
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
 
virtual G4VTrajectoryCloneForMaster () const
 

Detailed Description

Definition at line 60 of file G4ClonedTrajectory.hh.

Constructor & Destructor Documentation

◆ G4ClonedTrajectory() [1/2]

G4ClonedTrajectory::G4ClonedTrajectory ( )
default

◆ G4ClonedTrajectory() [2/2]

G4ClonedTrajectory::G4ClonedTrajectory ( const G4Trajectory & right)

Definition at line 54 of file G4ClonedTrajectory.cc.

55{
56 ParticleName = right.ParticleName;
57 PDGCharge = right.PDGCharge;
58 PDGEncoding = right.PDGEncoding;
59 fTrackID = right.fTrackID;
60 fParentID = right.fParentID;
61 initialKineticEnergy = right.initialKineticEnergy;
62 initialMomentum = right.initialMomentum;
63 positionRecord = new G4ClonedTrajectoryPointContainer();
64
65 for (auto& i : *right.positionRecord) {
66 auto rightPoint = (G4ClonedTrajectoryPoint*)i;
67 positionRecord->push_back(new G4ClonedTrajectoryPoint(*rightPoint));
68 }
69}

◆ ~G4ClonedTrajectory()

G4ClonedTrajectory::~G4ClonedTrajectory ( )
override

Definition at line 71 of file G4ClonedTrajectory.cc.

72{
73 if (positionRecord != nullptr) {
74 for (auto& i : *positionRecord) {
75 delete i;
76 }
77 positionRecord->clear();
78 delete positionRecord;
79 }
80}

Member Function Documentation

◆ AppendStep()

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

Implements G4VTrajectory.

Definition at line 163 of file G4ClonedTrajectory.cc.

164{
165 positionRecord->push_back(new G4ClonedTrajectoryPoint(aStep->GetPostStepPoint()->GetPosition()));
166}
const G4ThreeVector & GetPosition() const
G4StepPoint * GetPostStepPoint() const

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Definition at line 134 of file G4ClonedTrajectory.cc.

135{
136 auto values = new std::vector<G4AttValue>;
137
138 values->push_back(G4AttValue("ID", G4UIcommand::ConvertToString(fTrackID), ""));
139
140 values->push_back(G4AttValue("PID", G4UIcommand::ConvertToString(fParentID), ""));
141
142 values->push_back(G4AttValue("PN", ParticleName, ""));
143
144 values->push_back(G4AttValue("Ch", G4UIcommand::ConvertToString(PDGCharge), ""));
145
146 values->push_back(G4AttValue("PDG", G4UIcommand::ConvertToString(PDGEncoding), ""));
147
148 values->push_back(G4AttValue("IKE", G4BestUnit(initialKineticEnergy, "Energy"), ""));
149
150 values->push_back(G4AttValue("IMom", G4BestUnit(initialMomentum, "Energy"), ""));
151
152 values->push_back(G4AttValue("IMag", G4BestUnit(initialMomentum.mag(), "Energy"), ""));
153
154 values->push_back(G4AttValue("NTP", G4UIcommand::ConvertToString(GetPointEntries()), ""));
155
156#ifdef G4ATTDEBUG
157 G4cout << G4AttCheck(values, GetAttDefs());
158#endif
159
160 return values;
161}
#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 G4ClonedTrajectory::DrawTrajectory ( ) const
overridevirtual

Reimplemented from G4VTrajectory.

Definition at line 90 of file G4ClonedTrajectory.cc.

91{
92 // Invoke the default implementation in G4VTrajectory...
94
95 // ... or override with your own code here.
96}
virtual void DrawTrajectory() const

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Definition at line 98 of file G4ClonedTrajectory.cc.

99{
100 G4bool isNew;
101 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("G4Trajectory", isNew);
102 if (isNew) {
103 G4String ID("ID");
104 (*store)[ID] = G4AttDef(ID, "Track ID", "Physics", "", "G4int");
105
106 G4String PID("PID");
107 (*store)[PID] = G4AttDef(PID, "Parent ID", "Physics", "", "G4int");
108
109 G4String PN("PN");
110 (*store)[PN] = G4AttDef(PN, "Particle Name", "Physics", "", "G4String");
111
112 G4String Ch("Ch");
113 (*store)[Ch] = G4AttDef(Ch, "Charge", "Physics", "e+", "G4double");
114
115 G4String PDG("PDG");
116 (*store)[PDG] = G4AttDef(PDG, "PDG Encoding", "Physics", "", "G4int");
117
118 G4String IKE("IKE");
119 (*store)[IKE] = G4AttDef(IKE, "Initial kinetic energy", "Physics", "G4BestUnit", "G4double");
120
121 G4String IMom("IMom");
122 (*store)[IMom] = G4AttDef(IMom, "Initial momentum", "Physics", "G4BestUnit", "G4ThreeVector");
123
124 G4String IMag("IMag");
125 (*store)[IMag] =
126 G4AttDef(IMag, "Initial momentum magnitude", "Physics", "G4BestUnit", "G4double");
127
128 G4String NTP("NTP");
129 (*store)[NTP] = G4AttDef(NTP, "No. of points", "Physics", "", "G4int");
130 }
131 return store;
132}
bool G4bool
Definition G4Types.hh:86
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)

Referenced by CreateAttValues().

◆ GetCharge()

G4double G4ClonedTrajectory::GetCharge ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 82 of file G4ClonedTrajectory.hh.

82{ return PDGCharge; }

◆ GetInitialKineticEnergy()

G4double G4ClonedTrajectory::GetInitialKineticEnergy ( ) const
inline

Definition at line 84 of file G4ClonedTrajectory.hh.

84{ return initialKineticEnergy; }

◆ GetInitialMomentum()

G4ThreeVector G4ClonedTrajectory::GetInitialMomentum ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 85 of file G4ClonedTrajectory.hh.

85{ return initialMomentum; }

◆ GetParentID()

G4int G4ClonedTrajectory::GetParentID ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 80 of file G4ClonedTrajectory.hh.

80{ return fParentID; }

◆ GetParticleDefinition()

G4ParticleDefinition * G4ClonedTrajectory::GetParticleDefinition ( )

Definition at line 168 of file G4ClonedTrajectory.cc.

169{
170 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
171}
static G4ParticleTable * GetParticleTable()

◆ GetParticleName()

G4String G4ClonedTrajectory::GetParticleName ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 81 of file G4ClonedTrajectory.hh.

81{ return ParticleName; }

◆ GetPDGEncoding()

G4int G4ClonedTrajectory::GetPDGEncoding ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 83 of file G4ClonedTrajectory.hh.

83{ return PDGEncoding; }

◆ GetPoint()

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

Implements G4VTrajectory.

Definition at line 93 of file G4ClonedTrajectory.hh.

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

◆ GetPointEntries()

G4int G4ClonedTrajectory::GetPointEntries ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 92 of file G4ClonedTrajectory.hh.

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

Referenced by CreateAttValues().

◆ GetTrackID()

G4int G4ClonedTrajectory::GetTrackID ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 79 of file G4ClonedTrajectory.hh.

79{ return fTrackID; }

◆ MergeTrajectory()

void G4ClonedTrajectory::MergeTrajectory ( G4VTrajectory * secondTrajectory)
overridevirtual

Implements G4VTrajectory.

Definition at line 173 of file G4ClonedTrajectory.cc.

174{
175 if (secondTrajectory == nullptr) return;
176
177 auto seco = (G4ClonedTrajectory*)secondTrajectory;
178 G4int ent = seco->GetPointEntries();
179 for (G4int i = 1; i < ent; ++i) // initial pt of 2nd trajectory shouldn't be merged
180 {
181 positionRecord->push_back((*(seco->positionRecord))[i]);
182 }
183 delete (*seco->positionRecord)[0];
184 seco->positionRecord->clear();
185}
G4ClonedTrajectory()=default

◆ operator delete()

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

Definition at line 122 of file G4ClonedTrajectory.hh.

123{
124 aClonedTrajectoryAllocator()->FreeSingle((G4ClonedTrajectory*)aTrajectory);
125}
G4TRACKING_DLL G4Allocator< G4ClonedTrajectory > *& aClonedTrajectoryAllocator()

◆ operator new()

void * G4ClonedTrajectory::operator new ( size_t )
inline

Definition at line 114 of file G4ClonedTrajectory.hh.

115{
116 if (aClonedTrajectoryAllocator() == nullptr) {
117 aClonedTrajectoryAllocator() = new G4Allocator<G4ClonedTrajectory>;
118 }
119 return (void*)aClonedTrajectoryAllocator()->MallocSingle();
120}

◆ operator==()

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

Definition at line 127 of file G4ClonedTrajectory.hh.

127{ return (this == &r); }

◆ ShowTrajectory()

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

Reimplemented from G4VTrajectory.

Definition at line 82 of file G4ClonedTrajectory.cc.

83{
84 // Invoke the default implementation in G4VTrajectory...
86
87 // ... or override with your own code here.
88}
virtual void ShowTrajectory(std::ostream &os=G4cout) const

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