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

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

import all except CVS

File size: 15.0 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.11 2007/05/30 08:56: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
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    G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
244    if (0 == pVVisManager) return;
245   
246    if (myContext.GetDrawLine()) {
247      G4VisAttributes trajectoryLineAttribs(myContext.GetLineColour());
248      trajectoryLineAttribs.SetVisibility(myContext.GetLineVisible());
249      trajectoryLine.SetVisAttributes(&trajectoryLineAttribs);
250
251      pVVisManager->Draw(trajectoryLine);
252    }
253 
254    if (myContext.GetDrawAuxPts() && (auxiliaryPoints.size() > 0)) {
255      auxiliaryPoints.SetMarkerType(myContext.GetAuxPtsType());
256      auxiliaryPoints.SetScreenSize(myContext.GetAuxPtsSize());
257      auxiliaryPoints.SetFillStyle(myContext.GetAuxPtsFillStyle());
258
259      G4VisAttributes auxiliaryPointsAttribs(myContext.GetAuxPtsColour()); 
260      auxiliaryPointsAttribs.SetVisibility(myContext.GetAuxPtsVisible());
261      auxiliaryPoints.SetVisAttributes(&auxiliaryPointsAttribs);
262
263      pVVisManager->Draw(auxiliaryPoints);
264    }
265 
266    if (myContext.GetDrawStepPts() && (stepPoints.size() > 0)) {
267      stepPoints.SetMarkerType(myContext.GetStepPtsType());
268      stepPoints.SetScreenSize(myContext.GetStepPtsSize());
269      stepPoints.SetFillStyle(myContext.GetStepPtsFillStyle());
270
271      G4VisAttributes stepPointsAttribs(myContext.GetStepPtsColour()); 
272      stepPointsAttribs.SetVisibility(myContext.GetStepPtsVisible());
273      stepPoints.SetVisAttributes(&stepPointsAttribs);
274
275      pVVisManager->Draw(stepPoints);
276    }
277  }
278
279  static void DrawWithTime(const G4VisTrajContext& myContext,
280                           G4Polyline& trajectoryLine,
281                           G4Polymarker& auxiliaryPoints,
282                           G4Polymarker& stepPoints,
283                           std::vector<G4double>& trajectoryLineTimes,
284                           std::vector<G4double>& auxiliaryPointTimes,
285                           std::vector<G4double>& stepPointTimes)
286  {
287    // Draw with time slice information
288
289    G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
290    if (0 == pVVisManager) return;
291
292    if (myContext.GetDrawLine()) {
293      G4VisAttributes trajectoryLineAttribs(myContext.GetLineColour());
294      trajectoryLineAttribs.SetVisibility(myContext.GetLineVisible());
295
296      for (size_t i = 1; i < trajectoryLine.size(); ++i ) {
297        G4Polyline slice;
298        slice.push_back(trajectoryLine[i -1]);
299        slice.push_back(trajectoryLine[i]);
300        trajectoryLineAttribs.SetStartTime(trajectoryLineTimes[i - 1]);
301        trajectoryLineAttribs.SetEndTime(trajectoryLineTimes[i]);
302        slice.SetVisAttributes(&trajectoryLineAttribs);
303        pVVisManager->Draw(slice);
304      }
305    }
306
307    if (myContext.GetDrawAuxPts() && (auxiliaryPoints.size() > 0)) {
308      G4VisAttributes auxiliaryPointsAttribs(myContext.GetAuxPtsColour()); 
309      auxiliaryPointsAttribs.SetVisibility(myContext.GetAuxPtsVisible());
310
311      for (size_t i = 0; i < auxiliaryPoints.size(); ++i ) {
312        G4Polymarker point;
313        point.push_back(auxiliaryPoints[i]);
314        point.SetMarkerType(myContext.GetAuxPtsType());
315        point.SetScreenSize(myContext.GetAuxPtsSize());
316        point.SetFillStyle(myContext.GetAuxPtsFillStyle());
317        auxiliaryPointsAttribs.SetStartTime(auxiliaryPointTimes[i]);
318        auxiliaryPointsAttribs.SetEndTime(auxiliaryPointTimes[i]);
319        point.SetVisAttributes(&auxiliaryPointsAttribs);
320        pVVisManager->Draw(point);
321      }
322    }
323
324    if (myContext.GetDrawStepPts() && (stepPoints.size() > 0)) {
325      G4VisAttributes stepPointsAttribs(myContext.GetStepPtsColour()); 
326      stepPointsAttribs.SetVisibility(myContext.GetStepPtsVisible());
327
328      for (size_t i = 0; i < stepPoints.size(); ++i ) {
329        G4Polymarker point;
330        point.push_back(stepPoints[i]);
331        point.SetMarkerType(myContext.GetStepPtsType());
332        point.SetScreenSize(myContext.GetStepPtsSize());
333        point.SetFillStyle(myContext.GetStepPtsFillStyle());
334        stepPointsAttribs.SetStartTime(stepPointTimes[i]);
335        stepPointsAttribs.SetEndTime(stepPointTimes[i]);
336        point.SetVisAttributes(&stepPointsAttribs);
337        pVVisManager->Draw(point);
338      }
339    }
340  }
341
342  void DrawLineAndPoints(const G4VTrajectory& traj, const G4VisTrajContext& context, const G4int& i_mode) 
343  {
344    // Extra copy while i_mode is still around
345    G4VisTrajContext myContext(context);
346   
347    if (i_mode != 0) {
348      const G4double markerSize = std::abs(i_mode)/1000;
349      G4bool lineRequired (i_mode >= 0);
350      G4bool markersRequired (markerSize > 0.);       
351     
352      myContext.SetDrawLine(lineRequired);
353      myContext.SetDrawAuxPts(markersRequired);
354      myContext.SetDrawStepPts(markersRequired);
355
356      myContext.SetAuxPtsSize(markerSize);
357      myContext.SetStepPtsSize(markerSize);
358
359      static G4bool warnedAboutIMode = false;
360
361      if (!warnedAboutIMode) {
362        G4cout<<"Trajectory drawing configuration will be based on imode value of "<<i_mode<<G4endl;
363        warnedAboutIMode = true;
364      }
365    }
366
367    // Return if don't need to do anything
368    if (!myContext.GetDrawLine() && !myContext.GetDrawAuxPts() && !myContext.GetDrawStepPts()) return;
369   
370    // Get points to draw
371    G4Polyline trajectoryLine;
372    G4Polymarker stepPoints;
373    G4Polymarker auxiliaryPoints;
374   
375    GetPoints(traj, trajectoryLine, auxiliaryPoints, stepPoints);
376   
377    if (myContext.GetTimeSliceInterval()) {
378
379      // Get corresponding track time information, if any
380      std::vector<G4double> trajectoryLineTimes;
381      std::vector<G4double> stepPointTimes;
382      std::vector<G4double> auxiliaryPointTimes;
383 
384      GetTimes(traj, trajectoryLineTimes, auxiliaryPointTimes, stepPointTimes);
385
386      // Check validity
387      if (trajectoryLineTimes.size() != trajectoryLine.size() ||
388          stepPointTimes.size() != stepPoints.size() ||
389          auxiliaryPointTimes.size() != auxiliaryPoints.size()) {
390
391        // Revert to drawing without time information...
392        DrawWithoutTime(myContext, trajectoryLine, auxiliaryPoints, stepPoints);
393      } else {
394
395        SliceLine(myContext.GetTimeSliceInterval(),
396                  trajectoryLine, trajectoryLineTimes);
397
398        DrawWithTime(myContext,
399                     trajectoryLine, auxiliaryPoints, stepPoints,
400                     trajectoryLineTimes, auxiliaryPointTimes, stepPointTimes);
401      }
402
403    } else {
404
405      DrawWithoutTime(myContext, trajectoryLine, auxiliaryPoints, stepPoints);
406
407    }
408  }
409}
Note: See TracBrowser for help on using the repository browser.