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

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

ajout de comments de debug

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