Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-02-23 09:21:15

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 /// \file field/BlineTracer/src/G4BlineEventAction.cc
0027 /// \brief Implementation of the G4BlineEventAction class
0028 //
0029 //
0030 //
0031 //
0032 // --------------------------------------------------------------------
0033 //
0034 // G4BlineEventAction implementation
0035 //
0036 // --------------------------------------------------------------------
0037 // Author: Laurent Desorgher (desorgher@phim.unibe.ch)
0038 //         Created - 2003-10-06
0039 // --------------------------------------------------------------------
0040 
0041 #include "G4BlineEventAction.hh"
0042 
0043 #include "G4BlineTracer.hh"
0044 #include "G4Event.hh"
0045 #include "G4EventManager.hh"
0046 #include "G4Polyline.hh"
0047 #include "G4Polymarker.hh"
0048 #include "G4Trajectory.hh"
0049 #include "G4UImanager.hh"
0050 #include "G4VisManager.hh"
0051 
0052 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0053 
0054 G4BlineEventAction::G4BlineEventAction(G4BlineTracer* aBlineTool)
0055 {
0056   fBlineTool = aBlineTool;
0057 }
0058 
0059 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0060 
0061 G4BlineEventAction::~G4BlineEventAction()
0062 {
0063   for (size_t i = 0; i < fTrajectoryVisAttributes.size(); i++)
0064     delete fTrajectoryVisAttributes[i];
0065 }
0066 
0067 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0068 
0069 void G4BlineEventAction::BeginOfEventAction(const G4Event*) {}
0070 
0071 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0072 
0073 void G4BlineEventAction::EndOfEventAction(const G4Event* evt)
0074 {
0075   G4TrajectoryContainer* trajectoryContainer = evt->GetTrajectoryContainer();
0076   if (trajectoryContainer) {
0077     // visualisation
0078     // -------------
0079 
0080     if (fDrawBline || fDrawPoints) {
0081       G4int n_point = (*(evt->GetTrajectoryContainer()))[0]->GetPointEntries();
0082 
0083       G4Polyline pPolyline;
0084       G4Polymarker stepPoints;
0085       fTrajectoryVisAttributes.push_back(new G4VisAttributes(fDrawColour));
0086       stepPoints.SetMarkerType(G4Polymarker::circles);
0087       stepPoints.SetScreenSize(fPointSize);
0088       stepPoints.SetFillStyle(G4VMarker::filled);
0089       stepPoints.SetVisAttributes(fTrajectoryVisAttributes.back());
0090 
0091       for (G4int i = 0; i < n_point; i++) {
0092         G4ThreeVector pos =
0093           ((G4TrajectoryPoint*)((*(evt->GetTrajectoryContainer()))[0]->GetPoint(i)))->GetPosition();
0094         if (fDrawBline) pPolyline.push_back(pos);
0095         if (fDrawPoints) stepPoints.push_back(pos);
0096       }
0097 
0098       pPolyline.SetVisAttributes(fTrajectoryVisAttributes.back());
0099 
0100       fTrajectoryPolyline.push_back(pPolyline);
0101       fTrajectoryPoints.push_back(stepPoints);
0102     }
0103   }
0104 }
0105 
0106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0107 
0108 void G4BlineEventAction::DrawFieldLines(G4double, G4double, G4double)
0109 {
0110   size_t nline = fTrajectoryPolyline.size();
0111   size_t npoints = fTrajectoryPoints.size();
0112 
0113   G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
0114   if (!pVVisManager) {
0115     G4Exception("G4BlineEventAction::DrawFieldLines()", "NullPointer", JustWarning,
0116                 "Missing visualisation driver for visualising magnetic field lines!");
0117     return;
0118   }
0119 
0120   if (nline == 0) {
0121     G4cout << "WARNING - G4BlineEventAction::DrawFieldLines()" << G4endl
0122            << "          There is nothing to visualise !" << G4endl;
0123     return;
0124   }
0125   ((G4VisManager*)pVVisManager)->GetCurrentSceneHandler()->ClearStore();
0126   G4UImanager::GetUIpointer()->ApplyCommand("/vis/drawVolume");
0127 
0128   for (size_t i = 0; i < nline; i++)
0129     pVVisManager->Draw(fTrajectoryPolyline[i]);
0130   for (size_t i = 0; i < npoints; i++)
0131     pVVisManager->Draw(fTrajectoryPoints[i]);
0132 
0133   // ((G4VisManager*)pVVisManager)->GetCurrentViewer()->DrawView();
0134   // ((G4VisManager*)pVVisManager)->GetCurrentViewer()->ShowView();
0135 }
0136 
0137 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0138 
0139 void G4BlineEventAction::ResetVectorObjectToBeDrawn()
0140 {
0141   fTrajectoryVisAttributes.clear();
0142   fTrajectoryPolyline.clear();
0143   fTrajectoryPoints.clear();
0144 }