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

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

bugged version

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