CMS 3D CMS Logo

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

#include <HFFibre.h>

Classes

struct  Params
 

Public Member Functions

double attLength (double lambda) const
 
 HFFibre (const HcalDDDSimConstants *hcons, const HcalSimulationParameters *hps, edm::ParameterSet const &p)
 
 HFFibre (Params iP)
 
double tShift (const G4ThreeVector &point, int depth, int fromEndAbs=0) const
 
double zShift (const G4ThreeVector &point, int depth, int fromEndAbs=0) const
 

Private Attributes

std::vector< double > attL_
 
double cFibre_
 
std::vector< double > gpar_
 
std::array< double, 2 > lambLim_
 
std::vector< double > longFL_
 
int nBinAtt_
 
int nBinR_
 
std::vector< double > radius_
 
std::vector< double > shortFL_
 

Detailed Description

Definition at line 19 of file HFFibre.h.

Constructor & Destructor Documentation

◆ HFFibre() [1/2]

HFFibre::HFFibre ( const HcalDDDSimConstants hcons,
const HcalSimulationParameters hps,
edm::ParameterSet const &  p 
)

Definition at line 26 of file HFFibre.cc.

27  : HFFibre(Params(p.getParameter<edm::ParameterSet>("HFShower")
28  .getParameter<edm::ParameterSet>("HFShowerBlock")
29  .getParameter<double>("CFibre"),
30  hcons,
31  hps)) {}
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
HFFibre(const HcalDDDSimConstants *hcons, const HcalSimulationParameters *hps, edm::ParameterSet const &p)
Definition: HFFibre.cc:26

◆ HFFibre() [2/2]

HFFibre::HFFibre ( Params  iP)

Definition at line 33 of file HFFibre.cc.

References attL_, cFibre_, gpar_, mps_fire::i, lambLim_, longFL_, visualization-live-secondInstance_cfg::m, nBinAtt_, nBinR_, radius_, and shortFL_.

34  : cFibre_(c_light * iP.fractionOfSpeedOfLightInFibre_),
35  gpar_(std::move(iP.gParHF_)),
36  radius_(std::move(iP.rTableHF_)),
37  shortFL_(std::move(iP.shortFibreLength_)),
38  longFL_(std::move(iP.longFibreLength_)),
39  attL_(std::move(iP.attenuationLength_)),
40  lambLim_(iP.lambdaLimits_) {
41  edm::LogVerbatim("HFShower") << "HFFibre:: Speed of light in fibre " << cFibre_ / (CLHEP::m / CLHEP::ns) << " m/ns";
42  // Attenuation length
43  nBinAtt_ = static_cast<int>(attL_.size());
44 
45  edm::LogVerbatim("HFShower").log([&](auto& logger) {
46  logger << "HFFibre: " << nBinAtt_ << " attL(1/cm): ";
47  for (int it = 0; it < nBinAtt_; it++) {
48  if (it / 10 * 10 == it) {
49  logger << "\n";
50  }
51  logger << " " << attL_[it] * CLHEP::cm;
52  }
53  });
54  edm::LogVerbatim("HFShower") << "HFFibre: Limits on lambda " << lambLim_[0] << " and " << lambLim_[1];
55  // Fibre Lengths
56  edm::LogVerbatim("HFShower").log([&](auto& logger) {
57  logger << "HFFibre: " << longFL_.size() << " Long Fibre Length(cm):";
58  for (unsigned int it = 0; it < longFL_.size(); it++) {
59  if (it / 10 * 10 == it) {
60  logger << "\n";
61  }
62  logger << " " << longFL_[it] / CLHEP::cm;
63  }
64  });
65  edm::LogVerbatim("HFShower").log([&](auto& logger) {
66  logger << "HFFibre: " << shortFL_.size() << " Short Fibre Length(cm):";
67  for (unsigned int it = 0; it < shortFL_.size(); it++) {
68  if (it / 10 * 10 == it) {
69  logger << "\n";
70  }
71  logger << " " << shortFL_[it] / CLHEP::cm;
72  }
73  });
74  // Now geometry parameters
75 
76  nBinR_ = static_cast<int>(radius_.size());
77  edm::LogVerbatim("HFShower").log([&](auto& logger) {
78  logger << "HFFibre: " << radius_.size() << " rTable(cm):";
79  for (int i = 0; i < nBinR_; ++i) {
80  if (i / 10 * 10 == i) {
81  logger << "\n";
82  }
83  logger << " " << radius_[i] / CLHEP::cm;
84  }
85  });
86 
87  edm::LogVerbatim("HFShower").log([&](auto& logger) {
88  logger << "HFFibre: " << gpar_.size() << " gParHF:";
89  for (std::size_t i = 0; i < gpar_.size(); ++i) {
90  if (i / 10 * 10 == i) {
91  logger << "\n";
92  }
93  logger << " " << gpar_[i];
94  }
95  });
96 }
std::vector< double > radius_
Definition: HFFibre.h:46
std::vector< double > shortFL_
Definition: HFFibre.h:47
Log< level::Info, true > LogVerbatim
std::vector< double > longFL_
Definition: HFFibre.h:47
int nBinR_
Definition: HFFibre.h:49
std::vector< double > gpar_
Definition: HFFibre.h:46
Definition: logger.py:1
double cFibre_
Definition: HFFibre.h:45
int nBinAtt_
Definition: HFFibre.h:49
std::vector< double > attL_
Definition: HFFibre.h:48
std::array< double, 2 > lambLim_
Definition: HFFibre.h:50
def move(src, dest)
Definition: eostools.py:511

