CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TotemTestGem.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Forward
4 // Class : TotemTestGem
5 //
6 // Implementation:
7 // <Notes on implementation>
8 //
9 // Original Author:
10 // Created: Tue May 16 10:14:34 CEST 2006
11 // $Id: TotemTestGem.cc,v 1.3 2009/05/25 06:48:53 fabiocos Exp $
12 //
13 
14 // system include files
15 #include <iostream>
16 #include <iomanip>
17 #include <cmath>
18 
19 // user include files
22 
26 
29 
30 #include "G4SDManager.hh"
31 #include "G4HCofThisEvent.hh"
32 #include "CLHEP/Units/GlobalSystemOfUnits.h"
33 #include "CLHEP/Units/GlobalPhysicalConstants.h"
34 
35 //
36 // constructors and destructor
37 //
38 
40 
41  edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("TotemTestGem");
42  names = m_Anal.getParameter<std::vector<std::string> >("Names");
43 
44  edm::LogInfo("ForwardSim") << "TotemTestGem:: Initialised as observer of "
45  << "begin of job, begin/end events and of G4step";
46 }
47 
49 }
50 
51 //
52 // member functions
53 //
54 
56 
57  std::auto_ptr<TotemTestHistoClass> product(new TotemTestHistoClass);
58  fillEvent(*product);
59  e.put(product);
60 }
61 
63 
64  int iev = (*evt)()->GetEventID();
65  LogDebug("ForwardSim") << "TotemTestGem: Begin of event = " << iev;
66  clear();
67 
68 }
69 
70 void TotemTestGem::update(const EndOfEvent * evt) {
71 
72  evtnum = (*evt)()->GetEventID();
73  LogDebug("ForwardSim") << "TotemTestGem:: Fill event " << evtnum;
74 
75  // access to the G4 hit collections
76  G4HCofThisEvent* allHC = (*evt)()->GetHCofThisEvent();
77 
78  int nhit = 0;
79  for (unsigned int in=0; in<names.size(); in++) {
80  int HCid = G4SDManager::GetSDMpointer()->GetCollectionID(names[in]);
81  TotemG4HitCollection* theHC = (TotemG4HitCollection*) allHC->GetHC(HCid);
82  LogDebug("ForwardSim") << "TotemTestGem :: Hit Collection for " <<names[in]
83  << " of ID " << HCid << " is obtained at " << theHC;
84 
85  if (HCid >= 0 && theHC > 0) {
86  int nentries = theHC->entries();
87  LogDebug("ForwardSim") << "TotemTestGem :: " << names[in] << " with "
88  << nentries << " entries";
89  for (int ihit = 0; ihit <nentries; ihit++) {
90  TotemG4Hit* aHit = (*theHC)[ihit];
91  hits.push_back(aHit);
92  }
93  nhit += nentries;
94  }
95  }
96 
97  // Writing the data to the Tree
98  LogDebug("ForwardSim") << "TotemTestGem:: --- after fillTree with " << nhit
99  << " Hits";
100 }
101 
103 
104  product.setEVT(evtnum);
105 
106  for (unsigned ihit = 0; ihit < hits.size(); ihit++) {
107  TotemG4Hit* aHit = hits[ihit];
108  int UID = aHit->getUnitID();
109  int Ptype = aHit->getParticleType();
110  int TID = aHit->getTrackID();
111  int PID = aHit->getParentId();
112  float ELoss = aHit->getEnergyLoss();
113  float PABS = aHit->getPabs();
114  float x = aHit->getEntry().x();
115  float y = aHit->getEntry().y();
116  float z = aHit->getEntry().z();
117  float vx = aHit->getVx();
118  float vy = aHit->getVy();
119  float vz = aHit->getVz();
120  product.fillHit(UID, Ptype, TID, PID, ELoss, PABS, vx, vy, vz, x, y,z);
121  }
122 }
123 
125 
126  evtnum = 0;
127  hits.clear();
128 }
#define LogDebug(id)
T getParameter(std::string const &) const
void fillEvent(TotemTestHistoClass &)
float getEnergyLoss() const
Definition: TotemG4Hit.cc:151
virtual void produce(edm::Event &, const edm::EventSetup &)
Definition: TotemTestGem.cc:55
math::XYZPoint getEntry() const
Definition: TotemG4Hit.cc:124
float getPabs() const
Definition: TotemG4Hit.cc:149
float getVy() const
Definition: TotemG4Hit.cc:180
double double double z
void update(const BeginOfEvent *evt)
This routine will be called when the appropriate signal arrives.
Definition: TotemTestGem.cc:62
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
int getTrackID() const
Definition: TotemG4Hit.cc:135
std::vector< std::string > names
Definition: TotemTestGem.h:63
float getVx() const
Definition: TotemG4Hit.cc:177
float getVz() const
Definition: TotemG4Hit.cc:183
int getParentId() const
Definition: TotemG4Hit.cc:174
uint32_t getUnitID() const
Definition: TotemG4Hit.cc:138
std::vector< TotemG4Hit * > hits
Definition: TotemTestGem.h:65
void fillHit(int uID, int pType, int tID, int pID, float eLoss, float pAbs, float vX, float vY, float vZ, float x, float y, float z)
int getParticleType() const
Definition: TotemG4Hit.cc:152
Definition: DDAxes.h:10
TotemTestGem(const edm::ParameterSet &p)
Definition: TotemTestGem.cc:39
virtual ~TotemTestGem()
Definition: TotemTestGem.cc:48