Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RayTrajectory Class Reference

#include <G4RayTrajectory.hh>

+ Inheritance diagram for G4RayTrajectory:

Public Member Functions

 G4RayTrajectory ()
 
 G4RayTrajectory (G4RayTrajectory &right)
 
virtual ~G4RayTrajectory ()
 
void * operator new (std::size_t)
 
void operator delete (void *)
 
virtual void AppendStep (const G4Step *)
 
virtual void ShowTrajectory (std::ostream &) const
 
virtual void DrawTrajectory () const
 
virtual G4int GetPointEntries () const
 
virtual G4VTrajectoryPointGetPoint (G4int i) const
 
G4RayTrajectoryPointGetPointC (G4int i) const
 
virtual void MergeTrajectory (G4VTrajectory *secondTrajectory)
 
G4int GetTrackID () const
 
G4int GetParentID () const
 
G4String GetParticleName () const
 
G4double GetCharge () const
 
G4int GetPDGEncoding () const
 
G4ThreeVector GetInitialMomentum () const
 
- Public Member Functions inherited from G4VTrajectory
 G4VTrajectory ()=default
 
virtual ~G4VTrajectory ()=default
 
G4bool operator== (const G4VTrajectory &right) const
 
virtual const std::map< G4String, G4AttDef > * GetAttDefs () const
 
virtual std::vector< G4AttValue > * CreateAttValues () const
 

Detailed Description

Definition at line 56 of file G4RayTrajectory.hh.

Constructor & Destructor Documentation

◆ G4RayTrajectory() [1/2]

G4RayTrajectory::G4RayTrajectory ( )

Definition at line 52 of file G4RayTrajectory.cc.

53{
54 positionRecord = new std::vector<G4RayTrajectoryPoint*>;
55}

◆ G4RayTrajectory() [2/2]

G4RayTrajectory::G4RayTrajectory ( G4RayTrajectory & right)

Definition at line 57 of file G4RayTrajectory.cc.

59{
60 positionRecord = new std::vector<G4RayTrajectoryPoint*>;
61 for(size_t i=0;i<right.positionRecord->size();i++)
62 {
64 ((*(right.positionRecord))[i]);
65 positionRecord->push_back(new G4RayTrajectoryPoint(*rightPoint));
66 }
67}
G4VTrajectory()=default

◆ ~G4RayTrajectory()

G4RayTrajectory::~G4RayTrajectory ( )
virtual

Definition at line 69 of file G4RayTrajectory.cc.

70{
71 //positionRecord->clearAndDestroy();
72 for(size_t i=0;i<positionRecord->size();i++)
73 { delete (*positionRecord)[i]; }
74 positionRecord->clear();
75 delete positionRecord;
76}

Member Function Documentation

◆ AppendStep()

void G4RayTrajectory::AppendStep ( const G4Step * theStep)
virtual

Implements G4VTrajectory.

Definition at line 78 of file G4RayTrajectory.cc.

79{
80 G4RayTrajectoryPoint* trajectoryPoint = new G4RayTrajectoryPoint();
81
82 const G4Step* aStep = theStep;
83 G4Navigator* theNavigator
85
86 // Take care of parallel world(s)
88 {
91 std::vector<G4Navigator*>::iterator iNav =
93 GetActiveNavigatorsIterator();
94 theNavigator = iNav[navID];
95 }
96
97 trajectoryPoint->SetStepLength(aStep->GetStepLength());
98
99 // Surface normal
100 G4bool valid;
101 G4ThreeVector theLocalNormal = theNavigator->GetLocalExitNormal(&valid);
102 if(valid) { theLocalNormal = -theLocalNormal; }
103 G4ThreeVector theGrobalNormal
104 = theNavigator->GetLocalToGlobalTransform().TransformAxis(theLocalNormal);
105 trajectoryPoint->SetSurfaceNormal(theGrobalNormal);
106
108 G4RayTracerSceneHandler* sceneHandler =
109 static_cast<G4RayTracerSceneHandler*>(visManager->GetCurrentSceneHandler());
110 const auto& sceneVisAttsMap = sceneHandler->GetSceneVisAttsMap();
111
112 // Make a path from the preStepPoint touchable
113 G4StepPoint* preStepPoint = aStep -> GetPreStepPoint();
114 const G4VTouchable* preTouchable = preStepPoint->GetTouchable();
115 G4int preDepth = preTouchable->GetHistoryDepth();
116 G4ModelingParameters::PVPointerCopyNoPath localPrePVPointerCopyNoPath;
117 for (G4int i = preDepth; i >= 0; --i) {
118 localPrePVPointerCopyNoPath.push_back
120 (preTouchable->GetVolume(i),preTouchable->GetCopyNumber(i)));
121 }
122
123 // Pick up the vis atts, if any, from the scene handler
124 auto preIterator = sceneVisAttsMap.find(localPrePVPointerCopyNoPath);
125 const G4VisAttributes* preVisAtts;
126 if (preIterator != sceneVisAttsMap.end()) {
127 preVisAtts = &preIterator->second;
128 } else {
129 preVisAtts = 0;
130 }
131 trajectoryPoint->SetPreStepAtt(preVisAtts);
132
133 // Make a path from the postStepPoint touchable
134 G4StepPoint* postStepPoint = aStep -> GetPostStepPoint();
135 const G4VTouchable* postTouchable = postStepPoint->GetTouchable();
136 G4int postDepth = postTouchable->GetHistoryDepth();
137 G4ModelingParameters::PVPointerCopyNoPath localPostPVPointerCopyNoPath;
138 for (G4int i = postDepth; i >= 0; --i) {
139 localPostPVPointerCopyNoPath.push_back
141 (postTouchable->GetVolume(i),postTouchable->GetCopyNumber(i)));
142 }
143
144 // Pick up the vis atts, if any, from the scene handler
145 auto postIterator = sceneVisAttsMap.find(localPostPVPointerCopyNoPath);
146 const G4VisAttributes* postVisAtts;
147 if (postIterator != sceneVisAttsMap.end()) {
148 postVisAtts = &postIterator->second;
149 } else {
150 postVisAtts = 0;
151 }
152 trajectoryPoint->SetPostStepAtt(postVisAtts);
153
154 positionRecord->push_back(trajectoryPoint);
155}
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
std::vector< PVPointerCopyNo > PVPointerCopyNoPath
virtual G4ThreeVector GetLocalExitNormal(G4bool *valid)
const G4AffineTransform GetLocalToGlobalTransform() const
static const G4Step * GetHyperStep()
const std::map< G4ModelingParameters::PVPointerCopyNoPath, G4VisAttributes, PathLessThan > & GetSceneVisAttsMap() const
void SetStepLength(G4double val)
void SetSurfaceNormal(const G4ThreeVector &val)
void SetPreStepAtt(const G4VisAttributes *val)
void SetPostStepAtt(const G4VisAttributes *val)
const G4VTouchable * GetTouchable() const
G4double GetStepLength() const
G4int GetCopyNumber(G4int depth=0) const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
virtual G4int GetHistoryDepth() const
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
G4VSceneHandler * GetCurrentSceneHandler() const
static G4VisManager * GetInstance()

