CMS 3D CMS Logo

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