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

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

update

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