Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
LXeTrajectory Class Reference

#include <LXeTrajectory.hh>

Inheritance diagram for LXeTrajectory:
G4Trajectory G4VTrajectory

Public Member Functions

 LXeTrajectory ()
 
 LXeTrajectory (const G4Track *aTrack)
 
 LXeTrajectory (LXeTrajectory &)
 
virtual ~LXeTrajectory ()
 
virtual void DrawTrajectory () const
 
voidoperator new (size_t)
 
void operator delete (void *)
 
void SetDrawTrajectory (G4bool b)
 
void WLS ()
 
void SetForceDrawTrajectory (G4bool b)
 
void SetForceNoDrawTrajectory (G4bool b)
 
- Public Member Functions inherited from G4Trajectory
 G4Trajectory ()
 
 G4Trajectory (const G4Track *aTrack)
 
 G4Trajectory (G4Trajectory &)
 
virtual ~G4Trajectory ()
 
voidoperator new (size_t)
 
void operator delete (void *)
 
int operator== (const G4Trajectory &right) const
 
G4int GetTrackID () const
 
G4int GetParentID () const
 
G4String GetParticleName () const
 
G4double GetCharge () const
 
G4int GetPDGEncoding () const
 
G4double GetInitialKineticEnergy () const
 
G4ThreeVector GetInitialMomentum () const
 
virtual void ShowTrajectory (std::ostream &os=G4cout) const
 
virtual void AppendStep (const G4Step *aStep)
 
virtual int GetPointEntries () const
 
virtual G4VTrajectoryPointGetPoint (G4int i) const
 
virtual void MergeTrajectory (G4VTrajectory *secondTrajectory)
 
G4ParticleDefinitionGetParticleDefinition ()
 
virtual const std::map
< G4String, G4AttDef > * 
GetAttDefs () const
 
virtual std::vector< G4AttValue > * CreateAttValues () const
 
- Public Member Functions inherited from G4VTrajectory
 G4VTrajectory ()
 
virtual ~G4VTrajectory ()
 
G4bool operator== (const G4VTrajectory &right) const
 

Detailed Description

Definition at line 45 of file LXeTrajectory.hh.

Constructor & Destructor Documentation

LXeTrajectory::LXeTrajectory ( )

Definition at line 49 of file LXeTrajectory.cc.

50  :G4Trajectory(),fWls(false),fDrawit(false),fForceNoDraw(false),fForceDraw(false)
51 {
52  fParticleDefinition=0;
53 }
LXeTrajectory::LXeTrajectory ( const G4Track aTrack)

Definition at line 57 of file LXeTrajectory.cc.

References G4Track::GetDefinition().

58  :G4Trajectory(aTrack),fWls(false),fDrawit(false)
59 {
60  fParticleDefinition=aTrack->GetDefinition();
61 }
G4ParticleDefinition * GetDefinition() const
LXeTrajectory::LXeTrajectory ( LXeTrajectory right)

Definition at line 65 of file LXeTrajectory.cc.

66  :G4Trajectory(right),fWls(right.fWls),fDrawit(right.fDrawit)
67 {
68  fParticleDefinition=right.fParticleDefinition;
69 }
LXeTrajectory::~LXeTrajectory ( )
virtual

Definition at line 73 of file LXeTrajectory.cc.

73 {}

Member Function Documentation

void LXeTrajectory::DrawTrajectory ( ) const
virtual

Reimplemented from G4Trajectory.

Definition at line 77 of file LXeTrajectory.cc.

References G4Polymarker::circles, G4VVisManager::Draw(), G4VMarker::filled, G4VTrajectoryPoint::GetAuxiliaryPoints(), G4VVisManager::GetConcreteInstance(), G4Trajectory::GetPoint(), G4Trajectory::GetPointEntries(), G4VTrajectoryPoint::GetPosition(), G4OpticalPhoton::OpticalPhotonDefinition(), G4VMarker::SetFillStyle(), G4Polymarker::SetMarkerType(), G4VMarker::SetScreenSize(), G4Visible::SetVisAttributes(), and G4Polymarker::squares.

Referenced by LXeEventAction::EndOfEventAction().

