CMS 3D CMS Logo

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

#include <HFFibre.h>

Public Member Functions

double attLength (double lambda)
 
 HFFibre (const std::string &name, const DDCompactView &cpv, edm::ParameterSet const &p)
 
void initRun (const HcalDDDSimConstants *)
 
double tShift (const G4ThreeVector &point, int depth, int fromEndAbs=0)
 
double zShift (const G4ThreeVector &point, int depth, int fromEndAbs=0)
 
 ~HFFibre ()=default
 

Protected Member Functions

std::vector< double > getDDDArray (const std::string &, const DDsvalues_type &, int &)
 

Private Attributes

std::vector< double > attL
 
double cFibre
 
std::vector< double > gpar
 
double lambLim [2]
 
std::vector< double > longFL
 
int nBinAtt
 
int nBinR
 
std::vector< double > radius
 
std::vector< double > shortFL
 

Detailed Description

Definition at line 20 of file HFFibre.h.

Constructor & Destructor Documentation

HFFibre::HFFibre ( const std::string &  name,
const DDCompactView cpv,
edm::ParameterSet const &  p 
)

Definition at line 20 of file HFFibre.cc.

References attL, cFibre, Exception, DDFilteredView::firstChild(), getDDDArray(), edm::ParameterSet::getParameter(), lambLim, longFL, DDFilteredView::mergedSpecifics(), dataset::name, nBinAtt, shortFL, AlCaHLTBitMon_QueryRunRegistry::string, and pfDeepBoostedJetPreprocessParams_cfi::sv.

20  {
21  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFShower");
22  cFibre = c_light * (m_HF.getParameter<double>("CFibre"));
23 
24  edm::LogVerbatim("HFShower") << "HFFibre:: Speed of light in fibre " << cFibre << " m/ns";
25 
26  std::string attribute = "Volume";
27  std::string value = "HF";
28  DDSpecificsMatchesValueFilter filter1{DDValue(attribute, value, 0)};
29  DDFilteredView fv1(cpv, filter1);
30  bool dodet = fv1.firstChild();
31 
32  if (dodet) {
33  DDsvalues_type sv(fv1.mergedSpecifics());
34 
35  // Attenuation length
36  nBinAtt = -1;
37  attL = getDDDArray("attl", sv, nBinAtt);
38  std::stringstream ss1;
39  for (int it = 0; it < nBinAtt; it++) {
40  if (it / 10 * 10 == it) {
41  ss1 << "\n";
42  }
43  ss1 << " " << attL[it] * cm;
44  }
45  edm::LogVerbatim("HFShower") << "HFFibre: " << nBinAtt << " attL(1/cm): " << ss1.str();
46 
47  // Limits on Lambda
48  int nb = 2;
49  std::vector<double> nvec = getDDDArray("lambLim", sv, nb);
50  lambLim[0] = static_cast<int>(nvec[0]);
51  lambLim[1] = static_cast<int>(nvec[1]);
52  edm::LogVerbatim("HFShower") << "HFFibre: Limits on lambda " << lambLim[0] << " and " << lambLim[1];
53 
54  // Fibre Lengths
55  nb = 0;
56  longFL = getDDDArray("LongFL", sv, nb);
57  std::stringstream ss2;
58  for (int it = 0; it < nb; it++) {
59  if (it / 10 * 10 == it) {
60  ss2 << "\n";
61  }
62  ss2 << " " << longFL[it] / cm;
63  }
64  edm::LogVerbatim("HFShower") << "HFFibre: " << nb << " Long Fibre Length(cm):" << ss2.str();
65  nb = 0;
66  shortFL = getDDDArray("ShortFL", sv, nb);
67  std::stringstream ss3;
68  for (int it = 0; it < nb; it++) {
69  if (it / 10 * 10 == it) {
70  ss3 << "\n";
71  }
72  ss3 << " " << shortFL[it] / cm;
73  }
74  edm::LogVerbatim("HFShower") << "HFFibre: " << nb << " Short Fibre Length(cm):" << ss3.str();
75  } else {
76  edm::LogError("HFShower") << "HFFibre: cannot get filtered view for " << attribute << " matching " << name;
77  throw cms::Exception("Unknown", "HFFibre") << "cannot match " << attribute << " to " << name << "\n";
78  }
79 }
T getParameter(std::string const &) const
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &, int &)
Definition: HFFibre.cc:163
double lambLim[2]
Definition: HFFibre.h:40
int nBinAtt
Definition: HFFibre.h:39
double cFibre
Definition: HFFibre.h:35
std::vector< double > shortFL
Definition: HFFibre.h:37
Definition: value.py:1
std::vector< double > attL
Definition: HFFibre.h:38
std::vector< std::pair< unsigned int, DDValue > > DDsvalues_type
Definition: DDsvalues.h:12
std::vector< double > longFL
Definition: HFFibre.h:37
HFFibre::~HFFibre ( )
default

