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

#include <G4ClonedRichTrajectory.hh>

+ Inheritance diagram for G4ClonedRichTrajectory:

Public Member Functions

 G4ClonedRichTrajectory ()=default
 
 ~G4ClonedRichTrajectory () override
 
 G4ClonedRichTrajectory (const G4RichTrajectory &)
 
G4ClonedRichTrajectoryoperator= (const G4ClonedRichTrajectory &)=delete
 
G4bool operator== (const G4ClonedRichTrajectory &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
 
void MergeTrajectory (G4VTrajectory *secondTrajectory) override
 
G4int GetPointEntries () const override
 
G4VTrajectoryPointGetPoint (G4int i) const 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 G4ClonedRichTrajectory.hh.

Constructor & Destructor Documentation

◆ G4ClonedRichTrajectory() [1/2]

G4ClonedRichTrajectory::G4ClonedRichTrajectory ( )
default

◆ ~G4ClonedRichTrajectory()

G4ClonedRichTrajectory::~G4ClonedRichTrajectory ( )
override

Definition at line 89 of file G4ClonedRichTrajectory.cc.

90{
91 if (fpRichPointContainer != nullptr) {
92 for (auto& i : *fpRichPointContainer) {
93 delete i;
94 }
95 fpRichPointContainer->clear();
96 delete fpRichPointContainer;
97 }
98}

◆ G4ClonedRichTrajectory() [2/2]

G4ClonedRichTrajectory::G4ClonedRichTrajectory ( const G4RichTrajectory & right)

Definition at line 57 of file G4ClonedRichTrajectory.cc.

58{
59 ParticleName = right.ParticleName;
60 PDGCharge = right.PDGCharge;
61 PDGEncoding = right.PDGEncoding;
62 fTrackID = right.fTrackID;
63 fParentID = right.fParentID;
64 initialKineticEnergy = right.initialKineticEnergy;
65 initialMomentum = right.initialMomentum;
66 //positionRecord = new G4ClonedTrajectoryPointContainer();
67 //for (auto& i : *right.positionRecord) {
68 // auto rightPoint = (G4ClonedTrajectoryPoint*)i;
69 // positionRecord->push_back(new G4ClonedTrajectoryPoint(*rightPoint));
70 //}
71 fpInitialVolume = right.fpInitialVolume;
72 fpInitialNextVolume = right.fpInitialNextVolume;
73 fpCreatorProcess = right.fpCreatorProcess;
74 fCreatorModelID = right.fCreatorModelID;
75 fpFinalVolume = right.fpFinalVolume;
76 fpFinalNextVolume = right.fpFinalNextVolume;
77 fpEndingProcess = right.fpEndingProcess;
78 fFinalKineticEnergy = right.fFinalKineticEnergy;
79 if(right.fpRichPointContainer!=nullptr && right.fpRichPointContainer->size()>0)
80 {
81 fpRichPointContainer = new G4TrajectoryPointContainer;
82 for (auto& i : *right.fpRichPointContainer) {
83 auto rightPoint = (G4RichTrajectoryPoint*)i;
84 fpRichPointContainer->push_back(new G4ClonedRichTrajectoryPoint(*rightPoint));
85 }
86 }
87}

Member Function Documentation

◆ AppendStep()

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

Implements G4VTrajectory.

Definition at line 100 of file G4ClonedRichTrajectory.cc.

101{
102 fpRichPointContainer->push_back(new G4ClonedRichTrajectoryPoint(aStep));
103
104 // Except for first step, which is a sort of virtual step to start
105 // the track, compute the final values...
106 //
107 const G4Track* track = aStep->GetTrack();
108 const G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
109 if (track->GetCurrentStepNumber() > 0) {
110 fpFinalVolume = track->GetTouchableHandle();
111 fpFinalNextVolume = track->GetNextTouchableHandle();
112 fpEndingProcess = postStepPoint->GetProcessDefinedStep();
113 fFinalKineticEnergy =
115 }
116}
const G4VProcess * GetProcessDefinedStep() const
G4double GetKineticEnergy() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetTotalEnergyDeposit() const
G4StepPoint * GetPostStepPoint() const
const G4TouchableHandle & GetNextTouchableHandle() const
G4int GetCurrentStepNumber() const
const G4TouchableHandle & GetTouchableHandle() const

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Definition at line 233 of file G4ClonedRichTrajectory.cc.

234{
235 // Create base class att values...
236 //std::vector<G4AttValue>* values = G4VTrajectory::CreateAttValues();
237 auto values = new std::vector<G4AttValue>;
238 values->push_back(G4AttValue("ID", G4UIcommand::ConvertToString(fTrackID), ""));
239 values->push_back(G4AttValue("PID", G4UIcommand::ConvertToString(fParentID), ""));
240 values->push_back(G4AttValue("PN", ParticleName, ""));
241 values->push_back(G4AttValue("Ch", G4UIcommand::ConvertToString(PDGCharge), ""));
242 values->push_back(G4AttValue("PDG", G4UIcommand::ConvertToString(PDGEncoding), ""));
243 values->push_back(G4AttValue("IKE", G4BestUnit(initialKineticEnergy, "Energy"), ""));
244 values->push_back(G4AttValue("IMom", G4BestUnit(initialMomentum, "Energy"), ""));
245 values->push_back(G4AttValue("IMag", G4BestUnit(initialMomentum.mag(), "Energy"), ""));
246 values->push_back(G4AttValue("NTP", G4UIcommand::ConvertToString(GetPointEntries()), ""));
247
248 if (fpInitialVolume && (fpInitialVolume->GetVolume() != nullptr)) {
249 values->push_back(G4AttValue("IVPath", Path(fpInitialVolume), ""));
250 }
251 else {
252 values->push_back(G4AttValue("IVPath", "None", ""));
253 }
254
255 if (fpInitialNextVolume && (fpInitialNextVolume->GetVolume() != nullptr)) {
256 values->push_back(G4AttValue("INVPath", Path(fpInitialNextVolume), ""));
257 }
258 else {
259 values->push_back(G4AttValue("INVPath", "None", ""));
260 }
261
262 if (fpCreatorProcess != nullptr) {
263 values->push_back(G4AttValue("CPN", fpCreatorProcess->GetProcessName(), ""));
264 G4ProcessType type = fpCreatorProcess->GetProcessType();
265 values->push_back(G4AttValue("CPTN", G4VProcess::GetProcessTypeName(type), ""));
266 values->push_back(G4AttValue("CMID", G4UIcommand::ConvertToString(fCreatorModelID), ""));
267 const G4String& creatorModelName = G4PhysicsModelCatalog::GetModelNameFromID(fCreatorModelID);
268 values->push_back(G4AttValue("CMN", creatorModelName, ""));
269 }
270 else {
271 values->push_back(G4AttValue("CPN", "None", ""));
272 values->push_back(G4AttValue("CPTN", "None", ""));
273 values->push_back(G4AttValue("CMID", "None", ""));
274 values->push_back(G4AttValue("CMN", "None", ""));
275 }
276
277 if (fpFinalVolume && (fpFinalVolume->GetVolume() != nullptr)) {
278 values->push_back(G4AttValue("FVPath", Path(fpFinalVolume), ""));
279 }
280 else {
281 values->push_back(G4AttValue("FVPath", "None", ""));
282 }
283
284 if (fpFinalNextVolume && (fpFinalNextVolume->GetVolume() != nullptr)) {
285 values->push_back(G4AttValue("FNVPath", Path(fpFinalNextVolume), ""));
286 }
287 else {
288 values->push_back(G4AttValue("FNVPath", "None", ""));
289 }
290
291 if (fpEndingProcess != nullptr) {
292 values->push_back(G4AttValue("EPN", fpEndingProcess->GetProcessName(), ""));
293 G4ProcessType type = fpEndingProcess->GetProcessType();
294 values->push_back(G4AttValue("EPTN", G4VProcess::GetProcessTypeName(type), ""));
295 }
296 else {
297 values->push_back(G4AttValue("EPN", "None", ""));
298 values->push_back(G4AttValue("EPTN", "None", ""));
299 }
300
301 values->push_back(G4AttValue("FKE", G4BestUnit(fFinalKineticEnergy, "Energy"), ""));
302
303#ifdef G4ATTDEBUG
304 G4cout << G4AttCheck(values, GetAttDefs());
305#endif
306
307 return values;
308}
G4ProcessType
#define G4BestUnit(a, b)
G4GLOB_DLL std::ostream G4cout
const std::map< G4String, G4AttDef > * GetAttDefs() const override
G4int GetPointEntries() const override
static const G4String GetModelNameFromID(const G4int modelID)
static G4String ConvertToString(G4bool boolVal)
static const G4String & GetProcessTypeName(G4ProcessType)

◆ DrawTrajectory()

void G4ClonedRichTrajectory::DrawTrajectory ( ) const
overridevirtual

Reimplemented from G4VTrajectory.

Definition at line 142 of file G4ClonedRichTrajectory.cc.

143{
144 // Invoke the default implementation in G4VTrajectory...
145 //
147
148 // ... or override with your own code here.
149}
virtual void DrawTrajectory() const

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Definition at line 151 of file G4ClonedRichTrajectory.cc.

152{
153 G4bool isNew;
154 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("G4ClonedRichTrajectory", isNew);
155 if (isNew) {
156 G4String ID;
157
158 ID = "ID";
159 (*store)[ID] = G4AttDef(ID, "Track ID", "Physics", "", "G4int");
160
161 ID = "PID";
162 (*store)[ID] = G4AttDef(ID, "Parent ID", "Physics", "", "G4int");
163
164 ID = "PN";
165 (*store)[ID] = G4AttDef(ID, "Particle Name", "Physics", "", "G4String");
166
167 ID = "Ch";
168 (*store)[ID] = G4AttDef(ID, "Charge", "Physics", "e+", "G4double");
169
170 ID = "PDG";
171 (*store)[ID] = G4AttDef(ID, "PDG Encoding", "Physics", "", "G4int");
172
173 ID = "IKE";
174 (*store)[ID] = G4AttDef(ID, "Initial kinetic energy", "Physics", "G4BestUnit", "G4double");
175
176 ID = "IMom";
177 (*store)[ID] = G4AttDef(ID, "Initial momentum", "Physics", "G4BestUnit", "G4ThreeVector");
178
179 ID = "IMag";
180 (*store)[ID] = G4AttDef(ID, "Initial momentum magnitude", "Physics", "G4BestUnit", "G4double");
181
182 ID = "NTP";
183 (*store)[ID] = G4AttDef(ID, "No. of points", "Physics", "", "G4int");
184
185 ID = "IVPath";
186 (*store)[ID] = G4AttDef(ID, "Initial Volume Path", "Physics", "", "G4String");
187
188 ID = "INVPath";
189 (*store)[ID] = G4AttDef(ID, "Initial Next Volume Path", "Physics", "", "G4String");
190
191 ID = "CPN";
192 (*store)[ID] = G4AttDef(ID, "Creator Process Name", "Physics", "", "G4String");
193
194 ID = "CPTN";
195 (*store)[ID] = G4AttDef(ID, "Creator Process Type Name", "Physics", "", "G4String");
196
197 ID = "CMID";
198 (*store)[ID] = G4AttDef(ID, "Creator Model ID", "Physics", "", "G4int");
199
200 ID = "CMN";
201 (*store)[ID] = G4AttDef(ID, "Creator Model Name", "Physics", "", "G4String");
202
203 ID = "FVPath";
204 (*store)[ID] = G4AttDef(ID, "Final Volume Path", "Physics", "", "G4String");
205
206 ID = "FNVPath";
207 (*store)[ID] = G4AttDef(ID, "Final Next Volume Path", "Physics", "", "G4String");
208
209 ID = "EPN";
210 (*store)[ID] = G4AttDef(ID, "Ending Process Name", "Physics", "", "G4String");
211
212 ID = "EPTN";
213 (*store)[ID] = G4AttDef(ID, "Ending Process Type Name", "Physics", "", "G4String");
214
215 ID = "FKE";
216 (*store)[ID] = G4AttDef(ID, "Final kinetic energy", "Physics", "G4BestUnit", "G4double");
217 }
218
219 return store;
220}
bool G4bool
Definition G4Types.hh:86
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)

