CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HFFibre.cc
Go to the documentation of this file.
1 // File: HFFibre.cc
3 // Description: Loads the table for attenuation length and calculates it
5 
8 
9 #include "CLHEP/Units/GlobalSystemOfUnits.h"
10 #include "CLHEP/Units/GlobalPhysicalConstants.h"
11 #include <iostream>
12 #include <sstream>
13 
14 //#define EDM_ML_DEBUG
15 
17  const HcalDDDSimConstants* hcons,
18  const HcalSimulationParameters* hps,
19  edm::ParameterSet const& p)
20  : hcalConstant_(hcons), hcalsimpar_(hps) {
21  edm::ParameterSet m_HF =
22  (p.getParameter<edm::ParameterSet>("HFShower")).getParameter<edm::ParameterSet>("HFShowerBlock");
23  cFibre = c_light * (m_HF.getParameter<double>("CFibre"));
24  edm::LogVerbatim("HFShower") << "HFFibre:: Speed of light in fibre " << cFibre / (CLHEP::m / CLHEP::ns) << " m/ns";
25  // Attenuation length
27  nBinAtt = static_cast<int>(attL.size());
28  std::stringstream ss1;
29  for (int it = 0; it < nBinAtt; it++) {
30  if (it / 10 * 10 == it) {
31  ss1 << "\n";
32  }
33  ss1 << " " << attL[it] * CLHEP::cm;
34  }
35  edm::LogVerbatim("HFShower") << "HFFibre: " << nBinAtt << " attL(1/cm): " << ss1.str();
36  // Limits on Lambda
37  std::vector<int> nvec = hcalsimpar_->lambdaLimits_;
38  lambLim[0] = nvec[0];
39  lambLim[1] = nvec[1];
40  edm::LogVerbatim("HFShower") << "HFFibre: Limits on lambda " << lambLim[0] << " and " << lambLim[1];
41  // Fibre Lengths
43  std::stringstream ss2;
44  for (unsigned int it = 0; it < longFL.size(); it++) {
45  if (it / 10 * 10 == it) {
46  ss2 << "\n";
47  }
48  ss2 << " " << longFL[it] / CLHEP::cm;
49  }
50  edm::LogVerbatim("HFShower") << "HFFibre: " << longFL.size() << " Long Fibre Length(cm):" << ss2.str();
52  std::stringstream ss3;
53  for (unsigned int it = 0; it < shortFL.size(); it++) {
54  if (it / 10 * 10 == it) {
55  ss3 << "\n";
56  }
57  ss3 << " " << shortFL[it] / CLHEP::cm;
58  }
59  edm::LogVerbatim("HFShower") << "HFFibre: " << shortFL.size() << " Short Fibre Length(cm):" << ss3.str();
60 
61  // Now geometry parameters
64 
65  nBinR = static_cast<int>(radius.size());
66  std::stringstream sss;
67  for (int i = 0; i < nBinR; ++i) {
68  if (i / 10 * 10 == i) {
69  sss << "\n";
70  }
71  sss << " " << radius[i] / CLHEP::cm;
72  }
73  edm::LogVerbatim("HFShower") << "HFFibre: " << radius.size() << " rTable(cm):" << sss.str();
74 }
75 
76 double HFFibre::attLength(double lambda) {
77  int i = int(nBinAtt * (lambda - lambLim[0]) / (lambLim[1] - lambLim[0]));
78 
79  int j = i;
80  if (i >= nBinAtt)
81  j = nBinAtt - 1;
82  else if (i < 0)
83  j = 0;
84  double att = attL[j];
85 #ifdef EDM_ML_DEBUG
86  edm::LogVerbatim("HFShower") << "HFFibre::attLength for Lambda " << lambda << " index " << i << " " << j
87  << " Att. Length " << att;
88 #endif
89  return att;
90 }
91 
92 double HFFibre::tShift(const G4ThreeVector& point, int depth, int fromEndAbs) {
93  double zFibre = zShift(point, depth, fromEndAbs);
94  double time = zFibre / cFibre;
95 #ifdef EDM_ML_DEBUG
96  edm::LogVerbatim("HFShower") << "HFFibre::tShift for point " << point << " ( depth = " << depth
97  << ", traversed length = " << zFibre / CLHEP::cm << " cm) = " << time / CLHEP::ns
98  << " ns";
99 #endif
100  return time;
101 }
102 
103 double HFFibre::zShift(const G4ThreeVector& point, int depth, int fromEndAbs) { // point is z-local
104 
105  double zFibre = 0;
106  int ieta = 0;
107  double length = 250 * CLHEP::cm;
108  double hR = sqrt((point.x()) * (point.x()) + (point.y()) * (point.y()));
109 
110  // Defines the Radius bin by radial subdivision
111  if (fromEndAbs >= 0) {
112  for (int i = nBinR - 1; i > 0; --i)
113  if (hR < radius[i])
114  ieta = nBinR - i - 1;
115  }
116 
117  // Defines the full length of the fibre
118  if (depth == 2) {
119  if (static_cast<int>(shortFL.size()) > ieta)
120  length = shortFL[ieta] + gpar[0];
121  } else {
122  if (static_cast<int>(longFL.size()) > ieta)
123  length = longFL[ieta];
124  }
125  zFibre = length; // from beginning of abs (full length)
126 
127  if (fromEndAbs > 0) {
128  zFibre -= gpar[1]; // length from end of HF
129  } else {
130  double zz = 0.5 * gpar[1] + point.z(); // depth of point of photon emission (from beginning of HF)
131  zFibre -= zz; // length of fiber from point of photon emission
132  }
133 
134 #ifdef EDM_ML_DEBUG
135  edm::LogVerbatim("HFShower") << "HFFibre::zShift for point " << point << " (R = " << hR / CLHEP::cm
136  << " cm, Index = " << ieta << ", depth = " << depth
137  << ", Fibre Length = " << length / CLHEP::cm << " cm = " << zFibre / CLHEP::cm << " cm)";
138 #endif
139  return zFibre;
140 }
double tShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
Definition: HFFibre.cc:92
std::vector< double > gpar
Definition: HFFibre.h:35
Log< level::Info, true > LogVerbatim
double lambLim[2]
Definition: HFFibre.h:39
HFFibre(const std::string &name, const HcalDDDSimConstants *hcons, const HcalSimulationParameters *hps, edm::ParameterSet const &p)
Definition: HFFibre.cc:16
const HcalDDDSimConstants * hcalConstant_
Definition: HFFibre.h:32
std::vector< double > shortFiberLength_
std::vector< double > longFiberLength_
double attLength(double lambda)
Definition: HFFibre.cc:76
const HcalSimulationParameters * hcalsimpar_
Definition: HFFibre.h:33
int nBinAtt
Definition: HFFibre.h:38
double zShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
Definition: HFFibre.cc:103
std::vector< double > radius
Definition: HFFibre.h:35
double cFibre
Definition: HFFibre.h:34
int nBinR
Definition: HFFibre.h:38
const std::vector< double > & getRTableHF() const
T sqrt(T t)
Definition: SSEVec.h:19
std::vector< double > shortFL
Definition: HFFibre.h:36
std::vector< double > attL
Definition: HFFibre.h:37
const std::vector< double > & getGparHF() const
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< double > attenuationLength_
*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
std::vector< double > longFL
Definition: HFFibre.h:36