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