Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-07-04 01:50:58

0001 ///////////////////////////////////////////////////////////////////////////////
0002 // File: DDHGCalMixRotatedLayer.cc
0003 // Description: Geometry factory class for HGCal (Mix) adopted for DD4hep
0004 ///////////////////////////////////////////////////////////////////////////////
0005 
0006 #include <cmath>
0007 #include <memory>
0008 #include <string>
0009 #include <unordered_set>
0010 #include <vector>
0011 
0012 #include "Geometry/HGCalCommonData/interface/HGCalCell.h"
0013 #include "Geometry/HGCalCommonData/interface/HGCalCassette.h"
0014 #include "Geometry/HGCalCommonData/interface/HGCalGeomTools.h"
0015 #include "Geometry/HGCalCommonData/interface/HGCalParameters.h"
0016 #include "Geometry/HGCalCommonData/interface/HGCalProperty.h"
0017 #include "Geometry/HGCalCommonData/interface/HGCalTileIndex.h"
0018 #include "Geometry/HGCalCommonData/interface/HGCalTypes.h"
0019 #include "Geometry/HGCalCommonData/interface/HGCalWaferIndex.h"
0020 #include "Geometry/HGCalCommonData/interface/HGCalWaferType.h"
0021 #include "DD4hep/DetFactoryHelper.h"
0022 #include "DataFormats/Math/interface/angle_units.h"
0023 #include "DetectorDescription/DDCMS/interface/DDPlugins.h"
0024 #include "DetectorDescription/DDCMS/interface/DDutils.h"
0025 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0026 
0027 //#define EDM_ML_DEBUG
0028 using namespace angle_units::operators;
0029 
0030 struct HGCalMixRotatedLayer {
0031   HGCalMixRotatedLayer() { throw cms::Exception("HGCalGeom") << "Wrong initialization to HGCalMixRotatedLayer"; }
0032   HGCalMixRotatedLayer(cms::DDParsingContext& ctxt, xml_h e) {
0033     cms::DDNamespace ns(ctxt, e, true);
0034     cms::DDAlgoArguments args(ctxt, e);
0035 
0036 #ifdef EDM_ML_DEBUG
0037     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: Creating an instance";
0038 #endif
0039 
0040     static constexpr double tol1 = 0.01 * dd4hep::mm;
0041     dd4hep::Volume mother = ns.volume(args.parentName());
0042 
0043     waferTypes_ = args.value<int>("WaferTypes");
0044     facingTypes_ = args.value<int>("FacingTypes");
0045     orientationTypes_ = args.value<int>("OrientationTypes");
0046     placeOffset_ = args.value<int>("PlaceOffset");
0047     phiBinsScint_ = args.value<int>("NPhiBinScint");
0048 #ifdef EDM_ML_DEBUG
0049     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer::Number of types of wafers: " << waferTypes_
0050                                   << " facings: " << facingTypes_ << " Orientations: " << orientationTypes_
0051                                   << " PlaceOffset: " << placeOffset_ << "; number of cells along phi "
0052                                   << phiBinsScint_;
0053 #endif
0054     firstLayer_ = args.value<int>("FirstLayer");
0055     absorbMode_ = args.value<int>("AbsorberMode");
0056     sensitiveMode_ = args.value<int>("SensitiveMode");
0057 #ifdef EDM_ML_DEBUG
0058     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer::First Layer " << firstLayer_ << " and "
0059                                   << "Absober:Sensitive mode " << absorbMode_ << ":" << sensitiveMode_;
0060 #endif
0061     zMinBlock_ = args.value<double>("zMinBlock");
0062     waferSize_ = args.value<double>("waferSize");
0063     waferSepar_ = args.value<double>("SensorSeparation");
0064     sectors_ = args.value<int>("Sectors");
0065     cassettes_ = args.value<int>("Cassettes");
0066     alpha_ = (1._pi) / sectors_;
0067     cosAlpha_ = cos(alpha_);
0068 #ifdef EDM_ML_DEBUG
0069     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: zStart " << cms::convert2mm(zMinBlock_) << " wafer width "
0070                                   << cms::convert2mm(waferSize_) << " separations " << cms::convert2mm(waferSepar_)
0071                                   << " sectors " << sectors_ << ":" << convertRadToDeg(alpha_) << ":" << cosAlpha_
0072                                   << " with " << cassettes_ << " cassettes";
0073 #endif
0074     slopeB_ = args.value<std::vector<double>>("SlopeBottom");
0075     zFrontB_ = args.value<std::vector<double>>("ZFrontBottom");
0076     rMinFront_ = args.value<std::vector<double>>("RMinFront");
0077     slopeT_ = args.value<std::vector<double>>("SlopeTop");
0078     zFrontT_ = args.value<std::vector<double>>("ZFrontTop");
0079     rMaxFront_ = args.value<std::vector<double>>("RMaxFront");
0080 #ifdef EDM_ML_DEBUG
0081     for (unsigned int i = 0; i < slopeB_.size(); ++i)
0082       edm::LogVerbatim("HGCalGeom") << "Bottom Block [" << i << "] Zmin " << cms::convert2mm(zFrontB_[i]) << " Rmin "
0083                                     << cms::convert2mm(rMinFront_[i]) << " Slope " << slopeB_[i];
0084     for (unsigned int i = 0; i < slopeT_.size(); ++i)
0085       edm::LogVerbatim("HGCalGeom") << "Top Block [" << i << "] Zmin " << cms::convert2mm(zFrontT_[i]) << " Rmax "
0086                                     << cms::convert2mm(rMaxFront_[i]) << " Slope " << slopeT_[i];
0087 #endif
0088 
0089     waferFull_ = args.value<std::vector<std::string>>("WaferNamesFull");
0090     waferPart_ = args.value<std::vector<std::string>>("WaferNamesPartial");
0091 #ifdef EDM_ML_DEBUG
0092     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << waferFull_.size() << " full and "
0093                                   << waferPart_.size() << " partial modules\nDDHGCalMixRotatedLayer:Full Modules:";
0094     unsigned int i1max = static_cast<unsigned int>(waferFull_.size());
0095     for (unsigned int i1 = 0; i1 < i1max; i1 += 2) {
0096       std::ostringstream st1;
0097       unsigned int i2 = std::min((i1 + 2), i1max);
0098       for (unsigned int i = i1; i < i2; ++i)
0099         st1 << " [" << i << "] " << waferFull_[i];
0100       edm::LogVerbatim("HGCalGeom") << st1.str();
0101     }
0102     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: Partial Modules:";
0103     i1max = static_cast<unsigned int>(waferPart_.size());
0104     for (unsigned int i1 = 0; i1 < i1max; i1 += 2) {
0105       std::ostringstream st1;
0106       unsigned int i2 = std::min((i1 + 2), i1max);
0107       for (unsigned int i = i1; i < i2; ++i)
0108         st1 << " [" << i << "] " << waferPart_[i];
0109       edm::LogVerbatim("HGCalGeom") << st1.str();
0110     }
0111 #endif
0112 
0113     materials_ = args.value<std::vector<std::string>>("MaterialNames");
0114     names_ = args.value<std::vector<std::string>>("VolumeNames");
0115     thick_ = args.value<std::vector<double>>("Thickness");
0116     copyNumber_.resize(materials_.size(), 1);
0117 #ifdef EDM_ML_DEBUG
0118     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << materials_.size() << " types of volumes";
0119     for (unsigned int i = 0; i < names_.size(); ++i)
0120       edm::LogVerbatim("HGCalGeom") << "Volume [" << i << "] " << names_[i] << " of thickness "
0121                                     << cms::convert2mm(thick_[i]) << " filled with " << materials_[i]
0122                                     << " first copy number " << copyNumber_[i];
0123 #endif
0124     layers_ = args.value<std::vector<int>>("Layers");
0125     layerThick_ = args.value<std::vector<double>>("LayerThick");
0126 #ifdef EDM_ML_DEBUG
0127     edm::LogVerbatim("HGCalGeom") << "There are " << layers_.size() << " blocks";
0128     for (unsigned int i = 0; i < layers_.size(); ++i)
0129       edm::LogVerbatim("HGCalGeom") << "Block [" << i << "] of thickness " << cms::convert2mm(layerThick_[i])
0130                                     << " with " << layers_[i] << " layers";
0131 #endif
0132     layerType_ = args.value<std::vector<int>>("LayerType");
0133     layerSense_ = args.value<std::vector<int>>("LayerSense");
0134     layerOrient_ = args.value<std::vector<int>>("LayerTypes");
0135     for (unsigned int k = 0; k < layerOrient_.size(); ++k)
0136       layerOrient_[k] = HGCalTypes::layerType(layerOrient_[k]);
0137 #ifdef EDM_ML_DEBUG
0138     for (unsigned int i = 0; i < layerOrient_.size(); ++i)
0139       edm::LogVerbatim("HGCalGeom") << "LayerTypes [" << i << "] " << layerOrient_[i];
0140 #endif
0141     if (firstLayer_ > 0) {
0142       for (unsigned int i = 0; i < layerType_.size(); ++i) {
0143         if (layerSense_[i] > 0) {
0144           int ii = layerType_[i];
0145           copyNumber_[ii] = firstLayer_;
0146 #ifdef EDM_ML_DEBUG
0147           edm::LogVerbatim("HGCalGeom") << "First copy number for layer type " << i << ":" << ii << " with "
0148                                         << materials_[ii] << " changed to " << copyNumber_[ii];
0149 #endif
0150           break;
0151         }
0152       }
0153     } else {
0154       firstLayer_ = 1;
0155     }
0156 #ifdef EDM_ML_DEBUG
0157     edm::LogVerbatim("HGCalGeom") << "There are " << layerType_.size() << " layers";
0158     for (unsigned int i = 0; i < layerType_.size(); ++i)
0159       edm::LogVerbatim("HGCalGeom") << "Layer [" << i << "] with material type " << layerType_[i] << " sensitive class "
0160                                     << layerSense_[i];
0161 #endif
0162     materialTop_ = args.value<std::vector<std::string>>("TopMaterialNames");
0163     namesTop_ = args.value<std::vector<std::string>>("TopVolumeNames");
0164     layerThickTop_ = args.value<std::vector<double>>("TopLayerThickness");
0165     layerTypeTop_ = args.value<std::vector<int>>("TopLayerType");
0166     copyNumberTop_.resize(materialTop_.size(), firstLayer_);
0167 #ifdef EDM_ML_DEBUG
0168     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << materialTop_.size()
0169                                   << " types of volumes in the top part";
0170     for (unsigned int i = 0; i < materialTop_.size(); ++i)
0171       edm::LogVerbatim("HGCalGeom") << "Volume [" << i << "] " << namesTop_[i] << " of thickness "
0172                                     << cms::convert2mm(layerThickTop_[i]) << " filled with " << materialTop_[i]
0173                                     << " first copy number " << copyNumberTop_[i];
0174     edm::LogVerbatim("HGCalGeom") << "There are " << layerTypeTop_.size() << " layers in the top part";
0175     for (unsigned int i = 0; i < layerTypeTop_.size(); ++i)
0176       edm::LogVerbatim("HGCalGeom") << "Layer [" << i << "] with material type " << layerTypeTop_[i];
0177 #endif
0178     waferIndex_ = args.value<std::vector<int>>("WaferIndex");
0179     waferProperty_ = args.value<std::vector<int>>("WaferProperties");
0180     waferLayerStart_ = args.value<std::vector<int>>("WaferLayerStart");
0181     cassetteShift_ = args.value<std::vector<double>>("CassetteShift");
0182 #ifdef EDM_ML_DEBUG
0183     edm::LogVerbatim("HGCalGeom") << "waferProperties with " << waferIndex_.size() << " entries in "
0184                                   << waferLayerStart_.size() << " layers";
0185     for (unsigned int k = 0; k < waferLayerStart_.size(); ++k)
0186       edm::LogVerbatim("HGCalGeom") << "LayerStart[" << k << "] " << waferLayerStart_[k];
0187     for (unsigned int k = 0; k < waferIndex_.size(); ++k)
0188       edm::LogVerbatim("HGCalGeom") << "[" << k << "] " << waferIndex_[k] << " ("
0189                                     << HGCalWaferIndex::waferLayer(waferIndex_[k]) << ", "
0190                                     << HGCalWaferIndex::waferU(waferIndex_[k]) << ", "
0191                                     << HGCalWaferIndex::waferV(waferIndex_[k]) << ") : ("
0192                                     << HGCalProperty::waferThick(waferProperty_[k]) << ":"
0193                                     << HGCalProperty::waferPartial(waferProperty_[k]) << ":"
0194                                     << HGCalProperty::waferOrient(waferProperty_[k]) << ")";
0195     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << cassetteShift_.size()
0196                                   << " elements for cassette shifts";
0197     unsigned int j1max = cassetteShift_.size();
0198     for (unsigned int j1 = 0; j1 < j1max; j1 += 6) {
0199       std::ostringstream st1;
0200       unsigned int j2 = std::min((j1 + 6), j1max);
0201       for (unsigned int j = j1; j < j2; ++j)
0202         st1 << " [" << j << "] " << std::setw(9) << cassetteShift_[j];
0203       edm::LogVerbatim("HGCalGeom") << st1.str();
0204     }
0205 #endif
0206     tileRMin_ = args.value<std::vector<double>>("TileRMin");
0207     tileRMax_ = args.value<std::vector<double>>("TileRMax");
0208     tileIndex_ = args.value<std::vector<int>>("TileLayerRings");
0209     tilePhis_ = args.value<std::vector<int>>("TilePhiRange");
0210     tileLayerStart_ = args.value<std::vector<int>>("TileLayerStart");
0211 #ifdef EDM_ML_DEBUG
0212     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer:: with " << tileRMin_.size() << " rings";
0213     for (unsigned int k = 0; k < tileRMin_.size(); ++k)
0214       edm::LogVerbatim("HGCalGeom") << "Ring[" << k << "] " << cms::convert2mm(tileRMin_[k]) << " : "
0215                                     << cms::convert2mm(tileRMax_[k]);
0216     edm::LogVerbatim("HGCalGeom") << "TileProperties with " << tileIndex_.size() << " entries in "
0217                                   << tileLayerStart_.size() << " layers";
0218     for (unsigned int k = 0; k < tileLayerStart_.size(); ++k)
0219       edm::LogVerbatim("HGCalGeom") << "LayerStart[" << k << "] " << tileLayerStart_[k];
0220     for (unsigned int k = 0; k < tileIndex_.size(); ++k)
0221       edm::LogVerbatim("HGCalGeom") << "[" << k << "] " << tileIndex_[k] << " ("
0222                                     << "Layer " << std::get<0>(HGCalTileIndex::tileUnpack(tileIndex_[k])) << " Ring "
0223                                     << std::get<1>(HGCalTileIndex::tileUnpack(tileIndex_[k])) << ":"
0224                                     << std::get<2>(HGCalTileIndex::tileUnpack(tileIndex_[k])) << ") Phi "
0225                                     << std::get<1>(HGCalTileIndex::tileUnpack(tilePhis_[k])) << ":"
0226                                     << std::get<2>(HGCalTileIndex::tileUnpack(tilePhis_[k]));
0227 #endif
0228     cassette_.setParameter(cassettes_, cassetteShift_);
0229 
0230 #ifdef EDM_ML_DEBUG
0231     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: NameSpace " << ns.name();
0232 
0233     edm::LogVerbatim("HGCalGeom") << "==>> Constructing DDHGCalMixRotatedLayer...";
0234     copies_.clear();
0235 #endif
0236 
0237     double zi(zMinBlock_);
0238     int laymin(0);
0239     for (unsigned int i = 0; i < layers_.size(); i++) {
0240       double zo = zi + layerThick_[i];
0241       double routF = HGCalGeomTools::radius(zi, zFrontT_, rMaxFront_, slopeT_);
0242       int laymax = laymin + layers_[i];
0243       double zz = zi;
0244       double thickTot(0);
0245       for (int ly = laymin; ly < laymax; ++ly) {
0246         int ii = layerType_[ly];
0247         int copy = copyNumber_[ii];
0248         double hthick = 0.5 * thick_[ii];
0249         double rinB = HGCalGeomTools::radius(zo, zFrontB_, rMinFront_, slopeB_);
0250         zz += hthick;
0251         thickTot += thick_[ii];
0252 
0253         std::string name = names_[ii] + std::to_string(copy);
0254 #ifdef EDM_ML_DEBUG
0255         edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: Layer " << ly << ":" << ii << " Front "
0256                                       << cms::convert2mm(zi) << ", " << cms::convert2mm(routF) << " Back "
0257                                       << cms::convert2mm(zo) << ", " << cms::convert2mm(rinB)
0258                                       << " superlayer thickness " << cms::convert2mm(layerThick_[i]);
0259 #endif
0260 
0261         dd4hep::Material matter = ns.material(materials_[ii]);
0262         dd4hep::Volume glog;
0263 
0264         if (layerSense_[ly] < 1) {
0265           std::vector<double> pgonZ, pgonRin, pgonRout;
0266           double rmax =
0267               (std::min(routF, HGCalGeomTools::radius(zz + hthick, zFrontT_, rMaxFront_, slopeT_)) * cosAlpha_) - tol1;
0268           HGCalGeomTools::radius(zz - hthick,
0269                                  zz + hthick,
0270                                  zFrontB_,
0271                                  rMinFront_,
0272                                  slopeB_,
0273                                  zFrontT_,
0274                                  rMaxFront_,
0275                                  slopeT_,
0276                                  -layerSense_[ly],
0277                                  pgonZ,
0278                                  pgonRin,
0279                                  pgonRout);
0280           for (unsigned int isec = 0; isec < pgonZ.size(); ++isec) {
0281             pgonZ[isec] -= zz;
0282             if (layerSense_[ly] == 0 || absorbMode_ == 0)
0283               pgonRout[isec] = rmax;
0284             else
0285               pgonRout[isec] = pgonRout[isec] * cosAlpha_ - tol1;
0286           }
0287 
0288           dd4hep::Solid solid = dd4hep::Polyhedra(sectors_, -alpha_, 2._pi, pgonZ, pgonRin, pgonRout);
0289           ns.addSolidNS(ns.prepend(name), solid);
0290           glog = dd4hep::Volume(solid.name(), solid, matter);
0291           ns.addVolumeNS(glog);
0292 #ifdef EDM_ML_DEBUG
0293           edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << solid.name() << " polyhedra of " << sectors_
0294                                         << " sectors covering " << convertRadToDeg(-alpha_) << ":"
0295                                         << convertRadToDeg(-alpha_ + 2._pi) << " with " << pgonZ.size() << " sections";
0296           for (unsigned int k = 0; k < pgonZ.size(); ++k)
0297             edm::LogVerbatim("HGCalGeom") << "[" << k << "] z " << cms::convert2mm(pgonZ[k]) << " R "
0298                                           << cms::convert2mm(pgonRin[k]) << ":" << cms::convert2mm(pgonRout[k]);
0299 #endif
0300         } else {
0301           double rins =
0302               (sensitiveMode_ < 1) ? rinB : HGCalGeomTools::radius(zz + hthick, zFrontB_, rMinFront_, slopeB_);
0303           double routs =
0304               (sensitiveMode_ < 1) ? routF : HGCalGeomTools::radius(zz - hthick, zFrontT_, rMaxFront_, slopeT_);
0305           dd4hep::Solid solid = dd4hep::Tube(rins, routs, hthick, 0.0, 2._pi);
0306           ns.addSolidNS(ns.prepend(name), solid);
0307           glog = dd4hep::Volume(solid.name(), solid, matter);
0308           ns.addVolumeNS(glog);
0309 
0310 #ifdef EDM_ML_DEBUG
0311           edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << solid.name() << " Tubs made of "
0312                                         << matter.name() << " of dimensions " << cms::convert2mm(rinB) << ":"
0313                                         << cms::convert2mm(rins) << ", " << cms::convert2mm(routF) << ":"
0314                                         << cms::convert2mm(routs) << ", " << cms::convert2mm(hthick)
0315                                         << ", 0.0, 360.0 and positioned in: " << glog.name() << " number " << copy;
0316 #endif
0317           positionMix(ctxt, e, glog, name, copy, thick_[ii], matter);
0318         }
0319 
0320         dd4hep::Position r1(0, 0, zz);
0321         mother.placeVolume(glog, copy, r1);
0322         ++copyNumber_[ii];
0323 #ifdef EDM_ML_DEBUG
0324         edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << glog.name() << " number " << copy
0325                                       << " positioned in " << mother.name() << " at (0,0," << cms::convert2mm(zz)
0326                                       << ") with no rotation";
0327 #endif
0328         zz += hthick;
0329       }  // End of loop over layers in a block
0330       zi = zo;
0331       laymin = laymax;
0332       if (std::abs(thickTot - layerThick_[i]) > tol2_) {
0333         if (thickTot > layerThick_[i]) {
0334           edm::LogError("HGCalGeom") << "Thickness of the partition " << cms::convert2mm(layerThick_[i])
0335                                      << " is smaller than " << cms::convert2mm(thickTot)
0336                                      << ": thickness of all its components **** ERROR ****";
0337         } else {
0338           edm::LogWarning("HGCalGeom") << "Thickness of the partition " << cms::convert2mm(layerThick_[i])
0339                                        << " does not match with " << cms::convert2mm(thickTot) << " of the components";
0340         }
0341       }
0342     }  // End of loop over blocks
0343 
0344 #ifdef EDM_ML_DEBUG
0345     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << copies_.size() << " different wafer copy numbers";
0346     int k(0);
0347     for (std::unordered_set<int>::const_iterator itr = copies_.begin(); itr != copies_.end(); ++itr, ++k) {
0348       edm::LogVerbatim("HGCalGeom") << "Copy [" << k << "] : " << (*itr);
0349     }
0350     copies_.clear();
0351     edm::LogVerbatim("HGCalGeom") << "<<== End of DDHGCalMixRotatedLayer construction...";
0352 #endif
0353   }
0354 
0355   void positionMix(cms::DDParsingContext& ctxt,
0356                    xml_h e,
0357                    const dd4hep::Volume& glog,
0358                    const std::string& nameM,
0359                    int copyM,
0360                    double thick,
0361                    const dd4hep::Material& matter) {
0362     cms::DDNamespace ns(ctxt, e, true);
0363 
0364     // Make the top part first
0365     for (unsigned int ly = 0; ly < layerTypeTop_.size(); ++ly) {
0366       int ii = layerTypeTop_[ly];
0367       copyNumberTop_[ii] = copyM;
0368     }
0369     double hthick = 0.5 * thick;
0370     double dphi = (2._pi) / phiBinsScint_;
0371     double thickTot(0), zpos(-hthick);
0372     for (unsigned int ly = 0; ly < layerTypeTop_.size(); ++ly) {
0373       int ii = layerTypeTop_[ly];
0374       int copy = copyNumberTop_[ii];
0375       int layer = copy - firstLayer_;
0376       double hthickl = 0.5 * layerThickTop_[ii];
0377       thickTot += layerThickTop_[ii];
0378       zpos += hthickl;
0379       dd4hep::Material matter1 = ns.material(materialTop_[ii]);
0380       unsigned int k = 0;
0381       int firstTile = tileLayerStart_[layer];
0382       int lastTile = ((layer + 1 < static_cast<int>(tileLayerStart_.size())) ? tileLayerStart_[layer + 1]
0383                                                                              : static_cast<int>(tileIndex_.size()));
0384 #ifdef EDM_ML_DEBUG
0385       edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: Layer " << ly << ":" << ii << " Copy " << copy
0386                                     << " Tiles " << firstTile << ":" << lastTile;
0387 #endif
0388       for (int ti = firstTile; ti < lastTile; ++ti) {
0389         double r1 = tileRMin_[std::get<1>(HGCalTileIndex::tileUnpack(tileIndex_[ti])) - 1];
0390         double r2 = tileRMax_[std::get<2>(HGCalTileIndex::tileUnpack(tileIndex_[ti])) - 1];
0391         int cassette = std::get<0>(HGCalTileIndex::tileUnpack(tilePhis_[ti]));
0392         int fimin = std::get<1>(HGCalTileIndex::tileUnpack(tilePhis_[ti]));
0393         int fimax = std::get<2>(HGCalTileIndex::tileUnpack(tilePhis_[ti]));
0394         double phi1 = dphi * (fimin - 1);
0395         double phi2 = dphi * (fimax - fimin + 1);
0396         auto cshift = cassette_.getShift(layer + 1, 1, cassette);
0397 #ifdef EDM_ML_DEBUG
0398         edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: Layer " << copy << " iR "
0399                                       << std::get<1>(HGCalTileIndex::tileUnpack(tileIndex_[ly])) << ":"
0400                                       << std::get<2>(HGCalTileIndex::tileUnpack(tileIndex_[ly])) << " R "
0401                                       << cms::convert2mm(r1) << ":" << cms::convert2mm(r2) << " Thick "
0402                                       << cms::convert2mm((2.0 * hthickl)) << " phi " << fimin << ":" << fimax << ":"
0403                                       << convertRadToDeg(phi1) << ":" << convertRadToDeg(phi2) << " cassette "
0404                                       << cassette << " Shift " << cshift.first << ":" << cshift.second;
0405 #endif
0406         std::string name = namesTop_[ii] + "L" + std::to_string(copy) + "F" + std::to_string(k);
0407         ++k;
0408         dd4hep::Solid solid = dd4hep::Tube(r1, r2, hthickl, phi1, phi2);
0409         ns.addSolidNS(ns.prepend(name), solid);
0410         dd4hep::Volume glog1 = dd4hep::Volume(solid.name(), solid, matter1);
0411         ns.addVolumeNS(glog1);
0412 #ifdef EDM_ML_DEBUG
0413         edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << glog1.name() << " Tubs made of "
0414                                       << materialTop_[ii] << " of dimensions " << cms::convert2mm(r1) << ", "
0415                                       << cms::convert2mm(r2) << ", " << cms::convert2mm(hthickl) << ", "
0416                                       << convertRadToDeg(phi1) << ", " << convertRadToDeg(phi2);
0417 #endif
0418         dd4hep::Position tran(cshift.first, cshift.second, zpos);
0419         glog.placeVolume(glog1, copy, tran);
0420 #ifdef EDM_ML_DEBUG
0421         edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: Position " << glog1.name() << " number " << copy
0422                                       << " in " << glog.name() << " at (" << cms::convert2mm(cshift.first) << ", "
0423                                       << cms::convert2mm(cshift.second) << ", " << cms::convert2mm(zpos)
0424                                       << ") with no rotation";
0425 #endif
0426       }
0427       ++copyNumberTop_[ii];
0428       zpos += hthickl;
0429     }
0430     if (std::abs(thickTot - thick) > tol2_) {
0431       if (thickTot > thick) {
0432         edm::LogError("HGCalGeom") << "Thickness of the partition " << cms::convert2mm(thick) << " is smaller than "
0433                                    << cms::convert2mm(thickTot)
0434                                    << ": thickness of all its components in the top part **** ERROR ****";
0435       } else {
0436         edm::LogWarning("HGCalGeom") << "Thickness of the partition " << cms::convert2mm(thick)
0437                                      << " does not match with " << cms::convert2mm(thickTot)
0438                                      << " of the components in top part";
0439       }
0440     }
0441 
0442     // Make the bottom part next
0443     int layer = (copyM - firstLayer_);
0444     static const double sqrt3 = std::sqrt(3.0);
0445     int layercenter = layerOrient_[layer];
0446     int layertype = HGCalTypes::layerFrontBack(layerOrient_[layer]);
0447     int firstWafer = waferLayerStart_[layer];
0448     int lastWafer = ((layer + 1 < static_cast<int>(waferLayerStart_.size())) ? waferLayerStart_[layer + 1]
0449                                                                              : static_cast<int>(waferIndex_.size()));
0450     double delx = 0.5 * (waferSize_ + waferSepar_);
0451     double dely = 2.0 * delx / sqrt3;
0452     double dy = 0.75 * dely;
0453     const auto& xyoff = geomTools_.shiftXY(layercenter, (waferSize_ + waferSepar_));
0454 #ifdef EDM_ML_DEBUG
0455     int ium(0), ivm(0), kount(0);
0456     std::vector<int> ntype(3, 0);
0457     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << glog.name() << "  r " << cms::convert2mm(delx)
0458                                   << " R " << cms::convert2mm(dely) << " dy " << cms::convert2mm(dy) << " Shift "
0459                                   << cms::convert2mm(xyoff.first) << ":" << cms::convert2mm(xyoff.second)
0460                                   << " WaferSize " << cms::convert2mm((waferSize_ + waferSepar_)) << " index "
0461                                   << firstWafer << ":" << (lastWafer - 1);
0462 #endif
0463     for (int k = firstWafer; k < lastWafer; ++k) {
0464       int u = HGCalWaferIndex::waferU(waferIndex_[k]);
0465       int v = HGCalWaferIndex::waferV(waferIndex_[k]);
0466 #ifdef EDM_ML_DEBUG
0467       int iu = std::abs(u);
0468       int iv = std::abs(v);
0469 #endif
0470       int nr = 2 * v;
0471       int nc = -2 * u + v;
0472       int type = HGCalProperty::waferThick(waferProperty_[k]);
0473       int part = HGCalProperty::waferPartial(waferProperty_[k]);
0474       int orien = HGCalProperty::waferOrient(waferProperty_[k]);
0475       int cassette = HGCalProperty::waferCassette(waferProperty_[k]);
0476       int place = HGCalCell::cellPlacementIndex(1, layertype, orien);
0477 #ifdef EDM_ML_DEBUG
0478       edm::LogVerbatim("HGCalGeom") << " index:Property:layertype:type:part:orien:cassette:place:offsets:ind " << k
0479                                     << ":" << waferProperty_[k] << ":" << layertype << ":" << type << ":" << part << ":"
0480                                     << orien << ":" << cassette << ":" << place;
0481 #endif
0482       auto cshift = cassette_.getShift(layer + 1, 1, cassette);
0483       double xpos = xyoff.first + cshift.first + nc * delx;
0484       double ypos = xyoff.second + cshift.second + nr * dy;
0485       std::string wafer;
0486       int i(999);
0487       if (part == HGCalTypes::WaferFull) {
0488         i = type * facingTypes_ * orientationTypes_ + place - placeOffset_;
0489 #ifdef EDM_ML_DEBUG
0490         edm::LogVerbatim("HGCalGeom") << " FullWafer type:place:ind " << type << ":" << place << ":" << i << ":"
0491                                       << waferFull_.size();
0492 #endif
0493         wafer = waferFull_[i];
0494       } else {
0495         int partoffset =
0496             (part >= HGCalTypes::WaferHDTop) ? HGCalTypes::WaferPartHDOffset : HGCalTypes::WaferPartLDOffset;
0497         i = (part - partoffset) * facingTypes_ * orientationTypes_ +
0498             HGCalTypes::WaferTypeOffset[type] * facingTypes_ * orientationTypes_ + place - placeOffset_;
0499 #ifdef EDM_ML_DEBUG
0500         edm::LogVerbatim("HGCalGeom") << " layertype:type:part:orien:cassette:place:offsets:ind " << layertype << ":"
0501                                       << type << ":" << part << ":" << orien << ":" << cassette << ":" << place << ":"
0502                                       << partoffset << ":" << HGCalTypes::WaferTypeOffset[type] << ":" << i << ":"
0503                                       << waferPart_.size();
0504 #endif
0505         wafer = waferPart_[i];
0506       }
0507       int copy = HGCalTypes::packTypeUV(type, u, v);
0508 #ifdef EDM_ML_DEBUG
0509       edm::LogVerbatim("HGCalGeom") << " DDHGCalMixRotatedLayer: Layer " << HGCalWaferIndex::waferLayer(waferIndex_[k])
0510                                     << " Wafer " << wafer << " number " << copy << " type :part:orien:ind " << type
0511                                     << ":" << part << ":" << orien << ":" << i << " layer:u:v " << (layer + firstLayer_)
0512                                     << ":" << u << ":" << v;
0513       if (iu > ium)
0514         ium = iu;
0515       if (iv > ivm)
0516         ivm = iv;
0517       kount++;
0518       if (copies_.count(copy) == 0)
0519         copies_.insert(copy);
0520 #endif
0521       dd4hep::Position tran(xpos, ypos, 0.0);
0522       glog.placeVolume(ns.volume(wafer), copy, tran);
0523 #ifdef EDM_ML_DEBUG
0524       ++ntype[type];
0525       edm::LogVerbatim("HGCalGeom") << " DDHGCalMixRotatedLayer: " << wafer << " number " << copy << " type "
0526                                     << layertype << ":" << type << " positioned in " << glog.name() << " at ("
0527                                     << cms::convert2mm(xpos) << "," << cms::convert2mm(ypos) << ",0) with no rotation";
0528 #endif
0529     }
0530 
0531 #ifdef EDM_ML_DEBUG
0532     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: Maximum # of u " << ium << " # of v " << ivm << " and "
0533                                   << kount << " wafers (" << ntype[0] << ":" << ntype[1] << ":" << ntype[2] << ") for "
0534                                   << glog.name();
0535 #endif
0536   }
0537 
0538   //Required data members to cache the values from XML file
0539   HGCalGeomTools geomTools_;
0540   HGCalCassette cassette_;
0541 
0542   static constexpr double tol2_ = 0.00001 * dd4hep::mm;
0543 
0544   int waferTypes_;                        // Number of wafer types
0545   int facingTypes_;                       // Types of facings of modules toward IP
0546   int orientationTypes_;                  // Number of partial wafer orienations
0547   int placeOffset_;                       // Offset for placement
0548   int phiBinsScint_;                      // Maximum number of cells along phi
0549   int firstLayer_;                        // Copy # of the first sensitive layer
0550   int absorbMode_;                        // Absorber mode
0551   int sensitiveMode_;                     // Sensitive mode
0552   double zMinBlock_;                      // Starting z-value of the block
0553   double waferSize_;                      // Width of the wafer
0554   double waferSepar_;                     // Sensor separation
0555   int sectors_;                           // Sectors
0556   int cassettes_;                         // Cassettes
0557   std::vector<double> slopeB_;            // Slope at the lower R
0558   std::vector<double> zFrontB_;           // Starting Z values for the slopes
0559   std::vector<double> rMinFront_;         // Corresponding rMin's
0560   std::vector<double> slopeT_;            // Slopes at the larger R
0561   std::vector<double> zFrontT_;           // Starting Z values for the slopes
0562   std::vector<double> rMaxFront_;         // Corresponding rMax's
0563   std::vector<std::string> waferFull_;    // Names of full wafer modules
0564   std::vector<std::string> waferPart_;    // Names of partial wafer modules
0565   std::vector<std::string> materials_;    // Materials
0566   std::vector<std::string> names_;        // Names
0567   std::vector<double> thick_;             // Thickness of the material
0568   std::vector<int> copyNumber_;           // Initial copy numbers
0569   std::vector<int> layers_;               // Number of layers in a section
0570   std::vector<double> layerThick_;        // Thickness of each section
0571   std::vector<int> layerType_;            // Type of the layer
0572   std::vector<int> layerSense_;           // Content of a layer (sensitive?)
0573   std::vector<std::string> materialTop_;  // Materials of top layers
0574   std::vector<std::string> namesTop_;     // Names of top layers
0575   std::vector<double> layerThickTop_;     // Thickness of the top sections
0576   std::vector<int> layerTypeTop_;         // Type of the Top layer
0577   std::vector<int> copyNumberTop_;        // Initial copy numbers (top section)
0578   std::vector<int> layerOrient_;          // Layer orientation for the silicon component
0579   std::vector<int> waferIndex_;           // Wafer index for the types
0580   std::vector<int> waferProperty_;        // Wafer property
0581   std::vector<int> waferLayerStart_;      // Start index of wafers in each layer
0582   std::vector<double> cassetteShift_;     // Shifts of the cassetes
0583   std::vector<double> tileRMin_;          // Minimum radius of each ring
0584   std::vector<double> tileRMax_;          // Maximum radius of each ring
0585   std::vector<int> tileIndex_;            // Index of tile (layer/start|end ring)
0586   std::vector<int> tilePhis_;             // Tile phi range for each index
0587   std::vector<int> tileLayerStart_;       // Start index of tiles in each layer
0588   std::unordered_set<int> copies_;        // List of copy #'s
0589   double alpha_, cosAlpha_;
0590 };
0591 
0592 static long algorithm(dd4hep::Detector& /* description */, cms::DDParsingContext& ctxt, xml_h e) {
0593   HGCalMixRotatedLayer healgo(ctxt, e);
0594   return cms::s_executed;
0595 }
0596 
0597 DECLARE_DDCMS_DETELEMENT(DDCMS_hgcal_DDHGCalMixRotatedLayer, algorithm)