CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
HFShowerPMT Class Reference

#include <HFShowerPMT.h>

Public Member Functions

double getHits (const G4Step *aStep)
 
double getRadius ()
 
 HFShowerPMT (const std::string &name, const HcalDDDSimConstants *hcons, const HcalSimulationParameters *hps, edm::ParameterSet const &p)
 
virtual ~HFShowerPMT ()
 

Private Attributes

std::unique_ptr< HFCherenkovcherenkov_
 
const HcalDDDSimConstantshcalConstant_
 
const HcalSimulationParametershcalsimpar_
 
int indexF
 
int indexR
 
double pePerGeV
 
std::vector< int > pmtFib1
 
std::vector< int > pmtFib2
 
std::vector< int > pmtR1
 
std::vector< int > pmtR2
 
std::vector< double > rTable
 

Detailed Description

Definition at line 19 of file HFShowerPMT.h.

Constructor & Destructor Documentation

◆ HFShowerPMT()

HFShowerPMT::HFShowerPMT ( const std::string &  name,
const HcalDDDSimConstants hcons,
const HcalSimulationParameters hps,
edm::ParameterSet const &  p 
)

Definition at line 19 of file HFShowerPMT.cc.

References cherenkov_, edm::ParameterSet::getParameter(), HcalDDDSimConstants::getRTableHF(), hcalConstant_, hcalsimpar_, cuy::ii, AlCaHLTBitMon_ParallelJobs::p, pePerGeV, pmtFib1, pmtFib2, HcalSimulationParameters::pmtFiberLeft_, HcalSimulationParameters::pmtFiberRight_, HcalSimulationParameters::pmtLeft_, pmtR1, pmtR2, HcalSimulationParameters::pmtRight_, and rTable.

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] / CLHEP::cm;
51  }
52  edm::LogVerbatim("HFShowerPMT") << "HFShowerPMT: " << rTable.size() << " rTable(cm):" << sss.str();
53 #endif
54 }
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
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
const std::vector< double > & getRTableHF() const
double pePerGeV
Definition: HFShowerPMT.h:33
std::vector< int > pmtR2
Definition: HFShowerPMT.h:37
const HcalDDDSimConstants * hcalConstant_
Definition: HFShowerPMT.h:30
ii
Definition: cuy.py:589
const HcalSimulationParameters * hcalsimpar_
Definition: HFShowerPMT.h:31
std::vector< int > pmtFib2
Definition: HFShowerPMT.h:37
std::vector< int > pmtR1
Definition: HFShowerPMT.h:36

◆ ~HFShowerPMT()

HFShowerPMT::~HFShowerPMT ( )
virtual

Definition at line 56 of file HFShowerPMT.cc.

56 {}

Member Function Documentation

◆ getHits()

double HFShowerPMT::getHits ( const G4Step *  aStep)

Definition at line 58 of file HFShowerPMT.cc.

References HLT_2024v13_cff::beta, cherenkov_, indexF, indexR, pePerGeV, BPHMonitor_cfi::photons, pmtFib1, pmtFib2, pmtR1, and pmtR2.

58  {
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 / CLHEP::MeV << " MeV; PE "
77  << edep * pePerGeV / CLHEP::GeV;
78 #endif
79 
80  double photons = 0;
81  if (indexR >= 0 && indexF > 0) {
82  G4Track* aTrack = aStep->GetTrack();
83  G4ParticleDefinition* particleDef = aTrack->GetDefinition();
84  double stepl = aStep->GetStepLength();
85  double beta = preStepPoint->GetBeta();
86  G4ThreeVector pDir = aTrack->GetDynamicParticle()->GetMomentumDirection();
87  G4ThreeVector localMom = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformAxis(pDir);
88  photons = cherenkov_->computeNPEinPMT(particleDef, beta, localMom.x(), localMom.y(), localMom.z(), stepl);
89 #ifdef EDM_ML_DEBUG
90  edm::LogVerbatim("HFShower") << "HFShowerPMT::getHits: for particle " << particleDef->GetParticleName() << " Step "
91  << stepl << " Beta " << beta << " Direction " << pDir << " Local " << localMom
92  << " p.e. " << photons;
93 #endif
94  }
95  return photons;
96 }
Log< level::Info, true > LogVerbatim
std::vector< int > pmtFib1
Definition: HFShowerPMT.h:36
std::unique_ptr< HFCherenkov > cherenkov_
Definition: HFShowerPMT.h:32
double pePerGeV
Definition: HFShowerPMT.h:33
std::vector< int > pmtR2
Definition: HFShowerPMT.h:37
std::vector< int > pmtFib2
Definition: HFShowerPMT.h:37
std::vector< int > pmtR1
Definition: HFShowerPMT.h:36

◆ getRadius()

double HFShowerPMT::getRadius ( )

Definition at line 98 of file HFShowerPMT.cc.

References indexF, indexR, alignCSCRings::r, and rTable.

98  {
99  double r = 0.;
100  if (indexR >= 0 && indexR + 1 < (int)(rTable.size()))
101  r = 0.5 * (rTable[indexR] + rTable[indexR + 1]);
102 #ifdef EDM_ML_DEBUG
103  else
104  edm::LogVerbatim("HFShower") << "HFShowerPMT::getRadius: R " << indexR << " F " << indexF;
105 #endif
106  if (indexF == 2)
107  r = -r;
108 #ifdef EDM_ML_DEBUG
109  edm::LogVerbatim("HFShower") << "HFShowerPMT: Radius (" << indexR << "/" << indexF << ") " << r;
110 #endif
111  return r;
112 }
Log< level::Info, true > LogVerbatim
std::vector< double > rTable
Definition: HFShowerPMT.h:35

Member Data Documentation

◆ cherenkov_

std::unique_ptr<HFCherenkov> HFShowerPMT::cherenkov_
private

Definition at line 32 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

◆ hcalConstant_

const HcalDDDSimConstants* HFShowerPMT::hcalConstant_
private

Definition at line 30 of file HFShowerPMT.h.

Referenced by HFShowerPMT().

◆ hcalsimpar_

const HcalSimulationParameters* HFShowerPMT::hcalsimpar_
private

Definition at line 31 of file HFShowerPMT.h.

Referenced by HFShowerPMT().

◆ indexF

int HFShowerPMT::indexF
private

Definition at line 34 of file HFShowerPMT.h.

Referenced by getHits(), and getRadius().

◆ indexR

int HFShowerPMT::indexR
private

Definition at line 34 of file HFShowerPMT.h.

Referenced by getHits(), and getRadius().

◆ pePerGeV

double HFShowerPMT::pePerGeV
private

Definition at line 33 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

◆ pmtFib1

std::vector<int> HFShowerPMT::pmtFib1
private

Definition at line 36 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

◆ pmtFib2

std::vector<int> HFShowerPMT::pmtFib2
private

Definition at line 37 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

◆ pmtR1

std::vector<int> HFShowerPMT::pmtR1
private

Definition at line 36 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

◆ pmtR2

std::vector<int> HFShowerPMT::pmtR2
private

Definition at line 37 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

◆ rTable

std::vector<double> HFShowerPMT::rTable
private

Definition at line 35 of file HFShowerPMT.h.

Referenced by getRadius(), and HFShowerPMT().