Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VTrajectory Class Referenceabstract

#include <G4VTrajectory.hh>

+ Inheritance diagram for G4VTrajectory:

Public Member Functions

 G4VTrajectory ()=default
 
virtual ~G4VTrajectory ()=default
 
G4bool operator== (const G4VTrajectory &right) const
 
virtual G4VTrajectoryCloneForMaster () 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 G4VTrajectory.hh.

Constructor & Destructor Documentation

◆ G4VTrajectory()

◆ ~G4VTrajectory()

virtual G4VTrajectory::~G4VTrajectory ( )
virtualdefault

Member Function Documentation

◆ AppendStep()

virtual void G4VTrajectory::AppendStep ( const G4Step * aStep)
pure virtual

◆ CloneForMaster()

G4VTrajectory * G4VTrajectory::CloneForMaster ( ) const
virtual

Reimplemented in G4RichTrajectory, G4SmoothTrajectory, and G4Trajectory.

Definition at line 50 of file G4VTrajectory.cc.

51{
52 // Base-class method must not be invoked.
53 // Each concrete class should implement its own method.
54 G4Exception("G4VTrajectory::CloneForMaster()","traj0100",FatalException,
55 "CloneForMaster() of G4VTrajectory base-class is invoked. Each concrete class should implement its own method.");
56 return nullptr;
57}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)

◆ CreateAttValues()

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

◆ DrawTrajectory()

void G4VTrajectory::DrawTrajectory ( ) const
virtual

◆ GetAttDefs()

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

◆ GetCharge()

◆ GetInitialMomentum()

◆ GetParentID()

virtual G4int G4VTrajectory::GetParentID ( ) const
pure virtual

◆ GetParticleName()

◆ GetPDGEncoding()

virtual G4int G4VTrajectory::GetPDGEncoding ( ) const
pure virtual

◆ GetPoint()

◆ GetPointEntries()

◆ GetTrackID()

◆ MergeTrajectory()

virtual void G4VTrajectory::MergeTrajectory ( G4VTrajectory * secondTrajectory)
pure virtual

◆ operator==()

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

Definition at line 48 of file G4VTrajectory.cc.

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

◆ ShowTrajectory()

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

Reimplemented in G4ClonedRichTrajectory, G4ClonedSmoothTrajectory, G4ClonedTrajectory, G4RayTrajectory, G4RichTrajectory, G4SmoothTrajectory, and G4Trajectory.

Definition at line 59 of file G4VTrajectory.cc.

60{
61 // Makes use of attribute values implemented in the concrete class.
62 // Note: the user needs to follow with new-line or end-of-string,
63 // depending on the nature of os
64
65 std::vector<G4AttValue>* attValues = CreateAttValues();
66 const std::map<G4String, G4AttDef>* attDefs = GetAttDefs();
67
68 // Ensure validity...
69 //
70 if (G4AttCheck(attValues, attDefs).Check("G4VTrajectory::ShowTrajectory")) {
71 return;
72 }
73
74 os << "Trajectory:";
75
76 for (const auto& attValue : *attValues) {
77 auto iAttDef = attDefs->find(attValue.GetName());
78 os << "\n " << iAttDef->second.GetDesc() << " (" << attValue.GetName()
79 << "): " << attValue.GetValue();
80 }
81
82 delete attValues; // AttValues must be deleted after use.
83
84 // Now do trajectory points...
85
86 for (G4int i = 0; i < GetPointEntries(); ++i) {
87 G4VTrajectoryPoint* aTrajectoryPoint = GetPoint(i);
88 attValues = aTrajectoryPoint->CreateAttValues();
89 attDefs = aTrajectoryPoint->GetAttDefs();
90
91 // Ensure validity...
92 //
93 if (G4AttCheck(attValues, attDefs).Check("G4VTrajectory::ShowTrajectory")) {
94 return;
95 }
96
97 for (const auto& attValue : *attValues) {
98 const auto iAttDef = attDefs->find(attValue.GetName());
99 os << "\n " << iAttDef->second.GetDesc() << " (" << attValue.GetName()
100 << "): " << attValue.GetValue();
101 }
102
103 delete attValues; // AttValues must be deleted after use.
104 }
105 os << std::endl;
106}
int G4int
Definition G4Types.hh:85
virtual std::vector< G4AttValue > * CreateAttValues() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual G4VTrajectoryPoint * GetPoint(G4int i) const =0
virtual G4int GetPointEntries() const =0
virtual std::vector< G4AttValue > * CreateAttValues() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const

Referenced by G4ClonedRichTrajectory::ShowTrajectory(), G4ClonedSmoothTrajectory::ShowTrajectory(), G4ClonedTrajectory::ShowTrajectory(), G4RichTrajectory::ShowTrajectory(), G4SmoothTrajectory::ShowTrajectory(), and G4Trajectory::ShowTrajectory().


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