Geant4 9.6.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 ()
 
 G4RichTrajectory (const G4Track *aTrack)
 
 G4RichTrajectory (G4RichTrajectory &)
 
virtual ~G4RichTrajectory ()
 
void * operator new (size_t)
 
void operator delete (void *)
 
int operator== (const G4RichTrajectory &right) const
 
void AppendStep (const G4Step *aStep)
 
void MergeTrajectory (G4VTrajectory *secondTrajectory)
 
int GetPointEntries () const
 
G4VTrajectoryPointGetPoint (G4int i) const
 
virtual const std::map< G4String, G4AttDef > * GetAttDefs () const
 
virtual std::vector< G4AttValue > * CreateAttValues () const
 
- Public Member Functions inherited from G4Trajectory
 G4Trajectory ()
 
 G4Trajectory (const G4Track *aTrack)
 
 G4Trajectory (G4Trajectory &)
 
virtual ~G4Trajectory ()
 
void * operator new (size_t)
 
void operator delete (void *)
 
int operator== (const G4Trajectory &right) const
 
G4int GetTrackID () const
 
G4int GetParentID () const
 
G4String GetParticleName () const
 
G4double GetCharge () const
 
G4int GetPDGEncoding () const
 
G4double GetInitialKineticEnergy () const
 
G4ThreeVector GetInitialMomentum () const
 
virtual void ShowTrajectory (std::ostream &os=G4cout) const
 
virtual void DrawTrajectory (G4int i_mode=0) const
 
virtual void AppendStep (const G4Step *aStep)
 
virtual int GetPointEntries () const
 
virtual G4VTrajectoryPointGetPoint (G4int i) const
 
virtual void MergeTrajectory (G4VTrajectory *secondTrajectory)
 
G4ParticleDefinitionGetParticleDefinition ()
 
virtual const std::map< G4String, G4AttDef > * GetAttDefs () const
 
virtual std::vector< G4AttValue > * CreateAttValues () const
 
- Public Member Functions inherited from G4VTrajectory
 G4VTrajectory ()
 
virtual ~G4VTrajectory ()
 
G4bool operator== (const G4VTrajectory &right) const
 
virtual G4int GetTrackID () const =0
 
virtual G4int GetParentID () const =0
 
virtual G4String GetParticleName () const =0
 
virtual G4double GetCharge () const =0
 
virtual G4int GetPDGEncoding () const =0
 
virtual G4ThreeVector GetInitialMomentum () const =0
 
virtual int GetPointEntries () const =0
 
virtual G4VTrajectoryPointGetPoint (G4int i) const =0
 
virtual void ShowTrajectory (std::ostream &os=G4cout) const
 
virtual void DrawTrajectory (G4int i_mode=0) const
 
virtual const std::map< G4String, G4AttDef > * GetAttDefs () const
 
virtual std::vector< G4AttValue > * CreateAttValues () const
 
virtual void AppendStep (const G4Step *aStep)=0
 
virtual void MergeTrajectory (G4VTrajectory *secondTrajectory)=0
 

Detailed Description

Definition at line 65 of file G4RichTrajectory.hh.

Constructor & Destructor Documentation

◆ G4RichTrajectory() [1/3]

G4RichTrajectory::G4RichTrajectory ( )

Definition at line 62 of file G4RichTrajectory.cc.

62 :
63 fpRichPointsContainer(0),
64 fpCreatorProcess(0),
65 fpEndingProcess(0),
66 fFinalKineticEnergy(0.)
67{}

◆ G4RichTrajectory() [2/3]

G4RichTrajectory::G4RichTrajectory ( const G4Track aTrack)

Definition at line 69 of file G4RichTrajectory.cc.

69 :
70 G4Trajectory(aTrack) // Note: this initialises the base class data
71 // members and, unfortunately but never mind,
72 // creates a G4TrajectoryPoint in
73 // TrajectoryPointContainer that we cannot
74 // access because it's private. We store the
75 // same information (plus more) in a
76 // G4RichTrajectoryPoint in the
77 // RichTrajectoryPointsContainer
78{
79 fpInitialVolume = aTrack->GetTouchableHandle();
80 fpInitialNextVolume = aTrack->GetNextTouchableHandle();
81 fpCreatorProcess = aTrack->GetCreatorProcess();
82 // On construction, set final values to initial values.
83 // Final values are updated at the addition of every step - see AppendStep.
84 fpFinalVolume = aTrack->GetTouchableHandle();
85 fpFinalNextVolume = aTrack->GetNextTouchableHandle();
86 fpEndingProcess = aTrack->GetCreatorProcess();
87 fFinalKineticEnergy = aTrack->GetKineticEnergy();
88 // Insert the first rich trajectory point (see note above)...
89 fpRichPointsContainer = new RichTrajectoryPointsContainer;
90 fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(aTrack));
91}
std::vector< G4VTrajectoryPoint * > RichTrajectoryPointsContainer
const G4TouchableHandle & GetNextTouchableHandle() const
const G4VProcess * GetCreatorProcess() const
const G4TouchableHandle & GetTouchableHandle() const
G4double GetKineticEnergy() const

