CMS 3D CMS Logo

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