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

#include <G4ClonedSmoothTrajectory.hh>

+ Inheritance diagram for G4ClonedSmoothTrajectory:

Public Member Functions

 G4ClonedSmoothTrajectory ()=default
 
 G4ClonedSmoothTrajectory (const G4SmoothTrajectory &)
 
 ~G4ClonedSmoothTrajectory () override
 
G4bool operator== (const G4ClonedSmoothTrajectory &r) const
 
void * operator new (size_t)
 
void operator delete (void *)
 
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 58 of file G4ClonedSmoothTrajectory.hh.

Constructor & Destructor Documentation

◆ G4ClonedSmoothTrajectory() [1/2]

G4ClonedSmoothTrajectory::G4ClonedSmoothTrajectory ( )
default

◆ G4ClonedSmoothTrajectory() [2/2]

G4ClonedSmoothTrajectory::G4ClonedSmoothTrajectory ( const G4SmoothTrajectory & right)

Definition at line 54 of file G4ClonedSmoothTrajectory.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 G4TrajectoryPointContainer();
64
65 for (auto& i : *right.positionRecord) {
66 auto rightPoint = (G4ClonedSmoothTrajectoryPoint*)i;
67 positionRecord->push_back(new G4ClonedSmoothTrajectoryPoint(*rightPoint));
68 }
69}

◆ ~G4ClonedSmoothTrajectory()

G4ClonedSmoothTrajectory::~G4ClonedSmoothTrajectory ( )
override

Definition at line 71 of file G4ClonedSmoothTrajectory.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 G4ClonedSmoothTrajectory::AppendStep ( const G4Step * aStep)
overridevirtual

Implements G4VTrajectory.

Definition at line 165 of file G4ClonedSmoothTrajectory.cc.

166{
167 positionRecord->push_back(new G4ClonedSmoothTrajectoryPoint(
169}
const G4ThreeVector & GetPosition() const
std::vector< G4ThreeVector > * GetPointerToVectorOfAuxiliaryPoints() const
G4StepPoint * GetPostStepPoint() const

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Definition at line 136 of file G4ClonedSmoothTrajectory.cc.

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

Reimplemented from G4VTrajectory.

Definition at line 91 of file G4ClonedSmoothTrajectory.cc.

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

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Definition at line 100 of file G4ClonedSmoothTrajectory.cc.

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

Referenced by CreateAttValues().

◆ GetCharge()

G4double G4ClonedSmoothTrajectory::GetCharge ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 80 of file G4ClonedSmoothTrajectory.hh.

80{ return PDGCharge; }

◆ GetInitialKineticEnergy()

G4double G4ClonedSmoothTrajectory::GetInitialKineticEnergy ( ) const
inline

Definition at line 82 of file G4ClonedSmoothTrajectory.hh.

82{ return initialKineticEnergy; }

◆ GetInitialMomentum()

G4ThreeVector G4ClonedSmoothTrajectory::GetInitialMomentum ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 83 of file G4ClonedSmoothTrajectory.hh.

83{ return initialMomentum; }

◆ GetParentID()

G4int G4ClonedSmoothTrajectory::GetParentID ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 78 of file G4ClonedSmoothTrajectory.hh.

78{ return fParentID; }

◆ GetParticleDefinition()

G4ParticleDefinition * G4ClonedSmoothTrajectory::GetParticleDefinition ( )

Definition at line 171 of file G4ClonedSmoothTrajectory.cc.

172{
173 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
174}
static G4ParticleTable * GetParticleTable()

◆ GetParticleName()

G4String G4ClonedSmoothTrajectory::GetParticleName ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 79 of file G4ClonedSmoothTrajectory.hh.

79{ return ParticleName; }

◆ GetPDGEncoding()

G4int G4ClonedSmoothTrajectory::GetPDGEncoding ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 81 of file G4ClonedSmoothTrajectory.hh.

81{ return PDGEncoding; }

◆ GetPoint()

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

Implements G4VTrajectory.

Definition at line 91 of file G4ClonedSmoothTrajectory.hh.

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

◆ GetPointEntries()

G4int G4ClonedSmoothTrajectory::GetPointEntries ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 90 of file G4ClonedSmoothTrajectory.hh.

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

Referenced by CreateAttValues().

◆ GetTrackID()

G4int G4ClonedSmoothTrajectory::GetTrackID ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 77 of file G4ClonedSmoothTrajectory.hh.

77{ return fTrackID; }

◆ MergeTrajectory()

void G4ClonedSmoothTrajectory::MergeTrajectory ( G4VTrajectory * secondTrajectory)
overridevirtual

Implements G4VTrajectory.

Definition at line 176 of file G4ClonedSmoothTrajectory.cc.

177{
178 if (secondTrajectory == nullptr) return;
179
180 auto seco = (G4ClonedSmoothTrajectory*)secondTrajectory;
181 G4int ent = seco->GetPointEntries();
182
183 // initial point of the second trajectory should not be merged
184 //
185 for (G4int i = 1; i < ent; ++i) {
186 positionRecord->push_back((*(seco->positionRecord))[i]);
187 }
188 delete (*seco->positionRecord)[0];
189 seco->positionRecord->clear();
190}
G4ClonedSmoothTrajectory()=default

◆ operator delete()

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

Definition at line 122 of file G4ClonedSmoothTrajectory.hh.

123{
125}
G4TRACKING_DLL G4Allocator< G4ClonedSmoothTrajectory > *& aClonedSmoothTrajectoryAllocator()

◆ operator new()

void * G4ClonedSmoothTrajectory::operator new ( size_t )
inline

Definition at line 114 of file G4ClonedSmoothTrajectory.hh.

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

◆ operator==()

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

Definition at line 127 of file G4ClonedSmoothTrajectory.hh.

128{
129 return (this == &r);
130}

◆ ShowTrajectory()

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

Reimplemented from G4VTrajectory.

Definition at line 82 of file G4ClonedSmoothTrajectory.cc.

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

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