◆ G4RichTrajectory() [3/3]

G4RichTrajectory::G4RichTrajectory ( G4RichTrajectory right)

Definition at line 93 of file G4RichTrajectory.cc.

93 :
94 G4Trajectory(right)
95{
96 fpInitialVolume = right.fpInitialVolume;
97 fpInitialNextVolume = right.fpInitialNextVolume;
98 fpCreatorProcess = right.fpCreatorProcess;
99 fpFinalVolume = right.fpFinalVolume;
100 fpFinalNextVolume = right.fpFinalNextVolume;
101 fpEndingProcess = right.fpEndingProcess;
102 fFinalKineticEnergy = right.fFinalKineticEnergy;
103 fpRichPointsContainer = new RichTrajectoryPointsContainer;
104 for(size_t i=0;i<right.fpRichPointsContainer->size();i++)
105 {
106 G4RichTrajectoryPoint* rightPoint =
107 (G4RichTrajectoryPoint*)((*(right.fpRichPointsContainer))[i]);
108 fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(*rightPoint));
109 }
110}

◆ ~G4RichTrajectory()

G4RichTrajectory::~G4RichTrajectory ( )
virtual

Definition at line 112 of file G4RichTrajectory.cc.

113{
114 if (fpRichPointsContainer) {
115 // fpRichPointsContainer->clearAndDestroy();
116 size_t i;
117 for(i=0;i<fpRichPointsContainer->size();i++){
118 delete (*fpRichPointsContainer)[i];
119 }
120 fpRichPointsContainer->clear();
121 delete fpRichPointsContainer;
122 }
123}

Member Function Documentation

◆ AppendStep()

void G4RichTrajectory::AppendStep ( const G4Step aStep)
virtual

Reimplemented from G4Trajectory.

Definition at line 125 of file G4RichTrajectory.cc.

126{
127 fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(aStep));
128 // Except for first step, which is a sort of virtual step to start
129 // the track, compute the final values...
130 const G4Track* track = aStep->GetTrack();
131 const G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
132 if (track->GetCurrentStepNumber() > 0) {
133 fpFinalVolume = track->GetTouchableHandle();
134 fpFinalNextVolume = track->GetNextTouchableHandle();
135 fpEndingProcess = postStepPoint->GetProcessDefinedStep();
136 fFinalKineticEnergy =
138 aStep->GetTotalEnergyDeposit();
139 }
140}
const G4VProcess * GetProcessDefinedStep() const
G4double GetKineticEnergy() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetTotalEnergyDeposit() const
G4StepPoint * GetPostStepPoint() const
G4int GetCurrentStepNumber() const

◆ CreateAttValues()

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

Reimplemented from G4Trajectory.

Definition at line 222 of file G4RichTrajectory.cc.

223{
224 // Create base class att values...
225 std::vector<G4AttValue>* values = G4Trajectory::CreateAttValues();
226
227 if (fpInitialVolume && fpInitialVolume->GetVolume()) {
228 values->push_back(G4AttValue("IVPath",Path(fpInitialVolume),""));
229 } else {
230 values->push_back(G4AttValue("IVPath","None",""));
231 }
232
233 if (fpInitialNextVolume && fpInitialNextVolume->GetVolume()) {
234 values->push_back(G4AttValue("INVPath",Path(fpInitialNextVolume),""));
235 } else {
236 values->push_back(G4AttValue("INVPath","None",""));
237 }
238
239 if (fpCreatorProcess) {
240 values->push_back(G4AttValue("CPN",fpCreatorProcess->GetProcessName(),""));
241 G4ProcessType type = fpCreatorProcess->GetProcessType();
242 values->push_back(G4AttValue("CPTN",G4VProcess::GetProcessTypeName(type),""));
243 } else {
244 values->push_back(G4AttValue("CPN","None",""));
245 values->push_back(G4AttValue("CPTN","None",""));
246 }
247
248 if (fpFinalVolume && fpFinalVolume->GetVolume()) {
249 values->push_back(G4AttValue("FVPath",Path(fpFinalVolume),""));
250 } else {
251 values->push_back(G4AttValue("FVPath","None",""));
252 }
253
254 if (fpFinalNextVolume && fpFinalNextVolume->GetVolume()) {
255 values->push_back(G4AttValue("FNVPath",Path(fpFinalNextVolume),""));
256 } else {
257 values->push_back(G4AttValue("FNVPath","None",""));
258 }
259
260 if (fpEndingProcess) {
261 values->push_back(G4AttValue("EPN",fpEndingProcess->GetProcessName(),""));
262 G4ProcessType type = fpEndingProcess->GetProcessType();
263 values->push_back(G4AttValue("EPTN",G4VProcess::GetProcessTypeName(type),""));
264 } else {
265 values->push_back(G4AttValue("EPN","None",""));
266 values->push_back(G4AttValue("EPTN","None",""));
267 }
268
269 values->push_back
270 (G4AttValue("FKE",G4BestUnit(fFinalKineticEnergy,"Energy"),""));
271
272#ifdef G4ATTDEBUG
273 G4cout << G4AttCheck(values,GetAttDefs());
274#endif
275
276 return values;
277}
G4ProcessType
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4DLLIMPORT std::ostream G4cout
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual std::vector< G4AttValue > * CreateAttValues() const
static const G4String & GetProcessTypeName(G4ProcessType)
Definition: G4VProcess.cc:150
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:385
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44

