Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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