Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
LXeTrajectory.cc
이 파일의 문서화 페이지로 가기
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: LXeTrajectory.cc 110138 2018-05-16 07:31:43Z gcosmo $
27 //
30 //
31 //
32 #include "LXeTrajectory.hh"
33 #include "G4TrajectoryPoint.hh"
34 #include "G4Trajectory.hh"
35 #include "G4ParticleTable.hh"
36 #include "G4ParticleTypes.hh"
37 #include "G4ThreeVector.hh"
38 #include "G4Polyline.hh"
39 #include "G4Circle.hh"
40 #include "G4Colour.hh"
41 #include "G4VisAttributes.hh"
42 #include "G4VVisManager.hh"
43 #include "G4Polymarker.hh"
44 
46 
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48 
50  :G4Trajectory(),fWls(false),fDrawit(false),
51  fForceNoDraw(false),fForceDraw(false)
52 {
53  fParticleDefinition = nullptr;
54 }
55 
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57 
59  :G4Trajectory(aTrack),fWls(false),fDrawit(false)
60 {
62 }
63 
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65 
67  :G4Trajectory(right),fWls(right.fWls),fDrawit(right.fDrawit)
68 {
70 }
71 
72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73 
75 
76 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
77 
79 {
80  // i_mode is no longer available as an argument of G4VTrajectory.
81  // In this exampple it was always called with an argument of 50.
82  const G4int i_mode = 50;
83  // Consider using commands /vis/modeling/trajectories.
84 
85  //Taken from G4VTrajectory and modified to select colours based on particle
86  //type and to selectively eliminate drawing of certain trajectories.
87 
88  if(!fForceDraw && (!fDrawit || fForceNoDraw))
89  return;
90 
91  // If i_mode>=0, draws a trajectory as a polyline and, if i_mode!=0,
92  // adds markers - yellow circles for step points and magenta squares
93  // for auxiliary points, if any - whose screen size in pixels is
94  // given by std::abs(i_mode)/1000. E.g: i_mode = 5000 gives easily
95  // visible markers.
96 
98  if (!pVVisManager) return;
99 
100  const G4double markerSize = std::abs(i_mode)/1000;
101  G4bool lineRequired (i_mode >= 0);
102  G4bool markersRequired (markerSize > 0.);
103 
104  G4Polyline trajectoryLine;
105  G4Polymarker stepPoints;
106  G4Polymarker auxiliaryPoints;
107 
108  for (G4int i = 0; i < GetPointEntries() ; i++) {
109  G4VTrajectoryPoint* aTrajectoryPoint = GetPoint(i);
110  const std::vector<G4ThreeVector>* auxiliaries
111  = aTrajectoryPoint->GetAuxiliaryPoints();
112  if (auxiliaries) {
113  for (size_t iAux = 0; iAux < auxiliaries->size(); ++iAux) {
114  const G4ThreeVector pos((*auxiliaries)[iAux]);
115  if (lineRequired) {
116  trajectoryLine.push_back(pos);
117  }
118  if (markersRequired) {
119  auxiliaryPoints.push_back(pos);
120  }
121  }
122  }
123  const G4ThreeVector pos(aTrajectoryPoint->GetPosition());
124  if (lineRequired) {
125  trajectoryLine.push_back(pos);
126  }
127  if (markersRequired) {
128  stepPoints.push_back(pos);
129  }
130  }
131 
132  if (lineRequired) {
133  G4Colour colour;
134 
136  if(fWls) //WLS photons are red
137  colour = G4Colour(1.,0.,0.);
138  else{ //Scintillation and Cerenkov photons are green
139  colour = G4Colour(0.,1.,0.);
140  }
141  }
142  else //All other particles are blue
143  colour = G4Colour(0.,0.,1.);
144 
145  G4VisAttributes trajectoryLineAttribs(colour);
146  trajectoryLine.SetVisAttributes(&trajectoryLineAttribs);
147  pVVisManager->Draw(trajectoryLine);
148  }
149  if (markersRequired) {
150  auxiliaryPoints.SetMarkerType(G4Polymarker::squares);
151  auxiliaryPoints.SetScreenSize(markerSize);
152  auxiliaryPoints.SetFillStyle(G4VMarker::filled);
153  G4VisAttributes auxiliaryPointsAttribs(G4Colour(0.,1.,1.)); // Magenta
154  auxiliaryPoints.SetVisAttributes(&auxiliaryPointsAttribs);
155  pVVisManager->Draw(auxiliaryPoints);
156 
158  stepPoints.SetScreenSize(markerSize);
159  stepPoints.SetFillStyle(G4VMarker::filled);
160  G4VisAttributes stepPointsAttribs(G4Colour(1.,1.,0.)); // Yellow
161  stepPoints.SetVisAttributes(&stepPointsAttribs);
162  pVVisManager->Draw(stepPoints);
163  }
164 }
virtual const G4ThreeVector GetPosition() const =0
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0
static const G4double pos
G4ParticleDefinition * fParticleDefinition
virtual int GetPointEntries() const
virtual ~LXeTrajectory()
#define G4ThreadLocal
Definition: tls.hh:69
G4ParticleDefinition * GetDefinition() const
static G4VVisManager * GetConcreteInstance()
void SetScreenSize(G4double)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
static G4OpticalPhoton * OpticalPhotonDefinition()
G4bool fForceNoDraw
G4ThreadLocal G4Allocator< LXeTrajectory > * LXeTrajectoryAllocator
virtual void DrawTrajectory() const
int G4int
Definition: G4Types.hh:78
virtual const std::vector< G4ThreeVector > * GetAuxiliaryPoints() const
virtual G4VTrajectoryPoint * GetPoint(G4int i) const
void SetMarkerType(MarkerType)
Definition of the LXeTrajectory class.
void SetVisAttributes(const G4VisAttributes *)
Definition: G4Visible.cc:80
void SetFillStyle(FillStyle)