source: trunk/source/visualization/modeling/src/G4TrajectoryDrawerUtils.cc@ 1277

Last change on this file since 1277 was 1258, checked in by garnier, 16 years ago

cvs update

File size: 17.6 KB
Line 
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: G4TrajectoryDrawerUtils.cc,v 1.14 2010/05/11 11:58:18 allison Exp $
27// GEANT4 tag $Name: $
28//
29// Jane Tinslay, John Allison, Joseph Perl November 2005
30//
31#include "G4TrajectoryDrawerUtils.hh"
32#include "G4Colour.hh"
33#include "G4Polyline.hh"
34#include "G4Polymarker.hh"
35#include "G4VTrajectory.hh"
36#include "G4VTrajectoryPoint.hh"
37#include "G4VisAttributes.hh"
38#include "G4VisTrajContext.hh"
39#include "G4VVisManager.hh"
40#include "G4UIcommand.hh"
41#include "G4AttValue.hh"
42#include <sstream>
43
44namespace G4TrajectoryDrawerUtils {
45
46
47 void GetPoints(const G4VTrajectory& traj, G4Polyline& trajectoryLine,
48 G4Polymarker& auxiliaryPoints, G4Polymarker& stepPoints)
49 {
50 for (G4int i=0; i<traj.GetPointEntries(); i++) {
51 G4VTrajectoryPoint* aTrajectoryPoint = traj.GetPoint(i);
52
53 const std::vector<G4ThreeVector>* auxiliaries
54 = aTrajectoryPoint->GetAuxiliaryPoints();
55
56 if (0 != auxiliaries) {
57 for (size_t iAux=0; iAux<auxiliaries->size(); ++iAux) {
58 const G4ThreeVector pos((*auxiliaries)[iAux]);
59 if (trajectoryLine.size() >0) {
60 if (pos != trajectoryLine[trajectoryLine.size()-1]) {
61 trajectoryLine.push_back(pos);
62 auxiliaryPoints.push_back(pos);
63 }
64 } else {
65 trajectoryLine.push_back(pos);
66 auxiliaryPoints.push_back(pos);
67 }
68 }
69 }
70 const G4ThreeVector pos(aTrajectoryPoint->GetPosition());
71 if (trajectoryLine.size() >0) {
72 if (pos != trajectoryLine[trajectoryLine.size()-1]) {
73 trajectoryLine.push_back(pos);
74 stepPoints.push_back(pos);
75 }
76 } else {
77 trajectoryLine.push_back(pos);
78 stepPoints.push_back(pos);
79 }
80 }
81 }
82
83 /***
84 void DrawLineAndPoints(const G4VTrajectory& traj, const G4int& i_mode, const G4Colour& colour, const G4bool& visible) {
85 // If i_mode>=0, draws a trajectory as a polyline (default is blue for
86 // positive, red for negative, green for neutral) and, if i_mode!=0,
87 // adds markers - yellow circles for step points and magenta squares
88 // for auxiliary points, if any - whose screen size in pixels is
89 // given by std::abs(i_mode)/1000. E.g: i_mode = 5000 gives easily
90 // visible markers.
91
92 static G4bool warnedAboutIMode = false;
93 std::ostringstream oss;
94 oss << "WARNING: DEPRECATED use of i_mode (i_mode: " << i_mode
95 << "). Feature will be removed at a future major release.";
96 if (!warnedAboutIMode) {
97 G4Exception
98 ("G4TrajectoryDrawerUtils::DrawLineAndPoints(traj, i_mode, colour, visible)",
99 "",
100 JustWarning,
101 oss.str().c_str());
102 warnedAboutIMode = true;
103 }
104
105 G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
106 if (0 == pVVisManager) return;
107
108 const G4double markerSize = std::abs(i_mode)/1000;
109 G4bool lineRequired (i_mode >= 0);
110 G4bool markersRequired (markerSize > 0.);
111
112 // Return if don't need to do anything
113 if (!lineRequired && !markersRequired) return;
114
115 // Get points to draw
116 G4Polyline trajectoryLine;
117 G4Polymarker stepPoints;
118 G4Polymarker auxiliaryPoints;
119
120 GetPoints(traj, trajectoryLine, auxiliaryPoints, stepPoints);
121
122 if (lineRequired) {
123 G4VisAttributes trajectoryLineAttribs(colour);
124 trajectoryLineAttribs.SetVisibility(visible);
125 trajectoryLine.SetVisAttributes(&trajectoryLineAttribs);
126 pVVisManager->Draw(trajectoryLine);
127 }
128
129 if (markersRequired) {
130 auxiliaryPoints.SetMarkerType(G4Polymarker::squares);
131 auxiliaryPoints.SetScreenSize(markerSize);
132 auxiliaryPoints.SetFillStyle(G4VMarker::filled);
133 G4VisAttributes auxiliaryPointsAttribs(G4Colour(0.,1.,1.)); // Magenta
134 auxiliaryPointsAttribs.SetVisibility(visible);
135 auxiliaryPoints.SetVisAttributes(&auxiliaryPointsAttribs);
136 pVVisManager->Draw(auxiliaryPoints);
137
138 stepPoints.SetMarkerType(G4Polymarker::circles);
139 stepPoints.SetScreenSize(markerSize);
140 stepPoints.SetFillStyle(G4VMarker::filled);
141 G4VisAttributes stepPointsAttribs(G4Colour(1.,1.,0.)); // Yellow.
142 stepPoints.SetVisAttributes(&stepPointsAttribs);
143 pVVisManager->Draw(stepPoints);
144 }
145
146 }
147 ***/
148
149 static void GetTimes(const G4VTrajectory& traj,
150 std::vector<G4double>& trajectoryLineTimes,
151 std::vector<G4double>& auxiliaryPointTimes,
152 std::vector<G4double>& stepPointTimes)
153 {
154 // It is important check that the sizes of times vectors produced
155 // by this function matches those of points vectors from
156 // GetPoints. If not, assume that the time information is
157 // invalid.
158
159 G4ThreeVector lastTrajectoryPointPosition;
160 for (G4int i=0; i<traj.GetPointEntries(); i++) {
161 G4VTrajectoryPoint* aTrajectoryPoint = traj.GetPoint(i);
162
163 // Pre- and Post-Point times from the trajectory point...
164 G4double trajectoryPointPreTime = -std::numeric_limits<double>::max();
165 G4double trajectoryPointPostTime = std::numeric_limits<double>::max();
166 std::vector<G4AttValue>* trajectoryPointAttValues =
167 aTrajectoryPoint->CreateAttValues();
168 if (!trajectoryPointAttValues) {
169 G4cout << "G4TrajectoryDrawerUtils::GetTimes: no att values."
170 << G4endl;
171 return;
172 } else {
173 G4bool foundPreTime = false, foundPostTime = false;
174 for (std::vector<G4AttValue>::iterator i =
175 trajectoryPointAttValues->begin();
176 i != trajectoryPointAttValues->end(); ++i) {
177 if (i->GetName() == "PreT") {
178 trajectoryPointPreTime =
179 G4UIcommand::ConvertToDimensionedDouble(i->GetValue());
180 foundPreTime = true;
181 }
182 if (i->GetName() == "PostT") {
183 trajectoryPointPostTime =
184 G4UIcommand::ConvertToDimensionedDouble(i->GetValue());
185 foundPostTime = true;
186 }
187 }
188 if (!foundPreTime || !foundPostTime) {
189 static G4bool warnedTimesNotFound = false;
190 if (!warnedTimesNotFound) {
191 G4cout <<
192 "WARNING: G4TrajectoryDrawerUtils::GetTimes: times not found."
193 << G4endl;
194 warnedTimesNotFound = true;
195 }
196 return;
197 }
198 }
199
200 const G4ThreeVector& trajectoryPointPosition =
201 aTrajectoryPoint->GetPosition();
202
203 const std::vector<G4ThreeVector>* auxiliaries
204 = aTrajectoryPoint->GetAuxiliaryPoints();
205
206 if (0 != auxiliaries) {
207 for (size_t iAux=0; iAux<auxiliaries->size(); ++iAux) {
208 // Interpolate time for auxiliary points...
209 const G4ThreeVector pos((*auxiliaries)[iAux]);
210 G4double s1 = (pos - lastTrajectoryPointPosition).mag();
211 G4double s2 = (trajectoryPointPosition - pos).mag();
212 G4double t = trajectoryPointPreTime +
213 (trajectoryPointPostTime - trajectoryPointPreTime) *
214 (s1 / (s1 + s2));
215 trajectoryLineTimes.push_back(t);
216 auxiliaryPointTimes.push_back(t);
217 }
218 }
219
220 trajectoryLineTimes.push_back(trajectoryPointPostTime);
221 stepPointTimes.push_back(trajectoryPointPostTime);
222
223 lastTrajectoryPointPosition = trajectoryPointPosition;
224 }
225 }
226
227 static void SliceLine(G4double timeIncrement,
228 G4Polyline& trajectoryLine,
229 std::vector<G4double>& trajectoryLineTimes)
230 {
231 // Assumes valid arguments from GetPoints and GetTimes.
232
233 G4Polyline newTrajectoryLine;
234 std::vector<G4double> newTrajectoryLineTimes;
235
236 newTrajectoryLine.push_back(trajectoryLine[0]);
237 newTrajectoryLineTimes.push_back(trajectoryLineTimes[0]);
238 size_t lineSize = trajectoryLine.size();
239 if (lineSize > 1) {
240 for (size_t i = 1; i < trajectoryLine.size(); ++i) {
241 G4double deltaT = trajectoryLineTimes[i] - trajectoryLineTimes[i - 1];
242 if (deltaT > 0.) {
243 G4double practicalTimeIncrement =
244 std::max(timeIncrement, deltaT / 100.);
245 for (G4double t =
246 (int(trajectoryLineTimes[i - 1]/practicalTimeIncrement) + 1) *
247 practicalTimeIncrement;
248 t <= trajectoryLineTimes[i];
249 t += practicalTimeIncrement) {
250 G4ThreeVector pos = trajectoryLine[i - 1] +
251 (trajectoryLine[i] - trajectoryLine[i - 1]) *
252 ((t - trajectoryLineTimes[i - 1]) / deltaT);
253 newTrajectoryLine.push_back(pos);
254 newTrajectoryLineTimes.push_back(t);
255 }
256 }
257 newTrajectoryLine.push_back(trajectoryLine[i]);
258 newTrajectoryLineTimes.push_back(trajectoryLineTimes[i]);
259 }
260 }
261
262 trajectoryLine = newTrajectoryLine;
263 trajectoryLineTimes = newTrajectoryLineTimes;
264 }
265
266 static void DrawWithoutTime(const G4VisTrajContext& myContext,
267 G4Polyline& trajectoryLine,
268 G4Polymarker& auxiliaryPoints,
269 G4Polymarker& stepPoints)
270 {
271 // Draw without time slice information
272
273 G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
274 if (0 == pVVisManager) return;
275
276 if (myContext.GetDrawLine()) {
277 G4VisAttributes trajectoryLineAttribs(myContext.GetLineColour());
278 trajectoryLineAttribs.SetVisibility(myContext.GetLineVisible());
279 trajectoryLine.SetVisAttributes(&trajectoryLineAttribs);
280
281 pVVisManager->Draw(trajectoryLine);
282 }
283
284 if (myContext.GetDrawAuxPts() && (auxiliaryPoints.size() > 0)) {
285 auxiliaryPoints.SetMarkerType(myContext.GetAuxPtsType());
286 auxiliaryPoints.SetSize(myContext.GetAuxPtsSizeType(), myContext.GetAuxPtsSize());
287 auxiliaryPoints.SetFillStyle(myContext.GetAuxPtsFillStyle());
288
289 G4VisAttributes auxiliaryPointsAttribs(myContext.GetAuxPtsColour());
290 auxiliaryPointsAttribs.SetVisibility(myContext.GetAuxPtsVisible());
291 auxiliaryPoints.SetVisAttributes(&auxiliaryPointsAttribs);
292
293 pVVisManager->Draw(auxiliaryPoints);
294 }
295
296 if (myContext.GetDrawStepPts() && (stepPoints.size() > 0)) {
297 stepPoints.SetMarkerType(myContext.GetStepPtsType());
298 stepPoints.SetSize(myContext.GetStepPtsSizeType(), myContext.GetStepPtsSize());
299 stepPoints.SetFillStyle(myContext.GetStepPtsFillStyle());
300
301 G4VisAttributes stepPointsAttribs(myContext.GetStepPtsColour());
302 stepPointsAttribs.SetVisibility(myContext.GetStepPtsVisible());
303 stepPoints.SetVisAttributes(&stepPointsAttribs);
304
305 pVVisManager->Draw(stepPoints);
306 }
307 }
308
309 static void DrawWithTime(const G4VisTrajContext& myContext,
310 G4Polyline& trajectoryLine,
311 G4Polymarker& auxiliaryPoints,
312 G4Polymarker& stepPoints,
313 std::vector<G4double>& trajectoryLineTimes,
314 std::vector<G4double>& auxiliaryPointTimes,
315 std::vector<G4double>& stepPointTimes)
316 {
317 // Draw with time slice information
318
319 G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
320 if (0 == pVVisManager) return;
321
322 if (myContext.GetDrawLine()) {
323 G4VisAttributes trajectoryLineAttribs(myContext.GetLineColour());
324 trajectoryLineAttribs.SetVisibility(myContext.GetLineVisible());
325
326 for (size_t i = 1; i < trajectoryLine.size(); ++i ) {
327 G4Polyline slice;
328 slice.push_back(trajectoryLine[i -1]);
329 slice.push_back(trajectoryLine[i]);
330 trajectoryLineAttribs.SetStartTime(trajectoryLineTimes[i - 1]);
331 trajectoryLineAttribs.SetEndTime(trajectoryLineTimes[i]);
332 slice.SetVisAttributes(&trajectoryLineAttribs);
333 pVVisManager->Draw(slice);
334 }
335 }
336
337 if (myContext.GetDrawAuxPts() && (auxiliaryPoints.size() > 0)) {
338 G4VisAttributes auxiliaryPointsAttribs(myContext.GetAuxPtsColour());
339 auxiliaryPointsAttribs.SetVisibility(myContext.GetAuxPtsVisible());
340
341 for (size_t i = 0; i < auxiliaryPoints.size(); ++i ) {
342 G4Polymarker point;
343 point.push_back(auxiliaryPoints[i]);
344 point.SetMarkerType(myContext.GetAuxPtsType());
345 point.SetSize(myContext.GetAuxPtsSizeType(), myContext.GetAuxPtsSize());
346 point.SetFillStyle(myContext.GetAuxPtsFillStyle());
347 auxiliaryPointsAttribs.SetStartTime(auxiliaryPointTimes[i]);
348 auxiliaryPointsAttribs.SetEndTime(auxiliaryPointTimes[i]);
349 point.SetVisAttributes(&auxiliaryPointsAttribs);
350 pVVisManager->Draw(point);
351 }
352 }
353
354 if (myContext.GetDrawStepPts() && (stepPoints.size() > 0)) {
355 G4VisAttributes stepPointsAttribs(myContext.GetStepPtsColour());
356 stepPointsAttribs.SetVisibility(myContext.GetStepPtsVisible());
357
358 for (size_t i = 0; i < stepPoints.size(); ++i ) {
359 G4Polymarker point;
360 point.push_back(stepPoints[i]);
361 point.SetMarkerType(myContext.GetStepPtsType());
362 point.SetSize(myContext.GetStepPtsSizeType(), myContext.GetStepPtsSize());
363 point.SetFillStyle(myContext.GetStepPtsFillStyle());
364 stepPointsAttribs.SetStartTime(stepPointTimes[i]);
365 stepPointsAttribs.SetEndTime(stepPointTimes[i]);
366 point.SetVisAttributes(&stepPointsAttribs);
367 pVVisManager->Draw(point);
368 }
369 }
370 }
371
372 void DrawLineAndPoints(const G4VTrajectory& traj, const G4VisTrajContext& context, const G4int& i_mode)
373 {
374 static G4bool warnedAboutIMode = false;
375 std::ostringstream oss;
376 oss << "WARNING: DEPRECATED use of i_mode (i_mode: " << i_mode
377 << "). Feature will be removed at a future major release.";
378 if (!warnedAboutIMode) {
379 G4Exception
380 ("G4TrajectoryDrawerUtils::DrawLineAndPoints(traj, context, i_mode)",
381 "",
382 JustWarning,
383 oss.str().c_str());
384 warnedAboutIMode = true;
385 }
386
387 // Extra copy while i_mode is still around
388 G4VisTrajContext myContext(context);
389
390 if (i_mode != 0) {
391 const G4double markerSize = std::abs(i_mode)/1000;
392 G4bool lineRequired (i_mode >= 0);
393 G4bool markersRequired (markerSize > 0.);
394
395 myContext.SetDrawLine(lineRequired);
396 myContext.SetDrawAuxPts(markersRequired);
397 myContext.SetDrawStepPts(markersRequired);
398
399 myContext.SetAuxPtsSize(markerSize);
400 myContext.SetStepPtsSize(markerSize);
401 }
402
403 // Return if don't need to do anything
404 if (!myContext.GetDrawLine() && !myContext.GetDrawAuxPts() && !myContext.GetDrawStepPts()) return;
405
406 // Get points to draw
407 G4Polyline trajectoryLine;
408 G4Polymarker stepPoints;
409 G4Polymarker auxiliaryPoints;
410
411 GetPoints(traj, trajectoryLine, auxiliaryPoints, stepPoints);
412
413 if (myContext.GetTimeSliceInterval()) {
414
415 // Get corresponding track time information, if any
416 std::vector<G4double> trajectoryLineTimes;
417 std::vector<G4double> stepPointTimes;
418 std::vector<G4double> auxiliaryPointTimes;
419
420 GetTimes(traj, trajectoryLineTimes, auxiliaryPointTimes, stepPointTimes);
421
422 // Check validity
423 if (trajectoryLineTimes.size() != trajectoryLine.size() ||
424 stepPointTimes.size() != stepPoints.size() ||
425 auxiliaryPointTimes.size() != auxiliaryPoints.size()) {
426
427 // Revert to drawing without time information...
428 DrawWithoutTime(myContext, trajectoryLine, auxiliaryPoints, stepPoints);
429 } else {
430
431 SliceLine(myContext.GetTimeSliceInterval(),
432 trajectoryLine, trajectoryLineTimes);
433
434 DrawWithTime(myContext,
435 trajectoryLine, auxiliaryPoints, stepPoints,
436 trajectoryLineTimes, auxiliaryPointTimes, stepPointTimes);
437 }
438
439 } else {
440
441 DrawWithoutTime(myContext, trajectoryLine, auxiliaryPoints, stepPoints);
442
443 }
444 }
445
446 void DrawLineAndPoints(const G4VTrajectory& traj, const G4VisTrajContext& context)
447 {
448 // Return if don't need to do anything
449 if (!context.GetDrawLine() && !context.GetDrawAuxPts() && !context.GetDrawStepPts()) return;
450
451 // Get points to draw
452 G4Polyline trajectoryLine;
453 G4Polymarker stepPoints;
454 G4Polymarker auxiliaryPoints;
455
456 GetPoints(traj, trajectoryLine, auxiliaryPoints, stepPoints);
457
458 if (context.GetTimeSliceInterval()) {
459
460 // Get corresponding track time information, if any
461 std::vector<G4double> trajectoryLineTimes;
462 std::vector<G4double> stepPointTimes;
463 std::vector<G4double> auxiliaryPointTimes;
464
465 GetTimes(traj, trajectoryLineTimes, auxiliaryPointTimes, stepPointTimes);
466
467 // Check validity
468 if (trajectoryLineTimes.size() != trajectoryLine.size() ||
469 stepPointTimes.size() != stepPoints.size() ||
470 auxiliaryPointTimes.size() != auxiliaryPoints.size()) {
471
472 // Revert to drawing without time information...
473 DrawWithoutTime(context, trajectoryLine, auxiliaryPoints, stepPoints);
474 } else {
475
476 SliceLine(context.GetTimeSliceInterval(),
477 trajectoryLine, trajectoryLineTimes);
478
479 DrawWithTime(context,
480 trajectoryLine, auxiliaryPoints, stepPoints,
481 trajectoryLineTimes, auxiliaryPointTimes, stepPointTimes);
482 }
483
484 } else {
485
486 DrawWithoutTime(context, trajectoryLine, auxiliaryPoints, stepPoints);
487
488 }
489 }
490}
Note: See TracBrowser for help on using the repository browser.