Geant4 11.1.1
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 ()
 
G4RichTrajectoryoperator= (const G4RichTrajectory &)=delete
 
G4int operator== (const G4RichTrajectory &r) const
 
void * operator new (size_t)
 
void operator delete (void *)
 
void ShowTrajectory (std::ostream &os=G4cout) const
 
void DrawTrajectory () const
 
void AppendStep (const G4Step *aStep)
 
void MergeTrajectory (G4VTrajectory *secondTrajectory)
 
G4int 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 *)
 
G4int operator== (const G4Trajectory &r) 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 () const
 
virtual void AppendStep (const G4Step *aStep)
 
virtual G4int 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 G4int GetPointEntries () const =0
 
virtual G4VTrajectoryPointGetPoint (G4int i) const =0
 
virtual void ShowTrajectory (std::ostream &os=G4cout) const
 
virtual void DrawTrajectory () 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 57 of file G4RichTrajectory.hh.

Constructor & Destructor Documentation

◆ G4RichTrajectory() [1/3]

G4RichTrajectory::G4RichTrajectory ( )

Definition at line 62 of file G4RichTrajectory.cc.

63{
64}

◆ G4RichTrajectory() [2/3]

G4RichTrajectory::G4RichTrajectory ( const G4Track aTrack)

Definition at line 66 of file G4RichTrajectory.cc.

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

◆ G4RichTrajectory() [3/3]

G4RichTrajectory::G4RichTrajectory ( G4RichTrajectory right)

Definition at line 95 of file G4RichTrajectory.cc.

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

◆ ~G4RichTrajectory()

G4RichTrajectory::~G4RichTrajectory ( )
virtual

Definition at line 115 of file G4RichTrajectory.cc.

116{
117 if (fpRichPointsContainer)
118 {
119 for(std::size_t i=0; i<fpRichPointsContainer->size(); ++i)
120 {
121 delete (*fpRichPointsContainer)[i];
122 }
123 fpRichPointsContainer->clear();
124 delete fpRichPointsContainer;
125 }
126}

Member Function Documentation

◆ AppendStep()

void G4RichTrajectory::AppendStep ( const G4Step aStep)
virtual

Reimplemented from G4Trajectory.

Definition at line 128 of file G4RichTrajectory.cc.

129{
130 fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(aStep));
131
132 // Except for first step, which is a sort of virtual step to start
133 // the track, compute the final values...
134 //
135 const G4Track* track = aStep->GetTrack();
136 const G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
137 if (track->GetCurrentStepNumber() > 0)
138 {
139 fpFinalVolume = track->GetTouchableHandle();
140 fpFinalNextVolume = track->GetNextTouchableHandle();
141 fpEndingProcess = postStepPoint->GetProcessDefinedStep();
142 fFinalKineticEnergy = aStep->GetPreStepPoint()->GetKineticEnergy()
143 - aStep->GetTotalEnergyDeposit();
144 }
145}
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 255 of file G4RichTrajectory.cc.

256{
257 // Create base class att values...
258 std::vector<G4AttValue>* values = G4Trajectory::CreateAttValues();
259
260 if (fpInitialVolume && fpInitialVolume->GetVolume())
261 {
262 values->push_back(G4AttValue("IVPath",Path(fpInitialVolume),""));
263 }
264 else
265 {
266 values->push_back(G4AttValue("IVPath","None",""));
267 }
268
269 if (fpInitialNextVolume && fpInitialNextVolume->GetVolume())
270 {
271 values->push_back(G4AttValue("INVPath",Path(fpInitialNextVolume),""));
272 }
273 else
274 {
275 values->push_back(G4AttValue("INVPath","None",""));
276 }
277
278 if (fpCreatorProcess != nullptr)
279 {
280 values->push_back
281 (G4AttValue("CPN",fpCreatorProcess->GetProcessName(),""));
282 G4ProcessType type = fpCreatorProcess->GetProcessType();
283 values->push_back
285 values->push_back
286 (G4AttValue("CMID",G4UIcommand::ConvertToString(fCreatorModelID),""));
287 const G4String& creatorModelName =
289 values->push_back(G4AttValue("CMN",creatorModelName,""));
290 }
291 else
292 {
293 values->push_back(G4AttValue("CPN","None",""));
294 values->push_back(G4AttValue("CPTN","None",""));
295 values->push_back(G4AttValue("CMID","None",""));
296 values->push_back(G4AttValue("CMN","None",""));
297 }
298
299 if (fpFinalVolume && fpFinalVolume->GetVolume())
300 {
301 values->push_back(G4AttValue("FVPath",Path(fpFinalVolume),""));
302 }
303 else
304 {
305 values->push_back(G4AttValue("FVPath","None",""));
306 }
307
308 if (fpFinalNextVolume && fpFinalNextVolume->GetVolume())
309 {
310 values->push_back(G4AttValue("FNVPath",Path(fpFinalNextVolume),""));
311 }
312 else
313 {
314 values->push_back(G4AttValue("FNVPath","None",""));
315 }
316
317 if (fpEndingProcess != nullptr)
318 {
319 values->push_back(G4AttValue("EPN",fpEndingProcess->GetProcessName(),""));
320 G4ProcessType type = fpEndingProcess->GetProcessType();
321 values->push_back(G4AttValue("EPTN",G4VProcess::GetProcessTypeName(type),""));
322 }
323 else
324 {
325 values->push_back(G4AttValue("EPN","None",""));
326 values->push_back(G4AttValue("EPTN","None",""));
327 }
328
329 values->push_back
330 (G4AttValue("FKE",G4BestUnit(fFinalKineticEnergy,"Energy"),""));
331
332#ifdef G4ATTDEBUG
333 G4cout << G4AttCheck(values,GetAttDefs());
334#endif
335
336 return values;
337}
G4ProcessType
#define G4BestUnit(a, b)
G4GLOB_DLL std::ostream G4cout
static const G4String GetModelNameFromID(const G4int modelID)
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual std::vector< G4AttValue > * CreateAttValues() const
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:446
static const G4String & GetProcessTypeName(G4ProcessType)
Definition: G4VProcess.cc:134
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:392
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:34

