Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/field/BlineTracer/src/G4BlineEventAction.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 
 43 #include "G4BlineTracer.hh"
 44 #include "G4Event.hh"
 45 #include "G4EventManager.hh"
 46 #include "G4Polyline.hh"
 47 #include "G4Polymarker.hh"
 48 #include "G4Trajectory.hh"
 49 #include "G4UImanager.hh"
 50 #include "G4VisManager.hh"
 51 
 52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 53 
 54 G4BlineEventAction::G4BlineEventAction(G4BlineTracer* aBlineTool)
 55 {
 56   fBlineTool = aBlineTool;
 57 }
 58 
 59 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 60 
 61 G4BlineEventAction::~G4BlineEventAction()
 62 {
 63   for (size_t i = 0; i < fTrajectoryVisAttributes.size(); i++)
 64     delete fTrajectoryVisAttributes[i];
 65 }
 66 
 67 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 68 
 69 void G4BlineEventAction::BeginOfEventAction(const G4Event*) {}
 70 
 71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 72 
 73 void G4BlineEventAction::EndOfEventAction(const G4Event* evt)
 74 {
 75   G4TrajectoryContainer* trajectoryContainer = evt->GetTrajectoryContainer();
 76   if (trajectoryContainer) {
 77     // visualisation
 78     // -------------
 79 
 80     if (fDrawBline || fDrawPoints) {
 81       G4int n_point = (*(evt->GetTrajectoryContainer()))[0]->GetPointEntries();
 82 
 83       G4Polyline pPolyline;
 84       G4Polymarker stepPoints;
 85       fTrajectoryVisAttributes.push_back(new G4VisAttributes(fDrawColour));
 86       stepPoints.SetMarkerType(G4Polymarker::circles);
 87       stepPoints.SetScreenSize(fPointSize);
 88       stepPoints.SetFillStyle(G4VMarker::filled);
 89       stepPoints.SetVisAttributes(fTrajectoryVisAttributes.back());
 90 
 91       for (G4int i = 0; i < n_point; i++) {
 92         G4ThreeVector pos =
 93           ((G4TrajectoryPoint*)((*(evt->GetTrajectoryContainer()))[0]->GetPoint(i)))->GetPosition();
 94         if (fDrawBline) pPolyline.push_back(pos);
 95         if (fDrawPoints) stepPoints.push_back(pos);
 96       }
 97 
 98       pPolyline.SetVisAttributes(fTrajectoryVisAttributes.back());
 99 
100       fTrajectoryPolyline.push_back(pPolyline);
101       fTrajectoryPoints.push_back(stepPoints);
102     }
103   }
104 }
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107 
108 void G4BlineEventAction::DrawFieldLines(G4double, G4double, G4double)
109 {
110   size_t nline = fTrajectoryPolyline.size();
111   size_t npoints = fTrajectoryPoints.size();
112 
113   G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
114   if (!pVVisManager) {
115     G4Exception("G4BlineEventAction::DrawFieldLines()", "NullPointer", JustWarning,
116                 "Missing visualisation driver for visualising magnetic field lines!");
117     return;
118   }
119 
120   if (nline == 0) {
121     G4cout << "WARNING - G4BlineEventAction::DrawFieldLines()" << G4endl
122            << "          There is nothing to visualise !" << G4endl;
123     return;
124   }
125   ((G4VisManager*)pVVisManager)->GetCurrentSceneHandler()->ClearStore();
126   G4UImanager::GetUIpointer()->ApplyCommand("/vis/drawVolume");
127 
128   for (size_t i = 0; i < nline; i++)
129     pVVisManager->Draw(fTrajectoryPolyline[i]);
130   for (size_t i = 0; i < npoints; i++)
131     pVVisManager->Draw(fTrajectoryPoints[i]);
132 
133   // ((G4VisManager*)pVVisManager)->GetCurrentViewer()->DrawView();
134   // ((G4VisManager*)pVVisManager)->GetCurrentViewer()->ShowView();
135 }
136 
137 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
138 
139 void G4BlineEventAction::ResetVectorObjectToBeDrawn()
140 {
141   fTrajectoryVisAttributes.clear();
142   fTrajectoryPolyline.clear();
143   fTrajectoryPoints.clear();
144 }
145