Referenced by CreateAttValues().

◆ GetCharge()

G4double G4ClonedRichTrajectory::GetCharge ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 82 of file G4ClonedRichTrajectory.hh.

82{ return PDGCharge; }

◆ GetInitialKineticEnergy()

G4double G4ClonedRichTrajectory::GetInitialKineticEnergy ( ) const
inline

Definition at line 84 of file G4ClonedRichTrajectory.hh.

84{ return initialKineticEnergy; }

◆ GetInitialMomentum()

G4ThreeVector G4ClonedRichTrajectory::GetInitialMomentum ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 85 of file G4ClonedRichTrajectory.hh.

85{ return initialMomentum; }

◆ GetParentID()

G4int G4ClonedRichTrajectory::GetParentID ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 80 of file G4ClonedRichTrajectory.hh.

80{ return fParentID; }

◆ GetParticleDefinition()

G4ParticleDefinition * G4ClonedRichTrajectory::GetParticleDefinition ( )

Definition at line 310 of file G4ClonedRichTrajectory.cc.

311{
312 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
313}
static G4ParticleTable * GetParticleTable()

◆ GetParticleName()

G4String G4ClonedRichTrajectory::GetParticleName ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 81 of file G4ClonedRichTrajectory.hh.

81{ return ParticleName; }

