Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TrajectoryDrawerUtils.cc
Go to the documentation of this file.
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//
27// Jane Tinslay, John Allison, Joseph Perl November 2005
28//
30#include "G4Colour.hh"
31#include "G4Polyline.hh"
32#include "G4Polymarker.hh"
33#include "G4VTrajectory.hh"
34#include "G4VTrajectoryPoint.hh"
35#include "G4VisAttributes.hh"
36#include "G4VisTrajContext.hh"
37#include "G4VVisManager.hh"
38#include "G4UIcommand.hh"
39#include "G4AttValue.hh"
40
41#define G4warn G4cout
42
44
46
48 (const G4VTrajectory& traj,
49 const G4VisTrajContext& context,
50 G4Polyline& trajectoryLine,
51 G4Polymarker& auxiliaryPoints,
52 G4Polymarker& stepPoints,
53 std::vector<G4double>& trajectoryLineTimes,
54 std::vector<G4double>& auxiliaryPointTimes,
55 std::vector<G4double>& stepPointTimes)
56 {
57 TimesValidity validity = InvalidTimes;
58 if (context.GetTimeSliceInterval()) validity = ValidTimes;
59
60 // Memory for last trajectory point position for auxiliary point
61 // time interpolation algorithm. There are no auxiliary points
62 // for the first trajectory point, so its initial value is
63 // immaterial.
64 G4ThreeVector lastTrajectoryPointPosition;
65
66 // Keep positions. Don't store unless first or different.
67 std::vector<G4ThreeVector> positions;
68
69 for (G4int iPoint=0; iPoint<traj.GetPointEntries(); iPoint++) {
70
71 G4VTrajectoryPoint* aTrajectoryPoint = traj.GetPoint(iPoint);
72 const G4ThreeVector& trajectoryPointPosition =
73 aTrajectoryPoint->GetPosition();
74
75 // Only store if first or if different
76 if (positions.size() == 0 ||
77 trajectoryPointPosition != positions[positions.size()-1]) {
78
79 // Pre- and Post-Point times from the trajectory point...
80 G4double trajectoryPointPreTime = -std::numeric_limits<double>::max();
81 G4double trajectoryPointPostTime = std::numeric_limits<double>::max();
82
83 if (context.GetTimeSliceInterval() && validity == ValidTimes) {
84
85 std::vector<G4AttValue>* trajectoryPointAttValues =
86 aTrajectoryPoint->CreateAttValues();
87 if (!trajectoryPointAttValues) {
88 static G4bool warnedNoAttValues = false;
89 if (!warnedNoAttValues) {
90 G4warn <<
91 "*************************************************************************"
92 "\n* WARNING: G4TrajectoryDrawerUtils::GetPointsAndTimes: no att values."
93 "\n*************************************************************************"
94 << G4endl;
95 warnedNoAttValues = true;
96 }
97 validity = InvalidTimes;
98 } else {
99 G4bool foundPreTime = false, foundPostTime = false;
100 for (std::vector<G4AttValue>::iterator i =
101 trajectoryPointAttValues->begin();
102 i != trajectoryPointAttValues->end(); ++i) {
103 if (i->GetName() == "PreT") {
104 trajectoryPointPreTime =
106 foundPreTime = true;
107 }
108 if (i->GetName() == "PostT") {
109 trajectoryPointPostTime =
111 foundPostTime = true;
112 }
113 }
114 if (!foundPreTime || !foundPostTime) {
115 static G4bool warnedTimesNotFound = false;
116 if (!warnedTimesNotFound) {
117 G4warn <<
118 "*************************************************************************"
119 "\n* WARNING: G4TrajectoryDrawerUtils::GetPointsAndTimes: times not found."
120 "\n You need to specify \"/vis/scene/add/trajectories rich\""
121 "\n*************************************************************************"
122 << G4endl;
123 warnedTimesNotFound = true;
124 }
125 validity = InvalidTimes;
126 }
127 }
128 delete trajectoryPointAttValues; // (Must be deleted after use.)
129 }
130
131 const std::vector<G4ThreeVector>* auxiliaries
132 = aTrajectoryPoint->GetAuxiliaryPoints();
133 if (0 != auxiliaries) {
134 for (size_t iAux=0; iAux<auxiliaries->size(); ++iAux) {
135 const G4ThreeVector& auxPointPosition = (*auxiliaries)[iAux];
136 if (positions.size() == 0 ||
137 auxPointPosition != positions[positions.size()-1]) {
138 // Only store if first or if different
139 positions.push_back(trajectoryPointPosition);
140 trajectoryLine.push_back(auxPointPosition);
141 auxiliaryPoints.push_back(auxPointPosition);
142 if (validity == ValidTimes) {
143 // Interpolate time for auxiliary points...
144 G4double s1 =
145 (auxPointPosition - lastTrajectoryPointPosition).mag();
146 G4double s2 =
147 (trajectoryPointPosition - auxPointPosition).mag();
148 G4double t = trajectoryPointPreTime +
149 (trajectoryPointPostTime - trajectoryPointPreTime) *
150 (s1 / (s1 + s2));
151 trajectoryLineTimes.push_back(t);
152 auxiliaryPointTimes.push_back(t);
153 }
154 }
155 }
156 }
157
158 positions.push_back(trajectoryPointPosition);
159 trajectoryLine.push_back(trajectoryPointPosition);
160 stepPoints.push_back(trajectoryPointPosition);
161 if (validity == ValidTimes) {
162 trajectoryLineTimes.push_back(trajectoryPointPostTime);
163 stepPointTimes.push_back(trajectoryPointPostTime);
164 }
165 lastTrajectoryPointPosition = trajectoryPointPosition;
166 }
167 }
168 return validity;
169 }
170
171 static void SliceLine(G4double timeIncrement,
172 G4Polyline& trajectoryLine,
173 std::vector<G4double>& trajectoryLineTimes)
174 {
175 // Assumes valid arguments from GetPoints and GetTimes.
176
177 G4Polyline newTrajectoryLine;
178 std::vector<G4double> newTrajectoryLineTimes;
179
180 newTrajectoryLine.push_back(trajectoryLine[0]);
181 newTrajectoryLineTimes.push_back(trajectoryLineTimes[0]);
182 size_t lineSize = trajectoryLine.size();
183 if (lineSize > 1) {
184 for (size_t i = 1; i < trajectoryLine.size(); ++i) {
185 G4double deltaT = trajectoryLineTimes[i] - trajectoryLineTimes[i - 1];
186 if (deltaT > 0.) {
187 G4double practicalTimeIncrement =
188 std::max(timeIncrement, deltaT / 100.);
189 for (G4double t =
190 (int(trajectoryLineTimes[i - 1]/practicalTimeIncrement) + 1) *
191 practicalTimeIncrement;
192 t <= trajectoryLineTimes[i];
193 t += practicalTimeIncrement) {
194 G4ThreeVector pos = trajectoryLine[i - 1] +
195 (trajectoryLine[i] - trajectoryLine[i - 1]) *
196 ((t - trajectoryLineTimes[i - 1]) / deltaT);
197 newTrajectoryLine.push_back(pos);
198 newTrajectoryLineTimes.push_back(t);
199 }
200 }
201 newTrajectoryLine.push_back(trajectoryLine[i]);
202 newTrajectoryLineTimes.push_back(trajectoryLineTimes[i]);
203 }
204 }
205
206 trajectoryLine = newTrajectoryLine;
207 trajectoryLineTimes = newTrajectoryLineTimes;
208 }
209
210 static void DrawWithoutTime(const G4VisTrajContext& myContext,
211 G4Polyline& trajectoryLine,
212 G4Polymarker& auxiliaryPoints,
213 G4Polymarker& stepPoints)
214 {
215 // Draw without time slice information
216
218 if (0 == pVVisManager) return;
219
220 if (myContext.GetDrawLine() && myContext.GetLineVisible()) {
221 G4VisAttributes trajectoryLineAttribs(myContext.GetLineColour());
222 trajectoryLineAttribs.SetLineWidth(myContext.GetLineWidth());
223 trajectoryLine.SetVisAttributes(&trajectoryLineAttribs);
224
225 pVVisManager->Draw(trajectoryLine);
226 }
227
228 if (myContext.GetDrawAuxPts() && myContext.GetAuxPtsVisible()
229 && (auxiliaryPoints.size() > 0)) {
230 auxiliaryPoints.SetMarkerType(myContext.GetAuxPtsType());
231 auxiliaryPoints.SetSize(myContext.GetAuxPtsSizeType(), myContext.GetAuxPtsSize());
232 auxiliaryPoints.SetFillStyle(myContext.GetAuxPtsFillStyle());
233
234 G4VisAttributes auxiliaryPointsAttribs(myContext.GetAuxPtsColour());
235 auxiliaryPoints.SetVisAttributes(&auxiliaryPointsAttribs);
236
237 pVVisManager->Draw(auxiliaryPoints);
238 }
239
240 if (myContext.GetDrawStepPts() && myContext.GetStepPtsVisible()
241 && (stepPoints.size() > 0)) {
242 stepPoints.SetMarkerType(myContext.GetStepPtsType());
243 stepPoints.SetSize(myContext.GetStepPtsSizeType(), myContext.GetStepPtsSize());
244 stepPoints.SetFillStyle(myContext.GetStepPtsFillStyle());
245
246 G4VisAttributes stepPointsAttribs(myContext.GetStepPtsColour());
247 stepPoints.SetVisAttributes(&stepPointsAttribs);
248
249 pVVisManager->Draw(stepPoints);
250 }
251 }
252
253 static void DrawWithTime(const G4VisTrajContext& myContext,
254 G4Polyline& trajectoryLine,
255 G4Polymarker& auxiliaryPoints,
256 G4Polymarker& stepPoints,
257 std::vector<G4double>& trajectoryLineTimes,
258 std::vector<G4double>& auxiliaryPointTimes,
259 std::vector<G4double>& stepPointTimes)
260 {
261 // Draw with time slice information
262
264 if (0 == pVVisManager) return;
265
266 if (myContext.GetDrawLine() && myContext.GetLineVisible()) {
267 G4VisAttributes trajectoryLineAttribs(myContext.GetLineColour());
268 trajectoryLineAttribs.SetLineWidth(myContext.GetLineWidth());
269
270 for (size_t i = 1; i < trajectoryLine.size(); ++i ) {
271 G4Polyline slice;
272 slice.push_back(trajectoryLine[i -1]);
273 slice.push_back(trajectoryLine[i]);
274 trajectoryLineAttribs.SetStartTime(trajectoryLineTimes[i - 1]);
275 trajectoryLineAttribs.SetEndTime(trajectoryLineTimes[i]);
276 slice.SetVisAttributes(&trajectoryLineAttribs);
277 pVVisManager->Draw(slice);
278 }
279 }
280
281 if (myContext.GetDrawAuxPts() && myContext.GetAuxPtsVisible()
282 && (auxiliaryPoints.size() > 0)) {
283 G4VisAttributes auxiliaryPointsAttribs(myContext.GetAuxPtsColour());
284
285 for (size_t i = 0; i < auxiliaryPoints.size(); ++i ) {
286 G4Polymarker point;
287 point.push_back(auxiliaryPoints[i]);
288 point.SetMarkerType(myContext.GetAuxPtsType());
289 point.SetSize(myContext.GetAuxPtsSizeType(), myContext.GetAuxPtsSize());
290 point.SetFillStyle(myContext.GetAuxPtsFillStyle());
291 auxiliaryPointsAttribs.SetStartTime(auxiliaryPointTimes[i]);
292 auxiliaryPointsAttribs.SetEndTime(auxiliaryPointTimes[i]);
293 point.SetVisAttributes(&auxiliaryPointsAttribs);
294 pVVisManager->Draw(point);
295 }
296 }
297
298 if (myContext.GetDrawStepPts() && myContext.GetStepPtsVisible()
299 && (stepPoints.size() > 0)) {
300 G4VisAttributes stepPointsAttribs(myContext.GetStepPtsColour());
301
302 for (size_t i = 0; i < stepPoints.size(); ++i ) {
303 G4Polymarker point;
304 point.push_back(stepPoints[i]);
305 point.SetMarkerType(myContext.GetStepPtsType());
306 point.SetSize(myContext.GetStepPtsSizeType(), myContext.GetStepPtsSize());
307 point.SetFillStyle(myContext.GetStepPtsFillStyle());
308 stepPointsAttribs.SetStartTime(stepPointTimes[i]);
309 stepPointsAttribs.SetEndTime(stepPointTimes[i]);
310 point.SetVisAttributes(&stepPointsAttribs);
311 pVVisManager->Draw(point);
312 }
313 }
314 }
315
316 void DrawLineAndPoints(const G4VTrajectory& traj, const G4VisTrajContext& context)
317 {
318 // Return if don't need to do anything
319 if (!context.GetDrawLine() && !context.GetDrawAuxPts() && !context.GetDrawStepPts()) return;
320
321 // Get points and times (times are returned only if time-slicing
322 // is requested).
323 G4Polyline trajectoryLine;
324 G4Polymarker stepPoints;
325 G4Polymarker auxiliaryPoints;
326 std::vector<G4double> trajectoryLineTimes;
327 std::vector<G4double> stepPointTimes;
328 std::vector<G4double> auxiliaryPointTimes;
329
331 (traj, context,
332 trajectoryLine, auxiliaryPoints, stepPoints,
333 trajectoryLineTimes, auxiliaryPointTimes, stepPointTimes);
334
335 if (validity == ValidTimes) {
336
337 SliceLine(context.GetTimeSliceInterval(),
338 trajectoryLine, trajectoryLineTimes);
339
340 DrawWithTime(context,
341 trajectoryLine, auxiliaryPoints, stepPoints,
342 trajectoryLineTimes, auxiliaryPointTimes, stepPointTimes);
343
344 } else {
345
346 DrawWithoutTime(context, trajectoryLine, auxiliaryPoints, stepPoints);
347
348 }
349 }
350}
#define G4warn
Definition G4Scene.cc:41
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
void SetMarkerType(MarkerType)
static G4double ConvertToDimensionedDouble(const char *st)
void SetSize(SizeType, G4double)
Definition G4VMarker.cc:86
void SetFillStyle(FillStyle)
virtual std::vector< G4AttValue > * CreateAttValues() const
virtual const std::vector< G4ThreeVector > * GetAuxiliaryPoints() const
virtual const G4ThreeVector GetPosition() const =0
virtual G4VTrajectoryPoint * GetPoint(G4int i) const =0
virtual G4int GetPointEntries() const =0
static G4VVisManager * GetConcreteInstance()
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0
G4bool GetDrawAuxPts() const
G4Colour GetStepPtsColour() const
G4double GetLineWidth() const
G4double GetStepPtsSize() const
G4bool GetLineVisible() const
G4double GetTimeSliceInterval() const
G4bool GetDrawLine() const
G4VMarker::SizeType GetStepPtsSizeType() const
G4Polymarker::MarkerType GetAuxPtsType() const
G4double GetAuxPtsSize() const
G4VMarker::SizeType GetAuxPtsSizeType() const
G4Colour GetAuxPtsColour() const
G4VMarker::FillStyle GetStepPtsFillStyle() const
G4VMarker::FillStyle GetAuxPtsFillStyle() const
G4bool GetAuxPtsVisible() const
G4bool GetStepPtsVisible() const
G4Colour GetLineColour() const
G4Polymarker::MarkerType GetStepPtsType() const
G4bool GetDrawStepPts() const
void SetVisAttributes(const G4VisAttributes *)
Definition G4Visible.cc:98
TimesValidity GetPointsAndTimes(const G4VTrajectory &traj, const G4VisTrajContext &context, G4Polyline &trajectoryLine, G4Polymarker &auxiliaryPoints, G4Polymarker &stepPoints, std::vector< G4double > &trajectoryLineTimes, std::vector< G4double > &auxiliaryPointTimes, std::vector< G4double > &stepPointTimes)
void DrawLineAndPoints(const G4VTrajectory &traj, const G4VisTrajContext &)