00001
00002
00003
00005
00006 #include "SimG4CMS/Calo/interface/HFFibre.h"
00007
00008 #include "DetectorDescription/Core/interface/DDFilter.h"
00009 #include "DetectorDescription/Core/interface/DDFilteredView.h"
00010 #include "DetectorDescription/Core/interface/DDValue.h"
00011 #include "FWCore/Utilities/interface/Exception.h"
00012
00013 #include "CLHEP/Units/SystemOfUnits.h"
00014 #include <iostream>
00015
00016 HFFibre::HFFibre(std::string & name, const DDCompactView & cpv,
00017 edm::ParameterSet const & p) {
00018
00019 edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFShower");
00020 cFibre = c_light*(m_HF.getParameter<double>("CFibre"));
00021
00022 edm::LogInfo("HFShower") << "HFFibre:: Speed of light in fibre " << cFibre
00023 << " m/ns";
00024
00025 std::string attribute = "Volume";
00026 std::string value = "HF";
00027 DDSpecificsFilter filter1;
00028 DDValue ddv1(attribute,value,0);
00029 filter1.setCriteria(ddv1,DDSpecificsFilter::equals);
00030 DDFilteredView fv1(cpv);
00031 fv1.addFilter(filter1);
00032 bool dodet = fv1.firstChild();
00033
00034 if (dodet) {
00035 DDsvalues_type sv(fv1.mergedSpecifics());
00036
00037
00038 nBinAtt = -1;
00039 attL = getDDDArray("attl",sv,nBinAtt);
00040 edm::LogInfo("HFShower") << "HFFibre: " << nBinAtt << " attL ";
00041 for (int it=0; it<nBinAtt; it++)
00042 edm::LogInfo("HFShower") << "HFFibre: attL[" << it << "] = "
00043 << attL[it]*cm << "(1/cm)";
00044
00045
00046 int nb = 2;
00047 std::vector<double> nvec = getDDDArray("lambLim",sv,nb);
00048 lambLim[0] = static_cast<int>(nvec[0]);
00049 lambLim[1] = static_cast<int>(nvec[1]);
00050 edm::LogInfo("HFShower") << "HFFibre: Limits on lambda " << lambLim[0]
00051 << " and " << lambLim[1];
00052
00053
00054 nb = 0;
00055 longFL = getDDDArray("LongFL",sv,nb);
00056 edm::LogInfo("HFShower") << "HFFibre: " << nb << " Long Fibre Length";
00057 for (int it=0; it<nb; it++)
00058 edm::LogInfo("HFShower") << "HFFibre: longFL[" << it << "] = "
00059 << longFL[it]/cm << " cm";
00060 nb = 0;
00061 shortFL = getDDDArray("ShortFL",sv,nb);
00062 edm::LogInfo("HFShower") << "HFFibre: " << nb << " Short Fibre Length";
00063 for (int it=0; it<nb; it++)
00064 edm::LogInfo("HFShower") << "HFFibre: shortFL[" << it << "] = "
00065 << shortFL[it]/cm << " cm";
00066 } else {
00067 edm::LogError("HFShower") << "HFFibre: cannot get filtered "
00068 << " view for " << attribute << " matching "
00069 << name;
00070 throw cms::Exception("Unknown", "HFFibre")
00071 << "cannot match " << attribute << " to " << name <<"\n";
00072 }
00073
00074
00075 attribute = "ReadOutName";
00076 value = name;
00077 DDSpecificsFilter filter2;
00078 DDValue ddv2(attribute,value,0);
00079 filter2.setCriteria(ddv2,DDSpecificsFilter::equals);
00080 DDFilteredView fv2(cpv);
00081 fv2.addFilter(filter2);
00082 dodet = fv2.firstChild();
00083 if (dodet) {
00084 DDsvalues_type sv(fv2.mergedSpecifics());
00085
00086
00087 int nb = -1;
00088 gpar = getDDDArray("gparHF",sv,nb);
00089 edm::LogInfo("HFShower") << "HFFibre: " << nb <<" gpar (cm)";
00090 for (int i=0; i<nb; i++)
00091 edm::LogInfo("HFShower") << "HFFibre: gpar[" << i << "] = "
00092 << gpar[i]/cm << " cm";
00093
00094 nBinR = -1;
00095 radius = getDDDArray("rTable",sv,nBinR);
00096 edm::LogInfo("HFShower") << "HFFibre: " << nBinR <<" rTable (cm)";
00097 for (int i=0; i<nBinR; i++)
00098 edm::LogInfo("HFShower") << "HFFibre: radius[" << i << "] = "
00099 << radius[i]/cm << " cm";
00100 } else {
00101 edm::LogError("HFShower") << "HFFibre: cannot get filtered "
00102 << " view for " << attribute << " matching "
00103 << name;
00104 throw cms::Exception("Unknown", "HFFibre")
00105 << "cannot match " << attribute << " to " << name <<"\n";
00106 }
00107 }
00108
00109 HFFibre::~HFFibre() {}
00110
00111 double HFFibre::attLength(double lambda) {
00112
00113 int i = int(nBinAtt*(lambda - lambLim[0])/(lambLim[1]-lambLim[0]));
00114
00115 int j =i;
00116 if (i >= nBinAtt)
00117 j = nBinAtt-1;
00118 else if (i < 0)
00119 j = 0;
00120 double att = attL[j];
00121 LogDebug("HFShower") << "HFFibre::attLength for Lambda " << lambda
00122 << " index " << i << " " << j << " Att. Length "
00123 << att;
00124
00125 return att;
00126 }
00127
00128 double HFFibre::tShift(G4ThreeVector point, int depth, bool fromEndAbs) {
00129
00130 double hR = sqrt((point.x())*(point.x())+(point.y())*(point.y()));
00131 int ieta = 0;
00132 for (int i = nBinR-1; i > 0; i--)
00133 if (hR < radius[i]) ieta = nBinR - i - 1;
00134 double length = 250*cm;
00135 if (depth == 2) {
00136 if ((int)(shortFL.size()) > ieta) length = shortFL[ieta];
00137 } else {
00138 if ((int)(longFL.size()) > ieta) length = longFL[ieta];
00139 }
00140 double zFibre = length;
00141 if (fromEndAbs)
00142 zFibre -= gpar[1];
00143 else {
00144 double zint = point.z();
00145 double zz = std::abs(zint) - gpar[4];
00146 zFibre -= zz;
00147 }
00148 if (depth == 2) zFibre += gpar[0];
00149 double time = zFibre/cFibre;
00150 LogDebug("HFShower") << "HFFibre::tShift for point " << point
00151 << " (R = " << hR/cm << " cm, Index = " << ieta
00152 << ", depth = " << depth << ", Fibre Length = "
00153 << length/cm << " cm, traversed length = "
00154 << zFibre/cm << " cm) = " << time/ns << " ns";
00155 return time;
00156 }
00157
00158 std::vector<double> HFFibre::getDDDArray(const std::string & str,
00159 const DDsvalues_type & sv,
00160 int & nmin) {
00161
00162 LogDebug("HFShower") << "HFFibre:getDDDArray called for " << str
00163 << " with nMin " << nmin;
00164 DDValue value(str);
00165 if (DDfetch(&sv,value)) {
00166 LogDebug("HFShower") << value;
00167 const std::vector<double> & fvec = value.doubles();
00168 int nval = fvec.size();
00169 if (nmin > 0) {
00170 if (nval < nmin) {
00171 edm::LogError("HFShower") << "HFFibre : # of " << str << " bins "
00172 << nval << " < " << nmin << " ==> illegal";
00173 throw cms::Exception("Unknown", "HFFibre")
00174 << "nval < nmin for array " << str <<"\n";
00175 }
00176 } else {
00177 if (nval < 1 && nmin != 0) {
00178 edm::LogError("HFShower") << "HFFibre : # of " << str << " bins "
00179 << nval << " < 1 ==> illegal (nmin="
00180 << nmin << ")";
00181 throw cms::Exception("Unknown", "HFFibre")
00182 << "nval < 1 for array " << str <<"\n";
00183 }
00184 }
00185 nmin = nval;
00186 return fvec;
00187 } else {
00188 if (nmin != 0) {
00189 edm::LogError("HFShower") << "HFFibre : cannot get array " << str;
00190 throw cms::Exception("Unknown", "HFFibre")
00191 << "cannot get array " << str <<"\n";
00192 } else {
00193 std::vector<double> fvec;
00194 return fvec;
00195 }
00196 }
00197 }