Geant4-11
Public Member Functions | Private Types | Private Attributes
G4RichTrajectory Class Reference

#include <G4RichTrajectory.hh>

Inheritance diagram for G4RichTrajectory:
G4Trajectory G4VTrajectory

Public Member Functions

void AppendStep (const G4Step *aStep)
 
virtual std::vector< G4AttValue > * CreateAttValues () const
 
void DrawTrajectory () const
 
 G4RichTrajectory ()
 
 G4RichTrajectory (const G4Track *aTrack)
 
 G4RichTrajectory (G4RichTrajectory &)
 
virtual const std::map< G4String, G4AttDef > * GetAttDefs () const
 
G4double GetCharge () const
 
G4double GetInitialKineticEnergy () const
 
G4ThreeVector GetInitialMomentum () const
 
G4int GetParentID () const
 
G4ParticleDefinitionGetParticleDefinition ()
 
G4String GetParticleName () const
 
G4int GetPDGEncoding () const
 
G4VTrajectoryPointGetPoint (G4int i) const
 
G4int GetPointEntries () const
 
G4int GetTrackID () const
 
void MergeTrajectory (G4VTrajectory *secondTrajectory)
 
void operator delete (void *)
 
void * operator new (size_t)
 
G4RichTrajectoryoperator= (const G4RichTrajectory &)=delete
 
G4int operator== (const G4RichTrajectory &r) const
 
G4int operator== (const G4Trajectory &r) const
 
G4bool operator== (const G4VTrajectory &right) const
 
void ShowTrajectory (std::ostream &os=G4cout) const
 
virtual ~G4RichTrajectory ()
 

Private Types

using G4TrajectoryPointContainer = std::vector< G4VTrajectoryPoint * >
 
using RichTrajectoryPointsContainer = std::vector< G4VTrajectoryPoint * >
 

Private Attributes

G4int fCreatorModelID = 0
 
G4double fFinalKineticEnergy = 0.0
 
G4int fParentID = 0
 
const G4VProcessfpCreatorProcess = nullptr
 
const G4VProcessfpEndingProcess = nullptr
 
G4TouchableHandle fpFinalNextVolume
 
G4TouchableHandle fpFinalVolume
 
G4TouchableHandle fpInitialNextVolume
 
G4TouchableHandle fpInitialVolume
 
RichTrajectoryPointsContainerfpRichPointsContainer = nullptr
 
G4int fTrackID = 0
 
G4double initialKineticEnergy = 0.0
 
G4ThreeVector initialMomentum
 
G4String ParticleName = ""
 
G4double PDGCharge = 0.0
 
G4int PDGEncoding = 0
 
G4TrajectoryPointContainerpositionRecord = nullptr
 

Detailed Description

Definition at line 57 of file G4RichTrajectory.hh.

Member Typedef Documentation

◆ G4TrajectoryPointContainer

using G4Trajectory::G4TrajectoryPointContainer = std::vector<G4VTrajectoryPoint*>
privateinherited

Definition at line 64 of file G4Trajectory.hh.

◆ RichTrajectoryPointsContainer

Definition at line 59 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{
80
81 // On construction, set final values to initial values.
82 // Final values are updated at the addition of every step - see AppendStep.
83 //
88
89 // Insert the first rich trajectory point (see note above)...
90 //
92 fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(aTrack));
93}
G4TouchableHandle fpFinalVolume
G4TouchableHandle fpFinalNextVolume
G4TouchableHandle fpInitialNextVolume
G4double fFinalKineticEnergy
std::vector< G4VTrajectoryPoint * > RichTrajectoryPointsContainer
RichTrajectoryPointsContainer * fpRichPointsContainer
const G4VProcess * fpCreatorProcess
const G4VProcess * fpEndingProcess
G4TouchableHandle fpInitialVolume
const G4TouchableHandle & GetNextTouchableHandle() const
const G4VProcess * GetCreatorProcess() const
G4int GetCreatorModelID() const
const G4TouchableHandle & GetTouchableHandle() const
G4double GetKineticEnergy() const

References fCreatorModelID, fFinalKineticEnergy, fpCreatorProcess, fpEndingProcess, fpFinalNextVolume, fpFinalVolume, fpInitialNextVolume, fpInitialVolume, fpRichPointsContainer, G4Track::GetCreatorModelID(), G4Track::GetCreatorProcess(), G4Track::GetKineticEnergy(), G4Track::GetNextTouchableHandle(), and G4Track::GetTouchableHandle().