◆ GetPDGEncoding()

G4int G4ClonedRichTrajectory::GetPDGEncoding ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 83 of file G4ClonedRichTrajectory.hh.

83{ return PDGEncoding; }

◆ GetPoint()

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

Implements G4VTrajectory.

Definition at line 146 of file G4ClonedRichTrajectory.hh.

147{
148 return (*fpRichPointContainer)[i];
149}

◆ GetPointEntries()

G4int G4ClonedRichTrajectory::GetPointEntries ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 141 of file G4ClonedRichTrajectory.hh.

142{
143 return G4int(fpRichPointContainer->size());
144}
int G4int
Definition G4Types.hh:85

Referenced by CreateAttValues().

◆ GetTrackID()

G4int G4ClonedRichTrajectory::GetTrackID ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 79 of file G4ClonedRichTrajectory.hh.

79{ return fTrackID; }

◆ MergeTrajectory()

void G4ClonedRichTrajectory::MergeTrajectory ( G4VTrajectory * secondTrajectory)
overridevirtual

Implements G4VTrajectory.

Definition at line 118 of file G4ClonedRichTrajectory.cc.

119{
120 if (secondTrajectory == nullptr) return;
121
122 auto seco = (G4ClonedRichTrajectory*)secondTrajectory;
123 G4int ent = seco->GetPointEntries();
124 for (G4int i = 1; i < ent; ++i) {
125 // initial point of the second trajectory should not be merged
126 //
127 fpRichPointContainer->push_back((*(seco->fpRichPointContainer))[i]);
128 }
129 delete (*seco->fpRichPointContainer)[0];
130 seco->fpRichPointContainer->clear();
131}
G4ClonedRichTrajectory()=default

◆ operator delete()

void G4ClonedRichTrajectory::operator delete ( void * aRichTrajectory)
inline

Definition at line 136 of file G4ClonedRichTrajectory.hh.

137{
138 aClonedRichTrajectoryAllocator()->FreeSingle((G4ClonedRichTrajectory*)aRichTrajectory);
139}
G4TRACKING_DLL G4Allocator< G4ClonedRichTrajectory > *& aClonedRichTrajectoryAllocator()

◆ operator new()

void * G4ClonedRichTrajectory::operator new ( size_t )
inline

Definition at line 128 of file G4ClonedRichTrajectory.hh.

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

◆ operator=()

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

◆ operator==()

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

Definition at line 72 of file G4ClonedRichTrajectory.hh.

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

◆ ShowTrajectory()

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

Reimplemented from G4VTrajectory.

Definition at line 133 of file G4ClonedRichTrajectory.cc.

134{
135 // Invoke the default implementation in G4VTrajectory...
136 //
138
139 // ... or override with your own code here.
140}
virtual void ShowTrajectory(std::ostream &os=G4cout) const

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