Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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