◆ DrawTrajectory()

virtual void G4RayTrajectory::DrawTrajectory ( ) const
inlinevirtual

Reimplemented from G4VTrajectory.

Definition at line 76 of file G4RayTrajectory.hh.

76{;}

◆ GetCharge()

G4double G4RayTrajectory::GetCharge ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 88 of file G4RayTrajectory.hh.

88{ return 0.; }

◆ GetInitialMomentum()

G4ThreeVector G4RayTrajectory::GetInitialMomentum ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 90 of file G4RayTrajectory.hh.

90{ return G4ThreeVector(); }
CLHEP::Hep3Vector G4ThreeVector

◆ GetParentID()

G4int G4RayTrajectory::GetParentID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 86 of file G4RayTrajectory.hh.

86{ return 0; }

◆ GetParticleName()

G4String G4RayTrajectory::GetParticleName ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 87 of file G4RayTrajectory.hh.

87{ return ""; }

◆ GetPDGEncoding()

G4int G4RayTrajectory::GetPDGEncoding ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 89 of file G4RayTrajectory.hh.

89{ return 0; }

◆ GetPoint()

virtual G4VTrajectoryPoint * G4RayTrajectory::GetPoint ( G4int i) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 78 of file G4RayTrajectory.hh.

79 { return (*positionRecord)[i]; }

Referenced by MergeTrajectory().

◆ GetPointC()

G4RayTrajectoryPoint * G4RayTrajectory::GetPointC ( G4int i) const
inline

Definition at line 80 of file G4RayTrajectory.hh.

81 { return (*positionRecord)[i]; }

Referenced by G4TheRayTracer::GenerateColour(), and G4RTRun::RecordEvent().

◆ GetPointEntries()

virtual G4int G4RayTrajectory::GetPointEntries ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 77 of file G4RayTrajectory.hh.

77{return (G4int)positionRecord->size();}

Referenced by G4TheRayTracer::GenerateColour(), MergeTrajectory(), and G4RTRun::RecordEvent().

◆ GetTrackID()

G4int G4RayTrajectory::GetTrackID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 85 of file G4RayTrajectory.hh.

85{ return 0; }

◆ MergeTrajectory()

void G4RayTrajectory::MergeTrajectory ( G4VTrajectory * secondTrajectory)
virtual

Implements G4VTrajectory.

Definition at line 160 of file G4RayTrajectory.cc.

161{
162 if(!secondTrajectory) return;
163
164 G4RayTrajectory* seco = (G4RayTrajectory*)secondTrajectory;
165 G4int ent = seco->GetPointEntries();
166 for(G4int i=0;i<ent;i++)
167 { positionRecord->push_back((G4RayTrajectoryPoint*)seco->GetPoint(i)); }
168 seco->positionRecord->clear();
169}
virtual G4VTrajectoryPoint * GetPoint(G4int i) const
virtual G4int GetPointEntries() const

◆ operator delete()

void G4RayTrajectory::operator delete ( void * aTrajectory)
inline

Definition at line 110 of file G4RayTrajectory.hh.

111{
112 rayTrajectoryAllocator()->FreeSingle((G4RayTrajectory*)aTrajectory);
113}
G4DLLIMPORT G4Allocator< G4RayTrajectory > *& rayTrajectoryAllocator()

◆ operator new()

void * G4RayTrajectory::operator new ( std::size_t )
inline

Definition at line 103 of file G4RayTrajectory.hh.

104{
107 return (void*)rayTrajectoryAllocator()->MallocSingle();
108}

◆ ShowTrajectory()

void G4RayTrajectory::ShowTrajectory ( std::ostream & ) const
virtual

Reimplemented from G4VTrajectory.

Definition at line 157 of file G4RayTrajectory.cc.

158{ }

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