test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalTBH4Trigger.h
Go to the documentation of this file.
1 #ifndef HelpfulWatchers_EcalTBH4Trigger_h
2 #define HelpfulWatchers_EcalTBH4Trigger_h
3 // -*- C++ -*-
4 //
5 // Package: HelpfulWatchers
6 // Class : EcalTBH4Trigger
7 //
16 //
17 
18 // system include files
19 #include <iostream>
20 
21 // user include files
25 
27 
28 #include "G4Step.hh"
29 #include "G4VProcess.hh"
30 #include "G4VTouchable.hh"
31 
32 #include "CLHEP/Units/SystemOfUnits.h"
33 
34 // forward declarations
35 class DDDWorld;
36 class BeginOfJob;
37 class BeginOfRun;
38 class BeginOfEvent;
39 class BeginOfTrack;
40 
41 class EndOfRun;
42 class EndOfEvent;
43 class EndOfTrack;
44 
45 #define OBSERVES(type) public Observer<const type*>
46 #define UPDATEH4(type) void update(const type*) { }
47 class EcalTBH4Trigger : public SimWatcher,
48  //OBSERVES(DDDWorld),
49  //OBSERVES(BeginOfJob),
50  //OBSERVES(BeginOfRun),
51  OBSERVES(BeginOfEvent),
52  //OBSERVES(BeginOfTrack),
53  OBSERVES(G4Step),
54  //OBSERVES(EndOfRun),
55  OBSERVES(EndOfEvent)
56  //,
57  //OBSERVES(EndOfTrack)
58 {
59 
60  public:
62  m_verbose(pSet.getUntrackedParameter<bool>("verbose",false)), nTriggeredEvents_(0), trigEvents_(pSet.getUntrackedParameter<int>("trigEvents",-1)) {
63  }
64  //virtual ~EcalTBH4Trigger();
65 
66  // ---------- const member functions ---------------------
67 
68  // ---------- static member functions --------------------
69 
70  // ---------- member functions ---------------------------
71  // UPDATEH4(DDDWorld)
72  //UPDATEH4(BeginOfJob)
73  //UPDATEH4(BeginOfRun)
74  // UPDATEH4(BeginOfEvent)
75  void update(const BeginOfEvent* anEvent)
76  {
77  // std::cout <<"++ signal BeginOfEvent " ;
78  // m_enteringTBH4BeamLine=false;
79  // m_exitingTBH4BeamLine=false;
80  // m_passedTrigger=false;
81  m_passedTrg1=false;
82  m_passedTrg3=false;
83  m_passedTrg4=false;
84  m_passedTrg5=false;
85  m_passedTrg6=false;
86  }
87 
88  // UPDATEH4(BeginOfTrack)
89  void update(const G4Step* iStep)
90  {
92  throw SimG4Exception("Number of needed trigger events reached in ECALTBH4");
93 
94  const G4StepPoint* pre = iStep->GetPreStepPoint();
95  const G4StepPoint* post = iStep->GetPostStepPoint();
96  if(m_verbose) {
97  std::cout <<"++ signal G4Step" ;
98  const G4VTouchable* touch = iStep->GetPreStepPoint()->GetTouchable();
99  //Get name and copy numbers
100  if (touch->GetHistoryDepth() > 0) {
101  for (int ii = 0; ii <= touch->GetHistoryDepth() ; ii++) {
102  std::cout << "EcalTBH4::Level " << ii
103  << ": " << touch->GetVolume(ii)->GetName() << "["
104  << touch->GetReplicaNumber(ii) << "]";
105  }
106  }
107  std::cout <<std::endl;
108  const G4Track* theTrack = iStep->GetTrack();
109  const G4ThreeVector pos = post->GetPosition();
110  std::cout << "( "<<pos.x()<<","<<pos.y()<<","<<pos.z()<<") ";
111  std::cout << " released energy (MeV) "
112  << iStep->GetTotalEnergyDeposit()/CLHEP::MeV ;
113  if (theTrack)
114  {
115  const G4ThreeVector mom = theTrack->GetMomentum();
116  std::cout << " track length (cm) " << theTrack->GetTrackLength()/CLHEP::cm
117  << " particle type " << theTrack->GetDefinition()->GetParticleName()
118  << " momentum " << "( "<<mom.x()<<","<<mom.y()<<","<<mom.z()<<") ";
119  if (theTrack->GetCreatorProcess()) {
120  std::cout << " created by "
121  << theTrack->GetCreatorProcess()->GetProcessName();
122  }
123  }
124  if(post->GetPhysicalVolume()) {
125  std::cout << " " << pre->GetPhysicalVolume()->GetName()
126  << "->" << post->GetPhysicalVolume()->GetName();
127  }
128  std::cout <<std::endl;
129  }
130 
131  if (post && post->GetPhysicalVolume())
132  {
133 
134  if (!m_passedTrg1 && post->GetPhysicalVolume()->GetName() == "TRG1")
135  m_passedTrg1 = true;
136  if (!m_passedTrg3 && post->GetPhysicalVolume()->GetName() == "TRG3")
137  m_passedTrg3 = true;
138  if (!m_passedTrg4 && post->GetPhysicalVolume()->GetName() == "TRG4")
139  m_passedTrg4 = true;
140  if (!m_passedTrg5 && post->GetPhysicalVolume()->GetName() == "TRG5")
141  m_passedTrg5 = true;
142  if (!m_passedTrg6 && post->GetPhysicalVolume()->GetName() == "TRG6")
143  m_passedTrg6 = true;
144  if (post->GetPhysicalVolume()->GetName() == "CMSSE") //Exiting TBH4BeamLine
145  if (! (m_passedTrg1 && m_passedTrg6) ) // Trigger defined as Trg4 && Trg6
146  throw SimG4Exception("Event is not triggered by ECALTBH4");
147  }
148 
149 /* if (!m_enteringTBH4BeamLine && ( post->GetPhysicalVolume()->GetName() == */
150 
151 }
152 //UPDATEH4(G4Step)
153 //UPDATEH4(EndOfRun)
154 //UPDATEH4(EndOfEvent)
155  void update(const EndOfEvent* anEvent)
156  {
157  // std::cout <<"++ signal BeginOfEvent " ;
158  // m_enteringTBH4BeamLine=false;
159  // m_exitingTBH4BeamLine=false;
160  // m_passedTrigger=false;
162  }
163 //UPDATEH4(EndOfTrack)
164 
165  private:
166  //EcalTBH4Trigger(const EcalTBH4Trigger&); // stop default
167 
168  //const EcalTBH4Trigger& operator=(const EcalTBH4Trigger&); // stop default
169 
170  // ---------- member data --------------------------------
171  bool m_verbose;
172 // bool m_enteringTBH4BeamLine;
173 // bool m_exitingTBH4BeamLine;
181  // bool m_passedTrigger;
182 };
183 
184 
185 #endif
int ii
Definition: cuy.py:588
void update(const BeginOfEvent *anEvent)
This routine will be called when the appropriate signal arrives.
const double MeV
void update(const EndOfEvent *anEvent)
This routine will be called when the appropriate signal arrives.
void update(const G4Step *iStep)
This routine will be called when the appropriate signal arrives.
EcalTBH4Trigger(const edm::ParameterSet &pSet)
#define OBSERVES(type)
tuple cout
Definition: gather_cfg.py:121
volatile std::atomic< bool > shutdown_flag false