Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:14:58

0001 ///////////////////////////////////////////////////////////////////////////////
0002 // File: DDHGCalEEFileAlgo.cc
0003 // Description: Geometry factory class for HGCal (EE and HESil) using
0004 //              information from the file
0005 ///////////////////////////////////////////////////////////////////////////////
0006 
0007 #include "DataFormats/Math/interface/angle_units.h"
0008 #include "DetectorDescription/Core/interface/DDAlgorithm.h"
0009 #include "DetectorDescription/Core/interface/DDAlgorithmFactory.h"
0010 #include "DetectorDescription/Core/interface/DDCurrentNamespace.h"
0011 #include "DetectorDescription/Core/interface/DDLogicalPart.h"
0012 #include "DetectorDescription/Core/interface/DDMaterial.h"
0013 #include "DetectorDescription/Core/interface/DDSolid.h"
0014 #include "DetectorDescription/Core/interface/DDSplit.h"
0015 #include "DetectorDescription/Core/interface/DDTypes.h"
0016 #include "DetectorDescription/Core/interface/DDutils.h"
0017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0018 #include "FWCore/PluginManager/interface/PluginFactory.h"
0019 #include "Geometry/HGCalCommonData/interface/HGCalGeomTools.h"
0020 #include "Geometry/HGCalCommonData/interface/HGCalParameters.h"
0021 #include "Geometry/HGCalCommonData/interface/HGCalProperty.h"
0022 #include "Geometry/HGCalCommonData/interface/HGCalTypes.h"
0023 #include "Geometry/HGCalCommonData/interface/HGCalWaferIndex.h"
0024 #include "Geometry/HGCalCommonData/interface/HGCalWaferType.h"
0025 
0026 #include <cmath>
0027 #include <memory>
0028 #include <string>
0029 #include <unordered_set>
0030 #include <vector>
0031 
0032 //#define EDM_ML_DEBUG
0033 using namespace angle_units::operators;
0034 
0035 class DDHGCalEEFileAlgo : public DDAlgorithm {
0036 public:
0037   DDHGCalEEFileAlgo();
0038 
0039   void initialize(const DDNumericArguments& nArgs,
0040                   const DDVectorArguments& vArgs,
0041                   const DDMapArguments& mArgs,
0042                   const DDStringArguments& sArgs,
0043                   const DDStringVectorArguments& vsArgs) override;
0044   void execute(DDCompactView& cpv) override;
0045 
0046 protected:
0047   void constructLayers(const DDLogicalPart&, DDCompactView& cpv);
0048   void positionSensitive(
0049       const DDLogicalPart& glog, double rin, double rout, double zpos, int layertype, int layer, DDCompactView& cpv);
0050 
0051 private:
0052   HGCalGeomTools geomTools_;
0053 
0054   static constexpr double tol1_ = 0.01;
0055   static constexpr double tol2_ = 0.00001;
0056 
0057   std::vector<std::string> wafers_;     // Wafers
0058   std::vector<std::string> materials_;  // Materials
0059   std::vector<std::string> names_;      // Names
0060   std::vector<double> thick_;           // Thickness of the material
0061   std::vector<int> copyNumber_;         // Initial copy numbers
0062   std::vector<int> layers_;             // Number of layers in a section
0063   std::vector<double> layerThick_;      // Thickness of each section
0064   std::vector<int> layerType_;          // Type of the layer
0065   std::vector<int> layerSense_;         // Content of a layer (sensitive?)
0066   std::vector<int> layerCenter_;        // Centering of the wafers
0067   int firstLayer_;                      // Copy # of the first sensitive layer
0068   int absorbMode_;                      // Absorber mode
0069   int sensitiveMode_;                   // Sensitive mode
0070   double zMinBlock_;                    // Starting z-value of the block
0071   std::vector<int> waferIndex_;         // Wafer index for the types
0072   std::vector<int> waferProperty_;      // Wafer property
0073   double waferSize_;                    // Width of the wafer
0074   double waferSepar_;                   // Sensor separation
0075   int sectors_;                         // Sectors
0076   std::vector<double> slopeB_;          // Slope at the lower R
0077   std::vector<double> zFrontB_;         // Starting Z values for the slopes
0078   std::vector<double> rMinFront_;       // Corresponding rMin's
0079   std::vector<double> slopeT_;          // Slopes at the larger R
0080   std::vector<double> zFrontT_;         // Starting Z values for the slopes
0081   std::vector<double> rMaxFront_;       // Corresponding rMax's
0082   std::string nameSpace_;               // Namespace of this and ALL sub-parts
0083   std::unordered_set<int> copies_;      // List of copy #'s
0084   double alpha_, cosAlpha_;
0085 };
0086 
0087 DDHGCalEEFileAlgo::DDHGCalEEFileAlgo() {
0088 #ifdef EDM_ML_DEBUG
0089   edm::LogVerbatim("HGCalGeom") << "DDHGCalEEFileAlgo: Creating an instance";
0090 #endif
0091 }
0092 
0093 void DDHGCalEEFileAlgo::initialize(const DDNumericArguments& nArgs,
0094                                    const DDVectorArguments& vArgs,
0095                                    const DDMapArguments&,
0096                                    const DDStringArguments& sArgs,
0097                                    const DDStringVectorArguments& vsArgs) {
0098   wafers_ = vsArgs["WaferNames"];
0099 #ifdef EDM_ML_DEBUG
0100   edm::LogVerbatim("HGCalGeom") << "DDHGCalEEFileAlgo: " << wafers_.size() << " wafers";
0101   for (unsigned int i = 0; i < wafers_.size(); ++i)
0102     edm::LogVerbatim("HGCalGeom") << "Wafer[" << i << "] " << wafers_[i];
0103 #endif
0104   materials_ = vsArgs["MaterialNames"];
0105   names_ = vsArgs["VolumeNames"];
0106   thick_ = vArgs["Thickness"];
0107   copyNumber_.resize(materials_.size(), 1);
0108 #ifdef EDM_ML_DEBUG
0109   edm::LogVerbatim("HGCalGeom") << "DDHGCalEEFileAlgo: " << materials_.size() << " types of volumes";
0110   for (unsigned int i = 0; i < names_.size(); ++i)
0111     edm::LogVerbatim("HGCalGeom") << "Volume [" << i << "] " << names_[i] << " of thickness " << thick_[i]
0112                                   << " filled with " << materials_[i] << " first copy number " << copyNumber_[i];
0113 #endif
0114   layers_ = dbl_to_int(vArgs["Layers"]);
0115   layerThick_ = vArgs["LayerThick"];
0116 #ifdef EDM_ML_DEBUG
0117   edm::LogVerbatim("HGCalGeom") << "There are " << layers_.size() << " blocks";
0118   for (unsigned int i = 0; i < layers_.size(); ++i)
0119     edm::LogVerbatim("HGCalGeom") << "Block [" << i << "] of thickness " << layerThick_[i] << " with " << layers_[i]
0120                                   << " layers";
0121 #endif
0122   layerType_ = dbl_to_int(vArgs["LayerType"]);
0123   layerSense_ = dbl_to_int(vArgs["LayerSense"]);
0124   firstLayer_ = (int)(nArgs["FirstLayer"]);
0125   absorbMode_ = (int)(nArgs["AbsorberMode"]);
0126   sensitiveMode_ = (int)(nArgs["SensitiveMode"]);
0127 #ifdef EDM_ML_DEBUG
0128   edm::LogVerbatim("HGCalGeom") << "First Layer " << firstLayer_ << " and "
0129                                 << "Absober:Sensitive mode " << absorbMode_ << ":" << sensitiveMode_;
0130 #endif
0131   layerCenter_ = dbl_to_int(vArgs["LayerCenter"]);
0132 #ifdef EDM_ML_DEBUG
0133   for (unsigned int i = 0; i < layerCenter_.size(); ++i)
0134     edm::LogVerbatim("HGCalGeom") << "LayerCenter [" << i << "] " << layerCenter_[i];
0135 #endif
0136   if (firstLayer_ > 0) {
0137     for (unsigned int i = 0; i < layerType_.size(); ++i) {
0138       if (layerSense_[i] > 0) {
0139         int ii = layerType_[i];
0140         copyNumber_[ii] = firstLayer_;
0141 #ifdef EDM_ML_DEBUG
0142         edm::LogVerbatim("HGCalGeom") << "First copy number for layer type " << i << ":" << ii << " with "
0143                                       << materials_[ii] << " changed to " << copyNumber_[ii];
0144 #endif
0145         break;
0146       }
0147     }
0148   } else {
0149     firstLayer_ = 1;
0150   }
0151 #ifdef EDM_ML_DEBUG
0152   edm::LogVerbatim("HGCalGeom") << "There are " << layerType_.size() << " layers";
0153   for (unsigned int i = 0; i < layerType_.size(); ++i)
0154     edm::LogVerbatim("HGCalGeom") << "Layer [" << i << "] with material type " << layerType_[i] << " sensitive class "
0155                                   << layerSense_[i];
0156 #endif
0157   zMinBlock_ = nArgs["zMinBlock"];
0158   waferSize_ = nArgs["waferSize"];
0159   waferSepar_ = nArgs["SensorSeparation"];
0160   sectors_ = (int)(nArgs["Sectors"]);
0161   alpha_ = (1._pi) / sectors_;
0162   cosAlpha_ = cos(alpha_);
0163 #ifdef EDM_ML_DEBUG
0164   edm::LogVerbatim("HGCalGeom") << "zStart " << zMinBlock_ << " wafer width " << waferSize_ << " separations "
0165                                 << waferSepar_ << " sectors " << sectors_ << ":" << convertRadToDeg(alpha_) << ":"
0166                                 << cosAlpha_;
0167 #endif
0168   waferIndex_ = dbl_to_int(vArgs["WaferIndex"]);
0169   waferProperty_ = dbl_to_int(vArgs["WaferProperties"]);
0170 #ifdef EDM_ML_DEBUG
0171   edm::LogVerbatim("HGCalGeom") << "waferProperties with " << waferIndex_.size() << " entries";
0172   for (unsigned int k = 0; k < waferIndex_.size(); ++k)
0173     edm::LogVerbatim("HGCalGeom") << "[" << k << "] " << waferIndex_[k] << " ("
0174                                   << HGCalWaferIndex::waferLayer(waferIndex_[k]) << ", "
0175                                   << HGCalWaferIndex::waferU(waferIndex_[k]) << ", "
0176                                   << HGCalWaferIndex::waferV(waferIndex_[k]) << ") : ("
0177                                   << HGCalProperty::waferThick(waferProperty_[k]) << ":"
0178                                   << HGCalProperty::waferPartial(waferProperty_[k]) << ":"
0179                                   << HGCalProperty::waferOrient(waferProperty_[k]) << ")";
0180 #endif
0181   slopeB_ = vArgs["SlopeBottom"];
0182   zFrontB_ = vArgs["ZFrontBottom"];
0183   rMinFront_ = vArgs["RMinFront"];
0184   slopeT_ = vArgs["SlopeTop"];
0185   zFrontT_ = vArgs["ZFrontTop"];
0186   rMaxFront_ = vArgs["RMaxFront"];
0187 #ifdef EDM_ML_DEBUG
0188   for (unsigned int i = 0; i < slopeB_.size(); ++i)
0189     edm::LogVerbatim("HGCalGeom") << "Bottom Block [" << i << "] Zmin " << zFrontB_[i] << " Rmin " << rMinFront_[i]
0190                                   << " Slope " << slopeB_[i];
0191   for (unsigned int i = 0; i < slopeT_.size(); ++i)
0192     edm::LogVerbatim("HGCalGeom") << "Top Block [" << i << "] Zmin " << zFrontT_[i] << " Rmax " << rMaxFront_[i]
0193                                   << " Slope " << slopeT_[i];
0194 #endif
0195   nameSpace_ = DDCurrentNamespace::ns();
0196 #ifdef EDM_ML_DEBUG
0197   edm::LogVerbatim("HGCalGeom") << "DDHGCalEEFileAlgo: NameSpace " << nameSpace_ << ":";
0198 #endif
0199 }
0200 
0201 ////////////////////////////////////////////////////////////////////
0202 // DDHGCalEEFileAlgo methods...
0203 ////////////////////////////////////////////////////////////////////
0204 
0205 void DDHGCalEEFileAlgo::execute(DDCompactView& cpv) {
0206 #ifdef EDM_ML_DEBUG
0207   edm::LogVerbatim("HGCalGeom") << "==>> Constructing DDHGCalEEFileAlgo...";
0208   copies_.clear();
0209 #endif
0210   constructLayers(parent(), cpv);
0211 #ifdef EDM_ML_DEBUG
0212   edm::LogVerbatim("HGCalGeom") << "DDHGCalEEFileAlgo: " << copies_.size() << " different wafer copy numbers";
0213   int k(0);
0214   for (std::unordered_set<int>::const_iterator itr = copies_.begin(); itr != copies_.end(); ++itr, ++k) {
0215     edm::LogVerbatim("HGCalGeom") << "Copy [" << k << "] : " << (*itr);
0216   }
0217   copies_.clear();
0218   edm::LogVerbatim("HGCalGeom") << "<<== End of DDHGCalEEFileAlgo construction...";
0219 #endif
0220 }
0221 
0222 void DDHGCalEEFileAlgo::constructLayers(const DDLogicalPart& module, DDCompactView& cpv) {
0223   double zi(zMinBlock_);
0224   int laymin(0);
0225   for (unsigned int i = 0; i < layers_.size(); i++) {
0226     double zo = zi + layerThick_[i];
0227     double routF = HGCalGeomTools::radius(zi, zFrontT_, rMaxFront_, slopeT_);
0228     int laymax = laymin + layers_[i];
0229     double zz = zi;
0230     double thickTot(0);
0231     for (int ly = laymin; ly < laymax; ++ly) {
0232       int ii = layerType_[ly];
0233       int copy = copyNumber_[ii];
0234       double hthick = 0.5 * thick_[ii];
0235       double rinB = HGCalGeomTools::radius(zo - tol1_, zFrontB_, rMinFront_, slopeB_);
0236       zz += hthick;
0237       thickTot += thick_[ii];
0238 
0239       std::string name = names_[ii] + std::to_string(copy);
0240 #ifdef EDM_ML_DEBUG
0241       edm::LogVerbatim("HGCalGeom") << "DDHGCalEEFileAlgo: Layer " << ly << ":" << ii << " Front " << zi << ", "
0242                                     << routF << " Back " << zo << ", " << rinB << " superlayer thickness "
0243                                     << layerThick_[i];
0244 #endif
0245       DDName matName(DDSplit(materials_[ii]).first, DDSplit(materials_[ii]).second);
0246       DDMaterial matter(matName);
0247       DDLogicalPart glog;
0248       if (layerSense_[ly] < 1) {
0249         std::vector<double> pgonZ, pgonRin, pgonRout;
0250         double rmax = routF * cosAlpha_ - tol1_;
0251         HGCalGeomTools::radius(zz - hthick,
0252                                zz + hthick,
0253                                zFrontB_,
0254                                rMinFront_,
0255                                slopeB_,
0256                                zFrontT_,
0257                                rMaxFront_,
0258                                slopeT_,
0259                                -layerSense_[ly],
0260                                pgonZ,
0261                                pgonRin,
0262                                pgonRout);
0263         for (unsigned int isec = 0; isec < pgonZ.size(); ++isec) {
0264           pgonZ[isec] -= zz;
0265           if (layerSense_[ly] == 0 || absorbMode_ == 0)
0266             pgonRout[isec] = rmax;
0267           else
0268             pgonRout[isec] = pgonRout[isec] * cosAlpha_ - tol1_;
0269         }
0270         DDSolid solid =
0271             DDSolidFactory::polyhedra(DDName(name, nameSpace_), sectors_, -alpha_, 2._pi, pgonZ, pgonRin, pgonRout);
0272         glog = DDLogicalPart(solid.ddname(), matter, solid);
0273 #ifdef EDM_ML_DEBUG
0274         edm::LogVerbatim("HGCalGeom") << "DDHGCalEEFileAlgo: " << solid.name() << " polyhedra of " << sectors_
0275                                       << " sectors covering " << convertRadToDeg(-alpha_) << ":"
0276                                       << convertRadToDeg(-alpha_ + 2._pi) << " with " << pgonZ.size()
0277                                       << " sections and filled with " << matName;
0278         for (unsigned int k = 0; k < pgonZ.size(); ++k)
0279           edm::LogVerbatim("HGCalGeom") << "[" << k << "] z " << pgonZ[k] << " R " << pgonRin[k] << ":" << pgonRout[k];
0280 #endif
0281       } else {
0282         double rins =
0283             (sensitiveMode_ < 1) ? rinB : HGCalGeomTools::radius(zz + hthick - tol1_, zFrontB_, rMinFront_, slopeB_);
0284         double routs =
0285             (sensitiveMode_ < 1) ? routF : HGCalGeomTools::radius(zz - hthick, zFrontT_, rMaxFront_, slopeT_);
0286         DDSolid solid = DDSolidFactory::tubs(DDName(name, nameSpace_), hthick, rins, routs, 0.0, 2._pi);
0287         glog = DDLogicalPart(solid.ddname(), matter, solid);
0288 #ifdef EDM_ML_DEBUG
0289         edm::LogVerbatim("HGCalGeom") << "DDHGCalEEFileAlgo: " << solid.name() << " Tubs made of " << matName
0290                                       << " of dimensions " << rinB << ":" << rins << ", " << routF << ":" << routs
0291                                       << ", " << hthick << ", 0.0, 360.0 and position " << glog.name() << " number "
0292                                       << copy << ":" << layerCenter_[copy - firstLayer_];
0293 #endif
0294         positionSensitive(glog, rins, routs, zz, layerSense_[ly], (copy - firstLayer_), cpv);
0295       }
0296       DDTranslation r1(0, 0, zz);
0297       DDRotation rot;
0298       cpv.position(glog, module, copy, r1, rot);
0299       ++copyNumber_[ii];
0300 #ifdef EDM_ML_DEBUG
0301       edm::LogVerbatim("HGCalGeom") << "DDHGCalEEFileAlgo: " << glog.name() << " number " << copy << " positioned in "
0302                                     << module.name() << " at " << r1 << " with no rotation";
0303 #endif
0304       zz += hthick;
0305     }  // End of loop over layers in a block
0306     zi = zo;
0307     laymin = laymax;
0308     // Make consistency check of all the partitions of the block
0309     if (std::abs(thickTot - layerThick_[i]) >= tol2_) {
0310       if (thickTot > layerThick_[i]) {
0311         edm::LogError("HGCalGeom") << "Thickness of the partition " << layerThick_[i] << " is smaller than " << thickTot
0312                                    << ": thickness of all its components **** ERROR ****";
0313       } else {
0314         edm::LogWarning("HGCalGeom") << "Thickness of the partition " << layerThick_[i] << " does not match with "
0315                                      << thickTot << " of the components";
0316       }
0317     }
0318   }  // End of loop over blocks
0319 }
0320 
0321 void DDHGCalEEFileAlgo::positionSensitive(
0322     const DDLogicalPart& glog, double rin, double rout, double zpos, int layertype, int layer, DDCompactView& cpv) {
0323   static const double sqrt3 = std::sqrt(3.0);
0324   int layercenter = layerCenter_[layer];
0325   double r = 0.5 * (waferSize_ + waferSepar_);
0326   double R = 2.0 * r / sqrt3;
0327   double dy = 0.75 * R;
0328   int N = (int)(0.5 * rout / r) + 2;
0329   const auto& xyoff = geomTools_.shiftXY(layercenter, (waferSize_ + waferSepar_));
0330 #ifdef EDM_ML_DEBUG
0331   int ium(0), ivm(0), iumAll(0), ivmAll(0), kount(0), ntot(0), nin(0);
0332   std::vector<int> ntype(6, 0);
0333   edm::LogVerbatim("HGCalGeom") << "DDHGCalEEFileAlgo: " << glog.ddname() << " rin:rout " << rin << ":" << rout
0334                                 << " zpos " << zpos << " N " << N << " for maximum u, v;  r " << r << " R " << R
0335                                 << " dy " << dy << " Shift " << xyoff.first << ":" << xyoff.second << " WaferSize "
0336                                 << (waferSize_ + waferSepar_);
0337 #endif
0338   for (int u = -N; u <= N; ++u) {
0339     for (int v = -N; v <= N; ++v) {
0340 #ifdef EDM_ML_DEBUG
0341       int iu = std::abs(u);
0342       int iv = std::abs(v);
0343 #endif
0344       int nr = 2 * v;
0345       int nc = -2 * u + v;
0346       double xpos = xyoff.first + nc * r;
0347       double ypos = xyoff.second + nr * dy;
0348       const auto& corner = HGCalGeomTools::waferCorner(xpos, ypos, r, R, rin, rout, false);
0349 #ifdef EDM_ML_DEBUG
0350       ++ntot;
0351       if (((corner.first <= 0) && std::abs(u) < 5 && std::abs(v) < 5) || (std::abs(u) < 2 && std::abs(v) < 2)) {
0352         edm::LogVerbatim("HGCalGeom") << "DDHGCalEEFileAlgo: " << glog.ddname() << " R " << rin << ":" << rout
0353                                       << "\n Z " << zpos << " LayerType " << layertype << " u " << u << " v " << v
0354                                       << " with " << corner.first << " corners";
0355       }
0356 #endif
0357       int indx = HGCalWaferIndex::waferIndex((layer + firstLayer_), u, v, false);
0358       int type = HGCalWaferType::getType(indx, waferIndex_, waferProperty_);
0359       if (corner.first > 0 && type >= 0) {
0360         int copy = HGCalTypes::packTypeUV(type, u, v);
0361         if (layertype > 1)
0362           type += 3;
0363 #ifdef EDM_ML_DEBUG
0364         edm::LogVerbatim("HGCalGeom") << " DDHGCalEEFileAlgo: " << wafers_[type] << " number " << copy << " type "
0365                                       << type << " layer:u:v:indx " << (layer + firstLayer_) << ":" << u << ":" << v
0366                                       << ":" << indx;
0367         if (iu > ium)
0368           ium = iu;
0369         if (iv > ivm)
0370           ivm = iv;
0371         kount++;
0372         if (copies_.count(copy) == 0)
0373           copies_.insert(copy);
0374 #endif
0375         if (corner.first == (int)(HGCalParameters::k_CornerSize)) {
0376 #ifdef EDM_ML_DEBUG
0377           if (iu > iumAll)
0378             iumAll = iu;
0379           if (iv > ivmAll)
0380             ivmAll = iv;
0381           ++nin;
0382 #endif
0383           DDTranslation tran(xpos, ypos, 0.0);
0384           DDRotation rotation;
0385           DDName name = DDName(DDSplit(wafers_[type]).first, DDSplit(wafers_[type]).second);
0386           cpv.position(name, glog.ddname(), copy, tran, rotation);
0387 #ifdef EDM_ML_DEBUG
0388           ++ntype[type];
0389           edm::LogVerbatim("HGCalGeom") << " DDHGCalEEFileAlgo: " << name << " number " << copy << " type " << layertype
0390                                         << ":" << type << " positioned in " << glog.ddname() << " at " << tran
0391                                         << " with no rotation";
0392 #endif
0393         }
0394       }
0395     }
0396   }
0397 #ifdef EDM_ML_DEBUG
0398   edm::LogVerbatim("HGCalGeom") << "DDHGCalEEFileAlgo: Maximum # of u " << ium << ":" << iumAll << " # of v " << ivm
0399                                 << ":" << ivmAll << " and " << nin << ":" << kount << ":" << ntot << " wafers ("
0400                                 << ntype[0] << ":" << ntype[1] << ":" << ntype[2] << ":" << ntype[3] << ":" << ntype[4]
0401                                 << ":" << ntype[5] << ") for " << glog.ddname() << " R " << rin << ":" << rout;
0402 #endif
0403 }
0404 
0405 DEFINE_EDM_PLUGIN(DDAlgorithmFactory, DDHGCalEEFileAlgo, "hgcal:DDHGCalEEFileAlgo");