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

Last change on this file since 1108 was 959, checked in by garnier, 17 years ago

after maj on CVS good-HEAD

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