Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 13:02:55

0001 #include "CondFormats/GeometryObjects/interface/EcalSimulationParameters.h"
0002 #include "DetectorDescription/Core/interface/DDFilteredView.h"
0003 #include "DetectorDescription/Core/interface/DDFilter.h"
0004 #include "DetectorDescription/Core/interface/DDValue.h"
0005 #include "DetectorDescription/Core/interface/DDutils.h"
0006 #include "DetectorDescription/DDCMS/interface/DDFilteredView.h"
0007 #include "Geometry/EcalCommonData/interface/EcalSimParametersFromDD.h"
0008 #include <iostream>
0009 #include <iomanip>
0010 
0011 //#define EDM_ML_DEBUG
0012 
0013 template <typename T>
0014 void myPrint(std::string value, const std::vector<T>& vec) {
0015   edm::LogVerbatim("EcalGeom") << "EcalSimParametersFromDD: " << vec.size() << " entries for " << value << ":";
0016   unsigned int i(0);
0017   for (const auto& e : vec) {
0018     edm::LogVerbatim("EcalGeom") << " (" << i << ") " << e;
0019     ++i;
0020   }
0021 }
0022 
0023 bool EcalSimParametersFromDD::build(const DDCompactView* cpv, const std::string& name, EcalSimulationParameters& php) {
0024 #ifdef EDM_ML_DEBUG
0025   edm::LogVerbatim("EcalGeom")
0026       << "Inside EcalSimParametersFromDD::build(const DDCompactView*, const std::string&, EcalSimulationParameters&)";
0027 #endif
0028   // Get the filtered view
0029   std::string attribute = "ReadOutName";
0030   DDSpecificsMatchesValueFilter filter{DDValue(attribute, name, 0)};
0031   DDFilteredView fv(*cpv, filter);
0032   bool dodet = fv.firstChild();
0033   DDsvalues_type sv(fv.mergedSpecifics());
0034 
0035   //First the specpars
0036   php.useWeight_ = true;
0037   std::vector<double> tempD = getDDDArray("EnergyWeight", sv);
0038   if (!tempD.empty()) {
0039     if (tempD[0] < 0.1)
0040       php.useWeight_ = false;
0041   }
0042   tempD = getDDDArray("nxtalEta", sv);
0043   if (tempD.empty())
0044     php.nxtalEta_ = 0;
0045   else
0046     php.nxtalEta_ = static_cast<int>(tempD[0]);
0047   tempD = getDDDArray("nxtalPhi", sv);
0048   if (tempD.empty())
0049     php.nxtalPhi_ = 0;
0050   else
0051     php.nxtalPhi_ = static_cast<int>(tempD[0]);
0052   tempD = getDDDArray("PhiBaskets", sv);
0053   if (tempD.empty())
0054     php.phiBaskets_ = 0;
0055   else
0056     php.phiBaskets_ = static_cast<int>(tempD[0]);
0057   php.etaBaskets_ = dbl_to_int(getDDDArray("EtaBaskets", sv));
0058   tempD = getDDDArray("ncrys", sv);
0059   if (tempD.empty())
0060     php.ncrys_ = 0;
0061   else
0062     php.ncrys_ = static_cast<int>(tempD[0]);
0063   tempD = getDDDArray("nmods", sv);
0064   if (tempD.empty())
0065     php.nmods_ = 0;
0066   else
0067     php.nmods_ = static_cast<int>(tempD[0]);
0068 
0069   std::vector<std::string> tempS = getStringArray("Depth1Name", sv);
0070   if (!tempS.empty())
0071     php.depth1Name_ = tempS[0];
0072   else
0073     php.depth1Name_ = " ";
0074   tempS = getStringArray("Depth2Name", sv);
0075   if (!tempS.empty())
0076     php.depth2Name_ = tempS[0];
0077   else
0078     php.depth2Name_ = " ";
0079 
0080   //Then the logical volumes
0081   while (dodet) {
0082     if (std::find(php.lvNames_.begin(), php.lvNames_.end(), fv.logicalPart().name().name()) == php.lvNames_.end()) {
0083       php.matNames_.emplace_back(fv.logicalPart().material().name().name());
0084       php.lvNames_.emplace_back(fv.logicalPart().name().name());
0085       const DDSolid& sol = fv.logicalPart().solid();
0086       const std::vector<double>& paras = sol.parameters();
0087       double dz = (sol.shape() == DDSolidShape::ddtrap) ? (2 * paras[0]) : 0.0;
0088       php.dzs_.emplace_back(dz);
0089     }
0090     dodet = fv.next();
0091   }
0092   return this->buildParameters(php);
0093 }
0094 
0095 bool EcalSimParametersFromDD::build(const cms::DDCompactView* cpv,
0096                                     const std::string& name,
0097                                     EcalSimulationParameters& php) {
0098 #ifdef EDM_ML_DEBUG
0099   edm::LogVerbatim("EcalGeom") << "Inside EcalSimParametersFromDD::build(const cms::DDCompactView*, const std::string, "
0100                                   "EcalSimulationParameters&)";
0101 #endif
0102   // Get the filtered view
0103   std::string attribute = "ReadOutName";
0104   cms::DDFilteredView fv(cpv->detector(), cpv->detector()->worldVolume());
0105 
0106   //First the specpars
0107   std::string specName = ((name == "EcalHitsEE") ? "ecal_ee" : ((name == "EcalHitsES") ? "ecal_sf" : "ecal_eb"));
0108 
0109   php.useWeight_ = true;
0110   std::vector<double> tempD = fv.get<std::vector<double> >(specName, "EnergyWeight");
0111   if (!tempD.empty()) {
0112     if (tempD[0] < 0.1)
0113       php.useWeight_ = false;
0114   }
0115   tempD = fv.get<std::vector<double> >(specName, "nxtalEta");
0116   if (tempD.empty())
0117     php.nxtalEta_ = 0;
0118   else
0119     php.nxtalEta_ = static_cast<int>(tempD[0]);
0120   tempD = fv.get<std::vector<double> >(specName, "nxtalPhi");
0121   if (tempD.empty())
0122     php.nxtalPhi_ = 0;
0123   else
0124     php.nxtalPhi_ = static_cast<int>(tempD[0]);
0125   tempD = fv.get<std::vector<double> >(specName, "PhiBaskets");
0126   if (tempD.empty())
0127     php.phiBaskets_ = 0;
0128   else
0129     php.phiBaskets_ = static_cast<int>(tempD[0]);
0130   php.etaBaskets_ = dbl_to_int(fv.get<std::vector<double> >(specName, "EtaBaskets"));
0131   tempD = fv.get<std::vector<double> >(specName, "ncrys");
0132   if (tempD.empty())
0133     php.ncrys_ = 0;
0134   else
0135     php.ncrys_ = static_cast<int>(tempD[0]);
0136   tempD = fv.get<std::vector<double> >(specName, "nmods");
0137   if (tempD.empty())
0138     php.nmods_ = 0;
0139   else
0140     php.nmods_ = static_cast<int>(tempD[0]);
0141 
0142   std::vector<std::string> tempS = fv.get<std::vector<std::string> >(specName, "Depth1Name");
0143   if (!tempS.empty())
0144     php.depth1Name_ = tempS[0];
0145   else
0146     php.depth1Name_ = " ";
0147   tempS = fv.get<std::vector<std::string> >(specName, "Depth2Name");
0148   if (!tempS.empty())
0149     php.depth2Name_ = tempS[0];
0150   else
0151     php.depth2Name_ = " ";
0152 
0153   //Then the logical volumes
0154   cms::DDSpecParRefs refs;
0155   const cms::DDSpecParRegistry& mypar = cpv->specpars();
0156   mypar.filter(refs, attribute, name);
0157   fv.mergedSpecifics(refs);
0158   while (fv.firstChild()) {
0159     const std::string name{dd4hep::dd::noNamespace(fv.name()).data(), dd4hep::dd::noNamespace(fv.name()).size()};
0160     const std::string matName{dd4hep::dd::noNamespace(fv.materialName()).data(),
0161                               dd4hep::dd::noNamespace(fv.materialName()).size()};
0162     if (std::find(php.lvNames_.begin(), php.lvNames_.end(), name) == php.lvNames_.end()) {
0163       php.matNames_.emplace_back(matName);
0164       php.lvNames_.emplace_back(name);
0165       const std::vector<double>& paras = fv.parameters();
0166       double dz = (dd4hep::isA<dd4hep::Trap>(fv.solid())) ? ((2.0 * paras[0]) / dd4hep::mm) : 0.0;
0167       php.dzs_.emplace_back(dz);
0168     }
0169   };
0170 
0171   return this->buildParameters(php);
0172 }
0173 
0174 bool EcalSimParametersFromDD::buildParameters(const EcalSimulationParameters& php) {
0175 #ifdef EDM_ML_DEBUG
0176   edm::LogVerbatim("EcalGeom") << "EcalSimParametersFromDD:: nxtalEta:" << php.nxtalEta_
0177                                << " nxtalPhi:" << php.nxtalPhi_ << " phiBaskets:" << php.phiBaskets_
0178                                << " ncrys:" << php.ncrys_ << " nmods: " << php.nmods_ << " useWeight:" << php.useWeight_
0179                                << " DeothNames:" << php.depth1Name_ << ":" << php.depth2Name_;
0180   myPrint("etaBaskets", php.etaBaskets_);
0181   edm::LogVerbatim("EcalGeom") << "EcalSimParametersFromDD:: " << php.lvNames_.size() << " lvNames, "
0182                                << php.matNames_.size() << " matNames and " << php.dzs_.size() << "dzs";
0183 #endif
0184 
0185   return true;
0186 }
0187 
0188 std::vector<double> EcalSimParametersFromDD::getDDDArray(const std::string& str, const DDsvalues_type& sv) {
0189 #ifdef EDM_ML_DEBUG
0190   edm::LogVerbatim("EcalGeom") << "EcalSimParametersFromDD:getDDDArray called for " << str;
0191 #endif
0192   DDValue value(str);
0193   if (DDfetch(&sv, value)) {
0194 #ifdef EDM_ML_DEBUG
0195     edm::LogVerbatim("EcalGeom") << value;
0196 #endif
0197     const std::vector<double>& fvec = value.doubles();
0198     return fvec;
0199   } else {
0200     std::vector<double> fvec;
0201     return fvec;
0202   }
0203 }
0204 
0205 std::vector<std::string> EcalSimParametersFromDD::getStringArray(const std::string& str, const DDsvalues_type& sv) {
0206 #ifdef EDM_ML_DEBUG
0207   edm::LogVerbatim("EcalGeom") << "EcalSimParametersFromDD:getStringArray called for " << str;
0208 #endif
0209   DDValue value(str);
0210   if (DDfetch(&sv, value)) {
0211 #ifdef EDM_ML_DEBUG
0212     edm::LogVerbatim("EcalGeom") << value;
0213 #endif
0214     const std::vector<std::string>& fvec = value.strings();
0215     return fvec;
0216   } else {
0217     std::vector<std::string> fvec;
0218     return fvec;
0219   }
0220 }