Member Function Documentation

double HFFibre::attLength ( double  lambda)

Definition at line 97 of file HFFibre.cc.

References attL, mps_fire::i, createfilelist::int, lambLim, and nBinAtt.

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

97  {
98  int i = int(nBinAtt * (lambda - lambLim[0]) / (lambLim[1] - lambLim[0]));
99 
100  int j = i;
101  if (i >= nBinAtt)
102  j = nBinAtt - 1;
103  else if (i < 0)
104  j = 0;
105  double att = attL[j];
106 #ifdef EDM_ML_DEBUG
107  edm::LogVerbatim("HFShower") << "HFFibre::attLength for Lambda " << lambda << " index " << i << " " << j
108  << " Att. Length " << att;
109 #endif
110  return att;
111 }
double lambLim[2]
Definition: HFFibre.h:40
int nBinAtt
Definition: HFFibre.h:39
std::vector< double > attL
Definition: HFFibre.h:38
std::vector< double > HFFibre::getDDDArray ( const std::string &  str,
const DDsvalues_type sv,
int &  nmin 
)
protected

Definition at line 163 of file HFFibre.cc.

References DDfetch(), DDValue::doubles(), Exception, str, and relativeConstraints::value.

Referenced by HFFibre().

163  {
164 #ifdef EDM_ML_DEBUG
165  edm::LogVerbatim("HFShower") << "HFFibre:getDDDArray called for " << str << " with nMin " << nmin;
166 #endif
167  DDValue value(str);
168  if (DDfetch(&sv, value)) {
169 #ifdef EDM_ML_DEBUG
170  edm::LogVerbatim("HFShower") << value;
171 #endif
172  const std::vector<double>& fvec = value.doubles();
173  int nval = fvec.size();
174  if (nmin > 0) {
175  if (nval < nmin) {
176  edm::LogError("HFShower") << "HFFibre : # of " << str << " bins " << nval << " < " << nmin << " ==> illegal";
177  throw cms::Exception("Unknown", "HFFibre") << "nval < nmin for array " << str << "\n";
178  }
179  } else {
180  if (nval < 1 && nmin != 0) {
181  edm::LogError("HFShower") << "HFFibre : # of " << str << " bins " << nval << " < 1 ==> illegal (nmin=" << nmin
182  << ")";
183  throw cms::Exception("Unknown", "HFFibre") << "nval < 1 for array " << str << "\n";
184  }
185  }
186  nmin = nval;
187  return fvec;
188  } else {
189  if (nmin != 0) {
190  edm::LogError("HFShower") << "HFFibre : cannot get array " << str;
191  throw cms::Exception("Unknown", "HFFibre") << "cannot get array " << str << "\n";
192  } else {
193  std::vector<double> fvec;
194  return fvec;
195  }
196  }
197 }
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:81
Definition: value.py:1
#define str(s)
void HFFibre::initRun ( const HcalDDDSimConstants hcons)

Definition at line 81 of file HFFibre.cc.

References HcalDDDSimConstants::getGparHF(), HcalDDDSimConstants::getRTableHF(), gpar, mps_fire::i, createfilelist::int, nBinR, and radius.

Referenced by HFShower::initRun(), HFShowerParam::initRun(), and HFShowerLibrary::initRun().

81  {
82  // Now geometry parameters
83  gpar = hcons->getGparHF();
84  radius = hcons->getRTableHF();
85 
86  nBinR = (int)(radius.size());
87  std::stringstream sss;
88  for (int i = 0; i < nBinR; ++i) {
89  if (i / 10 * 10 == i) {
90  sss << "\n";
91  }
92  sss << " " << radius[i] / cm;
93  }
94  edm::LogVerbatim("HFShower") << "HFFibre: " << radius.size() << " rTable(cm):" << sss.str();
95 }
std::vector< double > gpar
Definition: HFFibre.h:36
std::vector< double > radius
Definition: HFFibre.h:36
int nBinR
Definition: HFFibre.h:39
const std::vector< double > & getRTableHF() const
const std::vector< double > & getGparHF() const
double HFFibre::tShift ( const G4ThreeVector &  point,
int  depth,
int  fromEndAbs = 0 
)

