ATLAS Offline Software
Loading...
Searching...
No Matches
AtlasTrajectory.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3*/
4
10
11// Visualization stuff
12#include "G4VVisManager.hh"
13#include "G4Polyline.hh"
14#include "G4Polymarker.hh"
15#include "G4VisAttributes.hh"
16#include "G4Colour.hh"
17
18
19AtlasTrajectory::AtlasTrajectory(const G4Track* track, int subDetVolLevel)
20 : G4Trajectory(track), m_subDetVolLevel(subDetVolLevel)
21{}
22
23
24void AtlasTrajectory::AppendStep(const G4Step* aStep)
25{
26 // only use truth service if there are new any secondaries
27 const int numSecondaries = aStep->GetSecondaryInCurrentStep()->size();
28 // This method not available until G4 10.2
29 //const int numSecondaries = aStep->GetNumberOfSecondariesInCurrentStep();
30
31 if (numSecondaries) {
32 // OK, there was an interation. look at the track, if it
33 // is not a secondary (i.e. we have a connected tree) we
34 // apply the MC truth machinery...
35 TrackHelper tHelper(aStep->GetTrack());
36 if (!tHelper.IsSecondary())
37 {
38 const TruthStrategyManager& sManager =
41 }
42 }
43
44 // Call the base class
45 G4Trajectory::AppendStep(aStep);
46}
47
48#if G4VERSION_NUMBER < 1010
50#else
51void AtlasTrajectory::DrawTrajectory(G4int i_mode) const
52#endif
53{
54 // If i_mode>=0, draws a trajectory as a polyline (blue for
55 // positive, red for negative, green for neutral) and, if i_mode!=0,
56 // adds markers - yellow circles for step points and magenta squares
57 // for auxiliary points, if any - whose screen size in pixels is
58 // given by abs(i_mode)/1000. E.g: i_mode = 5000 gives easily
59 // visible markers.
60
61 G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
62 if (!pVVisManager) return;
63
64 const G4double markerSize = abs(i_mode)*0.001;
65 G4bool lineRequired (i_mode >= 0);
66 G4bool markersRequired (markerSize > 0.);
67
68 G4Polyline trajectoryLine;
69 G4Polymarker stepPoints;
70 G4Polymarker auxiliaryPoints;
71
72 for (G4int i = 0; i < GetPointEntries() ; i++)
73 {
74 G4VTrajectoryPoint* aTrajectoryPoint = GetPoint(i);
75 const std::vector<G4ThreeVector>* auxiliaries
76 = aTrajectoryPoint->GetAuxiliaryPoints();
77 if (auxiliaries)
78 {
79 for (size_t iAux = 0; iAux < auxiliaries->size(); ++iAux)
80 {
81 const G4ThreeVector pos((*auxiliaries)[iAux]);
82 if (lineRequired)
83 {
84 trajectoryLine.push_back(pos);
85 }
86 if (markersRequired)
87 {
88 auxiliaryPoints.push_back(pos);
89 }
90 }
91 }
92 const G4ThreeVector pos(aTrajectoryPoint->GetPosition());
93 if (lineRequired)
94 {
95 trajectoryLine.push_back(pos);
96 }
97 if (markersRequired)
98 {
99 stepPoints.push_back(pos);
100 }
101 }
102
103 if (lineRequired)
104 {
106 getVisualizationHelper()->trackVisualizationScheme();
107 G4Colour colour;
108 const G4double charge = GetCharge();
109 if (visScheme==1)
110 {
111 // Geant3-like color code - neutrals are dashed
112 const G4String pname=GetParticleName();
113 if (pname=="e+" || pname=="e-") colour = G4Colour(1.,1.,0.);
114 else if (pname=="mu+" || pname=="mu-") colour = G4Colour(0.,1.,0.);
115 else if (pname=="gamma") colour = G4Colour(0.,0.,1.);
116 else if (pname=="neutron") colour = G4Colour(0.,0.,0.);
117 else if (charge) colour = G4Colour(1.,0.,0.);
118 }
119 else if (visScheme==2)
120 {
121 if(charge>0.) colour = G4Colour(0.,0.,1.); // Blue = positive.
122 else if(charge<0.) colour = G4Colour(1.,0.,0.); // Red = negative.
123 else colour = G4Colour(0.,1.,0.); // Green = neutral.
124 }
125 else if (visScheme==3)
126 {
127 const G4String pname=GetParticleName();
128 if (pname=="e+" || pname=="e-") colour = G4Colour(1.,1.,0.);
129 else if (pname=="mu+" || pname=="mu-") colour = G4Colour(0.,0.,1.);
130 else if (pname=="gamma") colour = G4Colour(60./256.,79./256.,113./256.);
131 else if (pname=="neutron") colour = G4Colour(1.,1.,1.);
132 else if (pname=="pi-"|| pname=="pi+")
133 colour = G4Colour(250/256.,128/256.,114/256.);
134 else colour = G4Colour(176/256.,48/256.,96/256.);
135 }
136 G4VisAttributes trajectoryLineAttribs(colour);
137 if(visScheme==1 && charge==0)
138 {
139 G4VisAttributes::LineStyle style=G4VisAttributes::dotted;
140 trajectoryLineAttribs.SetLineStyle(style);
141 }
142 trajectoryLine.SetVisAttributes(&trajectoryLineAttribs);
143 pVVisManager->Draw(trajectoryLine);
144 }
145 if (markersRequired)
146 {
147 auxiliaryPoints.SetMarkerType(G4Polymarker::squares);
148 auxiliaryPoints.SetScreenSize(markerSize);
149 auxiliaryPoints.SetFillStyle(G4VMarker::filled);
150 G4VisAttributes auxiliaryPointsAttribs(G4Colour(0.,1.,1.)); // Magenta
151 auxiliaryPoints.SetVisAttributes(&auxiliaryPointsAttribs);
152 pVVisManager->Draw(auxiliaryPoints);
153
154 stepPoints.SetMarkerType(G4Polymarker::circles);
155 stepPoints.SetScreenSize(markerSize);
156 stepPoints.SetFillStyle(G4VMarker::filled);
157 G4VisAttributes stepPointsAttribs(G4Colour(1.,1.,0.)); // Yellow.
158 stepPoints.SetVisAttributes(&stepPointsAttribs);
159 pVVisManager->Draw(stepPoints);
160 }
161}
double charge(const T &p)
Definition AtlasPID.h:997
AtlasTrajectory(const G4Track *aTrack, int subDetVolLevel)
Constructor.
void AppendStep(const G4Step *aStep)
Overriden from G4 in order to do vertex analysis.
int m_subDetVolLevel
The level in the G4 volume hierarchy at which can we find the sub-detector name.
void DrawTrajectory(G4int) const
Visualization stuff.
bool IsSecondary() const
static const TruthController * getTruthController()
Singleton class for creating truth incidents.
bool CreateTruthIncident(const G4Step *, int subDetVolLevel) const
Returns true if any of the truth strategies return true.
static const TruthStrategyManager & GetStrategyManager()
Retrieve the (const) singleton instance.