source: trunk/examples/extended/field/BlineTracer/src/G4BlineEventAction.cc @ 1309

Last change on this file since 1309 was 1230, checked in by garnier, 15 years ago

update to geant4.9.3

File size: 5.5 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//
27// $Id: G4BlineEventAction.cc,v 1.2 2006/06/29 17:15:09 gunter Exp $
28// GEANT4 tag $Name: geant4-09-03-cand-01 $
29//
30//
31// --------------------------------------------------------------------
32//
33// G4BlineEventAction implementation
34//
35// --------------------------------------------------------------------
36// Author: Laurent Desorgher (desorgher@phim.unibe.ch)
37//         Created - 2003-10-06
38// --------------------------------------------------------------------
39
40#include "G4BlineEventAction.hh"
41#include "G4BlineTracer.hh"
42
43#include "G4Event.hh"
44#include "G4Trajectory.hh"
45#include "G4EventManager.hh"
46#include "G4VisManager.hh"
47#include "G4UImanager.hh"
48#include "G4Polyline.hh"
49#include "G4Polymarker.hh"
50
51///////////////////////////////////////////////////////////////////////////
52
53G4BlineEventAction::G4BlineEventAction(G4BlineTracer* aBlineTool)
54{
55  fBlineTool=aBlineTool;
56}
57
58///////////////////////////////////////////////////////////////////////////
59
60G4BlineEventAction::~G4BlineEventAction()
61{
62  for (size_t i=0; i<TrajectoryVisAttributes.size(); i++)
63    delete TrajectoryVisAttributes[i];
64}
65
66///////////////////////////////////////////////////////////////////////////
67
68void G4BlineEventAction::BeginOfEventAction(const G4Event*)
69{
70}
71
72///////////////////////////////////////////////////////////////////////////
73
74void G4BlineEventAction::EndOfEventAction(const G4Event* evt)
75{
76  G4TrajectoryContainer * trajectoryContainer = evt->GetTrajectoryContainer();
77  G4int n_trajectories = 0;
78  if(trajectoryContainer) 
79  {
80    n_trajectories = trajectoryContainer->entries(); 
81   
82    // visualisation
83    // -------------
84 
85    if (DrawBline || DrawPoints)
86    {
87      G4int n_point = (*(evt->GetTrajectoryContainer()))[0]->GetPointEntries();
88       
89      G4Polyline pPolyline;
90      G4Polymarker stepPoints;
91      TrajectoryVisAttributes.push_back(new G4VisAttributes(DrawColour));
92      stepPoints.SetMarkerType(G4Polymarker::circles);
93      stepPoints.SetScreenSize(PointSize);
94      stepPoints.SetFillStyle(G4VMarker::filled);
95      stepPoints.SetVisAttributes(TrajectoryVisAttributes.back());
96
97      for(G4int i=0; i<n_point; i++)
98      {
99        G4ThreeVector pos = ((G4TrajectoryPoint*)
100          ((*(evt->GetTrajectoryContainer()))[0]->GetPoint(i)))->GetPosition();
101        if (DrawBline)  pPolyline.push_back( pos);
102        if (DrawPoints)  stepPoints.push_back(pos);
103      }
104
105      pPolyline.SetVisAttributes(TrajectoryVisAttributes.back());
106
107      TrajectoryPolyline.push_back(pPolyline);
108      TrajectoryPoints.push_back(stepPoints);
109    }
110  }
111}
112
113///////////////////////////////////////////////////////////////////////////
114
115void G4BlineEventAction::
116DrawFieldLines( G4double, G4double, G4double )
117{
118  size_t nline = TrajectoryPolyline.size();
119  size_t npoints =TrajectoryPoints.size();
120
121  G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
122  if (!pVVisManager) 
123  {
124    G4Exception("G4BlineEventAction::DrawFieldLines()",
125        "NullPointer", JustWarning,
126        "Missing visualisation driver for visualising magnetic field lines!");
127    return;
128  }
129     
130  if (nline ==0)
131  {
132    G4cout << "WARNING - G4BlineEventAction::DrawFieldLines()" << G4endl
133           << "          There is nothing to visualise !" << G4endl;
134    return;
135  }
136  ((G4VisManager*)pVVisManager)->GetCurrentSceneHandler()-> ClearStore ();
137  G4UImanager::GetUIpointer () -> ApplyCommand ("/vis/drawVolume");
138
139  for (size_t i=0;i<nline;i++)
140    pVVisManager->Draw(TrajectoryPolyline[i]);
141  for (size_t i=0;i<npoints;i++)
142    pVVisManager->Draw(TrajectoryPoints[i]); 
143
144  // ((G4VisManager*)pVVisManager)->GetCurrentViewer()->DrawView();   
145  // ((G4VisManager*)pVVisManager)->GetCurrentViewer()->ShowView();
146}
147
148///////////////////////////////////////////////////////////////////////////
149
150void G4BlineEventAction::ResetVectorObjectToBeDrawn()
151{
152  TrajectoryVisAttributes.clear();
153  TrajectoryPolyline.clear();
154  TrajectoryPoints.clear();
155}
Note: See TracBrowser for help on using the repository browser.