Member Function Documentation

◆ attLength()

double HFFibre::attLength ( double  lambda) const

Definition at line 98 of file HFFibre.cc.

References attL_, mps_fire::i, createfilelist::int, dqmiolumiharvest::j, lambLim_, and nBinAtt_.

Referenced by HFShowerLibrary::fillHits(), HFShower::getHits(), and HFShowerParam::HFShowerParam().

98  {
99  int i = int(nBinAtt_ * (lambda - lambLim_[0]) / (lambLim_[1] - lambLim_[0]));
100 
101  int j = i;
102  if (i >= nBinAtt_)
103  j = nBinAtt_ - 1;
104  else if (i < 0)
105  j = 0;
106  double att = attL_[j];
107 #ifdef EDM_ML_DEBUG
108  edm::LogVerbatim("HFShower") << "HFFibre::attLength for Lambda " << lambda << " index " << i << " " << j
109  << " Att. Length " << att;
110 #endif
111  return att;
112 }
Log< level::Info, true > LogVerbatim
int nBinAtt_
Definition: HFFibre.h:49
std::vector< double > attL_
Definition: HFFibre.h:48
std::array< double, 2 > lambLim_
Definition: HFFibre.h:50

◆ tShift()

double HFFibre::tShift ( const G4ThreeVector &  point,
int  depth,
int  fromEndAbs = 0 
) const

Definition at line 114 of file HFFibre.cc.

References cFibre_, LEDCalibrationChannels::depth, point, protons_cff::time, and zShift().

Referenced by HFShowerLibrary::fillHits(), HFShowerParam::getHits(), and HFShower::getHits().

114  {
115  double zFibre = zShift(point, depth, fromEndAbs);
116  double time = zFibre / cFibre_;
117 #ifdef EDM_ML_DEBUG
118  edm::LogVerbatim("HFShower") << "HFFibre::tShift for point " << point << " ( depth = " << depth
119  << ", traversed length = " << zFibre / CLHEP::cm << " cm) = " << time / CLHEP::ns
120  << " ns";
121 #endif
122  return time;
123 }
Log< level::Info, true > LogVerbatim
double zShift(const G4ThreeVector &point, int depth, int fromEndAbs=0) const
Definition: HFFibre.cc:125
double cFibre_
Definition: HFFibre.h:45
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5

◆ zShift()

double HFFibre::zShift ( const G4ThreeVector &  point,
int  depth,
int  fromEndAbs = 0 
) const

