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

#include <G4RichTrajectory.hh>

+ Inheritance diagram for G4RichTrajectory:

Public Member Functions

 G4RichTrajectory ()=default
 
 G4RichTrajectory (const G4Track *aTrack)
 
 ~G4RichTrajectory () override
 
 G4RichTrajectory (G4RichTrajectory &)
 
G4RichTrajectoryoperator= (const G4RichTrajectory &)=delete
 
G4bool operator== (const G4RichTrajectory &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
 
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
 

Friends

class G4ClonedRichTrajectory
 

Detailed Description

Definition at line 70 of file G4RichTrajectory.hh.

Constructor & Destructor Documentation

◆ G4RichTrajectory() [1/3]

G4RichTrajectory::G4RichTrajectory ( )
default

◆ G4RichTrajectory() [2/3]

G4RichTrajectory::G4RichTrajectory ( const G4Track * aTrack)

Definition at line 69 of file G4RichTrajectory.cc.

70{
71 G4ParticleDefinition* fpParticleDefinition = aTrack->GetDefinition();
72 ParticleName = fpParticleDefinition->GetParticleName();
73 PDGCharge = fpParticleDefinition->GetPDGCharge();
74 PDGEncoding = fpParticleDefinition->GetPDGEncoding();
75 fTrackID = aTrack->GetTrackID();
76 fParentID = aTrack->GetParentID();
77 initialKineticEnergy = aTrack->GetKineticEnergy();
78 initialMomentum = aTrack->GetMomentum();
79 positionRecord = new G4TrajectoryPointContainer();
80
81 // Following is for the first trajectory point
82 positionRecord->push_back(new G4RichTrajectoryPoint(aTrack));
83
84 fpInitialVolume = aTrack->GetTouchableHandle();
85 fpInitialNextVolume = aTrack->GetNextTouchableHandle();
86 fpCreatorProcess = aTrack->GetCreatorProcess();
87 fCreatorModelID = aTrack->GetCreatorModelID();
88
89 // On construction, set final values to initial values.
90 // Final values are updated at the addition of every step - see AppendStep.
91 //
92 fpFinalVolume = aTrack->GetTouchableHandle();
93 fpFinalNextVolume = aTrack->GetNextTouchableHandle();
94 fpEndingProcess = aTrack->GetCreatorProcess();
95 fFinalKineticEnergy = aTrack->GetKineticEnergy();
96
97 // Insert the first rich trajectory point (see note above)...
98 //
99 fpRichPointContainer = new G4TrajectoryPointContainer;
100 fpRichPointContainer->push_back(new G4RichTrajectoryPoint(aTrack));
101}
const G4String & GetParticleName() const
G4int GetTrackID() const
const G4TouchableHandle & GetNextTouchableHandle() const
const G4VProcess * GetCreatorProcess() const
G4int GetCreatorModelID() const
G4ThreeVector GetMomentum() const
G4ParticleDefinition * GetDefinition() const
const G4TouchableHandle & GetTouchableHandle() const
G4double GetKineticEnergy() const
G4int GetParentID() const

◆ ~G4RichTrajectory()

G4RichTrajectory::~G4RichTrajectory ( )
override

Definition at line 134 of file G4RichTrajectory.cc.

135{
136 if (fpRichPointContainer != nullptr) {
137 for (auto& i : *fpRichPointContainer) {
138 delete i;
139 }
140 fpRichPointContainer->clear();
141 delete fpRichPointContainer;
142 }
143}

◆ G4RichTrajectory() [3/3]

G4RichTrajectory::G4RichTrajectory ( G4RichTrajectory & right)

Definition at line 103 of file G4RichTrajectory.cc.

104{
105 ParticleName = right.ParticleName;
106 PDGCharge = right.PDGCharge;
107 PDGEncoding = right.PDGEncoding;
108 fTrackID = right.fTrackID;
109 fParentID = right.fParentID;
110 initialKineticEnergy = right.initialKineticEnergy;
111 initialMomentum = right.initialMomentum;
112 positionRecord = new G4TrajectoryPointContainer();
113
114 for (auto& i : *right.positionRecord) {
115 auto rightPoint = (G4RichTrajectoryPoint*)i;
116 positionRecord->push_back(new G4RichTrajectoryPoint(*rightPoint));
117 }
118
119 fpInitialVolume = right.fpInitialVolume;
120 fpInitialNextVolume = right.fpInitialNextVolume;
121 fpCreatorProcess = right.fpCreatorProcess;
122 fCreatorModelID = right.fCreatorModelID;
123 fpFinalVolume = right.fpFinalVolume;
124 fpFinalNextVolume = right.fpFinalNextVolume;
125 fpEndingProcess = right.fpEndingProcess;
126 fFinalKineticEnergy = right.fFinalKineticEnergy;
127 fpRichPointContainer = new G4TrajectoryPointContainer;
128 for (auto& i : *right.fpRichPointContainer) {
129 auto rightPoint = (G4RichTrajectoryPoint*)i;
130 fpRichPointContainer->push_back(new G4RichTrajectoryPoint(*rightPoint));
131 }
132}

Member Function Documentation

◆ AppendStep()

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

Implements G4VTrajectory.

Definition at line 145 of file G4RichTrajectory.cc.

146{
147 fpRichPointContainer->push_back(new G4RichTrajectoryPoint(aStep));
148
149 // Except for first step, which is a sort of virtual step to start
150 // the track, compute the final values...
151 //
152 const G4Track* track = aStep->GetTrack();
153 const G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
154 if (track->GetCurrentStepNumber() > 0) {
155 fpFinalVolume = track->GetTouchableHandle();
156 fpFinalNextVolume = track->GetNextTouchableHandle();
157 fpEndingProcess = postStepPoint->GetProcessDefinedStep();
158 fFinalKineticEnergy =
160 }
161}
const G4VProcess * GetProcessDefinedStep() const
G4double GetKineticEnergy() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetTotalEnergyDeposit() const
G4StepPoint * GetPostStepPoint() const
G4int GetCurrentStepNumber() const

◆ CloneForMaster()

G4VTrajectory * G4RichTrajectory::CloneForMaster ( ) const
overridevirtual

Reimplemented from G4VTrajectory.

Definition at line 360 of file G4RichTrajectory.cc.

361{
362 G4AutoLock lock(&CloneRichTrajectoryMutex);
363 auto* cloned = new G4ClonedRichTrajectory(*this);
364 return cloned;
365}
G4TemplateAutoLock< G4Mutex > G4AutoLock
friend class G4ClonedRichTrajectory

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Definition at line 278 of file G4RichTrajectory.cc.

279{
280 // Create base class att values...
281 //std::vector<G4AttValue>* values = G4VTrajectory::CreateAttValues();
282 auto values = new std::vector<G4AttValue>;
283 values->push_back(G4AttValue("ID", G4UIcommand::ConvertToString(fTrackID), ""));
284 values->push_back(G4AttValue("PID", G4UIcommand::ConvertToString(fParentID), ""));
285 values->push_back(G4AttValue("PN", ParticleName, ""));
286 values->push_back(G4AttValue("Ch", G4UIcommand::ConvertToString(PDGCharge), ""));
287 values->push_back(G4AttValue("PDG", G4UIcommand::ConvertToString(PDGEncoding), ""));
288 values->push_back(G4AttValue("IKE", G4BestUnit(initialKineticEnergy, "Energy"), ""));
289 values->push_back(G4AttValue("IMom", G4BestUnit(initialMomentum, "Energy"), ""));
290 values->push_back(G4AttValue("IMag", G4BestUnit(initialMomentum.mag(), "Energy"), ""));
291 values->push_back(G4AttValue("NTP", G4UIcommand::ConvertToString(GetPointEntries()), ""));
292
293 if (fpInitialVolume && (fpInitialVolume->GetVolume() != nullptr)) {
294 values->push_back(G4AttValue("IVPath", Path(fpInitialVolume), ""));
295 }
296 else {
297 values->push_back(G4AttValue("IVPath", "None", ""));
298 }
299
300 if (fpInitialNextVolume && (fpInitialNextVolume->GetVolume() != nullptr)) {
301 values->push_back(G4AttValue("INVPath", Path(fpInitialNextVolume), ""));
302 }
303 else {
304 values->push_back(G4AttValue("INVPath", "None", ""));
305 }
306
307 if (fpCreatorProcess != nullptr) {
308 values->push_back(G4AttValue("CPN", fpCreatorProcess->GetProcessName(), ""));
309 G4ProcessType type = fpCreatorProcess->GetProcessType();
310 values->push_back(G4AttValue("CPTN", G4VProcess::GetProcessTypeName(type), ""));
311 values->push_back(G4AttValue("CMID", G4UIcommand::ConvertToString(fCreatorModelID), ""));
312 const G4String& creatorModelName = G4PhysicsModelCatalog::GetModelNameFromID(fCreatorModelID);
313 values->push_back(G4AttValue("CMN", creatorModelName, ""));
314 }
315 else {
316 values->push_back(G4AttValue("CPN", "None", ""));
317 values->push_back(G4AttValue("CPTN", "None", ""));
318 values->push_back(G4AttValue("CMID", "None", ""));
319 values->push_back(G4AttValue("CMN", "None", ""));
320 }
321
322 if (fpFinalVolume && (fpFinalVolume->GetVolume() != nullptr)) {
323 values->push_back(G4AttValue("FVPath", Path(fpFinalVolume), ""));
324 }
325 else {
326 values->push_back(G4AttValue("FVPath", "None", ""));
327 }
328
329 if (fpFinalNextVolume && (fpFinalNextVolume->GetVolume() != nullptr)) {
330 values->push_back(G4AttValue("FNVPath", Path(fpFinalNextVolume), ""));
331 }
332 else {
333 values->push_back(G4AttValue("FNVPath", "None", ""));
334 }
335
336 if (fpEndingProcess != nullptr) {
337 values->push_back(G4AttValue("EPN", fpEndingProcess->GetProcessName(), ""));
338 G4ProcessType type = fpEndingProcess->GetProcessType();
339 values->push_back(G4AttValue("EPTN", G4VProcess::GetProcessTypeName(type), ""));
340 }
341 else {
342 values->push_back(G4AttValue("EPN", "None", ""));
343 values->push_back(G4AttValue("EPTN", "None", ""));
344 }
345
346 values->push_back(G4AttValue("FKE", G4BestUnit(fFinalKineticEnergy, "Energy"), ""));
347
348#ifdef G4ATTDEBUG
349 G4cout << G4AttCheck(values, GetAttDefs());
350#endif
351
352 return values;
353}
G4ProcessType
#define G4BestUnit(a, b)
G4GLOB_DLL std::ostream G4cout
static const G4String GetModelNameFromID(const G4int modelID)
const std::map< G4String, G4AttDef > * GetAttDefs() const override
G4int GetPointEntries() const override
static G4String ConvertToString(G4bool boolVal)
static const G4String & GetProcessTypeName(G4ProcessType)

◆ DrawTrajectory()

void G4RichTrajectory::DrawTrajectory ( ) const
overridevirtual

Reimplemented from G4VTrajectory.

Definition at line 187 of file G4RichTrajectory.cc.

188{
189 // Invoke the default implementation in G4VTrajectory...
190 //
192
193 // ... or override with your own code here.
194}
virtual void DrawTrajectory() const

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Definition at line 196 of file G4RichTrajectory.cc.

197{
198 G4bool isNew;
199 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("G4RichTrajectory", isNew);
200 if (isNew) {
201 G4String ID;
202
203 ID = "ID";
204 (*store)[ID] = G4AttDef(ID, "Track ID", "Physics", "", "G4int");
205
206 ID = "PID";
207 (*store)[ID] = G4AttDef(ID, "Parent ID", "Physics", "", "G4int");
208
209 ID = "PN";
210 (*store)[ID] = G4AttDef(ID, "Particle Name", "Physics", "", "G4String");
211
212 ID = "Ch";
213 (*store)[ID] = G4AttDef(ID, "Charge", "Physics", "e+", "G4double");
214
215 ID = "PDG";
216 (*store)[ID] = G4AttDef(ID, "PDG Encoding", "Physics", "", "G4int");
217
218 ID = "IKE";
219 (*store)[ID] = G4AttDef(ID, "Initial kinetic energy", "Physics", "G4BestUnit", "G4double");
220
221 ID = "IMom";
222 (*store)[ID] = G4AttDef(ID, "Initial momentum", "Physics", "G4BestUnit", "G4ThreeVector");
223
224 ID = "IMag";
225 (*store)[ID] = G4AttDef(ID, "Initial momentum magnitude", "Physics", "G4BestUnit", "G4double");
226
227 ID = "NTP";
228 (*store)[ID] = G4AttDef(ID, "No. of points", "Physics", "", "G4int");
229
230 ID = "IVPath";
231 (*store)[ID] = G4AttDef(ID, "Initial Volume Path", "Physics", "", "G4String");
232
233 ID = "INVPath";
234 (*store)[ID] = G4AttDef(ID, "Initial Next Volume Path", "Physics", "", "G4String");
235
236 ID = "CPN";
237 (*store)[ID] = G4AttDef(ID, "Creator Process Name", "Physics", "", "G4String");
238
239 ID = "CPTN";
240 (*store)[ID] = G4AttDef(ID, "Creator Process Type Name", "Physics", "", "G4String");
241
242 ID = "CMID";
243 (*store)[ID] = G4AttDef(ID, "Creator Model ID", "Physics", "", "G4int");
244
245 ID = "CMN";
246 (*store)[ID] = G4AttDef(ID, "Creator Model Name", "Physics", "", "G4String");
247
248 ID = "FVPath";
249 (*store)[ID] = G4AttDef(ID, "Final Volume Path", "Physics", "", "G4String");
250
251 ID = "FNVPath";
252 (*store)[ID] = G4AttDef(ID, "Final Next Volume Path", "Physics", "", "G4String");
253
254 ID = "EPN";
255 (*store)[ID] = G4AttDef(ID, "Ending Process Name", "Physics", "", "G4String");
256
257 ID = "EPTN";
258 (*store)[ID] = G4AttDef(ID, "Ending Process Type Name", "Physics", "", "G4String");
259
260 ID = "FKE";
261 (*store)[ID] = G4AttDef(ID, "Final kinetic energy", "Physics", "G4BestUnit", "G4double");
262 }
263
264 return store;
265}
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 G4RichTrajectory::GetCharge ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 100 of file G4RichTrajectory.hh.

100{ return PDGCharge; }

◆ GetInitialKineticEnergy()

G4double G4RichTrajectory::GetInitialKineticEnergy ( ) const
inline

Definition at line 102 of file G4RichTrajectory.hh.

102{ return initialKineticEnergy; }

◆ GetInitialMomentum()

G4ThreeVector G4RichTrajectory::GetInitialMomentum ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 103 of file G4RichTrajectory.hh.

103{ return initialMomentum; }

◆ GetParentID()

G4int G4RichTrajectory::GetParentID ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 98 of file G4RichTrajectory.hh.

98{ return fParentID; }

◆ GetParticleDefinition()

G4ParticleDefinition * G4RichTrajectory::GetParticleDefinition ( )

Definition at line 355 of file G4RichTrajectory.cc.

356{
357 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
358}
static G4ParticleTable * GetParticleTable()

◆ GetParticleName()

G4String G4RichTrajectory::GetParticleName ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 99 of file G4RichTrajectory.hh.

99{ return ParticleName; }

◆ GetPDGEncoding()

G4int G4RichTrajectory::GetPDGEncoding ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 101 of file G4RichTrajectory.hh.

101{ return PDGEncoding; }

◆ GetPoint()

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

Implements G4VTrajectory.

Definition at line 164 of file G4RichTrajectory.hh.

165{
166 return (*fpRichPointContainer)[i];
167}

◆ GetPointEntries()

G4int G4RichTrajectory::GetPointEntries ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 159 of file G4RichTrajectory.hh.

160{
161 return G4int(fpRichPointContainer->size());
162}
int G4int
Definition G4Types.hh:85

Referenced by CreateAttValues().

◆ GetTrackID()

G4int G4RichTrajectory::GetTrackID ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 97 of file G4RichTrajectory.hh.

97{ return fTrackID; }

◆ MergeTrajectory()

void G4RichTrajectory::MergeTrajectory ( G4VTrajectory * secondTrajectory)
overridevirtual

Implements G4VTrajectory.

Definition at line 163 of file G4RichTrajectory.cc.

164{
165 if (secondTrajectory == nullptr) return;
166
167 auto seco = (G4RichTrajectory*)secondTrajectory;
168 G4int ent = seco->GetPointEntries();
169 for (G4int i = 1; i < ent; ++i) {
170 // initial point of the second trajectory should not be merged
171 //
172 fpRichPointContainer->push_back((*(seco->fpRichPointContainer))[i]);
173 }
174 delete (*seco->fpRichPointContainer)[0];
175 seco->fpRichPointContainer->clear();
176}
G4RichTrajectory()=default

◆ operator delete()

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

Definition at line 154 of file G4RichTrajectory.hh.

155{
156 aRichTrajectoryAllocator()->FreeSingle((G4RichTrajectory*)aRichTrajectory);
157}
G4TRACKING_DLL G4Allocator< G4RichTrajectory > *& aRichTrajectoryAllocator()

◆ operator new()

void * G4RichTrajectory::operator new ( size_t )
inline

Definition at line 146 of file G4RichTrajectory.hh.

147{
148 if (aRichTrajectoryAllocator() == nullptr) {
149 aRichTrajectoryAllocator() = new G4Allocator<G4RichTrajectory>;
150 }
151 return (void*)aRichTrajectoryAllocator()->MallocSingle();
152}

◆ operator=()

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

◆ operator==()

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

Definition at line 87 of file G4RichTrajectory.hh.

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

◆ ShowTrajectory()

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

Reimplemented from G4VTrajectory.

Definition at line 178 of file G4RichTrajectory.cc.

179{
180 // Invoke the default implementation in G4VTrajectory...
181 //
183
184 // ... or override with your own code here.
185}
virtual void ShowTrajectory(std::ostream &os=G4cout) const

Friends And Related Symbol Documentation

◆ G4ClonedRichTrajectory

friend class G4ClonedRichTrajectory
friend

Definition at line 74 of file G4RichTrajectory.hh.

Referenced by CloneForMaster(), and G4ClonedRichTrajectory.


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