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

Friends

class G4ClonedRichTrajectory

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 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 135 of file G4RichTrajectory.cc.

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

◆ G4RichTrajectory() [3/3]

G4RichTrajectory::G4RichTrajectory ( G4RichTrajectory & right)

Definition at line 103 of file G4RichTrajectory.cc.

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

Member Function Documentation

◆ AppendStep()

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

Implements G4VTrajectory.

Definition at line 146 of file G4RichTrajectory.cc.

147{
148 fpRichPointContainer->push_back(new G4RichTrajectoryPoint(aStep));
149
150 // Except for first step, which is a sort of virtual step to start
151 // the track, compute the final values...
152 //
153 const G4Track* track = aStep->GetTrack();
154 const G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
155 if (track->GetCurrentStepNumber() > 0) {
156 fpFinalVolume = track->GetTouchableHandle();
157 fpFinalNextVolume = track->GetNextTouchableHandle();
158 fpEndingProcess = postStepPoint->GetProcessDefinedStep();
159 fFinalKineticEnergy =
161 }
162}
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 361 of file G4RichTrajectory.cc.

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

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Definition at line 279 of file G4RichTrajectory.cc.

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

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

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Definition at line 197 of file G4RichTrajectory.cc.

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

357{
358 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
359}
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 164 of file G4RichTrajectory.cc.

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

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

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