◆ G4RichTrajectory() [3/3]

G4RichTrajectory::G4RichTrajectory ( G4RichTrajectory right)

◆ ~G4RichTrajectory()

G4RichTrajectory::~G4RichTrajectory ( )
virtual

Definition at line 115 of file G4RichTrajectory.cc.

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

References fpRichPointsContainer.

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 {
141 fpEndingProcess = postStepPoint->GetProcessDefinedStep();
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

References fFinalKineticEnergy, fpEndingProcess, fpFinalNextVolume, fpFinalVolume, fpRichPointsContainer, G4Track::GetCurrentStepNumber(), G4StepPoint::GetKineticEnergy(), G4Track::GetNextTouchableHandle(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4StepPoint::GetProcessDefinedStep(), G4Step::GetTotalEnergyDeposit(), G4Track::GetTouchableHandle(), and G4Step::GetTrack().

◆ 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
261 {
262 values->push_back(G4AttValue("IVPath",Path(fpInitialVolume),""));
263 }
264 else
265 {
266 values->push_back(G4AttValue("IVPath","None",""));
267 }
268
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
283 values->push_back
285 values->push_back
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
300 {
301 values->push_back(G4AttValue("FVPath",Path(fpFinalVolume),""));
302 }
303 else
304 {
305 values->push_back(G4AttValue("FVPath","None",""));
306 }
307
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(),""));
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
static G4String Path(const G4TouchableHandle &th)
#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:445
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:34

References G4UIcommand::ConvertToString(), G4Trajectory::CreateAttValues(), fCreatorModelID, fFinalKineticEnergy, fpCreatorProcess, fpEndingProcess, fpFinalNextVolume, fpFinalVolume, fpInitialNextVolume, fpInitialVolume, G4BestUnit, G4cout, GetAttDefs(), G4PhysicsModelCatalog::GetModelNameFromID(), G4VProcess::GetProcessName(), G4VProcess::GetProcessType(), G4VProcess::GetProcessTypeName(), G4VTouchable::GetVolume(), and Path().

◆ 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

References G4VTrajectory::DrawTrajectory().

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

References G4Trajectory::GetAttDefs(), and G4AttDefStore::GetInstance().

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

◆ GetCharge()

G4double G4Trajectory::GetCharge ( ) const
inlinevirtualinherited

Implements G4VTrajectory.

Definition at line 89 of file G4Trajectory.hh.

90 { return PDGCharge; }
G4double PDGCharge

References G4Trajectory::PDGCharge.

◆ GetInitialKineticEnergy()

G4double G4Trajectory::GetInitialKineticEnergy ( ) const
inlineinherited

Definition at line 93 of file G4Trajectory.hh.

94 { return initialKineticEnergy; }
G4double initialKineticEnergy

References G4Trajectory::initialKineticEnergy.

◆ GetInitialMomentum()

G4ThreeVector G4Trajectory::GetInitialMomentum ( ) const
inlinevirtualinherited

Implements G4VTrajectory.

Definition at line 95 of file G4Trajectory.hh.

96 { return initialMomentum; }
G4ThreeVector initialMomentum

References G4Trajectory::initialMomentum.

◆ GetParentID()

G4int G4Trajectory::GetParentID ( ) const
inlinevirtualinherited

Implements G4VTrajectory.

Definition at line 85 of file G4Trajectory.hh.

86 { return fParentID; }

References G4Trajectory::fParentID.

◆ GetParticleDefinition()

G4ParticleDefinition * G4Trajectory::GetParticleDefinition ( void  )
inherited

Definition at line 211 of file G4Trajectory.cc.

212{
213 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
214}
static G4ParticleTable * GetParticleTable()
G4String ParticleName

References G4ParticleTable::GetParticleTable(), and G4Trajectory::ParticleName.

◆ GetParticleName()

G4String G4Trajectory::GetParticleName ( ) const
inlinevirtualinherited

Implements G4VTrajectory.

Definition at line 87 of file G4Trajectory.hh.

88 { return ParticleName; }

References G4Trajectory::ParticleName.

◆ GetPDGEncoding()

