Geant4-11
Public Member Functions | Private Member Functions | Private Attributes
G4RayTrajectory Class Reference

#include <G4RayTrajectory.hh>

Inheritance diagram for G4RayTrajectory:
G4VTrajectory

Public Member Functions

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

Private Member Functions

G4RayTrajectoryoperator= (const G4RayTrajectory &)
 

Private Attributes

std::vector< G4RayTrajectoryPoint * > * positionRecord
 

Detailed Description

Definition at line 56 of file G4RayTrajectory.hh.

Constructor & Destructor Documentation

◆ G4RayTrajectory() [1/2]

G4RayTrajectory::G4RayTrajectory ( )

Definition at line 53 of file G4RayTrajectory.cc.

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

References positionRecord.

◆ G4RayTrajectory() [2/2]

G4RayTrajectory::G4RayTrajectory ( G4RayTrajectory right)

Definition at line 58 of file G4RayTrajectory.cc.

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

References positionRecord.

◆ ~G4RayTrajectory()

G4RayTrajectory::~G4RayTrajectory ( )
virtual

Definition at line 70 of file G4RayTrajectory.cc.

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

References positionRecord.

Member Function Documentation

◆ AppendStep()

void G4RayTrajectory::AppendStep ( const G4Step theStep)
virtual

Implements G4VTrajectory.

Definition at line 79 of file G4RayTrajectory.cc.

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

References G4VTouchable::GetCopyNumber(), G4VisManager::GetCurrentSceneHandler(), G4VTouchable::GetHistoryDepth(), G4ParallelWorldProcess::GetHyperStep(), G4ParallelWorldProcess::GetHypNavigatorID(), G4VisManager::GetInstance(), G4Navigator::GetLocalExitNormal(), G4Navigator::GetLocalToGlobalTransform(), G4TransportationManager::GetNavigatorForTracking(), G4RayTracerSceneHandler::GetSceneVisAttsMap(), G4Step::GetStepLength(), G4StepPoint::GetTouchable(), G4TransportationManager::GetTransportationManager(), G4VTouchable::GetVolume(), positionRecord, G4RayTrajectoryPoint::SetPostStepAtt(), G4RayTrajectoryPoint::SetPreStepAtt(), G4RayTrajectoryPoint::SetStepLength(), G4RayTrajectoryPoint::SetSurfaceNormal(), and G4AffineTransform::TransformAxis().

◆ CreateAttValues()

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

◆ DrawTrajectory()

virtual void G4RayTrajectory::DrawTrajectory ( ) const
inlinevirtual

Reimplemented from G4VTrajectory.

Definition at line 76 of file G4RayTrajectory.hh.

76{;}

◆ GetAttDefs()

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

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

References positionRecord.

Referenced by MergeTrajectory().

◆ GetPointC()

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

Definition at line 80 of file G4RayTrajectory.hh.

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

References positionRecord.

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

◆ GetPointEntries()

virtual int G4RayTrajectory::GetPointEntries ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 77 of file G4RayTrajectory.hh.

77{return positionRecord->size();}

References positionRecord.

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 161 of file G4RayTrajectory.cc.

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

References GetPoint(), GetPointEntries(), and positionRecord.

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

References rayTrajectoryAllocator().

◆ operator new()

void * G4RayTrajectory::operator new ( size_t  )
inline

Definition at line 103 of file G4RayTrajectory.hh.

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

References rayTrajectoryAllocator().

◆ operator=()

G4RayTrajectory & G4RayTrajectory::operator= ( const G4RayTrajectory )
private

◆ operator==()

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

Definition at line 55 of file G4VTrajectory.cc.

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

◆ ShowTrajectory()

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

Reimplemented from G4VTrajectory.

Definition at line 158 of file G4RayTrajectory.cc.

159{ }

Field Documentation

◆ positionRecord

std::vector<G4RayTrajectoryPoint*>* G4RayTrajectory::positionRecord
private

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