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.

21  {
22 
23  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFShower");
24  cFibre = c_light*(m_HF.getParameter<double>("CFibre"));
25 
26  edm::LogInfo("HFShower") << "HFFibre:: Speed of light in fibre " << cFibre
27  << " m/ns";
28 
29  std::string attribute = "Volume";
30  std::string value = "HF";
31  DDSpecificsMatchesValueFilter filter1{DDValue(attribute,value,0)};
32  DDFilteredView fv1(cpv,filter1);
33  bool dodet = fv1.firstChild();
34 
35  if (dodet) {
36  DDsvalues_type sv(fv1.mergedSpecifics());
37 
38  // Attenuation length
39  nBinAtt = -1;
40  attL = getDDDArray("attl",sv,nBinAtt);
41  std::stringstream ss1;
42  for (int it=0; it<nBinAtt; it++) {
43  if(it/10*10 == it) { ss1 << "\n"; }
44  ss1 << " " << attL[it]*cm;
45  }
46  edm::LogInfo("HFShower") << "HFFibre: " << nBinAtt << " attL(1/cm): " << ss1.str();
47 
48  // Limits on Lambda
49  int nb = 2;
50  std::vector<double> nvec = getDDDArray("lambLim",sv,nb);
51  lambLim[0] = static_cast<int>(nvec[0]);
52  lambLim[1] = static_cast<int>(nvec[1]);
53  edm::LogInfo("HFShower") << "HFFibre: Limits on lambda " << lambLim[0]
54  << " and " << lambLim[1];
55 
56  // Fibre Lengths
57  nb = 0;
58  longFL = getDDDArray("LongFL",sv,nb);
59  std::stringstream ss2;
60  for (int it=0; it<nb; it++) {
61  if(it/10*10 == it) { ss2 << "\n"; }
62  ss2 << " " << longFL[it]/cm;
63  }
64  edm::LogInfo("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) { ss3 << "\n"; }
70  ss3 << " " << shortFL[it]/cm;
71  }
72  edm::LogInfo("HFShower") << "HFFibre: " << nb << " Short Fibre Length(cm):" << ss3.str();
73  } else {
74  edm::LogError("HFShower") << "HFFibre: cannot get filtered "
75  << " view for " << attribute << " matching "
76  << name;
77  throw cms::Exception("Unknown", "HFFibre")
78  << "cannot match " << attribute << " to " << name <<"\n";
79  }
80 }
T getParameter(std::string const &) const
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &, int &)
Definition: HFFibre.cc:164
double lambLim[2]
Definition: HFFibre.h:48
int nBinAtt
Definition: HFFibre.h:47
double cFibre
Definition: HFFibre.h:43
std::vector< std::pair< unsigned int, DDValue > > DDsvalues_type
std::maps an index to a DDValue. The index corresponds to the index assigned to the name of the std::...
Definition: DDsvalues.h:20
std::vector< double > shortFL
Definition: HFFibre.h:45
Definition: value.py:1
std::vector< double > attL
Definition: HFFibre.h:46
std::vector< double > longFL
Definition: HFFibre.h:45
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 
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 DebugLog
108  edm::LogInfo("HFShower") << "HFFibre::attLength for Lambda " << lambda
109  << " index " << i << " " << j << " Att. Length "
110  << att;
111 #endif
112  return att;
113 }
double lambLim[2]
Definition: HFFibre.h:48
int nBinAtt
Definition: HFFibre.h:47
std::vector< double > attL
Definition: HFFibre.h:46
std::vector< double > HFFibre::getDDDArray ( const std::string &  str,
const DDsvalues_type sv,
int &  nmin 
)
protected

Definition at line 164 of file HFFibre.cc.

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

Referenced by HFFibre().

166  {
167 
168 #ifdef DebugLog
169  LogDebug("HFShower") << "HFFibre:getDDDArray called for " << str
170  << " with nMin " << nmin;
171 #endif
172  DDValue value(str);
173  if (DDfetch(&sv,value)) {
174 #ifdef DebugLog
175  LogDebug("HFShower") << value;
176 #endif
177  const std::vector<double> & fvec = value.doubles();
178  int nval = fvec.size();
179  if (nmin > 0) {
180  if (nval < nmin) {
181  edm::LogError("HFShower") << "HFFibre : # of " << str << " bins "
182  << nval << " < " << nmin << " ==> illegal";
183  throw cms::Exception("Unknown", "HFFibre")
184  << "nval < nmin for array " << str <<"\n";
185  }
186  } else {
187  if (nval < 1 && nmin != 0) {
188  edm::LogError("HFShower") << "HFFibre : # of " << str << " bins "
189  << nval << " < 1 ==> illegal (nmin="
190  << nmin << ")";
191  throw cms::Exception("Unknown", "HFFibre")
192  << "nval < 1 for array " << str <<"\n";
193  }
194  }
195  nmin = nval;
196  return fvec;
197  } else {
198  if (nmin != 0) {
199  edm::LogError("HFShower") << "HFFibre : cannot get array " << str;
200  throw cms::Exception("Unknown", "HFFibre")
201  << "cannot get array " << str <<"\n";
202  } else {
203  std::vector<double> fvec;
204  return fvec;
205  }
206  }
207 }
#define LogDebug(id)
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 82 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().

82  {
83 
84  // Now geometry parameters
85  gpar = hcons->getGparHF();
86  radius = hcons->getRTableHF();
87 
88  nBinR = (int)(radius.size());
89  std::stringstream sss;
90  for (int i=0; i<nBinR; ++i) {
91  if(i/10*10 == i) { sss << "\n"; }
92  sss << " " << radius[i]/cm;
93  }
94  edm::LogInfo("HFShower") << "HFFibre: " << radius.size() <<" rTable(cm):" << sss.str();
95 }
std::vector< double > gpar
Definition: HFFibre.h:44
std::vector< double > radius
Definition: HFFibre.h:44
int nBinR
Definition: HFFibre.h:47
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 115 of file HFFibre.cc.

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

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

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

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

Member Data Documentation

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

Definition at line 46 of file HFFibre.h.

Referenced by attLength(), and HFFibre().

double HFFibre::cFibre
private

Definition at line 43 of file HFFibre.h.

Referenced by HFFibre(), and tShift().

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

Definition at line 44 of file HFFibre.h.

Referenced by initRun(), and zShift().

double HFFibre::lambLim[2]
private

Definition at line 48 of file HFFibre.h.

Referenced by attLength(), and HFFibre().

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

Definition at line 45 of file HFFibre.h.

Referenced by HFFibre(), and zShift().

int HFFibre::nBinAtt
private

Definition at line 47 of file HFFibre.h.

Referenced by attLength(), and HFFibre().

int HFFibre::nBinR
private

Definition at line 47 of file HFFibre.h.

Referenced by initRun(), and zShift().

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

Definition at line 44 of file HFFibre.h.

Referenced by initRun(), and zShift().

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

Definition at line 45 of file HFFibre.h.

Referenced by HFFibre(), and zShift().