◆ GetAttDefs()

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

Reimplemented from G4Trajectory.

Definition at line 157 of file G4RichTrajectory.cc.

158{
159 G4bool isNew;
160 std::map<G4String,G4AttDef>* store
161 = G4AttDefStore::GetInstance("G4RichTrajectory",isNew);
162 if (isNew) {
163
164 // Get att defs from base class...
165 *store = *(G4Trajectory::GetAttDefs());
166
167 G4String ID;
168
169 ID = "IVPath";
170 (*store)[ID] = G4AttDef(ID,"Initial Volume Path",
171 "Physics","","G4String");
172
173 ID = "INVPath";
174 (*store)[ID] = G4AttDef(ID,"Initial Next Volume Path",
175 "Physics","","G4String");
176
177 ID = "CPN";
178 (*store)[ID] = G4AttDef(ID,"Creator Process Name",
179 "Physics","","G4String");
180
181 ID = "CPTN";
182 (*store)[ID] = G4AttDef(ID,"Creator Process Type Name",
183 "Physics","","G4String");
184
185 ID = "FVPath";
186 (*store)[ID] = G4AttDef(ID,"Final Volume Path",
187 "Physics","","G4String");
188
189 ID = "FNVPath";
190 (*store)[ID] = G4AttDef(ID,"Final Next Volume Path",
191 "Physics","","G4String");
192
193 ID = "EPN";
194 (*store)[ID] = G4AttDef(ID,"Ending Process Name",
195 "Physics","","G4String");
196
197 ID = "EPTN";
198 (*store)[ID] = G4AttDef(ID,"Ending Process Type Name",
199 "Physics","","G4String");
200
201 ID = "FKE";
202 (*store)[ID] = G4AttDef(ID,"Final kinetic energy",
203 "Physics","G4BestUnit","G4double");
204
205 }
206
207 return store;
208}
bool G4bool
Definition: G4Types.hh:67
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
std::map< G4String, G4AttDef > * GetInstance(G4String storeKey, G4bool &isNew)

Referenced by CreateAttValues(), G4VisCommandList::SetNewValue(), and G4VisCommandSceneAddTrajectories::SetNewValue().

◆ GetPoint()

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

Reimplemented from G4Trajectory.

Definition at line 89 of file G4RichTrajectory.hh.

90 { return (*fpRichPointsContainer)[i]; }

◆ GetPointEntries()

int G4RichTrajectory::GetPointEntries ( ) const
inlinevirtual

Reimplemented from G4Trajectory.

Definition at line 88 of file G4RichTrajectory.hh.

88{ return fpRichPointsContainer->size(); }

Referenced by MergeTrajectory().

◆ MergeTrajectory()

void G4RichTrajectory::MergeTrajectory ( G4VTrajectory secondTrajectory)
virtual

Reimplemented from G4Trajectory.

Definition at line 142 of file G4RichTrajectory.cc.

143{
144 if(!secondTrajectory) return;
145
146 G4RichTrajectory* seco = (G4RichTrajectory*)secondTrajectory;
147 G4int ent = seco->GetPointEntries();
148 for(G4int i=1;i<ent;i++) {
149 // initial point of the second trajectory should not be merged
150 fpRichPointsContainer->push_back((*(seco->fpRichPointsContainer))[i]);
151 // fpRichPointsContainer->push_back(seco->fpRichPointsContainer->removeAt(1));
152 }
153 delete (*seco->fpRichPointsContainer)[0];
154 seco->fpRichPointsContainer->clear();
155}
int G4int
Definition: G4Types.hh:66
int GetPointEntries() const

◆ operator delete()

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

Definition at line 125 of file G4RichTrajectory.hh.

126{
127 aRichTrajectoryAllocator.FreeSingle
128 ((G4RichTrajectory*)aRichTrajectory);
129}
G4DLLIMPORT G4Allocator< G4RichTrajectory > aRichTrajectoryAllocator

◆ operator new()

void * G4RichTrajectory::operator new ( size_t  )
inline

Definition at line 118 of file G4RichTrajectory.hh.

119{
120 void* aRichTrajectory;
121 aRichTrajectory = (void*)aRichTrajectoryAllocator.MallocSingle();
122 return aRichTrajectory;
123}

◆ operator==()

int G4RichTrajectory::operator== ( const G4RichTrajectory right) const
inline

Definition at line 82 of file G4RichTrajectory.hh.

83 {return (this==&right);}

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