Definition at line 125 of file HFFibre.cc.

References LEDCalibrationChannels::depth, gpar_, mps_fire::i, LEDCalibrationChannels::ieta, longFL_, nBinR_, point, radius_, shortFL_, mathSSE::sqrt(), and geometryCSVtoXML::zz.

Referenced by HFShowerLibrary::fillHits(), HFShowerParam::getHits(), and tShift().

125  { // point is z-local
126 
127  double zFibre = 0;
128  int ieta = 0;
129  double length = 250 * CLHEP::cm;
130  double hR = sqrt((point.x()) * (point.x()) + (point.y()) * (point.y()));
131 
132  // Defines the Radius bin by radial subdivision
133  if (fromEndAbs >= 0) {
134  for (int i = nBinR_ - 1; i > 0; --i)
135  if (hR < radius_[i])
136  ieta = nBinR_ - i - 1;
137  }
138 
139  // Defines the full length of the fibre
140  if (depth == 2) {
141  if (static_cast<int>(shortFL_.size()) > ieta)
142  length = shortFL_[ieta] + gpar_[0];
143  } else {
144  if (static_cast<int>(longFL_.size()) > ieta)
145  length = longFL_[ieta];
146  }
147  zFibre = length; // from beginning of abs (full length)
148 
149  if (fromEndAbs > 0) {
150  zFibre -= gpar_[1]; // length from end of HF
151  } else {
152  double zz = 0.5 * gpar_[1] + point.z(); // depth of point of photon emission (from beginning of HF)
153  zFibre -= zz; // length of fiber from point of photon emission
154  }
155 
156 #ifdef EDM_ML_DEBUG
157  edm::LogVerbatim("HFShower") << "HFFibre::zShift for point " << point << " (R = " << hR / CLHEP::cm
158  << " cm, Index = " << ieta << ", depth = " << depth
159  << ", Fibre Length = " << length / CLHEP::cm << " cm = " << zFibre / CLHEP::cm << " cm)";
160 #endif
161  return zFibre;
162 }
std::vector< double > radius_
Definition: HFFibre.h:46
std::vector< double > shortFL_
Definition: HFFibre.h:47
Log< level::Info, true > LogVerbatim
std::vector< double > longFL_
Definition: HFFibre.h:47
int nBinR_
Definition: HFFibre.h:49
T sqrt(T t)
Definition: SSEVec.h:19
std::vector< double > gpar_
Definition: HFFibre.h:46
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5

Member Data Documentation

◆ attL_

std::vector<double> HFFibre::attL_
private

Definition at line 48 of file HFFibre.h.

Referenced by attLength(), and HFFibre().

◆ cFibre_

double HFFibre::cFibre_
private

Definition at line 45 of file HFFibre.h.

Referenced by HFFibre(), and tShift().

◆ gpar_

std::vector<double> HFFibre::gpar_
private

Definition at line 46 of file HFFibre.h.

Referenced by HFFibre(), and zShift().

◆ lambLim_

std::array<double, 2> HFFibre::lambLim_
private

Definition at line 50 of file HFFibre.h.

Referenced by attLength(), and HFFibre().

◆ longFL_

std::vector<double> HFFibre::longFL_
private

Definition at line 47 of file HFFibre.h.

Referenced by HFFibre(), and zShift().

◆ nBinAtt_

int HFFibre::nBinAtt_
private

Definition at line 49 of file HFFibre.h.

Referenced by attLength(), and HFFibre().

◆ nBinR_

int HFFibre::nBinR_
private

Definition at line 49 of file HFFibre.h.

Referenced by HFFibre(), and zShift().

◆ radius_

std::vector<double> HFFibre::radius_
private

Definition at line 46 of file HFFibre.h.

Referenced by HFFibre(), and zShift().

◆ shortFL_

std::vector<double> HFFibre::shortFL_
private

Definition at line 47 of file HFFibre.h.

Referenced by HFFibre(), and zShift().