CMS 3D CMS Logo

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

#include <HFShowerFibreBundle.h>

Public Member Functions

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

Private Attributes

std::unique_ptr< HFCherenkovcherenkov1_
 
std::unique_ptr< HFCherenkovcherenkov2_
 
double facCone
 
double facTube
 
const HcalDDDSimConstantshcalConstant_
 
const HcalSimulationParametershcalsimpar_
 
int indexF
 
int indexR
 
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 HFShowerFibreBundle.h.

Constructor & Destructor Documentation

◆ HFShowerFibreBundle()

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

Definition at line 19 of file HFShowerFibreBundle.cc.

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

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 }
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::vector< int > pmtFib1
std::vector< double > rTable
const HcalDDDSimConstants * hcalConstant_
const std::vector< double > & getRTableHF() const
std::vector< int > pmtR1
std::unique_ptr< HFCherenkov > cherenkov2_
ii
Definition: cuy.py:589
const HcalSimulationParameters * hcalsimpar_
std::unique_ptr< HFCherenkov > cherenkov1_
std::vector< int > pmtR2
std::vector< int > pmtFib2

◆ ~HFShowerFibreBundle()

HFShowerFibreBundle::~HFShowerFibreBundle ( )
virtual

Definition at line 61 of file HFShowerFibreBundle.cc.

61 {}

Member Function Documentation

◆ getHits()

double HFShowerFibreBundle::getHits ( const G4Step *  aStep,
bool  type 
)

Definition at line 63 of file HFShowerFibreBundle.cc.

References HLT_2024v10_cff::beta, cherenkov1_, cherenkov2_, facCone, facTube, indexF, indexR, BPHMonitor_cfi::photons, pmtFib1, pmtFib2, pmtR1, and pmtR2.

63  {
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 }
Log< level::Info, true > LogVerbatim
std::vector< int > pmtFib1
std::vector< int > pmtR1
std::unique_ptr< HFCherenkov > cherenkov2_
std::unique_ptr< HFCherenkov > cherenkov1_
std::vector< int > pmtR2
std::vector< int > pmtFib2

◆ getRadius()

double HFShowerFibreBundle::getRadius ( )

Definition at line 108 of file HFShowerFibreBundle.cc.

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

108  {
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 }
Log< level::Info, true > LogVerbatim
std::vector< double > rTable

Member Data Documentation

◆ cherenkov1_

std::unique_ptr<HFCherenkov> HFShowerFibreBundle::cherenkov1_
private

Definition at line 32 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

◆ cherenkov2_

std::unique_ptr<HFCherenkov> HFShowerFibreBundle::cherenkov2_
private

Definition at line 32 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

◆ facCone

double HFShowerFibreBundle::facCone
private

Definition at line 33 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

◆ facTube

double HFShowerFibreBundle::facTube
private

Definition at line 33 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

◆ hcalConstant_

const HcalDDDSimConstants* HFShowerFibreBundle::hcalConstant_
private

Definition at line 30 of file HFShowerFibreBundle.h.

Referenced by HFShowerFibreBundle().

◆ hcalsimpar_

const HcalSimulationParameters* HFShowerFibreBundle::hcalsimpar_
private

Definition at line 31 of file HFShowerFibreBundle.h.

Referenced by HFShowerFibreBundle().

◆ indexF

int HFShowerFibreBundle::indexF
private

Definition at line 34 of file HFShowerFibreBundle.h.

Referenced by getHits(), and getRadius().

◆ indexR

int HFShowerFibreBundle::indexR
private

Definition at line 34 of file HFShowerFibreBundle.h.

Referenced by getHits(), and getRadius().

◆ pmtFib1

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

Definition at line 36 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

◆ pmtFib2

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

Definition at line 37 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

◆ pmtR1

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

Definition at line 36 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

◆ pmtR2

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

Definition at line 37 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

◆ rTable

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

Definition at line 35 of file HFShowerFibreBundle.h.

Referenced by getRadius(), and HFShowerFibreBundle().