Geant4 10.7.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 ()
 
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 61 of file G4RichTrajectory.cc.

62{
63}

◆ G4RichTrajectory() [2/3]

G4RichTrajectory::G4RichTrajectory ( const G4Track aTrack)

Definition at line 65 of file G4RichTrajectory.cc.

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

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

◆ ~G4RichTrajectory()

G4RichTrajectory::~G4RichTrajectory ( )
virtual

Definition at line 114 of file G4RichTrajectory.cc.

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

Member Function Documentation

◆ AppendStep()

void G4RichTrajectory::AppendStep ( const G4Step aStep)
virtual

Reimplemented from G4Trajectory.

Definition at line 127 of file G4RichTrajectory.cc.

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

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

◆ DrawTrajectory()

void G4RichTrajectory::DrawTrajectory ( ) const
virtual

Reimplemented from G4Trajectory.

Definition at line 171 of file G4RichTrajectory.cc.

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

◆ GetAttDefs()

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

Reimplemented from G4Trajectory.

Definition at line 180 of file G4RichTrajectory.cc.

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

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

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

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