CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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 
7 
12 
13 #include "CLHEP/Units/GlobalSystemOfUnits.h"
14 #include "CLHEP/Units/GlobalPhysicalConstants.h"
15 #include <iostream>
16 
17 //#define DebugLog
18 
20  edm::ParameterSet const & p) {
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  DDSpecificsFilter filter1;
31  DDValue ddv1(attribute,value,0);
32  filter1.setCriteria(ddv1,DDCompOp::equals);
33  DDFilteredView fv1(cpv);
34  fv1.addFilter(filter1);
35  bool dodet = fv1.firstChild();
36 
37  if (dodet) {
39 
40  // Attenuation length
41  nBinAtt = -1;
42  attL = getDDDArray("attl",sv,nBinAtt);
43  edm::LogInfo("HFShower") << "HFFibre: " << nBinAtt << " attL ";
44  for (int it=0; it<nBinAtt; it++)
45  edm::LogInfo("HFShower") << "HFFibre: attL[" << it << "] = "
46  << attL[it]*cm << "(1/cm)";
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  edm::LogInfo("HFShower") << "HFFibre: " << nb << " Long Fibre Length";
60  for (int it=0; it<nb; it++)
61  edm::LogInfo("HFShower") << "HFFibre: longFL[" << it << "] = "
62  << longFL[it]/cm << " cm";
63  nb = 0;
64  shortFL = getDDDArray("ShortFL",sv,nb);
65  edm::LogInfo("HFShower") << "HFFibre: " << nb << " Short Fibre Length";
66  for (int it=0; it<nb; it++)
67  edm::LogInfo("HFShower") << "HFFibre: shortFL[" << it << "] = "
68  << shortFL[it]/cm << " cm";
69  } else {
70  edm::LogError("HFShower") << "HFFibre: cannot get filtered "
71  << " view for " << attribute << " matching "
72  << name;
73  throw cms::Exception("Unknown", "HFFibre")
74  << "cannot match " << attribute << " to " << name <<"\n";
75  }
76 }
77 
79 
81 
82  // Now geometry parameters
83  gpar = hcons->getGparHF();
84  edm::LogInfo("HFShower") << "HFFibre: " << gpar.size() <<" gpar (cm)";
85  for (unsigned int i=0; i<gpar.size(); i++)
86  edm::LogInfo("HFShower") << "HFFibre: gpar[" << i << "] = "
87  << gpar[i]/cm << " cm";
88 
89  radius = hcons->getRTableHF();
90  nBinR = (int)(radius.size());
91  edm::LogInfo("HFShower") << "HFFibre: " << radius.size() <<" rTable (cm)";
92  for (unsigned int i=0; i<radius.size(); i++)
93  edm::LogInfo("HFShower") << "HFFibre: radius[" << i << "] = "
94  << radius[i]/cm << " cm";
95 }
96 
97 double HFFibre::attLength(double lambda) {
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 }
114 
115 double HFFibre::tShift(const G4ThreeVector& point, int depth, int fromEndAbs) {
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 }
126 
127 double HFFibre::zShift(const G4ThreeVector& point, int depth, int fromEndAbs) { // 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 }
163 
164 std::vector<double> HFFibre::getDDDArray(const std::string & str,
165  const DDsvalues_type & sv,
166  int & nmin) {
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)
double tShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
Definition: HFFibre.cc:115
std::vector< double > gpar
Definition: HFFibre.h:44
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &, int &)
Definition: HFFibre.cc:164
double lambLim[2]
Definition: HFFibre.h:48
const std::vector< double > & doubles() const
a reference to the double-valued values stored in the given instance of DDValue
Definition: DDValue.cc:139
void initRun(HcalDDDSimConstants *)
Definition: HFFibre.cc:80
void addFilter(const DDFilter &, DDLogOp op=DDLogOp::AND)
~HFFibre()
Definition: HFFibre.cc:78
double attLength(double lambda)
Definition: HFFibre.cc:97
int nBinAtt
Definition: HFFibre.h:47
double zShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
Definition: HFFibre.cc:127
std::vector< double > radius
Definition: HFFibre.h:44
type of data representation of DDCompactView
Definition: DDCompactView.h:77
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:80
double cFibre
Definition: HFFibre.h:43
int nBinR
Definition: HFFibre.h:47
const std::vector< double > & getRTableHF() const
T sqrt(T t)
Definition: SSEVec.h:18
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:19
int j
Definition: DBlmapReader.cc:9
HFFibre(std::string &name, const DDCompactView &cpv, edm::ParameterSet const &p)
Definition: HFFibre.cc:19
std::vector< double > shortFL
Definition: HFFibre.h:45
std::vector< double > attL
Definition: HFFibre.h:46
const std::vector< double > & getGparHF() const
DDsvalues_type mergedSpecifics() const
bool firstChild()
set the current node to the first child ...
void setCriteria(const DDValue &nameVal, DDCompOp, DDLogOp l=DDLogOp::AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:245
*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
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:32