Loading...
Searching...
No Matches
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//
26/// \file field/BlineTracer/src/G4BlineEventAction.cc
27/// \brief Implementation of the G4BlineEventAction class
28//
29//
30//
31//
32// --------------------------------------------------------------------
33//
34// G4BlineEventAction implementation
35//
36// --------------------------------------------------------------------
37// Author: Laurent Desorgher (desorgher@phim.unibe.ch)
38// Created - 2003-10-06
39// --------------------------------------------------------------------
40
41#include "G4BlineEventAction.hh"
42#include "G4BlineTracer.hh"
43
44#include "G4Event.hh"
45#include "G4Trajectory.hh"
46#include "G4EventManager.hh"
47#include "G4VisManager.hh"
48#include "G4UImanager.hh"
49#include "G4Polyline.hh"
50#include "G4Polymarker.hh"
51
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53
58
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60
62{
63 for (size_t i=0; i<fTrajectoryVisAttributes.size(); i++)
65}
66
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68
72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
74
76{
77 G4TrajectoryContainer * trajectoryContainer = evt->GetTrajectoryContainer();
78 if(trajectoryContainer)
79 {
80 // visualisation
81 // -------------
82
84 {
85 G4int n_point = (*(evt->GetTrajectoryContainer()))[0]->GetPointEntries();
86
87 G4Polyline pPolyline;
88 G4Polymarker stepPoints;
90 stepPoints.SetMarkerType(G4Polymarker::circles);
91 stepPoints.SetScreenSize(fPointSize);
92 stepPoints.SetFillStyle(G4VMarker::filled);
93 stepPoints.SetVisAttributes(fTrajectoryVisAttributes.back());
94
95 for(G4int i=0; i<n_point; i++)
96 {
97 G4ThreeVector pos = ((G4TrajectoryPoint*)
98 ((*(evt->GetTrajectoryContainer()))[0]->GetPoint(i)))->GetPosition();
99 if (fDrawBline) pPolyline.push_back( pos);
100 if (fDrawPoints) stepPoints.push_back(pos);
101 }
102
103 pPolyline.SetVisAttributes(fTrajectoryVisAttributes.back());
104
105 fTrajectoryPolyline.push_back(pPolyline);
106 fTrajectoryPoints.push_back(stepPoints);
107 }
108 }
109}
110
111//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
112
114DrawFieldLines( G4double, G4double, G4double )
115{
116 size_t nline = fTrajectoryPolyline.size();
117 size_t npoints =fTrajectoryPoints.size();
118
119 G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
120 if (!pVVisManager)
121 {
122 G4Exception("G4BlineEventAction::DrawFieldLines()",
123 "NullPointer", JustWarning,
124 "Missing visualisation driver for visualising magnetic field lines!");
125 return;
126 }
127
128 if (nline ==0)
129 {
130 G4cout << "WARNING - G4BlineEventAction::DrawFieldLines()" << G4endl
131 << " There is nothing to visualise !" << G4endl;
132 return;
133 }
134 ((G4VisManager*)pVVisManager)->GetCurrentSceneHandler()-> ClearStore ();
135 G4UImanager::GetUIpointer () -> ApplyCommand ("/vis/drawVolume");
136
137 for (size_t i=0;i<nline;i++)
138 pVVisManager->Draw(fTrajectoryPolyline[i]);
139 for (size_t i=0;i<npoints;i++)
140 pVVisManager->Draw(fTrajectoryPoints[i]);
141
142 // ((G4VisManager*)pVVisManager)->GetCurrentViewer()->DrawView();
143 // ((G4VisManager*)pVVisManager)->GetCurrentViewer()->ShowView();
144}
145
146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
147
Definition of the G4BlineEventAction class.
Definition of the G4BlineTracer class.
std::vector< G4Polymarker > fTrajectoryPoints
void BeginOfEventAction(const G4Event *) override
G4BlineTracer * fBlineTool
G4BlineEventAction(G4BlineTracer *aBlineTool)
void DrawFieldLines(G4double zoom, G4double theta, G4double phi)
std::vector< G4Polyline > fTrajectoryPolyline
std::vector< G4VisAttributes * > fTrajectoryVisAttributes
void EndOfEventAction(const G4Event *) override

Applications | User Support | Publications | Collaboration