CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HFShowerPMT.cc
Go to the documentation of this file.
1 // File: HFShowerPMT.cc
3 // Description: Parametrized version of HF hits
5 
8 
9 #include "G4NavigationHistory.hh"
10 #include "G4VPhysicalVolume.hh"
11 #include "G4Step.hh"
12 #include "G4Track.hh"
13 #include "CLHEP/Units/GlobalPhysicalConstants.h"
14 #include "CLHEP/Units/GlobalSystemOfUnits.h"
15 #include <sstream>
16 
17 //#define EDM_ML_DEBUG
18 
20  const HcalDDDSimConstants* hcons,
21  const HcalSimulationParameters* hps,
22  edm::ParameterSet const& p)
23  : hcalConstant_(hcons), hcalsimpar_(hps) {
24  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFShowerPMT");
25  pePerGeV = m_HF.getParameter<double>("PEPerGeVPMT");
26 
27  //Special Geometry parameters
32 #ifdef EDM_ML_DEBUG
33  edm::LogVerbatim("HFShower") << "HFShowerPMT: gets the Index matches for " << pmtR1.size() << " PMTs";
34  for (unsigned int ii = 0; ii < pmtR1.size(); ii++) {
35  edm::LogVerbatim("HFShower") << "HFShowerPMT: rIndexR[" << ii << "] = " << pmtR1[ii] << " fibreR[" << ii
36  << "] = " << pmtFib1[ii] << " rIndexL[" << ii << "] = " << pmtR2[ii] << " fibreL["
37  << ii << "] = " << pmtFib2[ii];
38  }
39 #endif
40  cherenkov_ = std::make_unique<HFCherenkov>(m_HF);
41 
42  // Special Geometry parameters
44 #ifdef EDM_ML_DEBUG
45  std::stringstream sss;
46  for (unsigned int ig = 0; ig < rTable.size(); ++ig) {
47  if (ig / 10 * 10 == ig) {
48  sss << "\n";
49  }
50  sss << " " << rTable[ig] / cm;
51  }
52  edm::LogVerbatim("HFShowerPMT") << "HFShowerPMT: " << rTable.size() << " rTable(cm):" << sss.str();
53 #endif
54 }
55 
57 
58 double HFShowerPMT::getHits(const G4Step* aStep) {
59  indexR = indexF = -1;
60 
61  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
62  const G4VTouchable* touch = preStepPoint->GetTouchable();
63  int boxNo = touch->GetReplicaNumber(2);
64  int pmtNo = touch->GetReplicaNumber(1);
65  if (boxNo <= 1) {
66  indexR = pmtR1[pmtNo - 1];
67  indexF = pmtFib1[pmtNo - 1];
68  } else {
69  indexR = pmtR2[pmtNo - 1];
70  indexF = pmtFib2[pmtNo - 1];
71  }
72 
73 #ifdef EDM_ML_DEBUG
74  double edep = aStep->GetTotalEnergyDeposit();
75  edm::LogVerbatim("HFShower") << "HFShowerPMT: Box " << boxNo << " PMT " << pmtNo << " Mapped Indices " << indexR
76  << ", " << indexF << " Edeposit " << edep / MeV << " MeV; PE " << edep * pePerGeV / GeV;
77 #endif
78 
79  double photons = 0;
80  if (indexR >= 0 && indexF > 0) {
81  G4Track* aTrack = aStep->GetTrack();
82  G4ParticleDefinition* particleDef = aTrack->GetDefinition();
83  double stepl = aStep->GetStepLength();
84  double beta = preStepPoint->GetBeta();
85  G4ThreeVector pDir = aTrack->GetDynamicParticle()->GetMomentumDirection();
86  G4ThreeVector localMom = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformAxis(pDir);
87  photons = cherenkov_->computeNPEinPMT(particleDef, beta, localMom.x(), localMom.y(), localMom.z(), stepl);
88 #ifdef EDM_ML_DEBUG
89  edm::LogVerbatim("HFShower") << "HFShowerPMT::getHits: for particle " << particleDef->GetParticleName() << " Step "
90  << stepl << " Beta " << beta << " Direction " << pDir << " Local " << localMom
91  << " p.e. " << photons;
92 #endif
93  }
94  return photons;
95 }
96 
98  double r = 0.;
99  if (indexR >= 0 && indexR + 1 < (int)(rTable.size()))
100  r = 0.5 * (rTable[indexR] + rTable[indexR + 1]);
101 #ifdef EDM_ML_DEBUG
102  else
103  edm::LogVerbatim("HFShower") << "HFShowerPMT::getRadius: R " << indexR << " F " << indexF;
104 #endif
105  if (indexF == 2)
106  r = -r;
107 #ifdef EDM_ML_DEBUG
108  edm::LogVerbatim("HFShower") << "HFShowerPMT: Radius (" << indexR << "/" << indexF << ") " << r;
109 #endif
110  return r;
111 }
Log< level::Info, true > LogVerbatim
double getRadius()
Definition: HFShowerPMT.cc:97
HFShowerPMT(const std::string &name, const HcalDDDSimConstants *hcons, const HcalSimulationParameters *hps, edm::ParameterSet const &p)
Definition: HFShowerPMT.cc:19
std::vector< double > rTable
Definition: HFShowerPMT.h:35
std::vector< int > pmtFib1
Definition: HFShowerPMT.h:36
std::unique_ptr< HFCherenkov > cherenkov_
Definition: HFShowerPMT.h:32
int ii
Definition: cuy.py:589
virtual ~HFShowerPMT()
Definition: HFShowerPMT.cc:56
double pePerGeV
Definition: HFShowerPMT.h:33
double getHits(const G4Step *aStep)
Definition: HFShowerPMT.cc:58
std::vector< int > pmtR2
Definition: HFShowerPMT.h:37
const std::vector< double > & getRTableHF() const
const HcalDDDSimConstants * hcalConstant_
Definition: HFShowerPMT.h:30
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const HcalSimulationParameters * hcalsimpar_
Definition: HFShowerPMT.h:31
std::vector< int > pmtFib2
Definition: HFShowerPMT.h:37
std::vector< int > pmtR1
Definition: HFShowerPMT.h:36