Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
LXeTrajectory Class Reference

#include <Doxymodules_optical.h>

Inheritance diagram for LXeTrajectory:
G4Trajectory G4VTrajectory

Public Member Functions

 LXeTrajectory ()=default
 
 LXeTrajectory (const G4Track *aTrack)
 
 LXeTrajectory (LXeTrajectory &)
 
 ~LXeTrajectory ()=default
 
void DrawTrajectory () const override
 
void * operator new (size_t)
 
void operator delete (void *)
 
void SetDrawTrajectory (G4bool b)
 
void WLS ()
 
void SetForceDrawTrajectory (G4bool b)
 
void SetForceNoDrawTrajectory (G4bool b)
 

Private Attributes

G4bool fWls = false
 
G4bool fDrawit = false
 
G4bool fForceNoDraw = false
 
G4bool fForceDraw = false
 
G4ParticleDefinitionfParticleDefinition = nullptr
 

Detailed Description

Definition at line 69 of file Doxymodules_optical.h.

Constructor & Destructor Documentation

◆ LXeTrajectory() [1/3]

LXeTrajectory::LXeTrajectory ( )
default

◆ LXeTrajectory() [2/3]

LXeTrajectory::LXeTrajectory ( const G4Track aTrack)

Definition at line 49 of file LXeTrajectory.cc.

50 : G4Trajectory(aTrack)
51{
52 fParticleDefinition = aTrack->GetDefinition();
53}
G4ParticleDefinition * fParticleDefinition

◆ LXeTrajectory() [3/3]

LXeTrajectory::LXeTrajectory ( LXeTrajectory right)

Definition at line 57 of file LXeTrajectory.cc.

59 , fWls(right.fWls)
60 , fDrawit(right.fDrawit)
61{
62 fParticleDefinition = right.fParticleDefinition;
63}

◆ ~LXeTrajectory()

LXeTrajectory::~LXeTrajectory ( )
default

Member Function Documentation

◆ DrawTrajectory()

void LXeTrajectory::DrawTrajectory ( ) const
override

Definition at line 67 of file LXeTrajectory.cc.

68{
69 // Taken from G4VTrajectory and modified to select colours based on particle
70 // type and to selectively eliminate drawing of certain trajectories.
71
72 if(!fForceDraw && (!fDrawit || fForceNoDraw))
73 return;
74
75 G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
76 if(!pVVisManager)
77 return;
78
79 const G4double markerSize = 0.05;
80 G4bool lineRequired = true;
81 G4bool markersRequired = true;
82
83 G4Polyline trajectoryLine;
84 G4Polymarker stepPoints;
85 G4Polymarker auxiliaryPoints;
86
87 for(G4int i = 0; i < GetPointEntries(); ++i)
88 {
89 G4VTrajectoryPoint* aTrajectoryPoint = GetPoint(i);
90 const std::vector<G4ThreeVector>* auxiliaries =
91 aTrajectoryPoint->GetAuxiliaryPoints();
92 if(auxiliaries)
93 {
94 for(size_t iAux = 0; iAux < auxiliaries->size(); ++iAux)
95 {
96 const G4ThreeVector pos((*auxiliaries)[iAux]);
97 if(lineRequired)
98 {
99 trajectoryLine.push_back(pos);
100 }
101 if(markersRequired)
102 {
103 auxiliaryPoints.push_back(pos);
104 }
105 }
106 }
107 const G4ThreeVector pos(aTrajectoryPoint->GetPosition());
108 if(lineRequired)
109 {
110 trajectoryLine.push_back(pos);
111 }
112 if(markersRequired)
113 {
114 stepPoints.push_back(pos);
115 }
116 }
117
118 if(lineRequired)
119 {
120 G4Colour colour;
121
122 if(fParticleDefinition == G4OpticalPhoton::OpticalPhotonDefinition())
123 {
124 if(fWls) // WLS photons are red
125 colour = G4Colour(1., 0., 0.);
126 else
127 { // Scintillation and Cerenkov photons are green
128 colour = G4Colour(0., 1., 0.);
129 }
130 }
131 else // All other particles are blue
132 colour = G4Colour(0., 0., 1.);
133
134 G4VisAttributes trajectoryLineAttribs(colour);
135 trajectoryLine.SetVisAttributes(&trajectoryLineAttribs);
136 pVVisManager->Draw(trajectoryLine);
137 }
138 if(markersRequired)
139 {
140 auxiliaryPoints.SetMarkerType(G4Polymarker::squares);
141 auxiliaryPoints.SetScreenSize(markerSize);
142 auxiliaryPoints.SetFillStyle(G4VMarker::filled);
143 G4VisAttributes auxiliaryPointsAttribs(G4Colour(0., 1., 1.)); // Magenta
144 auxiliaryPoints.SetVisAttributes(&auxiliaryPointsAttribs);
145 pVVisManager->Draw(auxiliaryPoints);
146
147 stepPoints.SetMarkerType(G4Polymarker::circles);
148 stepPoints.SetScreenSize(markerSize);
149 stepPoints.SetFillStyle(G4VMarker::filled);
150 G4VisAttributes stepPointsAttribs(G4Colour(1., 1., 0.)); // Yellow
151 stepPoints.SetVisAttributes(&stepPointsAttribs);
152 pVVisManager->Draw(stepPoints);
153 }
154}

◆ operator new()

void * LXeTrajectory::operator new ( size_t  )
inline

Definition at line 68 of file LXeTrajectory.hh.

69{
72 return (void*) LXeTrajectoryAllocator->MallocSingle();
73}
G4ThreadLocal G4Allocator< LXeTrajectory > * LXeTrajectoryAllocator

◆ operator delete()

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

Definition at line 75 of file LXeTrajectory.hh.

76{
77 LXeTrajectoryAllocator->FreeSingle((LXeTrajectory*) aTrajectory);
78}

◆ SetDrawTrajectory()

void LXeTrajectory::SetDrawTrajectory ( G4bool  b)
inline

Definition at line 53 of file LXeTrajectory.hh.

53{ fDrawit = b; }

◆ WLS()

void LXeTrajectory::WLS ( )
inline

Definition at line 54 of file LXeTrajectory.hh.

54{ fWls = true; }

◆ SetForceDrawTrajectory()

void LXeTrajectory::SetForceDrawTrajectory ( G4bool  b)
inline

Definition at line 55 of file LXeTrajectory.hh.

55{ fForceDraw = b; }

◆ SetForceNoDrawTrajectory()

void LXeTrajectory::SetForceNoDrawTrajectory ( G4bool  b)
inline

Definition at line 56 of file LXeTrajectory.hh.

56{ fForceNoDraw = b; }

Member Data Documentation

◆ fWls

G4bool LXeTrajectory::fWls = false
private

Definition at line 59 of file LXeTrajectory.hh.

◆ fDrawit

G4bool LXeTrajectory::fDrawit = false
private

Definition at line 60 of file LXeTrajectory.hh.

◆ fForceNoDraw

G4bool LXeTrajectory::fForceNoDraw = false
private

Definition at line 61 of file LXeTrajectory.hh.

◆ fForceDraw

G4bool LXeTrajectory::fForceDraw = false
private

Definition at line 62 of file LXeTrajectory.hh.

◆ fParticleDefinition

G4ParticleDefinition* LXeTrajectory::fParticleDefinition = nullptr
private

Definition at line 63 of file LXeTrajectory.hh.


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

Applications | User Support | Publications | Collaboration