Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4BlineEventAction.cc
Go to the documentation of this file.
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 //
28 //
29 //
30 // $Id$
31 //
32 //
33 // --------------------------------------------------------------------
34 //
35 // G4BlineEventAction implementation
36 //
37 // --------------------------------------------------------------------
38 // Author: Laurent Desorgher (desorgher@phim.unibe.ch)
39 // Created - 2003-10-06
40 // --------------------------------------------------------------------
41 
42 #include "G4BlineEventAction.hh"
43 #include "G4BlineTracer.hh"
44 
45 #include "G4Event.hh"
46 #include "G4Trajectory.hh"
47 #include "G4EventManager.hh"
48 #include "G4VisManager.hh"
49 #include "G4UImanager.hh"
50 #include "G4Polyline.hh"
51 #include "G4Polymarker.hh"
52 
54 
56 {
57  fBlineTool=aBlineTool;
58 }
59 
61 
63 {
64  for (size_t i=0; i<fTrajectoryVisAttributes.size(); i++)
65  delete fTrajectoryVisAttributes[i];
66 }
67 
69 
71 {
72 }
73 
75 
77 {
78  G4TrajectoryContainer * trajectoryContainer = evt->GetTrajectoryContainer();
79  if(trajectoryContainer)
80  {
81  // visualisation
82  // -------------
83 
84  if (fDrawBline || fDrawPoints)
85  {
86  G4int n_point = (*(evt->GetTrajectoryContainer()))[0]->GetPointEntries();
87 
88  G4Polyline pPolyline;
89  G4Polymarker stepPoints;
90  fTrajectoryVisAttributes.push_back(new G4VisAttributes(fDrawColour));
92  stepPoints.SetScreenSize(fPointSize);
93  stepPoints.SetFillStyle(G4VMarker::filled);
94  stepPoints.SetVisAttributes(fTrajectoryVisAttributes.back());
95 
96  for(G4int i=0; i<n_point; i++)
97  {
99  ((*(evt->GetTrajectoryContainer()))[0]->GetPoint(i)))->GetPosition();
100  if (fDrawBline) pPolyline.push_back( pos);
101  if (fDrawPoints) stepPoints.push_back(pos);
102  }
103 
104  pPolyline.SetVisAttributes(fTrajectoryVisAttributes.back());
105 
106  fTrajectoryPolyline.push_back(pPolyline);
107  fTrajectoryPoints.push_back(stepPoints);
108  }
109  }
110 }
111 
113 
116 {
117  size_t nline = fTrajectoryPolyline.size();
118  size_t npoints =fTrajectoryPoints.size();
119 
121  if (!pVVisManager)
122  {
123  G4Exception("G4BlineEventAction::DrawFieldLines()",
124  "NullPointer", JustWarning,
125  "Missing visualisation driver for visualising magnetic field lines!");
126  return;
127  }
128 
129  if (nline ==0)
130  {
131  G4cout << "WARNING - G4BlineEventAction::DrawFieldLines()" << G4endl
132  << " There is nothing to visualise !" << G4endl;
133  return;
134  }
135  ((G4VisManager*)pVVisManager)->GetCurrentSceneHandler()-> ClearStore ();
136  G4UImanager::GetUIpointer () -> ApplyCommand ("/vis/drawVolume");
137 
138  for (size_t i=0;i<nline;i++)
139  pVVisManager->Draw(fTrajectoryPolyline[i]);
140  for (size_t i=0;i<npoints;i++)
141  pVVisManager->Draw(fTrajectoryPoints[i]);
142 
143  // ((G4VisManager*)pVVisManager)->GetCurrentViewer()->DrawView();
144  // ((G4VisManager*)pVVisManager)->GetCurrentViewer()->ShowView();
145 }
146 
148 
150 {
151  fTrajectoryVisAttributes.clear();
152  fTrajectoryPolyline.clear();
153  fTrajectoryPoints.clear();
154 }