◆ DrawTrajectory()

void G4RichTrajectory::DrawTrajectory ( ) const
virtual

Reimplemented from G4Trajectory.

Definition at line 172 of file G4RichTrajectory.cc.

173{
174 // Invoke the default implementation in G4VTrajectory...
175 //
177
178 // ... or override with your own code here.
179}
virtual void DrawTrajectory() const

◆ GetAttDefs()

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

Reimplemented from G4Trajectory.

Definition at line 181 of file G4RichTrajectory.cc.

182{
183 G4bool isNew;
184 std::map<G4String,G4AttDef>* store
185 = G4AttDefStore::GetInstance("G4RichTrajectory",isNew);
186 if (isNew)
187 {
188 // Get att defs from base class...
189 //
190 *store = *(G4Trajectory::GetAttDefs());
191
192 G4String ID;
193
194 ID = "IVPath";
195 (*store)[ID] = G4AttDef(ID,"Initial Volume Path",
196 "Physics","","G4String");
197
198 ID = "INVPath";
199 (*store)[ID] = G4AttDef(ID,"Initial Next Volume Path",
200 "Physics","","G4String");
201
202 ID = "CPN";
203 (*store)[ID] = G4AttDef(ID,"Creator Process Name",
204 "Physics","","G4String");
205
206 ID = "CPTN";
207 (*store)[ID] = G4AttDef(ID,"Creator Process Type Name",
208 "Physics","","G4String");
209
210 ID = "CMID";
211 (*store)[ID] = G4AttDef(ID,"Creator Model ID",
212 "Physics","","G4int");
213
214 ID = "CMN";
215 (*store)[ID] = G4AttDef(ID,"Creator Model Name",
216 "Physics","","G4String");
217
218 ID = "FVPath";
219 (*store)[ID] = G4AttDef(ID,"Final Volume Path",
220 "Physics","","G4String");
221
222 ID = "FNVPath";
223 (*store)[ID] = G4AttDef(ID,"Final Next Volume Path",
224 "Physics","","G4String");
225
226 ID = "EPN";
227 (*store)[ID] = G4AttDef(ID,"Ending Process Name",
228 "Physics","","G4String");
229
230 ID = "EPTN";
231 (*store)[ID] = G4AttDef(ID,"Ending Process Type Name",
232 "Physics","","G4String");
233
234 ID = "FKE";
235 (*store)[ID] = G4AttDef(ID,"Final kinetic energy",
236 "Physics","G4BestUnit","G4double");
237 }
238
239 return store;
240}
bool G4bool
Definition: G4Types.hh:86
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
std::map< G4String, G4AttDef > * GetInstance(const 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 128 of file G4RichTrajectory.hh.

129{
130 return (*fpRichPointsContainer)[i];
131}

Referenced by G4TrajectoryDrawByEncounteredVolume::Draw(), and G4TrajectoryEncounteredVolumeFilter::Evaluate().

◆ GetPointEntries()

G4int G4RichTrajectory::GetPointEntries ( ) const
inlinevirtual

Reimplemented from G4Trajectory.

Definition at line 123 of file G4RichTrajectory.hh.

124{
125 return G4int(fpRichPointsContainer->size());
126}
int G4int
Definition: G4Types.hh:85

Referenced by G4TrajectoryDrawByEncounteredVolume::Draw(), G4TrajectoryEncounteredVolumeFilter::Evaluate(), and MergeTrajectory().

◆ MergeTrajectory()

void G4RichTrajectory::MergeTrajectory ( G4VTrajectory secondTrajectory)
virtual

Reimplemented from G4Trajectory.

Definition at line 147 of file G4RichTrajectory.cc.

148{
149 if(secondTrajectory == nullptr) return;
150
151 G4RichTrajectory* seco = (G4RichTrajectory*)secondTrajectory;
152 G4int ent = seco->GetPointEntries();
153 for(G4int i=1; i<ent; ++i)
154 {
155 // initial point of the second trajectory should not be merged
156 //
157 fpRichPointsContainer->push_back((*(seco->fpRichPointsContainer))[i]);
158 }
159 delete (*seco->fpRichPointsContainer)[0];
160 seco->fpRichPointsContainer->clear();
161}
G4int GetPointEntries() const

◆ operator delete()

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

Definition at line 118 of file G4RichTrajectory.hh.

119{
120 aRichTrajectoryAllocator()->FreeSingle((G4RichTrajectory*)aRichTrajectory);
121}
G4TRACKING_DLL G4Allocator< G4RichTrajectory > *& aRichTrajectoryAllocator()

◆ operator new()

void * G4RichTrajectory::operator new ( size_t  )
inline

Definition at line 109 of file G4RichTrajectory.hh.

110{
111 if (aRichTrajectoryAllocator() == nullptr)
112 {
114 }
115 return (void*)aRichTrajectoryAllocator()->MallocSingle();
116}

◆ operator=()

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

◆ operator==()

G4int G4RichTrajectory::operator== ( const G4RichTrajectory r) const
inline

Definition at line 73 of file G4RichTrajectory.hh.

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

◆ ShowTrajectory()

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

Reimplemented from G4Trajectory.

Definition at line 163 of file G4RichTrajectory.cc.

164{
165 // Invoke the default implementation in G4VTrajectory...
166 //
168
169 // ... or override with your own code here.
170}
virtual void ShowTrajectory(std::ostream &os=G4cout) const

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