CMS 3D CMS Logo

HFShowerFibreBundle.cc
Go to the documentation of this file.
1 // File: HFShowerFibreBundle.cc
3 // Description: Hits in the fibre bundles
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_HF1 = p.getParameter<edm::ParameterSet>("HFShowerStraightBundle");
25  facTube = m_HF1.getParameter<double>("FactorBundle");
26  cherenkov1_ = std::make_unique<HFCherenkov>(m_HF1);
27  edm::ParameterSet m_HF2 = p.getParameter<edm::ParameterSet>("HFShowerConicalBundle");
28  facCone = m_HF2.getParameter<double>("FactorBundle");
29  cherenkov2_ = std::make_unique<HFCherenkov>(m_HF2);
30  edm::LogVerbatim("HFShower") << "HFShowerFibreBundle intialized with factors: " << facTube
31  << " for the straight portion and " << facCone << " for the curved portion";
32 
33  //Special Geometry parameters
38 #ifdef EDM_ML_DEBUG
39  edm::LogVerbatim("HFShower") << "HFShowerPMT: gets the Index matches for " << pmtR1.size() << " PMTs";
40  for (unsigned int ii = 0; ii < pmtR1.size(); ii++) {
41  edm::LogVerbatim("HFShower") << "HFShowerPMT: rIndexR[" << ii << "] = " << pmtR1[ii] << " fibreR[" << ii
42  << "] = " << pmtFib1[ii] << " rIndexL[" << ii << "] = " << pmtR2[ii] << " fibreL["
43  << ii << "] = " << pmtFib2[ii];
44  }
45 #endif
46 
47  // Other Geometry parameters
49 #ifdef EDM_ML_DEBUG
50  std::stringstream sss;
51  for (unsigned int ig = 0; ig < rTable.size(); ig++) {
52  if (ig / 10 * 10 == ig) {
53  sss << "\n";
54  }
55  sss << " " << rTable[ig] / CLHEP::cm;
56  }
57  edm::LogVerbatim("HFShower") << "HFShowerFibreBundle: " << rTable.size() << " rTable(cm):" << sss.str();
58 #endif
59 }
60 
62 
63 double HFShowerFibreBundle::getHits(const G4Step* aStep, bool type) {
64  indexR = indexF = -1;
65 
66  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
67  const G4VTouchable* touch = preStepPoint->GetTouchable();
68  int boxNo = touch->GetReplicaNumber(1);
69  int pmtNo = touch->GetReplicaNumber(0);
70  if (boxNo <= 1) {
71  indexR = pmtR1[pmtNo - 1];
72  indexF = pmtFib1[pmtNo - 1];
73  } else {
74  indexR = pmtR2[pmtNo - 1];
75  indexF = pmtFib2[pmtNo - 1];
76  }
77 
78 #ifdef EDM_ML_DEBUG
79  double edep = aStep->GetTotalEnergyDeposit();
80  edm::LogVerbatim("HFShower") << "HFShowerFibreBundle: Box " << boxNo << " PMT " << pmtNo << " Mapped Indices "
81  << indexR << ", " << indexF << " Edeposit " << edep / MeV << " MeV";
82 #endif
83 
84  double photons = 0;
85  if (indexR >= 0 && indexF > 0) {
86  const G4Track* aTrack = aStep->GetTrack();
87  const G4ParticleDefinition* particleDef = aTrack->GetDefinition();
88  double stepl = aStep->GetStepLength();
89  double beta = preStepPoint->GetBeta();
90  G4ThreeVector pDir = aTrack->GetDynamicParticle()->GetMomentumDirection();
91  G4ThreeVector localMom = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformAxis(pDir);
92  if (type) {
93  photons =
94  facCone * cherenkov2_->computeNPEinPMT(particleDef, beta, localMom.x(), localMom.y(), localMom.z(), stepl);
95  } else {
96  photons =
97  facTube * cherenkov1_->computeNPEinPMT(particleDef, beta, localMom.x(), localMom.y(), localMom.z(), stepl);
98  }
99 #ifdef EDM_ML_DEBUG
100  edm::LogVerbatim("HFShower") << "HFShowerFibreBundle::getHits: for particle " << particleDef->GetParticleName()
101  << " Step " << stepl << " Beta " << beta << " Direction " << pDir << " Local "
102  << localMom << " p.e. " << photons;
103 #endif
104  }
105  return photons;
106 }
107 
109  double r = 0.;
110  if (indexR >= 0 && indexR + 1 < (int)(rTable.size()))
111  r = 0.5 * (rTable[indexR] + rTable[indexR + 1]);
112 #ifdef EDM_ML_DEBUG
113  else
114  edm::LogVerbatim("HFShower") << "HFShowerFibreBundle::getRadius: R " << indexR << " F " << indexF;
115 #endif
116  if (indexF == 2)
117  r = -r;
118 #ifdef EDM_ML_DEBUG
119  edm::LogVerbatim("HFShower") << "HFShowerFibreBundle: Radius (" << indexR << "/" << indexF << ") " << r;
120 #endif
121  return r;
122 }
HFShowerFibreBundle::cherenkov2_
std::unique_ptr< HFCherenkov > cherenkov2_
Definition: HFShowerFibreBundle.h:32
HFShowerFibreBundle::pmtFib1
std::vector< int > pmtFib1
Definition: HFShowerFibreBundle.h:36
HFShowerFibreBundle::hcalConstant_
const HcalDDDSimConstants * hcalConstant_
Definition: HFShowerFibreBundle.h:30
HcalDDDSimConstants::getRTableHF
const std::vector< double > & getRTableHF() const
Definition: HcalDDDSimConstants.h:61
zMuMuMuonUserData.beta
beta
Definition: zMuMuMuonUserData.py:10
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
HFShowerFibreBundle::~HFShowerFibreBundle
virtual ~HFShowerFibreBundle()
Definition: HFShowerFibreBundle.cc:61
HcalSimulationParameters::pmtLeft_
std::vector< int > pmtLeft_
Definition: HcalSimulationParameters.h:18
MeV
const double MeV
HcalSimulationParameters::pmtFiberRight_
std::vector< int > pmtFiberRight_
Definition: HcalSimulationParameters.h:17
HFShowerFibreBundle::facTube
double facTube
Definition: HFShowerFibreBundle.h:33
HcalDDDSimConstants
Definition: HcalDDDSimConstants.h:24
HFShowerFibreBundle::facCone
double facCone
Definition: HFShowerFibreBundle.h:33
HFShowerFibreBundle::HFShowerFibreBundle
HFShowerFibreBundle(const std::string &name, const HcalDDDSimConstants *hcons, const HcalSimulationParameters *hps, edm::ParameterSet const &p)
Definition: HFShowerFibreBundle.cc:19
HFShowerFibreBundle::hcalsimpar_
const HcalSimulationParameters * hcalsimpar_
Definition: HFShowerFibreBundle.h:31
HcalSimulationParameters::pmtFiberLeft_
std::vector< int > pmtFiberLeft_
Definition: HcalSimulationParameters.h:19
HcalSimulationParameters::pmtRight_
std::vector< int > pmtRight_
Definition: HcalSimulationParameters.h:16
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
edm::ParameterSet
Definition: ParameterSet.h:36
HFShowerFibreBundle::indexF
int indexF
Definition: HFShowerFibreBundle.h:34
HFShowerFibreBundle::pmtR2
std::vector< int > pmtR2
Definition: HFShowerFibreBundle.h:37
HFShowerFibreBundle::pmtFib2
std::vector< int > pmtFib2
Definition: HFShowerFibreBundle.h:37
HFShowerFibreBundle::getHits
double getHits(const G4Step *aStep, bool type)
Definition: HFShowerFibreBundle.cc:63
edm::LogVerbatim
Definition: MessageLogger.h:297
HFShowerFibreBundle::indexR
int indexR
Definition: HFShowerFibreBundle.h:34
BPHMonitor_cfi.photons
photons
Definition: BPHMonitor_cfi.py:91
HFShowerFibreBundle.h
HFShowerFibreBundle::rTable
std::vector< double > rTable
Definition: HFShowerFibreBundle.h:35
alignCSCRings.r
r
Definition: alignCSCRings.py:93
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
type
type
Definition: HCALResponse.h:21
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
Exception.h
HFShowerFibreBundle::cherenkov1_
std::unique_ptr< HFCherenkov > cherenkov1_
Definition: HFShowerFibreBundle.h:32
HcalSimulationParameters
Definition: HcalSimulationParameters.h:6
cuy.ii
ii
Definition: cuy.py:590
HFShowerFibreBundle::getRadius
double getRadius()
Definition: HFShowerFibreBundle.cc:108
HFShowerFibreBundle::pmtR1
std::vector< int > pmtR1
Definition: HFShowerFibreBundle.h:36