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

Last change on this file since 1255 was 1247, checked in by garnier, 16 years ago

update

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