test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HcalForwardAnalysis.cc
Go to the documentation of this file.
1 // system include files
2 #include <cmath>
3 #include <iostream>
4 #include <iomanip>
5 
6 // user include files
8 
12 
14 
20 
21 #include "G4SDManager.hh"
22 #include "G4VProcess.hh"
23 #include "G4HCofThisEvent.hh"
24 #include "CLHEP/Units/GlobalSystemOfUnits.h"
25 #include "CLHEP/Units/GlobalPhysicalConstants.h"
26 
28 
29  edm::ParameterSet m_SLP = p.getParameter<edm::ParameterSet> ("HFShowerLibraryProducer");
30  theNames = m_SLP.getParameter<std::vector<std::string> > ("Names");
31  //LibVer = m_HS.getParameter<std::string> ("LibVer");
32  //produces<HFShowerPhotonCollection> ();
33  init();
34  theEventCounter = 0;
35  nphot = 0;
36  for (int i = 0; i < 10000; ++i) {
37  x[i] = 0.;
38  y[i] = 0.;
39  z[i] = 0.;
40  t[i] = 0.;
41  lambda[i] = 0.;
42  fiberId[i] = 0;
43  }
44  primX = primY = primZ = primT = 0.;
45  primMomX = primMomY = primMomZ = 0.;
46 
47 }
48 
50 }
51 
52 //
53 // member functions
54 //
55 
57 
58  //std::auto_ptr<HFShowerPhotonCollection> product(new HFShowerPhotonCollection);
59  //edm::LogInfo("HcalForwardLib") << "HcalForwardAnalysis: =====> Filling event";
60  //fillEvent(*product);
61  //iEvent.put(product);
62  //std::auto_ptr<PHcalForwardLibInfo> product(new PHcalForwardLibInfo);
63  //fillEvent(*product);
64  fillEvent();
65  //iEvent.put(product);
66 }
67 
69 
70  theTree = theFile->make<TTree> ("CherenkovPhotons", "Cherenkov Photons");
71  theTree->Branch("nphot", &nphot, "nphot/I");
72  theTree->Branch("x", &x, "x[nphot]/F");
73  theTree->Branch("y", &y, "y[nphot]/F");
74  theTree->Branch("z", &z, "z[nphot]/F");
75  theTree->Branch("t", &t, "t[nphot]/F");
76  theTree->Branch("lambda", &lambda, "lambda[nphot]/F");
77  theTree->Branch("fiberId", &fiberId, "fiberId[nphot]/I");
78  theTree->Branch("primX", &primX, "primX/F");
79  theTree->Branch("primY", &primY, "primY/F");
80  theTree->Branch("primZ", &primZ, "primZ/F");
81  theTree->Branch("primMomX", &primMomX, "primMomX/F");
82  theTree->Branch("primMomY", &primMomY, "primMomY/F");
83  theTree->Branch("primMomZ", &primMomZ, "primMomZ/F");
84  theTree->Branch("primT", &primT, "primT/F");
85 
86  // counter
87  count = 0;
88  evNum = 0;
89  clear();
90 }
91 
93 
94  int irun = (*run)()->GetRunID();
95  edm::LogInfo("HcalForwardLib") << " =====> Begin of Run = " << irun;
96 
97 }
98 
100 
101  evNum = (*evt)()->GetEventID();
102  clear();
103  edm::LogInfo("HcalForwardLib") << "HcalForwardAnalysis: =====> Begin of event = " << evNum;
104 }
105 
106 void HcalForwardAnalysis::update(const G4Step * aStep) {
107 }
108 
110 
111  count++;
112 
113  //fill the buffer
114  LogDebug("HcalForwardLib")<< "HcalForwardAnalysis::Fill event "
115  << (*evt)()->GetEventID();
116  setPhotons(evt);
117 
118  int iEvt = (*evt)()->GetEventID();
119  if (iEvt < 10)
120  edm::LogInfo("HcalForwardLib") << "HcalForwardAnalysis:: Event " << iEvt;
121  else if ((iEvt < 100) && (iEvt%10 == 0))
122  edm::LogInfo("HcalForwardLib") << "HcalForwardAnalysis:: Event " << iEvt;
123  else if ((iEvt < 1000) && (iEvt%100 == 0))
124  edm::LogInfo("HcalForwardLib") << "HcalForwardAnalysis:: Event " << iEvt;
125  else if ((iEvt < 10000) && (iEvt%1000 == 0))
126  edm::LogInfo("HcalForwardLib") << "HcalForwardAnalysis:: Event " << iEvt;
127 }
128 
130 
131  int idHC, j;
132  FiberG4HitsCollection* theHC;
133  // Look for the Hit Collection of HCal
134  G4HCofThisEvent* allHC = (*evt)()->GetHCofThisEvent();
135  std::string sdName = theNames[0];//name for fiber hits
136  idHC = G4SDManager::GetSDMpointer()->GetCollectionID(sdName);
137  theHC = (FiberG4HitsCollection*) allHC->GetHC(idHC);
138  LogDebug("HcalForwardLib")<< "HcalForwardAnalysis::setPhotons() Hit Collection for " << sdName
139  << " of ID " << idHC << " is obtained at " << theHC;
140  std::vector<HFShowerPhoton> ShortFiberPhotons;
141  std::vector<HFShowerPhoton> LongFiberPhotons;
142  LongFiberPhotons.clear();
143  ShortFiberPhotons.clear();
144  if (idHC >= 0 && theHC> 0) {
145  std::cout << "FiberhitSize " << theHC->entries() << std::endl;
146  for (j = 0; j < theHC->entries(); j++) {
147  FiberG4Hit* aHit = (*theHC)[j];
148  std::vector<HFShowerPhoton> thePhotonsFromHit = aHit->photon();
149  std::cout << "Fiberhit " << j << " has " << thePhotonsFromHit.size() << " photons." << std::endl;
150  int fTowerId = -1;
151  int fCellId = -1;
152  int fFiberId = -1;
153  parseDetId(aHit->towerId(), fTowerId, fCellId, fFiberId);
154  for(unsigned int iph = 0; iph < thePhotonsFromHit.size(); ++iph){
155  if(aHit->depth() == 1)LongFiberPhotons.push_back(thePhotonsFromHit[iph]);
156  if(aHit->depth() == 2)ShortFiberPhotons.push_back(thePhotonsFromHit[iph]);
157  }
158  LogDebug("HcalForwardLib") << "HcalForwardAnalysis::setPhotons() NbPhotons " << thePhotonsFromHit.size()
159  << " towerId " << fTowerId << " cellId " << fCellId << " fiberId " << fFiberId << " depth " << aHit->depth();
160  }
161  }else{
162  LogDebug("HcalForwardLib") << "HcalForwardAnalysis::setPhotons(): No Photons!";
163  return;
164  }
165  LogDebug("HcalForwardLib") << "HcalForwardAnalysis::setPhotons() LongFibPhotons: "
166  << LongFiberPhotons.size() << " ShortFibPhotons: " << ShortFiberPhotons.size();
167  std::cout << "HcalForwardAnalysis::setPhotons() LongFibPhotons: "
168  << LongFiberPhotons.size() << " ShortFibPhotons: " << ShortFiberPhotons.size() << std::endl;
169 
170  //Chamber hits to find information about primary particle on surface
171  HFShowerG4HitsCollection* theChamberHC;
172  G4HCofThisEvent* allChamberHC = (*evt)()->GetHCofThisEvent();
173  sdName= theNames[1];
174  idHC = G4SDManager::GetSDMpointer()->GetCollectionID(sdName);
175  theChamberHC = (HFShowerG4HitsCollection*) allChamberHC->GetHC(idHC);
176  math::XYZPoint primPosOnSurf(0,0,0);
177  math::XYZPoint primMomDirOnSurf(0,0,0);
178  float primTimeOnSurf = 0;
179  // the chamber hit is for primary particle, but step size can be small
180  // (in newer Geant4 versions) and as a result primary particle may have
181  // multiple hits. We want to take last one which is close the HF absorber
182  if(idHC >= 0 && theChamberHC> 0) {
183  LogDebug("HcalForwardLib") << "HcalForwardAnalysis::setPhotons() Chamber Hits size: " << theChamberHC->entries();
184  for(j = 0; j < theChamberHC->entries();++j) {
185  HFShowerG4Hit* aHit = (*theChamberHC)[j];
186  LogDebug("HcalForwardLib") << "HcalForwardAnalysis::setPhotons() Chamber Hit id " << aHit->hitId()
187  << " track id " << aHit->trackId() << " prim. pos. " << aHit->globalPosition() << " prom mom. dir. "
188  << aHit->primaryMomDir() << " time " << aHit->time();
189  primPosOnSurf.SetXYZ(aHit->globalPosition().x(), aHit->globalPosition().y(), aHit->globalPosition().z());
190  primMomDirOnSurf.SetXYZ(aHit->primaryMomDir().x(), aHit->primaryMomDir().y(), aHit->primaryMomDir().z());
191  primTimeOnSurf = aHit->time();
192  }
193  }else{
194  LogDebug("HcalForwardLib") << "HcalForwardAnalysis::setPhotons(): No Chamber hits. Something is wrong!";
195  return;
196  }
197  primX = primPosOnSurf.x();
198  primY = primPosOnSurf.y();
199  primZ = primPosOnSurf.z();
200  primT = primTimeOnSurf;
201  primMomX = primMomDirOnSurf.x();
202  primMomY = primMomDirOnSurf.y();
203  primMomZ = primMomDirOnSurf.z();
204  //production point of Cherenkov photons in long fibers w.r.t local coordinates
205  float photonProdX = 0.;
206  float photonProdY = 0.;
207  float photonProdZ = 0.;
208  float photonProdTime = 0.;
209  //angles for rotation matrices
210  double theta = primMomDirOnSurf.theta();
211  double phi = primMomDirOnSurf.phi();
212 
213  for(unsigned int k = 0; k < LongFiberPhotons.size(); ++k){
214  HFShowerPhoton aPhoton = LongFiberPhotons[k];
215  photonProdX = aPhoton.x()*cos(theta)*cos(phi) + aPhoton.y()*cos(theta)*sin(phi) - aPhoton.z()*sin(theta) - primPosOnSurf.x();
216  photonProdY = -aPhoton.x()*sin(phi) + aPhoton.y()*cos(phi) - primPosOnSurf.y();
217  photonProdZ = aPhoton.x()*sin(theta)*cos(phi) + aPhoton.y()*sin(theta)*sin(phi) + aPhoton.z()*cos(theta) - primPosOnSurf.z();
218  photonProdTime = aPhoton.t()- primTimeOnSurf;
219  thePhotons.push_back(Photon(1, photonProdX, photonProdY, photonProdZ, photonProdTime, aPhoton.lambda()));
220  }
221  for(unsigned int k = 0; k < ShortFiberPhotons.size(); ++k){
222  HFShowerPhoton aPhoton = ShortFiberPhotons[k];
223  photonProdX = aPhoton.x()*cos(theta)*cos(phi) + aPhoton.y()*cos(theta)*sin(phi) - aPhoton.z()*sin(theta) - primPosOnSurf.x();
224  photonProdY = -aPhoton.x()*sin(phi) + aPhoton.y()*cos(phi) - primPosOnSurf.y();
225  photonProdZ = aPhoton.x()*sin(theta)*cos(phi) + aPhoton.y()*sin(theta)*sin(phi) + aPhoton.z()*cos(theta) - primPosOnSurf.z();
226  photonProdTime = aPhoton.t()- primTimeOnSurf;
227  thePhotons.push_back(Photon(2, photonProdX, photonProdY, photonProdZ, photonProdTime, aPhoton.lambda()));
228  }
229 }
231  /*
232  edm::LogInfo("HcalForwardLib") << "HcalForwardAnalysis: =====> filledEvent";
233  */
234  nphot = int(thePhotons.size());
235  for (int i = 0; i < nphot; ++i) {
236  x[i] = thePhotons[i].x;
237  y[i] = thePhotons[i].y;
238  z[i] = thePhotons[i].z;
239  t[i] = thePhotons[i].t;
240  lambda[i] = thePhotons[i].lambda;
241  fiberId[i] = thePhotons[i].fiberId;
242  }
243  theTree->Fill();
244 }
245 void HcalForwardAnalysis::parseDetId(int id, int& tower, int& cell, int& fiber) {
246  tower = id / 10000;
247  cell = id / 10 - tower * 10;
248  fiber = id - tower * 10000 - cell * 10;
249 }
250 
252  nphot = 0;
253  for (int i = 0; i < 10000; ++i) {
254  x[i] = 0.;
255  y[i] = 0.;
256  z[i] = 0.;
257  t[i] = 0.;
258  lambda[i] = 0.;
259  fiberId[i] = 0;
260  }
261  primX = primY = primZ = primT = 0.;
262  primMomX = primMomY = primMomZ = 0.;
263 
264  thePhotons.clear();
265 }
#define LogDebug(id)
T getParameter(std::string const &) const
G4ThreeVector primaryMomDir() const
Definition: HFShowerG4Hit.h:57
int i
Definition: DBlmapReader.cc:9
G4int depth() const
Definition: FiberG4Hit.h:50
G4int hitId() const
Definition: HFShowerG4Hit.h:51
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
void setPhotons(const EndOfEvent *evt)
void update(const BeginOfRun *run)
This routine will be called when the appropriate signal arrives.
int iEvent
Definition: GenABIO.cc:230
float t() const
std::vector< Photon > thePhotons
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
G4int towerId() const
Definition: FiberG4Hit.h:49
int j
Definition: DBlmapReader.cc:9
G4int trackId() const
Definition: HFShowerG4Hit.h:52
std::vector< HFShowerPhoton > photon() const
Definition: FiberG4Hit.h:55
float z() const
float lambda() const
edm::Service< TFileService > theFile
float x() const
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::vector< std::string > theNames
G4THitsCollection< FiberG4Hit > FiberG4HitsCollection
Definition: FiberG4Hit.h:59
G4ThreeVector globalPosition() const
Definition: HFShowerG4Hit.h:56
tuple cout
Definition: gather_cfg.py:121
HcalForwardAnalysis(const edm::ParameterSet &p)
G4double time() const
Definition: HFShowerG4Hit.h:54
G4THitsCollection< HFShowerG4Hit > HFShowerG4HitsCollection
Definition: HFShowerG4Hit.h:60
void parseDetId(int id, int &tower, int &cell, int &fiber)
virtual void produce(edm::Event &, const edm::EventSetup &)
float y() const