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

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

update

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