CMS 3D CMS Logo

EcalTBH4Trigger.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HelpfulWatchers
4 // Class : EcalTBH4Trigger
5 //
16 //
17 
18 // system include files
19 #include <sstream>
20 
21 // user include files
28 
33 
34 #include "G4Step.hh"
35 #include "G4VProcess.hh"
36 #include "G4VTouchable.hh"
37 
38 #include <CLHEP/Units/SystemOfUnits.h>
39 
40 class EcalTBH4Trigger : public SimWatcher,
41  public Observer<const BeginOfEvent *>,
42  public Observer<const G4Step *>,
43  public Observer<const EndOfEvent *> {
44 public:
46  : m_verbose(pSet.getUntrackedParameter<bool>("verbose", false)),
48  trigEvents_(pSet.getUntrackedParameter<int>("trigEvents", -1)) {}
49  // virtual ~EcalTBH4Trigger();
50 
51  // ---------- const member functions ---------------------
52 
53  // ---------- static member functions --------------------
54 
55  // ---------- member functions ---------------------------
56  void update(const BeginOfEvent *anEvent) override;
57  void update(const G4Step *iStep) override;
58  void update(const EndOfEvent *anEvent) override;
59 
60 private:
61  // ---------- member data --------------------------------
62  bool m_verbose;
70 };
71 
72 void EcalTBH4Trigger::update(const BeginOfEvent *anEvent) {
73  m_passedTrg1 = false;
74  m_passedTrg3 = false;
75  m_passedTrg4 = false;
76  m_passedTrg5 = false;
77  m_passedTrg6 = false;
78 }
79 
80 void EcalTBH4Trigger::update(const G4Step *iStep) {
82  throw SimG4Exception("Number of needed trigger events reached in ECALTBH4");
83 
84  const G4StepPoint *pre = iStep->GetPreStepPoint();
85  const G4StepPoint *post = iStep->GetPostStepPoint();
86  if (m_verbose) {
87  std::ostringstream st1;
88  st1 << "++ signal G4Step";
89  const G4VTouchable *touch = iStep->GetPreStepPoint()->GetTouchable();
90  // Get name and copy numbers
91  if (touch->GetHistoryDepth() > 0) {
92  for (int ii = 0; ii <= touch->GetHistoryDepth(); ii++) {
93  st1 << "EcalTBH4::Level " << ii << ": " << touch->GetVolume(ii)->GetName() << "[" << touch->GetReplicaNumber(ii)
94  << "]";
95  }
96  }
97  edm::LogVerbatim("EcalTBInfo") << st1.str();
98 
99  std::ostringstream st2;
100  const G4Track *theTrack = iStep->GetTrack();
101  const G4ThreeVector &pos = post->GetPosition();
102  st2 << "( " << pos.x() << "," << pos.y() << "," << pos.z() << ") ";
103  st2 << " released energy (MeV) " << iStep->GetTotalEnergyDeposit() / CLHEP::MeV;
104  if (theTrack) {
105  const G4ThreeVector mom = theTrack->GetMomentum();
106  st2 << " track length (cm) " << theTrack->GetTrackLength() / CLHEP::cm << " particle type "
107  << theTrack->GetDefinition()->GetParticleName() << " momentum "
108  << "( " << mom.x() << "," << mom.y() << "," << mom.z() << ") ";
109  if (theTrack->GetCreatorProcess()) {
110  st2 << " created by " << theTrack->GetCreatorProcess()->GetProcessName();
111  }
112  }
113  if (post->GetPhysicalVolume()) {
114  st2 << " " << pre->GetPhysicalVolume()->GetName() << "->" << post->GetPhysicalVolume()->GetName();
115  }
116  edm::LogVerbatim("EcalTBInfo") << st2.str();
117  }
118 
119  if (post && post->GetPhysicalVolume()) {
120  if (!m_passedTrg1 && post->GetPhysicalVolume()->GetName() == "TRG1")
121  m_passedTrg1 = true;
122  if (!m_passedTrg3 && post->GetPhysicalVolume()->GetName() == "TRG3")
123  m_passedTrg3 = true;
124  if (!m_passedTrg4 && post->GetPhysicalVolume()->GetName() == "TRG4")
125  m_passedTrg4 = true;
126  if (!m_passedTrg5 && post->GetPhysicalVolume()->GetName() == "TRG5")
127  m_passedTrg5 = true;
128  if (!m_passedTrg6 && post->GetPhysicalVolume()->GetName() == "TRG6")
129  m_passedTrg6 = true;
130  if (post->GetPhysicalVolume()->GetName() == "CMSSE") // Exiting TBH4BeamLine
131  if (!(m_passedTrg1 && m_passedTrg6)) // Trigger defined as Trg4 && Trg6
132  throw SimG4Exception("Event is not triggered by ECALTBH4");
133  }
134 
135  /* if (!m_enteringTBH4BeamLine && ( post->GetPhysicalVolume()->GetName()
136  * == */
137 }
138 
139 void EcalTBH4Trigger::update(const EndOfEvent *anEvent) {
140  edm::LogVerbatim("EcalTBInfo") << "++ signal BeginOfEvent ";
142 }
143 
Log< level::Info, true > LogVerbatim
#define DEFINE_SIMWATCHER(type)
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:589