File indexing completed on 2024-04-06 12:14:47
0001 #include "Geometry/HcalCommonData/interface/HcalParametersFromDD.h"
0002 #include "Geometry/HcalCommonData/interface/HcalGeomParameters.h"
0003 #include "Geometry/HcalCommonData/interface/HcalTopologyMode.h"
0004 #include "CondFormats/GeometryObjects/interface/HcalParameters.h"
0005 #include "DataFormats/Math/interface/GeantUnits.h"
0006 #include "DetectorDescription/Core/interface/DDFilteredView.h"
0007 #include "DetectorDescription/Core/interface/DDutils.h"
0008
0009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0010 #include <iostream>
0011 #include <iomanip>
0012
0013
0014 using namespace geant_units::operators;
0015
0016 namespace {
0017 int getTopologyMode(const char* s, const DDsvalues_type& sv, bool type) {
0018 DDValue val(s);
0019 if (DDfetch(&sv, val)) {
0020 const std::vector<std::string>& fvec = val.strings();
0021 if (fvec.empty()) {
0022 throw cms::Exception("HcalParametersFromDD") << "Failed to get " << s << " tag.";
0023 }
0024
0025 int result(-1);
0026 if (type) {
0027 StringToEnumParser<HcalTopologyMode::Mode> eparser;
0028 HcalTopologyMode::Mode mode = (HcalTopologyMode::Mode)eparser.parseString(fvec[0]);
0029 result = static_cast<int>(mode);
0030 } else {
0031 StringToEnumParser<HcalTopologyMode::TriggerMode> eparser;
0032 HcalTopologyMode::TriggerMode mode = (HcalTopologyMode::TriggerMode)eparser.parseString(fvec[0]);
0033 result = static_cast<int>(mode);
0034 }
0035 return result;
0036 } else {
0037 throw cms::Exception("HcalParametersFromDD") << "Failed to get " << s << " tag.";
0038 }
0039 }
0040
0041 int getTopologyMode(const std::string& s, bool type) {
0042 int result(-1);
0043 if (type) {
0044 StringToEnumParser<HcalTopologyMode::Mode> eparser;
0045 HcalTopologyMode::Mode mode = (HcalTopologyMode::Mode)eparser.parseString(s);
0046 result = static_cast<int>(mode);
0047 } else {
0048 StringToEnumParser<HcalTopologyMode::TriggerMode> eparser;
0049 HcalTopologyMode::TriggerMode mode = (HcalTopologyMode::TriggerMode)eparser.parseString(s);
0050 result = static_cast<int>(mode);
0051 }
0052 return result;
0053 }
0054 }
0055
0056 bool HcalParametersFromDD::build(const DDCompactView* cpv, HcalParameters& php) {
0057 #ifdef EDM_ML_DEBUG
0058 edm::LogVerbatim("HCalGeom") << "HcalParametersFromDD::build(const DDCompactView*, HcalParameters&) is called";
0059 #endif
0060
0061 std::string attribute = "OnlyForHcalSimNumbering";
0062 DDSpecificsHasNamedValueFilter filter1{attribute};
0063 DDFilteredView fv1(*cpv, filter1);
0064 bool ok = fv1.firstChild();
0065
0066 const int nEtaMax = 100;
0067
0068 if (ok) {
0069 std::unique_ptr<HcalGeomParameters> geom = std::make_unique<HcalGeomParameters>();
0070 geom->loadGeometry(fv1, php);
0071 php.modHB = geom->getModHalfHBHE(0);
0072 php.modHE = geom->getModHalfHBHE(1);
0073 php.dzVcal = geom->getConstDzHF();
0074 geom->getConstRHO(php.rHO);
0075
0076 php.phioff = cpv->vector("phioff");
0077 php.etaTable = cpv->vector("etaTable");
0078 php.rTable = cpv->vector("rTable");
0079 rescale(php.rTable, HcalGeomParameters::k_ScaleFromDDDToG4);
0080 php.phibin = cpv->vector("phibin");
0081 php.phitable = cpv->vector("phitable");
0082 for (unsigned int i = 1; i <= nEtaMax; ++i) {
0083 std::stringstream sstm;
0084 sstm << "layerGroupSimEta" << i;
0085 std::string tempName = sstm.str();
0086 auto const& v = cpv->vector(tempName);
0087 if (!v.empty()) {
0088 HcalParameters::LayerItem layerGroupEta;
0089 layerGroupEta.layer = i;
0090 layerGroupEta.layerGroup = dbl_to_int(v);
0091 php.layerGroupEtaSim.emplace_back(layerGroupEta);
0092 }
0093 }
0094 php.etaMin = dbl_to_int(cpv->vector("etaMin"));
0095 php.etaMax = dbl_to_int(cpv->vector("etaMax"));
0096 php.etaRange = cpv->vector("etaRange");
0097 php.gparHF = cpv->vector("gparHF");
0098 rescale(php.gparHF, HcalGeomParameters::k_ScaleFromDDDToG4);
0099 php.noff = dbl_to_int(cpv->vector("noff"));
0100 php.Layer0Wt = cpv->vector("Layer0Wt");
0101 php.HBGains = cpv->vector("HBGains");
0102 php.HBShift = dbl_to_int(cpv->vector("HBShift"));
0103 php.HEGains = cpv->vector("HEGains");
0104 php.HEShift = dbl_to_int(cpv->vector("HEShift"));
0105 php.HFGains = cpv->vector("HFGains");
0106 php.HFShift = dbl_to_int(cpv->vector("HFShift"));
0107 php.maxDepth = dbl_to_int(cpv->vector("MaxDepth"));
0108 } else {
0109 throw cms::Exception("HcalParametersFromDD") << "Not found " << attribute.c_str() << " but needed.";
0110 }
0111
0112 attribute = "OnlyForHcalRecNumbering";
0113 DDSpecificsHasNamedValueFilter filter2{attribute};
0114 DDFilteredView fv2(*cpv, filter2);
0115 ok = fv2.firstChild();
0116 if (ok) {
0117 DDsvalues_type sv(fv2.mergedSpecifics());
0118 int topoMode = getTopologyMode("TopologyMode", sv, true);
0119 int trigMode = getTopologyMode("TriggerMode", sv, false);
0120 php.topologyMode = ((trigMode & 0xFF) << 8) | (topoMode & 0xFF);
0121 php.etagroup = dbl_to_int(cpv->vector("etagroup"));
0122 php.phigroup = dbl_to_int(cpv->vector("phigroup"));
0123 for (unsigned int i = 1; i <= nEtaMax; ++i) {
0124 std::stringstream sstm;
0125 sstm << "layerGroupRecEta" << i;
0126 std::string tempName = sstm.str();
0127 auto const& v = cpv->vector(tempName);
0128 if (!v.empty()) {
0129 HcalParameters::LayerItem layerGroupEta;
0130 layerGroupEta.layer = i;
0131 layerGroupEta.layerGroup = dbl_to_int(v);
0132 php.layerGroupEtaRec.emplace_back(layerGroupEta);
0133 }
0134 }
0135 } else {
0136 throw cms::Exception("HcalParametersFromDD") << "Not found " << attribute.c_str() << " but needed.";
0137 }
0138
0139 return build(php);
0140 }
0141
0142 bool HcalParametersFromDD::build(const cms::DDCompactView& cpv, HcalParameters& php) {
0143 #ifdef EDM_ML_DEBUG
0144 edm::LogVerbatim("HCalGeom") << "HcalParametersFromDD::build(const cms::DDCompactView*, HcalParameters&) is called";
0145 #endif
0146 const int nEtaMax = 100;
0147
0148 std::unique_ptr<HcalGeomParameters> geom = std::make_unique<HcalGeomParameters>();
0149 geom->loadGeometry(cpv, php);
0150 php.modHB = geom->getModHalfHBHE(0);
0151 php.modHE = geom->getModHalfHBHE(1);
0152 php.dzVcal = geom->getConstDzHF();
0153 geom->getConstRHO(php.rHO);
0154
0155 php.phioff = cpv.get<std::vector<double>>("phioff");
0156 php.etaTable = cpv.get<std::vector<double>>("etaTable");
0157 php.rTable = cpv.get<std::vector<double>>("rTable");
0158 php.phibin = cpv.get<std::vector<double>>("phibin");
0159 php.phitable = cpv.get<std::vector<double>>("phitable");
0160 php.etaMin = cpv.getVector<int>("etaMin");
0161 php.etaMax = cpv.getVector<int>("etaMax");
0162 php.etaRange = cpv.get<std::vector<double>>("etaRange");
0163 php.gparHF = cpv.get<std::vector<double>>("gparHF");
0164 php.noff = cpv.getVector<int>("noff");
0165 php.Layer0Wt = cpv.get<std::vector<double>>("Layer0Wt");
0166 php.HBGains = cpv.get<std::vector<double>>("HBGains");
0167 php.HBShift = cpv.getVector<int>("HBShift");
0168 php.HEGains = cpv.get<std::vector<double>>("HEGains");
0169 php.HEShift = cpv.getVector<int>("HEShift");
0170 php.HFGains = cpv.get<std::vector<double>>("HFGains");
0171 php.HFShift = cpv.getVector<int>("HFShift");
0172 php.maxDepth = cpv.getVector<int>("MaxDepth");
0173
0174 rescale(php.rTable, HcalGeomParameters::k_ScaleFromDD4hepToG4);
0175 rescale(php.gparHF, HcalGeomParameters::k_ScaleFromDD4hepToG4);
0176 for (unsigned int i = 1; i <= nEtaMax; ++i) {
0177 HcalParameters::LayerItem layerGroupEta;
0178 layerGroupEta.layer = i;
0179 layerGroupEta.layerGroup = cpv.getVector<int>(std::string("layerGroupSimEta") + std::to_string(i));
0180 if (!layerGroupEta.layerGroup.empty()) {
0181 php.layerGroupEtaSim.emplace_back(layerGroupEta);
0182 }
0183 }
0184
0185
0186 const cms::DDFilter filter("OnlyForHcalRecNumbering", "HCAL");
0187 cms::DDFilteredView fv(cpv, filter);
0188
0189 std::vector<std::string> tempS = fv.get<std::vector<std::string>>("hcal", "TopologyMode");
0190 std::string sv = (!tempS.empty()) ? tempS[0] : "HcalTopologyMode::SLHC";
0191 int topoMode = getTopologyMode(sv, true);
0192 tempS = fv.get<std::vector<std::string>>("hcal", "TriggerMode");
0193 sv = (!tempS.empty()) ? tempS[0] : "HcalTopologyMode::TriggerMode_2021";
0194 int trigMode = getTopologyMode(sv, false);
0195 php.topologyMode = ((trigMode & 0xFF) << 8) | (topoMode & 0xFF);
0196
0197 php.etagroup = cpv.getVector<int>("etagroup");
0198 php.phigroup = cpv.getVector<int>("phigroup");
0199
0200 for (unsigned int i = 1; i <= nEtaMax; ++i) {
0201 HcalParameters::LayerItem layerGroupEta;
0202 layerGroupEta.layer = i;
0203 layerGroupEta.layerGroup = cpv.getVector<int>(std::string("layerGroupRecEta") + std::to_string(i));
0204 if (!layerGroupEta.layerGroup.empty()) {
0205 php.layerGroupEtaRec.emplace_back(layerGroupEta);
0206 }
0207 }
0208
0209 return build(php);
0210 }
0211
0212 bool HcalParametersFromDD::build(HcalParameters& php) {
0213 php.etaMin[0] = 1;
0214 if (php.etaMax[1] >= php.etaMin[1])
0215 php.etaMax[1] = static_cast<int>(php.etaTable.size()) - 1;
0216 php.etaMax[2] = php.etaMin[2] + static_cast<int>(php.rTable.size()) - 2;
0217
0218 for (unsigned int i = 0; i < php.rTable.size(); ++i) {
0219 unsigned int k = php.rTable.size() - i - 1;
0220 php.etaTableHF.emplace_back(-log(tan(0.5 * atan(php.rTable[k] / php.gparHF[4]))));
0221 }
0222
0223 #ifdef EDM_ML_DEBUG
0224 int i(0);
0225 std::stringstream ss0;
0226 ss0 << "HcalParametersFromDD: MaxDepth[" << php.maxDepth.size() << "]: ";
0227 for (const auto& it : php.maxDepth)
0228 ss0 << it << ", ";
0229 edm::LogVerbatim("HCalGeom") << ss0.str();
0230 std::stringstream ss1;
0231 ss1 << "HcalParametersFromDD: ModHB [" << php.modHB.size() << "]: ";
0232 for (const auto& it : php.modHB)
0233 ss1 << it << ", ";
0234 edm::LogVerbatim("HCalGeom") << ss1.str();
0235 std::stringstream ss2;
0236 ss2 << "HcalParametersFromDD: ModHE [" << php.modHE.size() << "]: ";
0237 for (const auto& it : php.modHE)
0238 ss2 << it << ", ";
0239 edm::LogVerbatim("HCalGeom") << ss2.str();
0240 std::stringstream ss3;
0241 ss3 << "HcalParametersFromDD: " << php.phioff.size() << " phioff values:";
0242 std::vector<double>::const_iterator it;
0243 for (it = php.phioff.begin(), i = 0; it != php.phioff.end(); ++it)
0244 ss3 << " [" << ++i << "] = " << convertRadToDeg(*it);
0245 edm::LogVerbatim("HCalGeom") << ss3.str();
0246 std::stringstream ss4;
0247 ss4 << "HcalParametersFromDD: " << php.etaTable.size() << " entries for etaTable:";
0248 for (it = php.etaTable.begin(), i = 0; it != php.etaTable.end(); ++it)
0249 ss4 << " [" << ++i << "] = " << (*it);
0250 edm::LogVerbatim("HCalGeom") << ss4.str();
0251 std::stringstream ss5;
0252 ss5 << "HcalParametersFromDD: " << php.rTable.size() << " entries for rTable:";
0253 for (it = php.rTable.begin(), i = 0; it != php.rTable.end(); ++it)
0254 ss5 << " [" << ++i << "] = " << convertMmToCm(*it);
0255 edm::LogVerbatim("HCalGeom") << ss5.str();
0256 std::stringstream ss6;
0257 ss6 << "HcalParametersFromDD: " << php.phibin.size() << " phibin values:";
0258 for (it = php.phibin.begin(), i = 0; it != php.phibin.end(); ++it)
0259 ss6 << " [" << ++i << "] = " << convertRadToDeg(*it);
0260 edm::LogVerbatim("HCalGeom") << ss6.str();
0261 std::stringstream ss7;
0262 ss7 << "HcalParametersFromDD: " << php.phitable.size() << " phitable values:";
0263 for (it = php.phitable.begin(), i = 0; it != php.phitable.end(); ++it)
0264 ss7 << " [" << ++i << "] = " << convertRadToDeg(*it);
0265 edm::LogVerbatim("HCalGeom") << ss7.str();
0266 edm::LogVerbatim("HCalGeom") << "HcalParametersFromDD: " << php.layerGroupEtaSim.size() << " layerGroupEtaSim blocks"
0267 << std::endl;
0268 std::vector<int>::const_iterator kt;
0269 for (unsigned int k = 0; k < php.layerGroupEtaSim.size(); ++k) {
0270 std::stringstream ss8;
0271 ss8 << "layerGroupEtaSim[" << k << "] Layer " << php.layerGroupEtaSim[k].layer;
0272 for (kt = php.layerGroupEtaSim[k].layerGroup.begin(), i = 0; kt != php.layerGroupEtaSim[k].layerGroup.end(); ++kt)
0273 ss8 << " " << ++i << ":" << (*kt);
0274 edm::LogVerbatim("HCalGeom") << ss8.str();
0275 }
0276 std::stringstream ss8;
0277 ss8 << "HcalParametersFromDD: " << php.etaMin.size() << " etaMin values:";
0278 for (kt = php.etaMin.begin(), i = 0; kt != php.etaMin.end(); ++kt)
0279 ss8 << " [" << ++i << "] = " << (*kt);
0280 edm::LogVerbatim("HCalGeom") << ss8.str();
0281 std::stringstream ss9;
0282 ss9 << "HcalParametersFromDD: " << php.etaMax.size() << " etaMax values:";
0283 for (kt = php.etaMax.begin(), i = 0; kt != php.etaMax.end(); ++kt)
0284 ss9 << " [" << ++i << "] = " << (*kt);
0285 edm::LogVerbatim("HCalGeom") << ss9.str();
0286 std::stringstream ss10;
0287 ss10 << "HcalParametersFromDD: " << php.etaRange.size() << " etaRange values:";
0288 for (it = php.etaRange.begin(), i = 0; it != php.etaRange.end(); ++it)
0289 ss10 << " [" << ++i << "] = " << (*it);
0290 edm::LogVerbatim("HCalGeom") << ss10.str();
0291 std::stringstream ss11;
0292 ss11 << "HcalParametersFromDD: " << php.gparHF.size() << " gparHF values:";
0293 for (it = php.gparHF.begin(), i = 0; it != php.gparHF.end(); ++it)
0294 ss11 << " [" << ++i << "] = " << convertMmToCm(*it);
0295 edm::LogVerbatim("HCalGeom") << ss11.str();
0296 std::stringstream ss12;
0297 ss12 << "HcalParametersFromDD: " << php.noff.size() << " noff values:";
0298 for (kt = php.noff.begin(), i = 0; kt != php.noff.end(); ++kt)
0299 ss12 << " [" << ++i << "] = " << (*kt);
0300 edm::LogVerbatim("HCalGeom") << ss12.str();
0301 std::stringstream ss13;
0302 ss13 << "HcalParametersFromDD: " << php.Layer0Wt.size() << " Layer0Wt values:";
0303 for (it = php.Layer0Wt.begin(), i = 0; it != php.Layer0Wt.end(); ++it)
0304 ss13 << " [" << ++i << "] = " << (*it);
0305 edm::LogVerbatim("HCalGeom") << ss13.str();
0306 std::stringstream ss14;
0307 ss14 << "HcalParametersFromDD: " << php.HBGains.size() << " Shift/Gains values for HB:";
0308 for (unsigned k = 0; k < php.HBGains.size(); ++k)
0309 ss14 << " [" << k << "] = " << php.HBShift[k] << ":" << php.HBGains[k];
0310 edm::LogVerbatim("HCalGeom") << ss14.str();
0311 std::stringstream ss15;
0312 ss15 << "HcalParametersFromDD: " << php.HEGains.size() << " Shift/Gains values for HE:";
0313 for (unsigned k = 0; k < php.HEGains.size(); ++k)
0314 ss15 << " [" << k << "] = " << php.HEShift[k] << ":" << php.HEGains[k];
0315 edm::LogVerbatim("HCalGeom") << ss15.str();
0316 std::stringstream ss16;
0317 ss16 << "HcalParametersFromDD: " << php.HFGains.size() << " Shift/Gains values for HF:";
0318 for (unsigned k = 0; k < php.HFGains.size(); ++k)
0319 ss16 << " [" << k << "] = " << php.HFShift[k] << ":" << php.HFGains[k];
0320 edm::LogVerbatim("HCalGeom") << ss16.str();
0321 std::stringstream ss17;
0322 ss17 << "HcalParametersFromDD: " << php.etagroup.size() << " etagroup values:";
0323 for (kt = php.etagroup.begin(), i = 0; kt != php.etagroup.end(); ++kt)
0324 ss17 << " [" << ++i << "] = " << (*kt);
0325 edm::LogVerbatim("HCalGeom") << ss17.str();
0326 std::stringstream ss18;
0327 ss18 << "HcalParametersFromDD: " << php.phigroup.size() << " phigroup values:";
0328 for (kt = php.phigroup.begin(), i = 0; kt != php.phigroup.end(); ++kt)
0329 ss18 << " [" << ++i << "] = " << (*kt);
0330 edm::LogVerbatim("HCalGeom") << ss18.str();
0331 edm::LogVerbatim("HCalGeom") << "HcalParametersFromDD: " << php.layerGroupEtaRec.size() << " layerGroupEtaRec blocks";
0332 for (unsigned int k = 0; k < php.layerGroupEtaRec.size(); ++k) {
0333 std::stringstream ss19;
0334 ss19 << "layerGroupEtaRec[" << k << "] Layer " << php.layerGroupEtaRec[k].layer;
0335 for (kt = php.layerGroupEtaRec[k].layerGroup.begin(), i = 0; kt != php.layerGroupEtaRec[k].layerGroup.end(); ++kt)
0336 ss19 << " " << ++i << ":" << (*kt);
0337 edm::LogVerbatim("HCalGeom") << ss19.str();
0338 }
0339 edm::LogVerbatim("HCalGeom") << "HcalParametersFromDD: (topology|trigger)Mode " << std::hex << php.topologyMode
0340 << std::dec;
0341 #endif
0342
0343 return true;
0344 }
0345
0346 void HcalParametersFromDD::rescale(std::vector<double>& v, const double s) {
0347 std::for_each(v.begin(), v.end(), [s](double& n) { n *= s; });
0348 }