source: trunk/examples/advanced/air_shower/src/UltraEventAction.cc @ 1303

Last change on this file since 1303 was 807, checked in by garnier, 16 years ago

update

File size: 5.3 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// --------------------------------------------------------------
28//                 GEANT 4 - ULTRA experiment example
29// --------------------------------------------------------------
30//
31// Code developed by:
32// B. Tome, M.C. Espirito-Santo, A. Trindade, P. Rodrigues
33//
34//    ****************************************************
35//    *      UltraEventAction.cc
36//    ****************************************************
37//
38//    Ultra EventAction class. The UltraAnalysisManager class is used for histogram
39//    filling if the G4ANALYSIS_USE environment variable is set.
40//
41#include "UltraEventAction.hh"
42#include "UltraRunAction.hh"
43#include "UltraPrimaryGeneratorAction.hh"
44#include "UltraOpticalHit.hh"
45
46#include "G4RunManager.hh"
47#include "G4Event.hh"
48#include "G4EventManager.hh"
49#include "G4SDManager.hh"
50#include "G4HCofThisEvent.hh"
51#include "G4VHitsCollection.hh"
52#include "G4GeneralParticleSource.hh"
53
54#include "G4TrajectoryContainer.hh"
55#include "G4Trajectory.hh"
56#include "G4VVisManager.hh"
57
58#ifdef G4ANALYSIS_USE
59#include "UltraAnalysisManager.hh"
60#endif
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
63UltraEventAction::UltraEventAction(UltraRunAction* run)
64  :UltraRun(run),OpticalHitsCollID(-1)
65{;}
66
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
68
69UltraEventAction::~UltraEventAction(){;}
70
71//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
72
73void UltraEventAction::BeginOfEventAction(const G4Event* evt)
74{
75  G4int printModulo = 100;
76
77  evtNb = evt->GetEventID();
78
79  G4SDManager * SDman = G4SDManager::GetSDMpointer(); 
80
81
82  if(OpticalHitsCollID==-1) {
83    OpticalHitsCollID = SDman->GetCollectionID("OpticalHitsCollection");
84  }
85
86
87  if (evtNb%printModulo == 0)
88    G4cout << "\n---> Begin of Event: " << evtNb << G4endl;
89
90}
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92
93void UltraEventAction::EndOfEventAction(const G4Event* evt)
94{
95
96G4HCofThisEvent* HCE = evt->GetHCofThisEvent();
97UltraOpticalHitsCollection* OpticalHitsColl = 0;
98 
99if(HCE){
100
101  if(OpticalHitsCollID != -1) OpticalHitsColl = 
102  (UltraOpticalHitsCollection*)(HCE->GetHC(OpticalHitsCollID));
103
104}
105  G4int nOptHits = 0 ; 
106
107  if(OpticalHitsColl){
108
109    nOptHits = OpticalHitsColl->entries();
110
111#ifdef ULTRA_VERBOSE
112    if (nOptHits > 0){
113     G4cout << " Optical Hit # " << " " << "Energy (eV)" <<  " " << "x,y,z (cm)" << G4endl ;
114    }
115#endif
116
117    for(G4int iHit=0; iHit<nOptHits; iHit++){
118
119      G4double HitEnergy ;
120      HitEnergy = (*OpticalHitsColl)[iHit]->GetEnergy() ;
121      G4ThreeVector HitPos = (*OpticalHitsColl)[iHit]->GetPosition()/cm;
122
123#ifdef G4ANALYSIS_USE
124      UltraAnalysisManager* analysis = UltraAnalysisManager::getInstance();
125      analysis->FillHistogram(1,HitEnergy/eV);
126#endif
127
128    }
129
130  }
131
132#ifdef G4ANALYSIS_USE
133  UltraAnalysisManager* analysis = UltraAnalysisManager::getInstance();
134  analysis->FillHistogram(2,nOptHits);
135#endif
136
137  // Define visualization atributes:
138
139  G4String drawFlag = "all";
140  G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
141
142  if(pVVisManager)
143  {
144   G4TrajectoryContainer* trajectoryContainer = evt->GetTrajectoryContainer();
145   G4int n_trajectories = 0;
146   if (trajectoryContainer) n_trajectories = trajectoryContainer->entries(); 
147   for(G4int i=0; i<n_trajectories; i++) 
148      { G4Trajectory* trj = (G4Trajectory *)((*(evt->GetTrajectoryContainer()))[i]);
149        if (drawFlag == "all") trj->DrawTrajectory();
150        else if ((drawFlag == "charged")&&(trj->GetCharge() != 0.))
151                               trj->DrawTrajectory(); 
152      }
153  }
154
155}
156
157//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
Note: See TracBrowser for help on using the repository browser.