Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
WLSEventAction.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 // $Id: WLSEventAction.cc 70603 2013-06-03 11:23:16Z gcosmo $
27 //
28 /// \file optical/wls/src/WLSEventAction.cc
29 /// \brief Implementation of the WLSEventAction class
30 //
31 //
32 #include "WLSEventAction.hh"
33 
34 #include "WLSRunAction.hh"
35 
37 
38 #include "WLSPhotonDetHit.hh"
39 
40 #include "G4Event.hh"
41 #include "G4EventManager.hh"
42 #include "WLSTrajectory.hh"
43 #include "G4TrajectoryContainer.hh"
44 #include "G4VVisManager.hh"
45 #include "G4SDManager.hh"
46 
47 #include "Randomize.hh"
48 
49 // Purpose: Invoke visualization at the end
50 // Also can accumulate statistics regarding hits
51 // in the PhotonDet detector
52 
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54 
56  : fRunAction(runaction), fVerboseLevel(0),
57  fPrintModulo(100), fDrawFlag("all")
58 {
59  fMPPCCollID = 0;
60 
61  fEventMessenger = new WLSEventActionMessenger(this);
62 
63  fForceDrawPhotons = false;
64  fForceNoPhotons = false;
65 
66 }
67 
68 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
69 
71 {
72  delete fEventMessenger;
73 }
74 
75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76 
78 {
79  G4int evtNb = evt->GetEventID();
80  if (evtNb%fPrintModulo == 0)
81  G4cout << "\n---> Begin of Event: " << evtNb << G4endl;
82 
83  if(fVerboseLevel>0)
84  G4cout << "<<< Event " << evtNb << " started." << G4endl;
85 }
86 
87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
88 
90 {
92 
93  // Visualization of Trajectory
94  if(pVVisManager)
95  {
96  G4TrajectoryContainer* trajectoryContainer = evt->GetTrajectoryContainer();
97 
98  G4int n_trajectories = 0;
99  if (trajectoryContainer) n_trajectories = trajectoryContainer->entries();
100  G4cout << "n_trajectories: " << n_trajectories << G4endl;
101  if (fDrawFlag == "all") G4cout << "draw all trajectories" << G4endl;
102  if (fDrawFlag == "charged") G4cout << "draw only charged" << G4endl;
103 
104  for(G4int i=0; i<n_trajectories; i++)
105  { WLSTrajectory* trj =
106  (WLSTrajectory *)((*(evt->GetTrajectoryContainer()))[i]);
107  if (fDrawFlag == "all") {
108  G4cout << "Now calling DrawTrajectory" << G4endl;
109  G4cout << "Particle Name: " << trj->GetParticleName() << G4endl;
110  trj->DrawTrajectory();
111  }
112  else if ((fDrawFlag == "charged")&&(trj->GetCharge() != 0.))
113  trj->DrawTrajectory();
114  else if (trj->GetParticleName()=="opticalphoton")
115  {
116  G4cout << "We should be drawing an opticalphoton" << G4endl;
117  trj->SetForceDrawTrajectory(fForceDrawPhotons);
118  trj->SetForceNoDrawTrajectory(fForceNoPhotons);
119  trj->DrawTrajectory();
120  }
121  }
122  }
123 
124  if (fVerboseLevel>0)
125  G4cout << "<<< Event " << evt->GetEventID() << " ended." << G4endl;
126 
127  // Save the Random Engine Status at the of event
128  if (fRunAction->GetRndmFreq() == 2)
129  {
130  G4Random::saveEngineStatus("endOfEvent.rndm");
131  G4int evtNb = evt->GetEventID();
132  if (evtNb%fPrintModulo == 0)
133  {
134  G4cout << "\n---> End of Event: " << evtNb << G4endl;
135  G4Random::showEngineStatus();
136  }
137  }
138 
139  // Get Hits from the detector if any
141  G4String colName = "PhotonDetHitCollection";
142  fMPPCCollID = SDman->GetCollectionID(colName);
143 
144  G4HCofThisEvent* HCE = evt->GetHCofThisEvent();
145  WLSPhotonDetHitsCollection* mppcHC = 0;
146 
147  // Get the hit collections
148  if (HCE)
149  {
150  if (fMPPCCollID>=0) mppcHC =
151  (WLSPhotonDetHitsCollection*)(HCE->GetHC(fMPPCCollID));
152  }
153 
154  // Get hit information about photons that reached the detector in this event
155  if (mppcHC)
156  {
157 // G4int n_hit = mppcHC->entries();
158  }
159 }
160 
161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
162 
164 {
166 }
167 
168 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
169 
171 {
172  fVerboseLevel = level;
173 }
void SetEventVerbose(G4int)
virtual ~WLSEventAction()
G4int GetCollectionID(G4String colName)
Definition: G4SDManager.cc:131
WLSEventAction(WLSRunAction *)
virtual void EndOfEventAction(const G4Event *)
static G4VVisManager * GetConcreteInstance()
virtual G4String GetParticleName() const
Definition of the WLSPhotonDetHit class.
G4EventManager * fpEventManager
int G4int
Definition: G4Types.hh:78
G4TrajectoryContainer * GetTrajectoryContainer() const
Definition: G4Event.hh:178
Definition of the WLSEventActionMessenger class.
G4int GetEventID() const
Definition: G4Event.hh:140
G4GLOB_DLL std::ostream G4cout
Definition of the WLSRunAction class.
virtual G4double GetCharge() const
virtual void BeginOfEventAction(const G4Event *)
void SetForceNoDrawTrajectory(G4bool b)
Definition of the WLSTrajectory class.
G4int GetRndmFreq()
Definition: WLSRunAction.hh:59
void SetForceDrawTrajectory(G4bool b)
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
#define G4endl
Definition: G4ios.hh:61
Definition of the WLSEventAction class.
G4HCofThisEvent * GetHCofThisEvent() const
Definition: G4Event.hh:174
const G4Event * GetConstCurrentEvent()
virtual void DrawTrajectory() const