Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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// $Id$
27//
28// Jane Tinslay, John Allison, Joseph Perl November 2005
29//
31#include "G4Colour.hh"
32#include "G4Polyline.hh"
33#include "G4Polymarker.hh"
34#include "G4VTrajectory.hh"
35#include "G4VTrajectoryPoint.hh"
36#include "G4VisAttributes.hh"
37#include "G4VisTrajContext.hh"
38#include "G4VVisManager.hh"
39#include "G4UIcommand.hh"
40#include "G4AttValue.hh"
41
43
45
47 (const G4VTrajectory& traj,
48 const G4VisTrajContext& context,
49 G4Polyline& trajectoryLine,
50 G4Polymarker& auxiliaryPoints,
51 G4Polymarker& stepPoints,
52 std::vector<G4double>& trajectoryLineTimes,
53 std::vector<G4double>& auxiliaryPointTimes,
54 std::vector<G4double>& stepPointTimes)
55 {
56 TimesValidity validity = InvalidTimes;
57 if (context.GetTimeSliceInterval()) validity = ValidTimes;
58
59 // Memory for last trajectory point position for auxiliary point
60 // time interpolation algorithm. There are no auxiliary points
61 // for the first trajectory point, so its initial value is
62 // immaterial.
63 G4ThreeVector lastTrajectoryPointPosition;
64
65 // Keep positions. Don't store unless first or different.
66 std::vector<G4ThreeVector> positions;
67
68 for (G4int iPoint=0; iPoint<traj.GetPointEntries(); iPoint++) {
69
70 G4VTrajectoryPoint* aTrajectoryPoint = traj.GetPoint(iPoint);
71 const G4ThreeVector& trajectoryPointPosition =
72 aTrajectoryPoint->GetPosition();
73
74 // Only store if first or if different
75 if (positions.size() == 0 ||
76 trajectoryPointPosition != positions[positions.size()-1]) {
77
78 // Pre- and Post-Point times from the trajectory point...
79 G4double trajectoryPointPreTime = -std::numeric_limits<double>::max();
80 G4double trajectoryPointPostTime = std::numeric_limits<double>::max();
81
82 if (context.GetTimeSliceInterval() && validity == ValidTimes) {
83
84 std::vector<G4AttValue>* trajectoryPointAttValues =
85 aTrajectoryPoint->CreateAttValues();
86 if (!trajectoryPointAttValues) {
87 static G4bool warnedNoAttValues = false;
88 if (!warnedNoAttValues) {
89 G4cout <<
90 "*************************************************************************"
91 "\n* WARNING: G4TrajectoryDrawerUtils::GetPointsAndTimes: no att values."
92 "\n*************************************************************************"
93 << G4endl;
94 warnedNoAttValues = true;
95 }
96 validity = InvalidTimes;
97 } else {
98 G4bool foundPreTime = false, foundPostTime = false;
99 for (std::vector<G4AttValue>::iterator i =
100 trajectoryPointAttValues->begin();
101 i != trajectoryPointAttValues->end(); ++i) {
102 if (i->GetName() == "PreT") {
103 trajectoryPointPreTime =
105 foundPreTime = true;
106 }
107 if (i->GetName() == "PostT") {
108 trajectoryPointPostTime =
110 foundPostTime = true;
111 }
112 }
113 if (!foundPreTime || !foundPostTime) {
114 static G4bool warnedTimesNotFound = false;
115 if (!warnedTimesNotFound) {
116 G4cout <<
117 "*************************************************************************"
118 "\n* WARNING: G4TrajectoryDrawerUtils::GetPointsAndTimes: times not found."
119 "\n*************************************************************************"
120 << G4endl;
121 warnedTimesNotFound = true;
122 }
123 validity = InvalidTimes;
124 }
125 }
126 }
127
128 const std::vector<G4ThreeVector>* auxiliaries
129 = aTrajectoryPoint->GetAuxiliaryPoints();
130 if (0 != auxiliaries) {
131 for (size_t iAux=0; iAux<auxiliaries->size(); ++iAux) {
132 const G4ThreeVector& auxPointPosition = (*auxiliaries)[iAux];
133 if (positions.size() == 0 ||
134 auxPointPosition != positions[positions.size()-1]) {
135 // Only store if first or if different
136 positions.push_back(trajectoryPointPosition);
137 trajectoryLine.push_back(auxPointPosition);
138 auxiliaryPoints.push_back(auxPointPosition);
139 if (validity == ValidTimes) {
140 // Interpolate time for auxiliary points...
141 G4double s1 =
142 (auxPointPosition - lastTrajectoryPointPosition).mag();
143 G4double s2 =
144 (trajectoryPointPosition - auxPointPosition).mag();
145 G4double t = trajectoryPointPreTime +
146 (trajectoryPointPostTime - trajectoryPointPreTime) *
147 (s1 / (s1 + s2));
148 trajectoryLineTimes.push_back(t);
149 auxiliaryPointTimes.push_back(t);
150 }
151 }
152 }
153 }
154
155 positions.push_back(trajectoryPointPosition);
156 trajectoryLine.push_back(trajectoryPointPosition);
157 stepPoints.push_back(trajectoryPointPosition);
158 if (validity == ValidTimes) {
159 trajectoryLineTimes.push_back(trajectoryPointPostTime);
160 stepPointTimes.push_back(trajectoryPointPostTime);
161 }
162 lastTrajectoryPointPosition = trajectoryPointPosition;
163 }
164 }
165 return validity;
166 }
167
168 static void SliceLine(G4double timeIncrement,
169 G4Polyline& trajectoryLine,
170 std::vector<G4double>& trajectoryLineTimes)
171 {
172 // Assumes valid arguments from GetPoints and GetTimes.
173
174 G4Polyline newTrajectoryLine;
175 std::vector<G4double> newTrajectoryLineTimes;
176
177 newTrajectoryLine.push_back(trajectoryLine[0]);
178 newTrajectoryLineTimes.push_back(trajectoryLineTimes[0]);
179 size_t lineSize = trajectoryLine.size();
180 if (lineSize > 1) {
181 for (size_t i = 1; i < trajectoryLine.size(); ++i) {
182 G4double deltaT = trajectoryLineTimes[i] - trajectoryLineTimes[i - 1];
183 if (deltaT > 0.) {
184 G4double practicalTimeIncrement =
185 std::max(timeIncrement, deltaT / 100.);
186 for (G4double t =
187 (int(trajectoryLineTimes[i - 1]/practicalTimeIncrement) + 1) *
188 practicalTimeIncrement;
189 t <= trajectoryLineTimes[i];
190 t += practicalTimeIncrement) {
191 G4ThreeVector pos = trajectoryLine[i - 1] +
192 (trajectoryLine[i] - trajectoryLine[i - 1]) *
193 ((t - trajectoryLineTimes[i - 1]) / deltaT);
194 newTrajectoryLine.push_back(pos);
195 newTrajectoryLineTimes.push_back(t);
196 }
197 }
198 newTrajectoryLine.push_back(trajectoryLine[i]);
199 newTrajectoryLineTimes.push_back(trajectoryLineTimes[i]);
200 }
201 }
202
203 trajectoryLine = newTrajectoryLine;
204 trajectoryLineTimes = newTrajectoryLineTimes;
205 }
206
207 static void DrawWithoutTime(const G4VisTrajContext& myContext,
208 G4Polyline& trajectoryLine,
209 G4Polymarker& auxiliaryPoints,
210 G4Polymarker& stepPoints)
211 {
212 // Draw without time slice information
213
215 if (0 == pVVisManager) return;
216
217 if (myContext.GetDrawLine()) {
218 G4VisAttributes trajectoryLineAttribs(myContext.GetLineColour());
219 trajectoryLineAttribs.SetVisibility(myContext.GetLineVisible());
220 trajectoryLine.SetVisAttributes(&trajectoryLineAttribs);
221
222 pVVisManager->Draw(trajectoryLine);
223 }
224
225 if (myContext.GetDrawAuxPts() && (auxiliaryPoints.size() > 0)) {
226 auxiliaryPoints.SetMarkerType(myContext.GetAuxPtsType());
227 auxiliaryPoints.SetSize(myContext.GetAuxPtsSizeType(), myContext.GetAuxPtsSize());
228 auxiliaryPoints.SetFillStyle(myContext.GetAuxPtsFillStyle());
229
230 G4VisAttributes auxiliaryPointsAttribs(myContext.GetAuxPtsColour());
231 auxiliaryPointsAttribs.SetVisibility(myContext.GetAuxPtsVisible());
232 auxiliaryPoints.SetVisAttributes(&auxiliaryPointsAttribs);
233
234 pVVisManager->Draw(auxiliaryPoints);
235 }
236
237 if (myContext.GetDrawStepPts() && (stepPoints.size() > 0)) {
238 stepPoints.SetMarkerType(myContext.GetStepPtsType());
239 stepPoints.SetSize(myContext.GetStepPtsSizeType(), myContext.GetStepPtsSize());
240 stepPoints.SetFillStyle(myContext.GetStepPtsFillStyle());
241
242 G4VisAttributes stepPointsAttribs(myContext.GetStepPtsColour());
243 stepPointsAttribs.SetVisibility(myContext.GetStepPtsVisible());
244 stepPoints.SetVisAttributes(&stepPointsAttribs);
245
246 pVVisManager->Draw(stepPoints);
247 }
248 }
249
250 static void DrawWithTime(const G4VisTrajContext& myContext,
251 G4Polyline& trajectoryLine,
252 G4Polymarker& auxiliaryPoints,
253 G4Polymarker& stepPoints,
254 std::vector<G4double>& trajectoryLineTimes,
255 std::vector<G4double>& auxiliaryPointTimes,
256 std::vector<G4double>& stepPointTimes)
257 {
258 // Draw with time slice information
259
261 if (0 == pVVisManager) return;
262
263 if (myContext.GetDrawLine()) {
264 G4VisAttributes trajectoryLineAttribs(myContext.GetLineColour());
265 trajectoryLineAttribs.SetVisibility(myContext.GetLineVisible());
266
267 for (size_t i = 1; i < trajectoryLine.size(); ++i ) {
268 G4Polyline slice;
269 slice.push_back(trajectoryLine[i -1]);
270 slice.push_back(trajectoryLine[i]);
271 trajectoryLineAttribs.SetStartTime(trajectoryLineTimes[i - 1]);
272 trajectoryLineAttribs.SetEndTime(trajectoryLineTimes[i]);
273 slice.SetVisAttributes(&trajectoryLineAttribs);
274 pVVisManager->Draw(slice);
275 }
276 }
277
278 if (myContext.GetDrawAuxPts() && (auxiliaryPoints.size() > 0)) {
279 G4VisAttributes auxiliaryPointsAttribs(myContext.GetAuxPtsColour());
280 auxiliaryPointsAttribs.SetVisibility(myContext.GetAuxPtsVisible());
281
282 for (size_t i = 0; i < auxiliaryPoints.size(); ++i ) {
283 G4Polymarker point;
284 point.push_back(auxiliaryPoints[i]);
285 point.SetMarkerType(myContext.GetAuxPtsType());
286 point.SetSize(myContext.GetAuxPtsSizeType(), myContext.GetAuxPtsSize());
287 point.SetFillStyle(myContext.GetAuxPtsFillStyle());
288 auxiliaryPointsAttribs.SetStartTime(auxiliaryPointTimes[i]);
289 auxiliaryPointsAttribs.SetEndTime(auxiliaryPointTimes[i]);
290 point.SetVisAttributes(&auxiliaryPointsAttribs);
291 pVVisManager->Draw(point);
292 }
293 }
294
295 if (myContext.GetDrawStepPts() && (stepPoints.size() > 0)) {
296 G4VisAttributes stepPointsAttribs(myContext.GetStepPtsColour());
297 stepPointsAttribs.SetVisibility(myContext.GetStepPtsVisible());
298
299 for (size_t i = 0; i < stepPoints.size(); ++i ) {
300 G4Polymarker point;
301 point.push_back(stepPoints[i]);
302 point.SetMarkerType(myContext.GetStepPtsType());
303 point.SetSize(myContext.GetStepPtsSizeType(), myContext.GetStepPtsSize());
304 point.SetFillStyle(myContext.GetStepPtsFillStyle());
305 stepPointsAttribs.SetStartTime(stepPointTimes[i]);
306 stepPointsAttribs.SetEndTime(stepPointTimes[i]);
307 point.SetVisAttributes(&stepPointsAttribs);
308 pVVisManager->Draw(point);
309 }
310 }
311 }
312
313 void DrawLineAndPoints(const G4VTrajectory& traj, const G4VisTrajContext& context, const G4int& i_mode)
314 {
315 static G4bool warnedAboutIMode = false;
317 ed << "WARNING: DEPRECATED use of i_mode (i_mode: " << i_mode
318 << "). Feature will be removed at a future major release.";
319 if (!warnedAboutIMode) {
321 ("G4TrajectoryDrawerUtils::DrawLineAndPoints(traj, context, i_mode)",
322 "modeling0125", JustWarning, ed);
323 warnedAboutIMode = true;
324 }
325
326 // Extra copy while i_mode is still around
327 G4VisTrajContext myContext(context);
328
329 if (i_mode != 0) {
330 const G4double markerSize = std::abs(i_mode)/1000;
331 G4bool lineRequired (i_mode >= 0);
332 G4bool markersRequired (markerSize > 0.);
333
334 myContext.SetDrawLine(lineRequired);
335 myContext.SetDrawAuxPts(markersRequired);
336 myContext.SetDrawStepPts(markersRequired);
337
338 myContext.SetAuxPtsSize(markerSize);
339 myContext.SetStepPtsSize(markerSize);
340 }
341
342 // Return if don't need to do anything
343 if (!myContext.GetDrawLine() && !myContext.GetDrawAuxPts() && !myContext.GetDrawStepPts()) return;
344
345 // Get points and times (times are returned only if time-slicing
346 // is requested).
347 G4Polyline trajectoryLine;
348 G4Polymarker stepPoints;
349 G4Polymarker auxiliaryPoints;
350 std::vector<G4double> trajectoryLineTimes;
351 std::vector<G4double> stepPointTimes;
352 std::vector<G4double> auxiliaryPointTimes;
353
355 (traj, context,
356 trajectoryLine, auxiliaryPoints, stepPoints,
357 trajectoryLineTimes, auxiliaryPointTimes, stepPointTimes);
358
359 if (validity == ValidTimes) {
360
361 SliceLine(context.GetTimeSliceInterval(),
362 trajectoryLine, trajectoryLineTimes);
363
364 DrawWithTime(context,
365 trajectoryLine, auxiliaryPoints, stepPoints,
366 trajectoryLineTimes, auxiliaryPointTimes, stepPointTimes);
367
368 } else {
369
370 DrawWithoutTime(context, trajectoryLine, auxiliaryPoints, stepPoints);
371
372 }
373 }
374
375 void DrawLineAndPoints(const G4VTrajectory& traj, const G4VisTrajContext& context)
376 {
377 // Return if don't need to do anything
378 if (!context.GetDrawLine() && !context.GetDrawAuxPts() && !context.GetDrawStepPts()) return;
379
380 // Get points and times (times are returned only if time-slicing
381 // is requested).
382 G4Polyline trajectoryLine;
383 G4Polymarker stepPoints;
384 G4Polymarker auxiliaryPoints;
385 std::vector<G4double> trajectoryLineTimes;
386 std::vector<G4double> stepPointTimes;
387 std::vector<G4double> auxiliaryPointTimes;
388
390 (traj, context,
391 trajectoryLine, auxiliaryPoints, stepPoints,
392 trajectoryLineTimes, auxiliaryPointTimes, stepPointTimes);
393
394 if (validity == ValidTimes) {
395
396 SliceLine(context.GetTimeSliceInterval(),
397 trajectoryLine, trajectoryLineTimes);
398
399 DrawWithTime(context,
400 trajectoryLine, auxiliaryPoints, stepPoints,
401 trajectoryLineTimes, auxiliaryPointTimes, stepPointTimes);
402
403 } else {
404
405 DrawWithoutTime(context, trajectoryLine, auxiliaryPoints, stepPoints);
406
407 }
408 }
409}
@ JustWarning
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void SetMarkerType(MarkerType)
static G4double ConvertToDimensionedDouble(const char *st)
Definition: G4UIcommand.cc:429
void SetSize(SizeType, G4double)
Definition: G4VMarker.cc:119
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 int GetPointEntries() const =0
static G4VVisManager * GetConcreteInstance()
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0
void SetStepPtsSize(const G4double &size)
G4bool GetDrawAuxPts() const
G4Colour GetStepPtsColour() const
G4double GetStepPtsSize() const
G4bool GetLineVisible() const
void SetDrawAuxPts(const G4bool &draw)
G4double GetTimeSliceInterval() const
void SetDrawLine(const G4bool &draw)
G4bool GetDrawLine() const
G4VMarker::SizeType GetStepPtsSizeType() const
G4Polymarker::MarkerType GetAuxPtsType() const
void SetDrawStepPts(const G4bool &draw)
G4double GetAuxPtsSize() const
G4VMarker::SizeType GetAuxPtsSizeType() const
G4Colour GetAuxPtsColour() const
G4VMarker::FillStyle GetStepPtsFillStyle() const
void SetAuxPtsSize(const G4double &size)
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:80
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
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 &, const G4int &i_mode)