Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-08-15 04:54:02

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