source: trunk/source/visualization/modeling/src/G4TrajectoryDrawerUtils.cc~ @ 1350

Last change on this file since 1350 was 1350, checked in by garnier, 13 years ago

update to last version 4.9.4

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