Definition at line 113 of file HFFibre.cc.

References cFibre, protons_cff::time, and zShift().

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

113  {
114  double zFibre = zShift(point, depth, fromEndAbs);
115  double time = zFibre / cFibre;
116 #ifdef EDM_ML_DEBUG
117  edm::LogVerbatim("HFShower") << "HFFibre::tShift for point " << point << " ( depth = " << depth
118  << ", traversed length = " << zFibre / cm << " cm) = " << time / ns << " ns";
119 #endif
120  return time;
121 }
double zShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
Definition: HFFibre.cc:123
double cFibre
Definition: HFFibre.h:35
*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
double HFFibre::zShift ( const G4ThreeVector &  point,
int  depth,
int  fromEndAbs = 0 
)

Definition at line 123 of file HFFibre.cc.

References gpar, mps_fire::i, longFL, nBinR, radius, shortFL, mathSSE::sqrt(), and geometryCSVtoXML::zz.

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

123  { // point is z-local
124 
125  double zFibre = 0;
126  double hR = sqrt((point.x()) * (point.x()) + (point.y()) * (point.y()));
127  int ieta = 0;
128  double length = 250 * cm;
129  if (fromEndAbs < 0) {
130  zFibre = 0.5 * gpar[1] - point.z(); // Never, as fromEndAbs=0 (?)
131  } else {
132  // Defines the Radius bin by radial subdivision
133  for (int i = nBinR - 1; i > 0; --i)
134  if (hR < radius[i])
135  ieta = nBinR - i - 1;
136  // define the length of the fibre
137  if (depth == 2) {
138  if ((int)(shortFL.size()) > ieta)
139  length = shortFL[ieta];
140  } else {
141  if ((int)(longFL.size()) > ieta)
142  length = longFL[ieta];
143  }
144  zFibre = length;
145  if (fromEndAbs > 0) {
146  zFibre -= gpar[1]; // Never, as fromEndAbs=0 (M.K. ?)
147  } else {
148  double zz = 0.5 * gpar[1] + point.z();
149  zFibre -= zz;
150  }
151  if (depth == 2)
152  zFibre += gpar[0]; // here zFibre is reduced for Short
153  }
154 
155 #ifdef EDM_ML_DEBUG
156  edm::LogVerbatim("HFShower") << "HFFibre::zShift for point " << point << " (R = " << hR / cm
157  << " cm, Index = " << ieta << ", depth = " << depth << ", Fibre Length = " << length / cm
158  << " cm = " << zFibre / cm << " cm)";
159 #endif
160  return zFibre;
161 }
std::vector< double > gpar
Definition: HFFibre.h:36
std::vector< double > radius
Definition: HFFibre.h:36
int nBinR
Definition: HFFibre.h:39
T sqrt(T t)
Definition: SSEVec.h:18
std::vector< double > shortFL
Definition: HFFibre.h:37
*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:37

Member Data Documentation

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

Definition at line 38 of file HFFibre.h.

Referenced by attLength(), and HFFibre().

double HFFibre::cFibre
private

Definition at line 35 of file HFFibre.h.

Referenced by HFFibre(), and tShift().

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

Definition at line 36 of file HFFibre.h.

Referenced by initRun(), and zShift().

double HFFibre::lambLim[2]
private

Definition at line 40 of file HFFibre.h.

Referenced by attLength(), and HFFibre().

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

Definition at line 37 of file HFFibre.h.

Referenced by HFFibre(), and zShift().

int HFFibre::nBinAtt
private

Definition at line 39 of file HFFibre.h.

Referenced by attLength(), and HFFibre().

int HFFibre::nBinR
private

Definition at line 39 of file HFFibre.h.

Referenced by initRun(), and zShift().

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

Definition at line 36 of file HFFibre.h.

Referenced by initRun(), and zShift().

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

Definition at line 37 of file HFFibre.h.

Referenced by HFFibre(), and zShift().