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 (std::string &name, const DDCompactView &cpv, edm::ParameterSet const &p)
 
void initRun (HcalDDDSimConstants *)
 
double tShift (const G4ThreeVector &point, int depth, int fromEndAbs=0)
 
double zShift (const G4ThreeVector &point, int depth, int fromEndAbs=0)
 
 ~HFFibre ()
 

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 ( std::string &  name,
const DDCompactView cpv,
edm::ParameterSet const &  p 
)

Definition at line 19 of file HFFibre.cc.

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

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

Definition at line 75 of file HFFibre.cc.

75 {}

Member Function Documentation

double HFFibre::attLength ( double  lambda)

Definition at line 94 of file HFFibre.cc.

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

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

94  {
95 
96  int i = int(nBinAtt*(lambda - lambLim[0])/(lambLim[1]-lambLim[0]));
97 
98  int j =i;
99  if (i >= nBinAtt)
100  j = nBinAtt-1;
101  else if (i < 0)
102  j = 0;
103  double att = attL[j];
104 #ifdef DebugLog
105  edm::LogInfo("HFShower") << "HFFibre::attLength for Lambda " << lambda
106  << " index " << i << " " << j << " Att. Length "
107  << att;
108 #endif
109  return att;
110 }
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 161 of file HFFibre.cc.

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

Referenced by HFFibre().

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

Definition at line 77 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().

77  {
78 
79  // Now geometry parameters
80  gpar = hcons->getGparHF();
81  edm::LogInfo("HFShower") << "HFFibre: " << gpar.size() <<" gpar (cm)";
82  for (unsigned int i=0; i<gpar.size(); i++)
83  edm::LogInfo("HFShower") << "HFFibre: gpar[" << i << "] = "
84  << gpar[i]/cm << " cm";
85 
86  radius = hcons->getRTableHF();
87  nBinR = (int)(radius.size());
88  edm::LogInfo("HFShower") << "HFFibre: " << radius.size() <<" rTable (cm)";
89  for (unsigned int i=0; i<radius.size(); i++)
90  edm::LogInfo("HFShower") << "HFFibre: radius[" << i << "] = "
91  << radius[i]/cm << " cm";
92 }
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 112 of file HFFibre.cc.

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

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

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

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