Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-08-22 04:57:42

0001 ///////////////////////////////////////////////////////////////////////////////
0002 // File: DDHGCalMixRotatedFineFineCassette.cc
0003 // Description: Geometry factory class for HGCal (Mix)
0004 ///////////////////////////////////////////////////////////////////////////////
0005 
0006 #include "DataFormats/Math/interface/angle_units.h"
0007 #include "DetectorDescription/DDCMS/interface/DDPlugins.h"
0008 #include "DetectorDescription/DDCMS/interface/DDutils.h"
0009 #include "DD4hep/DetFactoryHelper.h"
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011 #include "Geometry/HGCalCommonData/interface/HGCalCell.h"
0012 #include "Geometry/HGCalCommonData/interface/HGCalCassette.h"
0013 #include "Geometry/HGCalCommonData/interface/HGCalGeomTools.h"
0014 #include "Geometry/HGCalCommonData/interface/HGCalParameters.h"
0015 #include "Geometry/HGCalCommonData/interface/HGCalProperty.h"
0016 #include "Geometry/HGCalCommonData/interface/HGCalTileIndex.h"
0017 #include "Geometry/HGCalCommonData/interface/HGCalTypes.h"
0018 #include "Geometry/HGCalCommonData/interface/HGCalWaferIndex.h"
0019 #include "Geometry/HGCalCommonData/interface/HGCalWaferType.h"
0020 
0021 #include <cmath>
0022 #include <memory>
0023 #include <string>
0024 #include <unordered_set>
0025 #include <vector>
0026 
0027 //#define EDM_ML_DEBUG
0028 using namespace angle_units::operators;
0029 
0030 struct HGCalMixRotatedFineCassette {
0031   HGCalMixRotatedFineCassette() {
0032     throw cms::Exception("HGCalGeom") << "Wrong initialization to HGCalMixRotatedFineCassette";
0033   }
0034 
0035   HGCalMixRotatedFineCassette(cms::DDParsingContext& ctxt, xml_h e) {
0036     cms::DDNamespace ns(ctxt, e, true);
0037     cms::DDAlgoArguments args(ctxt, e);
0038 
0039 #ifdef EDM_ML_DEBUG
0040     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: Creating an instance";
0041 #endif
0042     static constexpr double tol1 = 0.01 * dd4hep::mm;
0043 
0044     dd4hep::Volume mother = ns.volume(args.parentName());
0045 
0046     waferTypes_ = args.value<int>("WaferTypes");
0047     passiveTypes_ = args.value<int>("PassiveTypes");
0048     facingTypes_ = args.value<int>("FacingTypes");
0049     orientationTypes_ = args.value<int>("OrientationTypes");
0050     partialTypes_ = args.value<int>("PartialTypes");
0051     placeOffset_ = args.value<int>("PlaceOffset");
0052     phiBinsScint_ = args.value<int>("NPhiBinScint");
0053     phiBinsFineScint_ = args.value<int>("NPhiBinFineScint");
0054 #ifdef EDM_ML_DEBUG
0055     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette::Number of types of wafers: " << waferTypes_
0056                                   << " passives: " << passiveTypes_ << " facings: " << facingTypes_
0057                                   << " Orientations: " << orientationTypes_ << " PartialTypes: " << partialTypes_
0058                                   << " PlaceOffset: " << placeOffset_ << "; number of cells along phi "
0059                                   << phiBinsFineScint_ << ":" << phiBinsScint_;
0060 #endif
0061     firstFineLayer_ = args.value<int>("FirstFineLayer");
0062     firstCoarseLayer_ = args.value<int>("FirstCoarseLayer");
0063     absorbMode_ = args.value<int>("AbsorberMode");
0064     sensitiveMode_ = args.value<int>("SensitiveMode");
0065     passiveMode_ = args.value<int>("PassiveMode");
0066 #ifdef EDM_ML_DEBUG
0067     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette::First Layers " << firstFineLayer_ << ":"
0068                                   << firstCoarseLayer_ << " and Absober:Sensitive mode " << absorbMode_ << ":"
0069                                   << sensitiveMode_ << ":" << passiveMode_;
0070 #endif
0071     zMinBlock_ = args.value<double>("zMinBlock");
0072     waferSize_ = args.value<double>("waferSize");
0073     waferSepar_ = args.value<double>("SensorSeparation");
0074     sectors_ = args.value<int>("Sectors");
0075     cassettes_ = args.value<int>("Cassettes");
0076     alpha_ = (1._pi) / sectors_;
0077     cosAlpha_ = cos(alpha_);
0078 #ifdef EDM_ML_DEBUG
0079     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: zStart " << zMinBlock_ << " wafer width "
0080                                   << waferSize_ << " separations " << waferSepar_ << " sectors " << sectors_ << ":"
0081                                   << convertRadToDeg(alpha_) << ":" << cosAlpha_ << " with " << cassettes_
0082                                   << " cassettes";
0083 #endif
0084     slopeB_ = args.value<std::vector<double>>("SlopeBottom");
0085     zFrontB_ = args.value<std::vector<double>>("ZFrontBottom");
0086     rMinFront_ = args.value<std::vector<double>>("RMinFront");
0087     slopeT_ = args.value<std::vector<double>>("SlopeTop");
0088     zFrontT_ = args.value<std::vector<double>>("ZFrontTop");
0089     rMaxFront_ = args.value<std::vector<double>>("RMaxFront");
0090 #ifdef EDM_ML_DEBUG
0091     for (unsigned int i = 0; i < slopeB_.size(); ++i)
0092       edm::LogVerbatim("HGCalGeom") << "Bottom Block [" << i << "] Zmin " << zFrontB_[i] << " Rmin " << rMinFront_[i]
0093                                     << " Slope " << slopeB_[i];
0094     for (unsigned int i = 0; i < slopeT_.size(); ++i)
0095       edm::LogVerbatim("HGCalGeom") << "Top Block [" << i << "] Zmin " << zFrontT_[i] << " Rmax " << rMaxFront_[i]
0096                                     << " Slope " << slopeT_[i];
0097 #endif
0098     waferFull_ = args.value<std::vector<std::string>>("WaferNamesFull");
0099     waferPart_ = args.value<std::vector<std::string>>("WaferNamesPartial");
0100 #ifdef EDM_ML_DEBUG
0101     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: " << waferFull_.size() << " full and "
0102                                   << waferPart_.size()
0103                                   << " partial modules\nDDHGCalMixRotatedFineCassette:Full Modules:";
0104     unsigned int i1max = static_cast<unsigned int>(waferFull_.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 << "] " << waferFull_[i];
0110       edm::LogVerbatim("HGCalGeom") << st1.str();
0111     }
0112     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: Partial Modules:";
0113     i1max = static_cast<unsigned int>(waferPart_.size());
0114     for (unsigned int i1 = 0; i1 < i1max; i1 += 2) {
0115       std::ostringstream st1;
0116       unsigned int i2 = std::min((i1 + 2), i1max);
0117       for (unsigned int i = i1; i < i2; ++i)
0118         st1 << " [" << i << "] " << waferPart_[i];
0119       edm::LogVerbatim("HGCalGeom") << st1.str();
0120     }
0121 #endif
0122     passiveFull_ = args.value<std::vector<std::string>>("PassiveNamesFull");
0123     passivePart_ = args.value<std::vector<std::string>>("PassiveNamesPartial");
0124 #ifdef EDM_ML_DEBUG
0125     edm::LogVerbatim("HGCalGeom") << "DDHGCalSiliconRotatedCassette: " << passiveFull_.size() << " full and "
0126                                   << passivePart_.size() << " partial passive modules";
0127     i1max = static_cast<unsigned int>(passiveFull_.size());
0128     for (unsigned int i1 = 0; i1 < i1max; i1 += 2) {
0129       std::ostringstream st1;
0130       unsigned int i2 = std::min((i1 + 2), i1max);
0131       for (unsigned int i = i1; i < i2; ++i)
0132         st1 << " [" << i << "] " << passiveFull_[i];
0133       edm::LogVerbatim("HGCalGeom") << st1.str();
0134     }
0135     edm::LogVerbatim("HGCalGeom") << "DDHGCalSiliconRotatedCassette: Partial Modules:";
0136     i1max = static_cast<unsigned int>(passivePart_.size());
0137     for (unsigned int i1 = 0; i1 < i1max; i1 += 2) {
0138       std::ostringstream st1;
0139       unsigned int i2 = std::min((i1 + 2), i1max);
0140       for (unsigned int i = i1; i < i2; ++i)
0141         st1 << " [" << i << "] " << passivePart_[i];
0142       edm::LogVerbatim("HGCalGeom") << st1.str();
0143     }
0144 #endif
0145     materials_ = args.value<std::vector<std::string>>("MaterialNames");
0146     names_ = args.value<std::vector<std::string>>("VolumeNames");
0147     thick_ = args.value<std::vector<double>>("Thickness");
0148     copyNumber_.resize(materials_.size(), 1);
0149 #ifdef EDM_ML_DEBUG
0150     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: " << materials_.size() << " types of volumes";
0151     for (unsigned int i = 0; i < names_.size(); ++i)
0152       edm::LogVerbatim("HGCalGeom") << "Volume [" << i << "] " << names_[i] << " of thickness " << thick_[i]
0153                                     << " filled with " << materials_[i] << " first copy number " << copyNumber_[i];
0154 #endif
0155     layers_ = args.value<std::vector<int>>("Layers");
0156     layerThick_ = args.value<std::vector<double>>("LayerThick");
0157 #ifdef EDM_ML_DEBUG
0158     edm::LogVerbatim("HGCalGeom") << "There are " << layers_.size() << " blocks";
0159     for (unsigned int i = 0; i < layers_.size(); ++i)
0160       edm::LogVerbatim("HGCalGeom") << "Block [" << i << "] of thickness " << layerThick_[i] << " with " << layers_[i]
0161                                     << " layers";
0162 #endif
0163     layerType_ = args.value<std::vector<int>>("LayerType");
0164     layerSense_ = args.value<std::vector<int>>("LayerSense");
0165     layerOrient_ = args.value<std::vector<int>>("LayerTypes");
0166     for (unsigned int k = 0; k < layerOrient_.size(); ++k)
0167       layerOrient_[k] = HGCalTypes::layerType(layerOrient_[k]);
0168 #ifdef EDM_ML_DEBUG
0169     for (unsigned int i = 0; i < layerOrient_.size(); ++i)
0170       edm::LogVerbatim("HGCalGeom") << "LayerOrient [" << i << "] " << layerOrient_[i];
0171 #endif
0172     if (firstFineLayer_ > 0) {
0173       for (unsigned int i = 0; i < layerType_.size(); ++i) {
0174         if (layerSense_[i] != 0) {
0175           int ii = layerType_[i];
0176           copyNumber_[ii] = firstFineLayer_;
0177 #ifdef EDM_ML_DEBUG
0178           edm::LogVerbatim("HGCalGeom") << "First copy number for layer type " << i << ":" << ii << " with "
0179                                         << materials_[ii] << " changed to " << copyNumber_[ii];
0180 #endif
0181         }
0182       }
0183     } else {
0184       firstFineLayer_ = 1;
0185     }
0186 #ifdef EDM_ML_DEBUG
0187     edm::LogVerbatim("HGCalGeom") << "There are " << layerType_.size() << " layers";
0188     for (unsigned int i = 0; i < layerType_.size(); ++i)
0189       edm::LogVerbatim("HGCalGeom") << "Layer [" << i << "] with material type " << layerType_[i] << " sensitive class "
0190                                     << layerSense_[i];
0191 #endif
0192     materialTop_ = args.value<std::vector<std::string>>("TopMaterialNames");
0193     namesTop_ = args.value<std::vector<std::string>>("TopVolumeNames");
0194     layerThickTop_ = args.value<std::vector<double>>("TopLayerThickness");
0195     layerTypeTop_ = args.value<std::vector<int>>("TopLayerType");
0196     copyNumberTop_.resize(materialTop_.size(), firstFineLayer_);
0197     coverTypeTop_ = args.value<int>("TopCoverLayerType");
0198     coverTopLayers_ = args.value<int>("TopCoverLayers");
0199     copyNumberCoverTop_.resize(coverTopLayers_, firstFineLayer_);
0200 #ifdef EDM_ML_DEBUG
0201     std::ostringstream st0;
0202     for (int k = 0; k < coverTopLayers_; ++k)
0203       st0 << " " << copyNumberCoverTop_[k];
0204     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: " << materialTop_.size()
0205                                   << " types of volumes in the top part; " << coverTopLayers_ << " covers of Type "
0206                                   << coverTypeTop_ << " with initial copy numbers: " << st0.str();
0207     for (unsigned int i = 0; i < materialTop_.size(); ++i)
0208       edm::LogVerbatim("HGCalGeom") << "Volume [" << i << "] " << namesTop_[i] << " of thickness " << layerThickTop_[i]
0209                                     << " filled with " << materialTop_[i] << " first copy number " << copyNumberTop_[i];
0210     edm::LogVerbatim("HGCalGeom") << "There are " << layerTypeTop_.size() << " layers in the top part";
0211     for (unsigned int i = 0; i < layerTypeTop_.size(); ++i)
0212       edm::LogVerbatim("HGCalGeom") << "Layer [" << i << "] with material type " << layerTypeTop_[i];
0213 #endif
0214     waferIndex_ = args.value<std::vector<int>>("WaferIndex");
0215     waferProperty_ = args.value<std::vector<int>>("WaferProperties");
0216     waferLayerStart_ = args.value<std::vector<int>>("WaferLayerStart");
0217     cassetteShift_ = args.value<std::vector<double>>("CassetteShift");
0218 #ifdef EDM_ML_DEBUG
0219     edm::LogVerbatim("HGCalGeom") << "waferProperties with " << waferIndex_.size() << " entries in "
0220                                   << waferLayerStart_.size() << " layers";
0221     for (unsigned int k = 0; k < waferLayerStart_.size(); ++k)
0222       edm::LogVerbatim("HGCalGeom") << "LayerStart[" << k << "] " << waferLayerStart_[k];
0223     for (unsigned int k = 0; k < waferIndex_.size(); ++k)
0224       edm::LogVerbatim("HGCalGeom") << "[" << k << "] " << waferIndex_[k] << " ("
0225                                     << HGCalWaferIndex::waferLayer(waferIndex_[k]) << ", "
0226                                     << HGCalWaferIndex::waferU(waferIndex_[k]) << ", "
0227                                     << HGCalWaferIndex::waferV(waferIndex_[k]) << ") : ("
0228                                     << HGCalProperty::waferThick(waferProperty_[k]) << ":"
0229                                     << HGCalProperty::waferPartial(waferProperty_[k]) << ":"
0230                                     << HGCalProperty::waferOrient(waferProperty_[k]) << ")";
0231     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: " << cassetteShift_.size()
0232                                   << " elements for cassette shifts";
0233     unsigned int j1max = cassetteShift_.size();
0234     for (unsigned int j1 = 0; j1 < j1max; j1 += 6) {
0235       std::ostringstream st1;
0236       unsigned int j2 = std::min((j1 + 6), j1max);
0237       for (unsigned int j = j1; j < j2; ++j)
0238         st1 << " [" << j << "] " << std::setw(9) << cassetteShift_[j];
0239       edm::LogVerbatim("HGCalGeom") << st1.str();
0240     }
0241 #endif
0242     tileFineRMin_ = args.value<std::vector<double>>("Tile6RMin");
0243     tileFineRMax_ = args.value<std::vector<double>>("Tile6RMax");
0244     tileFineIndex_ = args.value<std::vector<int>>("Tile6LayerRings");
0245     tileFinePhis_ = args.value<std::vector<int>>("Tile6PhiRange");
0246     tileFineLayerStart_ = args.value<std::vector<int>>("Tile6LayerStart");
0247     tileCoarseRMin_ = args.value<std::vector<double>>("TileRMin");
0248     tileCoarseRMax_ = args.value<std::vector<double>>("TileRMax");
0249     tileCoarseIndex_ = args.value<std::vector<int>>("TileLayerRings");
0250     tileCoarsePhis_ = args.value<std::vector<int>>("TilePhiRange");
0251     tileCoarseLayerStart_ = args.value<std::vector<int>>("TileLayerStart");
0252 #ifdef EDM_ML_DEBUG
0253     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette:: with " << tileFineRMin_.size() << ":"
0254                                   << tileCoarseRMin_.size() << " rings";
0255     for (unsigned int k = 0; k < tileFineRMin_.size(); ++k)
0256       edm::LogVerbatim("HGCalGeom") << "Fine Ring[" << k << "] " << tileFineRMin_[k] << " : " << tileFineRMax_[k];
0257     for (unsigned int k = 0; k < tileCoarseRMin_.size(); ++k)
0258       edm::LogVerbatim("HGCalGeom") << "Coarse Ring[" << k << "] " << tileCoarseRMin_[k] << " : " << tileCoarseRMax_[k];
0259     edm::LogVerbatim("HGCalGeom") << "TileProperties with " << tileFineIndex_.size() << ":" << tileCoarseIndex_.size()
0260                                   << " entries in " << tileFineLayerStart_.size() << ":" << tileCoarseLayerStart_.size()
0261                                   << " layers";
0262     for (unsigned int k = 0; k < tileFineLayerStart_.size(); ++k)
0263       edm::LogVerbatim("HGCalGeom") << "FineLayerStart[" << k << "] " << tileFineLayerStart_[k];
0264     for (unsigned int k = 0; k < tileCoarseLayerStart_.size(); ++k)
0265       edm::LogVerbatim("HGCalGeom") << "CoarseLayerStart[" << k << "] " << tileCoarseLayerStart_[k];
0266     for (unsigned int k = 0; k < tileFineIndex_.size(); ++k)
0267       edm::LogVerbatim("HGCalGeom") << "[" << k << "] " << tileFineIndex_[k] << " ("
0268                                     << "Layer " << std::get<0>(HGCalTileIndex::tileUnpack(tileFineIndex_[k]))
0269                                     << " Ring " << std::get<1>(HGCalTileIndex::tileUnpack(tileFineIndex_[k])) << ":"
0270                                     << std::get<2>(HGCalTileIndex::tileUnpack(tileFineIndex_[k])) << ") Phi "
0271                                     << std::get<1>(HGCalTileIndex::tileUnpack(tileFinePhis_[k])) << ":"
0272                                     << std::get<2>(HGCalTileIndex::tileUnpack(tileFinePhis_[k]));
0273     for (unsigned int k = 0; k < tileCoarseIndex_.size(); ++k)
0274       edm::LogVerbatim("HGCalGeom") << "[" << k << "] " << tileCoarseIndex_[k] << " ("
0275                                     << "Layer " << std::get<0>(HGCalTileIndex::tileUnpack(tileCoarseIndex_[k]))
0276                                     << " Ring " << std::get<1>(HGCalTileIndex::tileUnpack(tileCoarseIndex_[k])) << ":"
0277                                     << std::get<2>(HGCalTileIndex::tileUnpack(tileCoarseIndex_[k])) << ") Phi "
0278                                     << std::get<1>(HGCalTileIndex::tileUnpack(tileCoarsePhis_[k])) << ":"
0279                                     << std::get<2>(HGCalTileIndex::tileUnpack(tileCoarsePhis_[k]));
0280 #endif
0281     std::vector<double> retract = args.value<std::vector<double>>("ScintRetract");
0282     double dphi = M_PI / cassettes_;
0283     for (int k = 0; k < cassettes_; ++k) {
0284       double phi = (2 * k + 1) * dphi;
0285       cassetteShiftScnt_.emplace_back(retract[k] * cos(phi));
0286       cassetteShiftScnt_.emplace_back(retract[k] * sin(phi));
0287     }
0288 #ifdef EDM_ML_DEBUG
0289     unsigned int j2max = cassetteShiftScnt_.size();
0290     for (unsigned int j1 = 0; j1 < j2max; j1 += 6) {
0291       std::ostringstream st1;
0292       unsigned int j2 = std::min((j1 + 6), j2max);
0293       for (unsigned int j = j1; j < j2; ++j)
0294         st1 << " [" << j << "] " << std::setw(9) << cassetteShiftScnt_[j];
0295       edm::LogVerbatim("HGCalGeom") << st1.str();
0296     }
0297 #endif
0298     cassette_.setParameter(cassettes_, cassetteShift_, false);
0299     cassette_.setParameterScint(cassetteShiftScnt_);
0300 
0301     ////////////////////////////////////////////////////////////////////
0302     // DDHGCalMixRotatedFineCassette methods...
0303     ////////////////////////////////////////////////////////////////////
0304 
0305 #ifdef EDM_ML_DEBUG
0306     edm::LogVerbatim("HGCalGeom") << "==>> Constructing DDHGCalMixRotatedFineCassette...";
0307     copies_.clear();
0308 #endif
0309 
0310     double zi(zMinBlock_);
0311     int laymin(0);
0312     unsigned int fineLayers = static_cast<unsigned int>(firstCoarseLayer_ - firstFineLayer_);
0313     for (unsigned int i = 0; i < layers_.size(); i++) {
0314       double zo = zi + layerThick_[i];
0315       double routF = HGCalGeomTools::radius(zi, zFrontT_, rMaxFront_, slopeT_);
0316       int laymax = laymin + layers_[i];
0317       double zz = zi;
0318       double thickTot(0);
0319       bool fine = (i < fineLayers) ? true : false;
0320       for (int ly = laymin; ly < laymax; ++ly) {
0321         int ii = layerType_[ly];
0322         int copy = copyNumber_[ii];
0323         double hthick = 0.5 * thick_[ii];
0324         double rinB = HGCalGeomTools::radius(zo, zFrontB_, rMinFront_, slopeB_);
0325         zz += hthick;
0326         thickTot += thick_[ii];
0327 
0328         std::string name = names_[ii] + std::to_string(copy);
0329 #ifdef EDM_ML_DEBUG
0330         edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: Layer " << ly << ":" << ii << " Front "
0331                                       << cms::convert2mm(zi) << ", " << cms::convert2mm(routF) << " Back "
0332                                       << cms::convert2mm(zo) << ", " << cms::convert2mm(rinB)
0333                                       << " superlayer thickness " << layerThick_[i];
0334 #endif
0335         dd4hep::Material matter = ns.material(materials_[ii]);
0336         dd4hep::Volume glog;
0337         if (layerSense_[ly] == 0) {
0338           std::vector<double> pgonZ, pgonRin, pgonRout;
0339           double rmax =
0340               (std::min(routF, HGCalGeomTools::radius(zz + hthick, zFrontT_, rMaxFront_, slopeT_)) * cosAlpha_) - tol1;
0341           HGCalGeomTools::radius(zz - hthick,
0342                                  zz + hthick,
0343                                  zFrontB_,
0344                                  rMinFront_,
0345                                  slopeB_,
0346                                  zFrontT_,
0347                                  rMaxFront_,
0348                                  slopeT_,
0349                                  -layerSense_[ly],
0350                                  pgonZ,
0351                                  pgonRin,
0352                                  pgonRout);
0353           for (unsigned int isec = 0; isec < pgonZ.size(); ++isec) {
0354             pgonZ[isec] -= zz;
0355             if (layerSense_[ly] == 0 || absorbMode_ == 0)
0356               pgonRout[isec] = rmax;
0357             else
0358               pgonRout[isec] = pgonRout[isec] * cosAlpha_ - tol1;
0359           }
0360           dd4hep::Solid solid = dd4hep::Polyhedra(sectors_, -alpha_, 2._pi, pgonZ, pgonRin, pgonRout);
0361           ns.addSolidNS(ns.prepend(name), solid);
0362           glog = dd4hep::Volume(solid.name(), solid, matter);
0363           ns.addVolumeNS(glog);
0364 #ifdef EDM_ML_DEBUG
0365           edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: " << solid.name() << " polyhedra of "
0366                                         << sectors_ << " sectors covering " << convertRadToDeg(-alpha_) << ":"
0367                                         << convertRadToDeg(-alpha_ + 2._pi) << " with " << pgonZ.size() << " sections";
0368           for (unsigned int k = 0; k < pgonZ.size(); ++k)
0369             edm::LogVerbatim("HGCalGeom") << "[" << k << "] z " << cms::convert2mm(pgonZ[k]) << " R "
0370                                           << cms::convert2mm(pgonRin[k]) << ":" << cms::convert2mm(pgonRout[k]);
0371 #endif
0372         } else {
0373           int mode = (layerSense_[ly] > 0) ? sensitiveMode_ : absorbMode_;
0374           double rins = (mode < 1) ? rinB : HGCalGeomTools::radius(zz + hthick, zFrontB_, rMinFront_, slopeB_);
0375           double routs = (mode < 1) ? routF : HGCalGeomTools::radius(zz - hthick, zFrontT_, rMaxFront_, slopeT_);
0376           dd4hep::Solid solid = dd4hep::Tube(rins, routs, hthick, 0.0, 2._pi);
0377           ns.addSolidNS(ns.prepend(name), solid);
0378           glog = dd4hep::Volume(solid.name(), solid, matter);
0379           ns.addVolumeNS(glog);
0380 #ifdef EDM_ML_DEBUG
0381           edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: " << solid.name() << " Tubs made of "
0382                                         << matter.name() << " of dimensions " << cms::convert2mm(rinB) << ":"
0383                                         << cms::convert2mm(rins) << ", " << cms::convert2mm(routF) << ":"
0384                                         << cms::convert2mm(routs) << ", " << cms::convert2mm(hthick)
0385                                         << ", 0.0, 360.0 and positioned in: " << glog.name() << " number " << copy;
0386 #endif
0387           positionMix(ctxt, e, glog, name, copy, thick_[ii], matter, -layerSense_[ly], fine);
0388         }
0389         dd4hep::Position r1(0, 0, zz);
0390         mother.placeVolume(glog, copy, r1);
0391         ++copyNumber_[ii];
0392 #ifdef EDM_ML_DEBUG
0393         edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: " << glog.name() << " number " << copy
0394                                       << " positioned in " << mother.name() << " at (0,0," << cms::convert2mm(zz)
0395                                       << ") with no rotation";
0396 #endif
0397         zz += hthick;
0398       }  // End of loop over layers in a block
0399       zi = zo;
0400       laymin = laymax;
0401       // Make consistency check of all the partitions of the block
0402       if (std::abs(thickTot - layerThick_[i]) >= tol2_) {
0403         if (thickTot > layerThick_[i]) {
0404           edm::LogError("HGCalGeom") << "Thickness of the partition " << cms::convert2mm(layerThick_[i])
0405                                      << " is smaller than " << cms::convert2mm(thickTot)
0406                                      << ": thickness of all its components **** ERROR ****";
0407         } else {
0408           edm::LogWarning("HGCalGeom") << "Thickness of the partition " << cms::convert2mm(layerThick_[i])
0409                                        << " does not match with " << cms::convert2mm(thickTot) << " of the components";
0410         }
0411       }
0412     }  // End of loop over blocks
0413 #ifdef EDM_ML_DEBUG
0414     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: " << copies_.size()
0415                                   << " different wafer copy numbers";
0416     int k(0);
0417     for (std::unordered_set<int>::const_iterator itr = copies_.begin(); itr != copies_.end(); ++itr, ++k) {
0418       edm::LogVerbatim("HGCalGeom") << "Copy [" << k << "] : " << (*itr);
0419     }
0420     copies_.clear();
0421     edm::LogVerbatim("HGCalGeom") << "<<== End of DDHGCalMixRotatedFineCassette construction...";
0422 #endif
0423   }
0424 
0425   void positionMix(cms::DDParsingContext& ctxt,
0426                    xml_h e,
0427                    const dd4hep::Volume& glog,
0428                    const std::string& nameM,
0429                    int copyM,
0430                    double thick,
0431                    const dd4hep::Material& matter,
0432                    int absType,
0433                    bool fine) {
0434     cms::DDNamespace ns(ctxt, e, true);
0435 
0436     // Make the top part first
0437     for (unsigned int ly = 0; ly < layerTypeTop_.size(); ++ly) {
0438       int ii = layerTypeTop_[ly];
0439       copyNumberTop_[ii] = copyM;
0440     }
0441     double hthick = 0.5 * thick;
0442     double dphi = (fine) ? ((2._pi) / phiBinsFineScint_) : ((2._pi) / phiBinsScint_);
0443     double thickTot(0), zpos(-hthick);
0444 #ifdef EDM_ML_DEBUG
0445     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: Entry to positionMix with Name " << nameM
0446                                   << " copy " << copyM << " Thick " << thick << " AbsType " << absType << " Fine "
0447                                   << fine << " dphi " << convertRadToDeg(dphi);
0448 #endif
0449     if (absType < 0) {
0450       for (unsigned int ly = 0; ly < layerTypeTop_.size(); ++ly) {
0451         int ii = layerTypeTop_[ly];
0452         int copy = copyNumberTop_[ii];
0453         int layer = (fine) ? (copy - firstFineLayer_) : (copy - firstCoarseLayer_);
0454         double hthickl = 0.5 * layerThickTop_[ii];
0455         thickTot += layerThickTop_[ii];
0456         zpos += hthickl;
0457         dd4hep::Material matter1 = ns.material(materialTop_[ii]);
0458         unsigned int k = 0;
0459         int firstTile = (fine) ? tileFineLayerStart_[layer] : tileCoarseLayerStart_[layer];
0460         int lastTile = (fine) ? ((layer + 1 < static_cast<int>(tileFineLayerStart_.size()))
0461                                      ? tileFineLayerStart_[layer + 1]
0462                                      : static_cast<int>(tileFineIndex_.size()))
0463                               : ((layer + 1 < static_cast<int>(tileCoarseLayerStart_.size()))
0464                                      ? tileCoarseLayerStart_[layer + 1]
0465                                      : static_cast<int>(tileCoarseIndex_.size()));
0466 #ifdef EDM_ML_DEBUG
0467         edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: Layer " << ly << ":" << ii << ":" << layer
0468                                       << " Copy " << copy << " Tiles " << firstTile << ":" << lastTile << " Size "
0469                                       << tileFineIndex_.size() << ":" << tileCoarseIndex_.size() << " Fine " << fine
0470                                       << " absType " << absType;
0471 #endif
0472         for (int ti = firstTile; ti < lastTile; ++ti) {
0473           double r1, r2;
0474           int cassette, fimin, fimax;
0475 #ifdef EDM_ML_DEBUG
0476           edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: ti " << ti << ":" << fine << " index "
0477                                         << tileFineIndex_.size() << ":" << tileCoarseIndex_.size() << " Phis "
0478                                         << tileFinePhis_.size() << ":" << tileCoarsePhis_.size();
0479 #endif
0480           if (fine) {
0481             r1 = tileFineRMin_[std::get<1>(HGCalTileIndex::tileUnpack(tileFineIndex_[ti])) - 1];
0482             r2 = tileFineRMax_[std::get<2>(HGCalTileIndex::tileUnpack(tileFineIndex_[ti])) - 1];
0483             cassette = std::get<0>(HGCalTileIndex::tileUnpack(tileFinePhis_[ti]));
0484             fimin = std::get<1>(HGCalTileIndex::tileUnpack(tileFinePhis_[ti]));
0485             fimax = std::get<2>(HGCalTileIndex::tileUnpack(tileFinePhis_[ti]));
0486           } else {
0487             r1 = tileCoarseRMin_[std::get<1>(HGCalTileIndex::tileUnpack(tileCoarseIndex_[ti])) - 1];
0488             r2 = tileCoarseRMax_[std::get<2>(HGCalTileIndex::tileUnpack(tileCoarseIndex_[ti])) - 1];
0489             cassette = std::get<0>(HGCalTileIndex::tileUnpack(tileCoarsePhis_[ti]));
0490             fimin = std::get<1>(HGCalTileIndex::tileUnpack(tileCoarsePhis_[ti]));
0491             fimax = std::get<2>(HGCalTileIndex::tileUnpack(tileCoarsePhis_[ti]));
0492           }
0493 #ifdef EDM_ML_DEBUG
0494           edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: Casstee|Fimin|Fimax " << cassette << ":"
0495                                         << fimin << ":" << fimax;
0496 #endif
0497           double phi1 = dphi * (fimin - 1);
0498           double phi2 = dphi * (fimax - fimin + 1);
0499           auto cshift = cassette_.getShift(layer + 1, 1, cassette, true);
0500 #ifdef EDM_ML_DEBUG
0501           int cassette0 = HGCalCassette::cassetteType(2, 1, cassette);  //
0502           int ir1 = (fine) ? std::get<1>(HGCalTileIndex::tileUnpack(tileFineIndex_[ti]))
0503                            : std::get<1>(HGCalTileIndex::tileUnpack(tileCoarseIndex_[ti]));
0504           int ir2 = (fine) ? std::get<2>(HGCalTileIndex::tileUnpack(tileFineIndex_[ti]))
0505                            : std::get<2>(HGCalTileIndex::tileUnpack(tileCoarseIndex_[ti]));
0506           edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: Layer " << copy << ":" << (layer + 1)
0507                                         << " iR " << ir1 << ":" << ir2 << " R " << cms::convert2mm(r1) << ":"
0508                                         << cms::convert2mm(r2) << " Thick " << cms::convert2mm(2.0 * hthickl) << " phi "
0509                                         << fimin << ":" << fimax << ":" << convertRadToDeg(phi1) << ":"
0510                                         << convertRadToDeg(phi2) << " cassette " << cassette << ":" << cassette0
0511                                         << " Shift " << cms::convert2mm(cshift.first) << ":"
0512                                         << cms::convert2mm(cshift.second);
0513 #endif
0514           std::string name = namesTop_[ii] + "L" + std::to_string(copy) + "F" + std::to_string(k);
0515           ++k;
0516           dd4hep::Solid solid = dd4hep::Tube(r1, r2, hthickl, phi1, phi2);
0517           ns.addSolidNS(ns.prepend(name), solid);
0518           dd4hep::Volume glog1 = dd4hep::Volume(solid.name(), solid, matter1);
0519           ns.addVolumeNS(glog1);
0520 #ifdef EDM_ML_DEBUG
0521           edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: " << glog1.name() << " Tubs made of "
0522                                         << matter1.name() << " of dimensions " << cms::convert2mm(r1) << ", "
0523                                         << cms::convert2mm(r2) << ", " << cms::convert2mm(hthickl) << ", "
0524                                         << convertRadToDeg(phi1) << ", " << convertRadToDeg(phi2);
0525 #endif
0526           dd4hep::Position tran(-cshift.first, cshift.second, zpos);
0527           glog.placeVolume(glog1, copy, tran);
0528 #ifdef EDM_ML_DEBUG
0529           edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: Position " << glog1.name() << " number "
0530                                         << copy << " in " << glog.name() << " at (" << cms::convert2mm(cshift.first)
0531                                         << ", " << cms::convert2mm(cshift.second) << ", " << cms::convert2mm(zpos)
0532                                         << ") with no rotation";
0533 #endif
0534         }
0535         ++copyNumberTop_[ii];
0536         zpos += hthickl;
0537       }
0538       if (std::abs(thickTot - thick) >= tol2_) {
0539         if (thickTot > thick) {
0540           edm::LogError("HGCalGeom") << "DDHGCalMixRotatedFineCassette: Thickness of the partition "
0541                                      << cms::convert2mm(thick) << " is smaller than " << cms::convert2mm(thickTot)
0542                                      << ": thickness of all its components in the top part **** ERROR ****";
0543         } else {
0544           edm::LogWarning("HGCalGeom") << "DDHGCalMixRotatedFineCassette: Thickness of the partition "
0545                                        << cms::convert2mm(thick) << " does not match with " << cms::convert2mm(thickTot)
0546                                        << " of the components in top part";
0547         }
0548       }
0549     } else {
0550       int ii = coverTypeTop_;
0551       int copy = copyNumberCoverTop_[absType - 1];
0552       int layer = (fine) ? (copy - firstFineLayer_) : (copy - firstCoarseLayer_);
0553       double hthickl = 0.5 * layerThickTop_[ii];
0554       zpos += hthickl;
0555       dd4hep::Material matter1 = ns.material(materialTop_[ii]);
0556       unsigned int k = 0;
0557       int firstTile = (fine) ? tileFineLayerStart_[layer] : tileCoarseLayerStart_[layer];
0558       int lastTile = (fine) ? (((layer + 1) < static_cast<int>(tileFineLayerStart_.size()))
0559                                    ? tileFineLayerStart_[layer + 1]
0560                                    : static_cast<int>(tileFineIndex_.size()))
0561                             : (((layer + 1) < static_cast<int>(tileCoarseLayerStart_.size()))
0562                                    ? tileCoarseLayerStart_[layer + 1]
0563                                    : static_cast<int>(tileCoarseIndex_.size()));
0564 #ifdef EDM_ML_DEBUG
0565       edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: TOP Layer " << ii << ":" << layer << " Copy "
0566                                     << copy << " Tiles " << firstTile << ":" << lastTile << " Size "
0567                                     << tileFineIndex_.size() << ":" << tileCoarseIndex_.size() << " Fine " << fine
0568                                     << " absType " << absType;
0569 #endif
0570       for (int ti = firstTile; ti < lastTile; ++ti) {
0571         double r1, r2;
0572         int cassette, fimin, fimax;
0573 #ifdef EDM_ML_DEBUG
0574         edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: ti " << ti << ":" << fine << " index "
0575                                       << tileFineIndex_.size() << ":" << tileCoarseIndex_.size() << " Phis "
0576                                       << tileFinePhis_.size() << ":" << tileCoarsePhis_.size();
0577 #endif
0578         if (fine) {
0579           r1 = tileFineRMin_[std::get<1>(HGCalTileIndex::tileUnpack(tileFineIndex_[ti])) - 1];
0580           r2 = tileFineRMax_[std::get<2>(HGCalTileIndex::tileUnpack(tileFineIndex_[ti])) - 1];
0581           cassette = std::get<0>(HGCalTileIndex::tileUnpack(tileFinePhis_[ti]));
0582           fimin = std::get<1>(HGCalTileIndex::tileUnpack(tileFinePhis_[ti]));
0583           fimax = std::get<2>(HGCalTileIndex::tileUnpack(tileFinePhis_[ti]));
0584         } else {
0585           r1 = tileCoarseRMin_[std::get<1>(HGCalTileIndex::tileUnpack(tileCoarseIndex_[ti])) - 1];
0586           r2 = tileCoarseRMax_[std::get<2>(HGCalTileIndex::tileUnpack(tileCoarseIndex_[ti])) - 1];
0587           cassette = std::get<0>(HGCalTileIndex::tileUnpack(tileCoarsePhis_[ti]));
0588           fimin = std::get<1>(HGCalTileIndex::tileUnpack(tileCoarsePhis_[ti]));
0589           fimax = std::get<2>(HGCalTileIndex::tileUnpack(tileCoarsePhis_[ti]));
0590         }
0591 #ifdef EDM_ML_DEBUG
0592         edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: Casstee|Fimin|Fimax " << cassette << ":"
0593                                       << fimin << ":" << fimax;
0594 #endif
0595         double phi1 = dphi * (fimin - 1);
0596         double phi2 = dphi * (fimax - fimin + 1);
0597         auto cshift = cassette_.getShift(layer + 1, 1, cassette, true);
0598 #ifdef EDM_ML_DEBUG
0599         int cassette0 = HGCalCassette::cassetteType(2, 1, cassette);  //
0600         int ir1 = (fine) ? std::get<1>(HGCalTileIndex::tileUnpack(tileFineIndex_[ti]))
0601                          : std::get<1>(HGCalTileIndex::tileUnpack(tileCoarseIndex_[ti]));
0602         int ir2 = (fine) ? std::get<2>(HGCalTileIndex::tileUnpack(tileFineIndex_[ti]))
0603                          : std::get<2>(HGCalTileIndex::tileUnpack(tileCoarseIndex_[ti]));
0604         edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: Layer " << copy << ":" << (layer + 1) << " iR "
0605                                       << ir1 << ":" << ir2 << " R " << cms::convert2mm(r1) << ":" << cms::convert2mm(r2)
0606                                       << " Thick " << cms::convert2mm(2.0 * hthickl) << " phi " << fimin << ":" << fimax
0607                                       << ":" << convertRadToDeg(phi1) << ":" << convertRadToDeg(phi2) << " cassette "
0608                                       << cassette << ":" << cassette0 << " Shift " << cms::convert2mm(cshift.first)
0609                                       << ":" << cms::convert2mm(cshift.second);
0610 #endif
0611         std::string name = namesTop_[ii] + "L" + std::to_string(copy) + "F" + std::to_string(k);
0612         ++k;
0613         dd4hep::Solid solid = dd4hep::Tube(r1, r2, hthickl, phi1, phi2);
0614         ns.addSolidNS(ns.prepend(name), solid);
0615         dd4hep::Volume glog1 = dd4hep::Volume(solid.name(), solid, matter1);
0616         ns.addVolumeNS(glog1);
0617 #ifdef EDM_ML_DEBUG
0618         edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: " << glog1.name() << " Tubs made of "
0619                                       << matter1.name() << " of dimensions " << cms::convert2mm(r1) << ", "
0620                                       << cms::convert2mm(r2) << ", " << cms::convert2mm(hthickl) << ", "
0621                                       << convertRadToDeg(phi1) << ", " << convertRadToDeg(phi2);
0622 #endif
0623         dd4hep::Position tran(-cshift.first, cshift.second, zpos);
0624         glog.placeVolume(glog1, copy, tran);
0625 #ifdef EDM_ML_DEBUG
0626         edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: Position " << glog1.name() << " number "
0627                                       << copy << " in " << glog.name() << " at (" << cms::convert2mm(-cshift.first)
0628                                       << ", " << cms::convert2mm(cshift.second) << ", " << cms::convert2mm(zpos)
0629                                       << ") with no rotation";
0630 #endif
0631       }
0632       ++copyNumberCoverTop_[absType - 1];
0633     }
0634 
0635     // Make the bottom part next
0636     int layer = (copyM - firstFineLayer_);
0637 #ifdef EDM_ML_DEBUG
0638     edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: Start bottom section for layer " << layer
0639                                   << " absType " << absType;
0640 #endif
0641     if (absType > 0) {
0642 #ifdef EDM_ML_DEBUG
0643       int kount(0);
0644 #endif
0645       for (int k = 0; k < cassettes_; ++k) {
0646         int cassette = k + 1;
0647         auto cshift = cassette_.getShift(layer + 1, -1, cassette);
0648         double xpos = -cshift.first;
0649         double ypos = cshift.second;
0650         int i = layer * cassettes_ + k;
0651 #ifdef EDM_ML_DEBUG
0652         edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette::Passive: layer " << layer + 1 << " cassette "
0653                                       << cassette << " Shift " << cms::convert2mm(-cshift.first) << ":"
0654                                       << cms::convert2mm(cshift.second) << " PassiveIndex " << i << ":"
0655                                       << passiveFull_.size() << ":" << passivePart_.size();
0656 #endif
0657         std::string passive = (absType <= waferTypes_) ? passiveFull_[i] : passivePart_[i];
0658 #ifdef EDM_ML_DEBUG
0659         edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: Passive " << passive << " number " << cassette
0660                                       << " pos " << cms::convert2mm(xpos) << ":" << cms::convert2mm(ypos);
0661         kount++;
0662 #endif
0663         dd4hep::Position tran(xpos, ypos, 0.0);
0664         glog.placeVolume(ns.volume(passive), cassette, tran);
0665 #ifdef EDM_ML_DEBUG
0666         edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: " << passive << " number " << cassette
0667                                       << " positioned in " << glog.name() << " at (" << cms::convert2mm(xpos) << ","
0668                                       << cms::convert2mm(ypos) << ",0) with no rotation";
0669 #endif
0670       }
0671 #ifdef EDM_ML_DEBUG
0672       edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: " << kount << " passives of type " << absType
0673                                     << " for " << glog.name();
0674 #endif
0675     } else {
0676       static const double sqrt3 = std::sqrt(3.0);
0677       int layercenter = layerOrient_[layer];
0678       int layertype = HGCalTypes::layerFrontBack(layerOrient_[layer]);
0679       int firstWafer = waferLayerStart_[layer];
0680       int lastWafer =
0681           (((layer + 1) < static_cast<int>(waferLayerStart_.size())) ? waferLayerStart_[layer + 1]
0682                                                                      : static_cast<int>(waferIndex_.size()));
0683       double delx = 0.5 * (waferSize_ + waferSepar_);
0684       double dely = 2.0 * delx / sqrt3;
0685       double dy = 0.75 * dely;
0686       const auto& xyoff = geomTools_.shiftXY(layercenter, (waferSize_ + waferSepar_));
0687 #ifdef EDM_ML_DEBUG
0688       int ium(0), ivm(0), kount(0);
0689       std::vector<int> ntype(3, 0);
0690       edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette::Bottom: " << glog.name() << "  r "
0691                                     << cms::convert2mm(delx) << " R " << cms::convert2mm(dely) << " dy "
0692                                     << cms::convert2mm(dy) << " Shift " << cms::convert2mm(xyoff.first) << ":"
0693                                     << cms::convert2mm(xyoff.second) << " WaferSize "
0694                                     << cms::convert2mm(waferSize_ + waferSepar_) << " index " << firstWafer << ":"
0695                                     << (lastWafer - 1) << " Copy " << copyM << ":" << layer;
0696 #endif
0697       for (int k = firstWafer; k < lastWafer; ++k) {
0698         int u = HGCalWaferIndex::waferU(waferIndex_[k]);
0699         int v = HGCalWaferIndex::waferV(waferIndex_[k]);
0700 #ifdef EDM_ML_DEBUG
0701         int iu = std::abs(u);
0702         int iv = std::abs(v);
0703 #endif
0704         int nr = 2 * v;
0705         int nc = -2 * u + v;
0706         int type = HGCalProperty::waferThick(waferProperty_[k]);
0707         int part = HGCalProperty::waferPartial(waferProperty_[k]);
0708         int orien = HGCalProperty::waferOrient(waferProperty_[k]);
0709         int cassette = HGCalProperty::waferCassette(waferProperty_[k]);
0710         int place = HGCalCell::cellPlacementIndex(1, layertype, orien);
0711 #ifdef EDM_ML_DEBUG
0712         edm::LogVerbatim("HGCalGeom")
0713             << "DDHGCalMixRotatedFineCassette::index:Property:layertype:type:part:orien:cassette:place:offsets:ind "
0714             << k << waferProperty_[k] << ":" << layertype << ":" << type << ":" << part << ":" << orien << ":"
0715             << cassette << ":" << place;
0716 #endif
0717         auto cshift = cassette_.getShift(layer + 1, -1, cassette, false);
0718         double xpos = xyoff.first - cshift.first + nc * delx;
0719         double ypos = xyoff.second + cshift.second + nr * dy;
0720 #ifdef EDM_ML_DEBUG
0721         double xorig = xyoff.first + nc * delx;
0722         double yorig = xyoff.second + nr * dy;
0723         double angle = std::atan2(yorig, xorig);
0724         edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette::Wafer: layer " << layer + 1 << " cassette "
0725                                       << cassette << " Shift " << cms::convert2mm(cshift.first) << ":"
0726                                       << cms::convert2mm(cshift.second) << " Original " << cms::convert2mm(xorig) << ":"
0727                                       << cms::convert2mm(yorig) << ":" << convertRadToDeg(angle) << " Final "
0728                                       << cms::convert2mm(xpos) << ":" << cms::convert2mm(ypos);
0729 #endif
0730         std::string wafer;
0731         int i(999);
0732         if (part == HGCalTypes::WaferFull) {
0733           i = type * facingTypes_ * orientationTypes_ + place - placeOffset_;
0734 #ifdef EDM_ML_DEBUG
0735           edm::LogVerbatim("HGCalGeom") << " FullWafer type:place:ind " << type << ":" << place << ":" << i << ":"
0736                                         << waferFull_.size();
0737 #endif
0738           wafer = waferFull_[i];
0739         } else {
0740           int partoffset =
0741               (part >= HGCalTypes::WaferHDTop) ? HGCalTypes::WaferPartHDOffset : HGCalTypes::WaferPartLDOffset;
0742           i = (part - partoffset) * facingTypes_ * orientationTypes_ +
0743               HGCalTypes::WaferTypeOffset[type] * facingTypes_ * orientationTypes_ + place - placeOffset_;
0744 #ifdef EDM_ML_DEBUG
0745           edm::LogVerbatim("HGCalGeom") << " layertype:type:part:orien:cassette:place:offsets:ind " << layertype << ":"
0746                                         << type << ":" << part << ":" << orien << ":" << cassette << ":" << place << ":"
0747                                         << partoffset << ":" << HGCalTypes::WaferTypeOffset[type] << ":" << i << ":"
0748                                         << waferPart_.size();
0749 #endif
0750           wafer = waferPart_[i];
0751         }
0752         int copy = HGCalTypes::packTypeUV(type, u, v);
0753 #ifdef EDM_ML_DEBUG
0754         edm::LogVerbatim("HGCalGeom") << " DDHGCalMixRotatedFineCassette: Layer "
0755                                       << HGCalWaferIndex::waferLayer(waferIndex_[k]) << " Wafer " << wafer << " number "
0756                                       << copy << " type :part:orien:ind " << type << ":" << part << ":" << orien << ":"
0757                                       << i << " layer:u:v " << (layer + firstFineLayer_) << ":" << u << ":" << v;
0758         if (iu > ium)
0759           ium = iu;
0760         if (iv > ivm)
0761           ivm = iv;
0762         kount++;
0763         if (copies_.count(copy) == 0)
0764           copies_.insert(copy);
0765 #endif
0766         dd4hep::Position tran(xpos, ypos, 0.0);
0767         glog.placeVolume(ns.volume(wafer), copy, tran);
0768 #ifdef EDM_ML_DEBUG
0769         ++ntype[type];
0770         edm::LogVerbatim("HGCalGeom") << " DDHGCalMixRotatedFineCassette: " << wafer << " number " << copy << " type "
0771                                       << layertype << ":" << type << " positioned in " << glog.name() << " at ("
0772                                       << cms::convert2mm(xpos) << "," << cms::convert2mm(ypos)
0773                                       << ",0) with no rotation";
0774 #endif
0775       }
0776 #ifdef EDM_ML_DEBUG
0777       edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedFineCassette: Maximum # of u " << ium << " # of v " << ivm
0778                                     << " and " << kount << " wafers (" << ntype[0] << ":" << ntype[1] << ":" << ntype[2]
0779                                     << ") for " << glog.name();
0780 #endif
0781     }
0782   }
0783 
0784   HGCalGeomTools geomTools_;
0785   HGCalCassette cassette_;
0786 
0787   int waferTypes_;                         // Number of wafer types
0788   int passiveTypes_;                       // Number of passive types
0789   int facingTypes_;                        // Types of facings of modules toward IP
0790   int orientationTypes_;                   // Number of partial wafer orienations
0791   int partialTypes_;                       // Number of partial types
0792   int placeOffset_;                        // Offset for placement
0793   int phiBinsScint_;                       // Maximum number of cells along phi (coarse)
0794   int phiBinsFineScint_;                   // Maximum number of cells along phi (fine)
0795   int firstFineLayer_;                     // Copy # of the first Fine sensitive layer
0796   int firstCoarseLayer_;                   // Copy # of the first Coarse sensitive layer
0797   int absorbMode_;                         // Absorber mode
0798   int sensitiveMode_;                      // Sensitive mode
0799   int passiveMode_;                        // Mode for passive components
0800   double zMinBlock_;                       // Starting z-value of the block
0801   double waferSize_;                       // Width of the wafer
0802   double waferSepar_;                      // Sensor separation
0803   int sectors_;                            // Sectors
0804   int cassettes_;                          // Cassettes
0805   std::vector<double> slopeB_;             // Slope at the lower R
0806   std::vector<double> zFrontB_;            // Starting Z values for the slopes
0807   std::vector<double> rMinFront_;          // Corresponding rMin's
0808   std::vector<double> slopeT_;             // Slopes at the larger R
0809   std::vector<double> zFrontT_;            // Starting Z values for the slopes
0810   std::vector<double> rMaxFront_;          // Corresponding rMax's
0811   std::vector<std::string> waferFull_;     // Names of full wafer modules
0812   std::vector<std::string> waferPart_;     // Names of partial wafer modules
0813   std::vector<std::string> passiveFull_;   // Names of full passive modules
0814   std::vector<std::string> passivePart_;   // Names of partial passive modules
0815   std::vector<std::string> materials_;     // Materials
0816   std::vector<std::string> names_;         // Names
0817   std::vector<double> thick_;              // Thickness of the material
0818   std::vector<int> copyNumber_;            // Initial copy numbers
0819   std::vector<int> layers_;                // Number of layers in a section
0820   std::vector<double> layerThick_;         // Thickness of each section
0821   std::vector<int> layerType_;             // Type of the layer
0822   std::vector<int> layerSense_;            // Content of a layer (sensitive?)
0823   std::vector<std::string> materialTop_;   // Materials of top layers
0824   std::vector<std::string> namesTop_;      // Names of top layers
0825   std::vector<double> layerThickTop_;      // Thickness of the top sections
0826   std::vector<int> layerTypeTop_;          // Type of the Top layer
0827   std::vector<int> copyNumberTop_;         // Initial copy numbers (top section)
0828   int coverTypeTop_;                       // Type of the Top layer cover
0829   int coverTopLayers_;                     // Number of cover layers in top section
0830   std::vector<int> copyNumberCoverTop_;    // Initial copy number of top cover
0831   std::vector<int> layerOrient_;           // Layer orientation for the silicon component
0832   std::vector<int> waferIndex_;            // Wafer index for the types
0833   std::vector<int> waferProperty_;         // Wafer property
0834   std::vector<int> waferLayerStart_;       // Start index of wafers in each layer
0835   std::vector<double> cassetteShift_;      // Shifts of the cassetes
0836   std::vector<double> tileFineRMin_;       // Minimum radius of each fine ring
0837   std::vector<double> tileFineRMax_;       // Maximum radius of each fine ring
0838   std::vector<int> tileFineIndex_;         // Index of tile (layer/start|end fine ring)
0839   std::vector<int> tileFinePhis_;          // Tile phi range for each index in fine ring
0840   std::vector<int> tileFineLayerStart_;    // Start index of tiles in each fine layer
0841   std::vector<double> tileCoarseRMin_;     // Minimum radius of each coarse ring
0842   std::vector<double> tileCoarseRMax_;     // Maximum radius of each coarse ring
0843   std::vector<int> tileCoarseIndex_;       // Index of tile (layer/start|end coarse ring)
0844   std::vector<int> tileCoarsePhis_;        // Tile phi range for each index in coarse ring
0845   std::vector<int> tileCoarseLayerStart_;  // Start index of tiles in each coarse layer
0846   std::vector<double> cassetteShiftScnt_;  // Shifts of the cassetes for scintillators
0847   std::string nameSpace_;                  // Namespace of this and ALL sub-parts
0848   std::unordered_set<int> copies_;         // List of copy #'s
0849   double alpha_, cosAlpha_;
0850   static constexpr double tol2_ = 0.00001 * dd4hep::mm;
0851 };
0852 
0853 static long algorithm(dd4hep::Detector& /* description */, cms::DDParsingContext& ctxt, xml_h e) {
0854   HGCalMixRotatedFineCassette siliconRotatedCassetteAlgo(ctxt, e);
0855   return cms::s_executed;
0856 }
0857 
0858 DECLARE_DDCMS_DETELEMENT(DDCMS_hgcal_DDHGCalMixRotatedFineCassette, algorithm)