G4int G4Trajectory::GetPDGEncoding ( ) const
inlinevirtualinherited

Implements G4VTrajectory.

Definition at line 91 of file G4Trajectory.hh.

92 { return PDGEncoding; }

References G4Trajectory::PDGEncoding.

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

References fpRichPointsContainer.

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

References fpRichPointsContainer.

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

◆ GetTrackID()

G4int G4Trajectory::GetTrackID ( ) const
inlinevirtualinherited

Implements G4VTrajectory.

Definition at line 83 of file G4Trajectory.hh.

84 { return fTrackID; }

References G4Trajectory::fTrackID.

◆ 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

References fpRichPointsContainer, and GetPointEntries().

◆ 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()

References 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}

References aRichTrajectoryAllocator().

◆ operator=()

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

◆ operator==() [1/3]

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

Definition at line 73 of file G4RichTrajectory.hh.

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

◆ operator==() [2/3]

G4int G4Trajectory::operator== ( const G4Trajectory r) const
inlineinherited

Definition at line 142 of file G4Trajectory.hh.

143{
144 return (this==&r);
145}

◆ operator==() [3/3]

G4bool G4VTrajectory::operator== ( const G4VTrajectory right) const
inherited

Definition at line 55 of file G4VTrajectory.cc.

56{
57 return (this==&right);
58}

◆ 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

References G4VTrajectory::ShowTrajectory().

Field Documentation

◆ fCreatorModelID

G4int G4RichTrajectory::fCreatorModelID = 0
private

Definition at line 100 of file G4RichTrajectory.hh.

Referenced by CreateAttValues(), and G4RichTrajectory().

◆ fFinalKineticEnergy

G4double G4RichTrajectory::fFinalKineticEnergy = 0.0
private

Definition at line 104 of file G4RichTrajectory.hh.

Referenced by AppendStep(), CreateAttValues(), and G4RichTrajectory().

◆ fParentID

G4int G4Trajectory::fParentID = 0
privateinherited

◆ fpCreatorProcess

const G4VProcess* G4RichTrajectory::fpCreatorProcess = nullptr
private

Definition at line 99 of file G4RichTrajectory.hh.

Referenced by CreateAttValues(), and G4RichTrajectory().

◆ fpEndingProcess

const G4VProcess* G4RichTrajectory::fpEndingProcess = nullptr
private

Definition at line 103 of file G4RichTrajectory.hh.

Referenced by AppendStep(), CreateAttValues(), and G4RichTrajectory().

◆ fpFinalNextVolume

G4TouchableHandle G4RichTrajectory::fpFinalNextVolume
private

Definition at line 102 of file G4RichTrajectory.hh.

Referenced by AppendStep(), CreateAttValues(), and G4RichTrajectory().

◆ fpFinalVolume

G4TouchableHandle G4RichTrajectory::fpFinalVolume
private

Definition at line 101 of file G4RichTrajectory.hh.

Referenced by AppendStep(), CreateAttValues(), and G4RichTrajectory().

◆ fpInitialNextVolume

G4TouchableHandle G4RichTrajectory::fpInitialNextVolume
private

Definition at line 98 of file G4RichTrajectory.hh.

Referenced by CreateAttValues(), and G4RichTrajectory().

◆ fpInitialVolume

G4TouchableHandle G4RichTrajectory::fpInitialVolume
private

Definition at line 97 of file G4RichTrajectory.hh.

Referenced by CreateAttValues(), and G4RichTrajectory().

◆ fpRichPointsContainer

RichTrajectoryPointsContainer* G4RichTrajectory::fpRichPointsContainer = nullptr
private

◆ fTrackID

G4int G4Trajectory::fTrackID = 0
privateinherited

◆ initialKineticEnergy

G4double G4Trajectory::initialKineticEnergy = 0.0
privateinherited

◆ initialMomentum

G4ThreeVector G4Trajectory::initialMomentum
privateinherited

◆ ParticleName

G4String G4Trajectory::ParticleName = ""
privateinherited

◆ PDGCharge

G4double G4Trajectory::PDGCharge = 0.0
privateinherited

◆ PDGEncoding

G4int G4Trajectory::PDGEncoding = 0
privateinherited

◆ positionRecord

G4TrajectoryPointContainer* G4Trajectory::positionRecord = nullptr
privateinherited

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