CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EventAction.cc
Go to the documentation of this file.
7 
9 
10 #include <fstream>
11 #include "Randomize.hh"
12 
13 //using std::cout;
14 //using std::endl;
15 
18  SimTrackManager* iManager)
19  : m_runInterface(rm),
20  m_trackManager(iManager),
21  m_stopFile(p.getParameter<std::string>("StopFile")),
22  m_printRandom(p.getParameter<bool>("PrintRandomSeed")),
23  m_debug(p.getUntrackedParameter<bool>("debug",false))
24 {
25  m_trackManager->setCollapsePrimaryVertices(p.getParameter<bool>("CollapsePrimaryVertices"));
26 }
27 
29 
30 void EventAction::BeginOfEventAction(const G4Event * anEvent)
31 {
32  if (std::ifstream(m_stopFile.c_str()))
33  {
34  edm::LogWarning("SimG4CoreApplication")
35  << "BeginOfEventAction: termination signal received at event "
36  << anEvent->GetEventID();
37  /*
38  cout << "BeginOfEventAction: termination signal received at event "
39  << anEvent->GetEventID() << endl;
40  */
41  m_runInterface->abortRun(true);
42  }
44  BeginOfEvent e(anEvent);
46 }
47 
48 void EventAction::EndOfEventAction(const G4Event * anEvent)
49 {
50  if(m_printRandom)
51  {
52  edm::LogInfo("SimG4CoreApplication") << " Event " << anEvent->GetEventID()
53  << " Random number: " << G4UniformRand();
54  //std::cout << " Event " << anEvent->GetEventID()
55  // << " Random number: " << G4UniformRand() << std::endl;
56  //CLHEP::HepRandom::showEngineStatus();
57  }
58  if (std::ifstream(m_stopFile.c_str()))
59  {
60  edm::LogWarning("SimG4CoreApplication")
61  << "EndOfEventAction: termination signal received at event "
62  << anEvent->GetEventID();
63  // soft abort run
64  m_runInterface->abortRun(true);
65  }
66  if (anEvent->GetNumberOfPrimaryVertex()==0)
67  {
68  edm::LogWarning("SimG4CoreApplication")
69  << "EndOfEventAction: event " << anEvent->GetEventID()
70  << " must have failed (no G4PrimaryVertices found) and will be skipped ";
71  return;
72  }
73 
75 
76  // dispatch now end of event, and only then delete tracks...
77  EndOfEvent e(anEvent);
79 
82 }
83 
84 void EventAction::addTrack(TrackWithHistory* iTrack, bool inHistory,
85  bool withAncestor)
86 {
87  m_trackManager->addTrack(iTrack, inHistory, withAncestor);
88 }
89 
90 void EventAction::addTkCaloStateInfo(uint32_t t,const std::pair< math::XYZVectorD,
92 {
94 }
95 
97 {
99 }
void addTrack(TrackWithHistory *iTrack, bool inHistory, bool withAncestor)
T getParameter(std::string const &) const
std::string m_stopFile
Definition: EventAction.h:52
tuple t
Definition: tree.py:139
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
SimRunInterface * m_runInterface
Definition: EventAction.h:50
EventAction(const edm::ParameterSet &ps, SimRunInterface *, SimTrackManager *)
Definition: EventAction.cc:16
SimActivityRegistry::EndOfEventSignal m_endOfEventSignal
Definition: EventAction.h:46
void EndOfEventAction(const G4Event *evt)
Definition: EventAction.cc:48
bool m_printRandom
Definition: EventAction.h:53
G4SimEvent * simEvent()
void cleanTkCaloStateInfoMap()
void addTkCaloStateInfo(uint32_t t, const std::pair< math::XYZVectorD, math::XYZTLorentzVectorD > &p)
void BeginOfEventAction(const G4Event *evt)
Definition: EventAction.cc:30
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
void abortRun(bool softAbort)
string rm
Definition: submit.py:76
void storeTracks(G4SimEvent *simEvent)
SimTrackManager * m_trackManager
Definition: EventAction.h:51
void addTkCaloStateInfo(uint32_t t, const std::pair< math::XYZVectorD, math::XYZTLorentzVectorD > &p)
Definition: EventAction.cc:90
void abortEvent()
Definition: EventAction.cc:96
volatile std::atomic< bool > shutdown_flag false
void addTrack(TrackWithHistory *iTrack, bool inHistory, bool withAncestor)
Definition: EventAction.cc:84
void setCollapsePrimaryVertices(bool iSet)
SimActivityRegistry::BeginOfEventSignal m_beginOfEventSignal
Definition: EventAction.h:45