78 {
79  // i_mode is no longer available as an argument of G4VTrajectory.
80  // In this exampple it was always called with an argument of 50.
81  const G4int i_mode = 50;
82  // Consider using commands /vis/modeling/trajectories.
83 
84  //Taken from G4VTrajectory and modified to select colours based on particle
85  //type and to selectively eliminate drawing of certain trajectories.
86 
87  if(!fForceDraw && (!fDrawit || fForceNoDraw))
88  return;
89 
90  // If i_mode>=0, draws a trajectory as a polyline and, if i_mode!=0,
91  // adds markers - yellow circles for step points and magenta squares
92  // for auxiliary points, if any - whose screen size in pixels is
93  // given by std::abs(i_mode)/1000. E.g: i_mode = 5000 gives easily
94  // visible markers.
95 
97  if (!pVVisManager) return;
98 
99  const G4double markerSize = std::abs(i_mode)/1000;
100  G4bool lineRequired (i_mode >= 0);
101  G4bool markersRequired (markerSize > 0.);
102 
103  G4Polyline trajectoryLine;
104  G4Polymarker stepPoints;
105  G4Polymarker auxiliaryPoints;
106 
107  for (G4int i = 0; i < GetPointEntries() ; i++) {
108  G4VTrajectoryPoint* aTrajectoryPoint = GetPoint(i);
109  const std::vector<G4ThreeVector>* auxiliaries
110  = aTrajectoryPoint->GetAuxiliaryPoints();
111  if (auxiliaries) {
112  for (size_t iAux = 0; iAux < auxiliaries->size(); ++iAux) {
113  const G4ThreeVector pos((*auxiliaries)[iAux]);
114  if (lineRequired) {
115  trajectoryLine.push_back(pos);
116  }
117  if (markersRequired) {
118  auxiliaryPoints.push_back(pos);
119  }
120  }
121  }
122  const G4ThreeVector pos(aTrajectoryPoint->GetPosition());
123  if (lineRequired) {
124  trajectoryLine.push_back(pos);
125  }
126  if (markersRequired) {
127  stepPoints.push_back(pos);
128  }
129  }
130 
131  if (lineRequired) {
132  G4Colour colour;
133 
134  if(fParticleDefinition==G4OpticalPhoton::OpticalPhotonDefinition()){
135  if(fWls) //WLS photons are red
136  colour = G4Colour(1.,0.,0.);
137  else{ //Scintillation and Cerenkov photons are green
138  colour = G4Colour(0.,1.,0.);
139  }
140  }
141  else //All other particles are blue
142  colour = G4Colour(0.,0.,1.);
143 
144  G4VisAttributes trajectoryLineAttribs(colour);
145  trajectoryLine.SetVisAttributes(&trajectoryLineAttribs);
146  pVVisManager->Draw(trajectoryLine);
147  }
148  if (markersRequired) {
149  auxiliaryPoints.SetMarkerType(G4Polymarker::squares);
150  auxiliaryPoints.SetScreenSize(markerSize);
151  auxiliaryPoints.SetFillStyle(G4VMarker::filled);
152  G4VisAttributes auxiliaryPointsAttribs(G4Colour(0.,1.,1.)); // Magenta
153  auxiliaryPoints.SetVisAttributes(&auxiliaryPointsAttribs);
154  pVVisManager->Draw(auxiliaryPoints);
155 
157  stepPoints.SetScreenSize(markerSize);
158  stepPoints.SetFillStyle(G4VMarker::filled);
159  G4VisAttributes stepPointsAttribs(G4Colour(1.,1.,0.)); // Yellow
160  stepPoints.SetVisAttributes(&stepPointsAttribs);
161  pVVisManager->Draw(stepPoints);
162  }
163 }
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0
virtual const std::vector< G4ThreeVector > * GetAuxiliaryPoints() const
static G4VVisManager * GetConcreteInstance()
void SetMarkerType(MarkerType)
void SetFillStyle(FillStyle)
int G4int
Definition: G4Types.hh:78
virtual int GetPointEntries() const
bool G4bool
Definition: G4Types.hh:79
void SetVisAttributes(const G4VisAttributes *)
Definition: G4Visible.cc:80
virtual const G4ThreeVector GetPosition() const =0
virtual G4VTrajectoryPoint * GetPoint(G4int i) const
static G4OpticalPhoton * OpticalPhotonDefinition()
double G4double
Definition: G4Types.hh:76
void SetScreenSize(G4double)
void LXeTrajectory::operator delete ( void aTrajectory)
inline

Definition at line 82 of file LXeTrajectory.hh.

References LXeTrajectoryAllocator.

83 {
84  LXeTrajectoryAllocator->FreeSingle((LXeTrajectory*)aTrajectory);
85 }
G4ThreadLocal G4Allocator< LXeTrajectory > * LXeTrajectoryAllocator
void * LXeTrajectory::operator new ( size_t  )
inline

Definition at line 75 of file LXeTrajectory.hh.

References LXeTrajectoryAllocator.

76 {
79  return (void*)LXeTrajectoryAllocator->MallocSingle();
80 }
G4ThreadLocal G4Allocator< LXeTrajectory > * LXeTrajectoryAllocator
void LXeTrajectory::SetDrawTrajectory ( G4bool  b)
inline

Definition at line 59 of file LXeTrajectory.hh.

References test::b.

Referenced by LXeTrackingAction::PostUserTrackingAction().

59 {fDrawit=b;}
void LXeTrajectory::SetForceDrawTrajectory ( G4bool  b)
inline

Definition at line 61 of file LXeTrajectory.hh.

References test::b.

Referenced by LXeEventAction::EndOfEventAction().

61 {fForceDraw=b;}
void LXeTrajectory::SetForceNoDrawTrajectory ( G4bool  b)
inline

Definition at line 62 of file LXeTrajectory.hh.

References test::b.

Referenced by LXeEventAction::EndOfEventAction().

62 {fForceNoDraw=b;}
void LXeTrajectory::WLS ( )
inline

Definition at line 60 of file LXeTrajectory.hh.

Referenced by LXeTrackingAction::PostUserTrackingAction().

60 {fWls=true;}

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