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 <iostream>
15 
16 //#define DebugLog
17 
19  edm::ParameterSet const & p) {
20 
21  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFShower");
22  cFibre = c_light*(m_HF.getParameter<double>("CFibre"));
23 
24  edm::LogInfo("HFShower") << "HFFibre:: Speed of light in fibre " << cFibre
25  << " m/ns";
26 
27  std::string attribute = "Volume";
28  std::string value = "HF";
29  DDSpecificsFilter filter1;
30  DDValue ddv1(attribute,value,0);
32  DDFilteredView fv1(cpv);
33  fv1.addFilter(filter1);
34  bool dodet = fv1.firstChild();
35 
36  if (dodet) {
38 
39  // Attenuation length
40  nBinAtt = -1;
41  attL = getDDDArray("attl",sv,nBinAtt);
42  edm::LogInfo("HFShower") << "HFFibre: " << nBinAtt << " attL ";
43  for (int it=0; it<nBinAtt; it++)
44  edm::LogInfo("HFShower") << "HFFibre: attL[" << it << "] = "
45  << attL[it]*cm << "(1/cm)";
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::LogInfo("HFShower") << "HFFibre: Limits on lambda " << lambLim[0]
53  << " and " << lambLim[1];
54 
55  // Fibre Lengths
56  nb = 0;
57  longFL = getDDDArray("LongFL",sv,nb);
58  edm::LogInfo("HFShower") << "HFFibre: " << nb << " Long Fibre Length";
59  for (int it=0; it<nb; it++)
60  edm::LogInfo("HFShower") << "HFFibre: longFL[" << it << "] = "
61  << longFL[it]/cm << " cm";
62  nb = 0;
63  shortFL = getDDDArray("ShortFL",sv,nb);
64  edm::LogInfo("HFShower") << "HFFibre: " << nb << " Short Fibre Length";
65  for (int it=0; it<nb; it++)
66  edm::LogInfo("HFShower") << "HFFibre: shortFL[" << it << "] = "
67  << shortFL[it]/cm << " cm";
68  } else {
69  edm::LogError("HFShower") << "HFFibre: cannot get filtered "
70  << " view for " << attribute << " matching "
71  << name;
72  throw cms::Exception("Unknown", "HFFibre")
73  << "cannot match " << attribute << " to " << name <<"\n";
74  }
75 
76  // Now geometry parameters
77  attribute = "ReadOutName";
78  value = name;
79  DDSpecificsFilter filter2;
80  DDValue ddv2(attribute,value,0);
82  DDFilteredView fv2(cpv);
83  fv2.addFilter(filter2);
84  dodet = fv2.firstChild();
85  if (dodet) {
87 
88  //Special Geometry parameters
89  int nb = -1;
90  gpar = getDDDArray("gparHF",sv,nb);
91  edm::LogInfo("HFShower") << "HFFibre: " << nb <<" gpar (cm)";
92  for (int i=0; i<nb; i++)
93  edm::LogInfo("HFShower") << "HFFibre: gpar[" << i << "] = "
94  << gpar[i]/cm << " cm";
95 
96  nBinR = -1;
97  radius = getDDDArray("rTable",sv,nBinR);
98  edm::LogInfo("HFShower") << "HFFibre: " << nBinR <<" rTable (cm)";
99  for (int i=0; i<nBinR; i++)
100  edm::LogInfo("HFShower") << "HFFibre: radius[" << i << "] = "
101  << radius[i]/cm << " cm";
102  } else {
103  edm::LogError("HFShower") << "HFFibre: cannot get filtered "
104  << " view for " << attribute << " matching "
105  << name;
106  throw cms::Exception("Unknown", "HFFibre")
107  << "cannot match " << attribute << " to " << name <<"\n";
108  }
109 }
110 
112 
113 double HFFibre::attLength(double lambda) {
114 
115  int i = int(nBinAtt*(lambda - lambLim[0])/(lambLim[1]-lambLim[0]));
116 
117  int j =i;
118  if (i >= nBinAtt)
119  j = nBinAtt-1;
120  else if (i < 0)
121  j = 0;
122  double att = attL[j];
123 #ifdef DebugLog
124  edm::LogInfo("HFShower") << "HFFibre::attLength for Lambda " << lambda
125  << " index " << i << " " << j << " Att. Length "
126  << att;
127 #endif
128  return att;
129 }
130 
131 double HFFibre::tShift(G4ThreeVector point, int depth, int fromEndAbs) {
132 
133  double zFibre = zShift(point, depth, fromEndAbs);
134  double time = zFibre/cFibre;
135 #ifdef DebugLog
136  edm::LogInfo("HFShower") << "HFFibre::tShift for point " << point
137  << " ( depth = " << depth <<", traversed length = "
138  << zFibre/cm << " cm) = " << time/ns << " ns";
139 #endif
140  return time;
141 }
142 
143 double HFFibre::zShift(G4ThreeVector point, int depth, int fromEndAbs) { // point is z-local
144 
145  double zFibre = 0;
146  double hR = sqrt((point.x())*(point.x())+(point.y())*(point.y()));
147  int ieta = 0;
148  double length = 250*cm;
149  if (fromEndAbs < 0) {
150  zFibre = 0.5*gpar[1] - point.z(); // Never, as fromEndAbs=0 (?)
151  } else {
152  // Defines the Radius bin by radial subdivision
153  for (int i = nBinR-1; i > 0; --i) if (hR < radius[i]) ieta = nBinR - i - 1;
154  // define the length of the fibre
155  if (depth == 2) {
156  if ((int)(shortFL.size()) > ieta) length = shortFL[ieta];
157  } else {
158  if ((int)(longFL.size()) > ieta) length = longFL[ieta];
159  }
160  zFibre = length;
161  if (fromEndAbs > 0) {
162  zFibre -= gpar[1]; // Never, as fromEndAbs=0 (M.K. ?)
163  } else {
164  double zz = 0.5*gpar[1] + point.z();
165  zFibre -= zz;
166  }
167  if (depth == 2) zFibre += gpar[0]; // here zFibre is reduced for Short
168  }
169 
170 #ifdef DebugLog
171  edm::LogInfo("HFShower") << "HFFibre::zShift for point " << point
172  << " (R = " << hR/cm << " cm, Index = " << ieta
173  << ", depth = " << depth << ", Fibre Length = "
174  << length/cm << " cm = " << zFibre/cm
175  << " cm)";
176 #endif
177  return zFibre;
178 }
179 
180 std::vector<double> HFFibre::getDDDArray(const std::string & str,
181  const DDsvalues_type & sv,
182  int & nmin) {
183 
184 #ifdef DebugLog
185  LogDebug("HFShower") << "HFFibre:getDDDArray called for " << str
186  << " with nMin " << nmin;
187 #endif
188  DDValue value(str);
189  if (DDfetch(&sv,value)) {
190 #ifdef DebugLog
191  LogDebug("HFShower") << value;
192 #endif
193  const std::vector<double> & fvec = value.doubles();
194  int nval = fvec.size();
195  if (nmin > 0) {
196  if (nval < nmin) {
197  edm::LogError("HFShower") << "HFFibre : # of " << str << " bins "
198  << nval << " < " << nmin << " ==> illegal";
199  throw cms::Exception("Unknown", "HFFibre")
200  << "nval < nmin for array " << str <<"\n";
201  }
202  } else {
203  if (nval < 1 && nmin != 0) {
204  edm::LogError("HFShower") << "HFFibre : # of " << str << " bins "
205  << nval << " < 1 ==> illegal (nmin="
206  << nmin << ")";
207  throw cms::Exception("Unknown", "HFFibre")
208  << "nval < 1 for array " << str <<"\n";
209  }
210  }
211  nmin = nval;
212  return fvec;
213  } else {
214  if (nmin != 0) {
215  edm::LogError("HFShower") << "HFFibre : cannot get array " << str;
216  throw cms::Exception("Unknown", "HFFibre")
217  << "cannot get array " << str <<"\n";
218  } else {
219  std::vector<double> fvec;
220  return fvec;
221  }
222  }
223 }
#define LogDebug(id)
std::vector< double > gpar
Definition: HFFibre.h:42
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:180
double lambLim[2]
Definition: HFFibre.h:46
const std::vector< double > & doubles() const
a reference to the double-valued values stored in the given instance of DDValue
Definition: DDValue.cc:134
void addFilter(const DDFilter &, log_op op=AND)
~HFFibre()
Definition: HFFibre.cc:111
double attLength(double lambda)
Definition: HFFibre.cc:113
int nBinAtt
Definition: HFFibre.h:45
std::vector< double > radius
Definition: HFFibre.h:42
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:102
double zShift(G4ThreeVector point, int depth, int fromEndAbs=0)
Definition: HFFibre.cc:143
double cFibre
Definition: HFFibre.h:41
int nBinR
Definition: HFFibre.h:45
T sqrt(T t)
Definition: SSEVec.h:48
double tShift(G4ThreeVector point, int depth, int fromEndAbs=0)
Definition: HFFibre.cc:131
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:18
std::vector< double > shortFL
Definition: HFFibre.h:43
std::vector< double > attL
Definition: HFFibre.h:44
DDsvalues_type mergedSpecifics() const
bool firstChild()
set the current node to the first child ...
void setCriteria(const DDValue &nameVal, comp_op, log_op l=AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:285
*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:43
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:37