Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:15:00

0001 ///////////////////////////////////////////////////////////////////////////////
0002 // File: DDHGCalSiliconModule.cc
0003 // Description: Geometry factory class for HGCal (EE and HESil) using
0004 //              information from the file
0005 ///////////////////////////////////////////////////////////////////////////////
0006 
0007 #include "DataFormats/Math/interface/angle_units.h"
0008 #include "DetectorDescription/Core/interface/DDAlgorithm.h"
0009 #include "DetectorDescription/Core/interface/DDAlgorithmFactory.h"
0010 #include "DetectorDescription/Core/interface/DDCurrentNamespace.h"
0011 #include "DetectorDescription/Core/interface/DDLogicalPart.h"
0012 #include "DetectorDescription/Core/interface/DDMaterial.h"
0013 #include "DetectorDescription/Core/interface/DDSolid.h"
0014 #include "DetectorDescription/Core/interface/DDSplit.h"
0015 #include "DetectorDescription/Core/interface/DDTypes.h"
0016 #include "DetectorDescription/Core/interface/DDutils.h"
0017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0018 #include "FWCore/PluginManager/interface/PluginFactory.h"
0019 #include "Geometry/HGCalCommonData/interface/HGCalGeomTools.h"
0020 #include "Geometry/HGCalCommonData/interface/HGCalParameters.h"
0021 #include "Geometry/HGCalCommonData/interface/HGCalProperty.h"
0022 #include "Geometry/HGCalCommonData/interface/HGCalTypes.h"
0023 #include "Geometry/HGCalCommonData/interface/HGCalWaferIndex.h"
0024 #include "Geometry/HGCalCommonData/interface/HGCalWaferType.h"
0025 
0026 #include <cmath>
0027 #include <memory>
0028 #include <sstream>
0029 #include <string>
0030 #include <unordered_set>
0031 #include <vector>
0032 
0033 //#define EDM_ML_DEBUG
0034 using namespace angle_units::operators;
0035 
0036 class DDHGCalSiliconModule : public DDAlgorithm {
0037 public:
0038   DDHGCalSiliconModule();
0039 
0040   void initialize(const DDNumericArguments& nArgs,
0041                   const DDVectorArguments& vArgs,
0042                   const DDMapArguments& mArgs,
0043                   const DDStringArguments& sArgs,
0044                   const DDStringVectorArguments& vsArgs) override;
0045   void execute(DDCompactView& cpv) override;
0046 
0047 protected:
0048   void constructLayers(const DDLogicalPart&, DDCompactView& cpv);
0049   void positionSensitive(const DDLogicalPart& glog, int layer, DDCompactView& cpv);
0050 
0051 private:
0052   HGCalGeomTools geomTools_;
0053 
0054   static constexpr double tol1_ = 0.01;
0055   static constexpr double tol2_ = 0.00001;
0056 
0057   int waferTypes_;                      // Number of wafer types
0058   int facingTypes_;                     // Types of facings of modules toward IP
0059   int partialTypes_;                    // Number of partial wafer types
0060   int orientationTypes_;                // Number of partial wafer orienations
0061   int firstLayer_;                      // Copy # of the first sensitive layer
0062   int absorbMode_;                      // Absorber mode
0063   int sensitiveMode_;                   // Sensitive mode
0064   double zMinBlock_;                    // Starting z-value of the block
0065   double waferSize_;                    // Width of the wafer
0066   double waferSepar_;                   // Sensor separation
0067   int sectors_;                         // Sectors
0068   std::string rotstr_;                  // Rotation matrix (if needed)
0069   std::vector<std::string> waferFull_;  // Names of full wafer modules
0070   std::vector<std::string> waferPart_;  // Names of partial wafer modules
0071   std::vector<std::string> materials_;  // names of materials
0072   std::vector<std::string> names_;      // Names of volumes
0073   std::vector<double> thick_;           // Thickness of the material
0074   std::vector<int> copyNumber_;         // Initial copy numbers
0075   std::vector<int> layers_;             // Number of layers in a section
0076   std::vector<double> layerThick_;      // Thickness of each section
0077   std::vector<int> layerType_;          // Type of the layer
0078   std::vector<int> layerSense_;         // Content of a layer (sensitive?)
0079   std::vector<double> slopeB_;          // Slope at the lower R
0080   std::vector<double> zFrontB_;         // Starting Z values for the slopes
0081   std::vector<double> rMinFront_;       // Corresponding rMin's
0082   std::vector<double> slopeT_;          // Slopes at the larger R
0083   std::vector<double> zFrontT_;         // Starting Z values for the slopes
0084   std::vector<double> rMaxFront_;       // Corresponding rMax's
0085   std::vector<int> layerOrient_;        // Layer orientation (Centering, rotations..)
0086   std::vector<int> waferIndex_;         // Wafer index for the types
0087   std::vector<int> waferProperty_;      // Wafer property
0088   std::vector<int> waferLayerStart_;    // Index of wafers in each layer
0089   std::string nameSpace_;               // Namespace of this and ALL sub-parts
0090   std::unordered_set<int> copies_;      // List of copy #'s
0091   double alpha_, cosAlpha_;
0092 };
0093 
0094 DDHGCalSiliconModule::DDHGCalSiliconModule() {
0095 #ifdef EDM_ML_DEBUG
0096   edm::LogVerbatim("HGCalGeom") << "DDHGCalSiliconModule: Creating an instance";
0097 #endif
0098 }
0099 
0100 void DDHGCalSiliconModule::initialize(const DDNumericArguments& nArgs,
0101                                       const DDVectorArguments& vArgs,
0102                                       const DDMapArguments&,
0103                                       const DDStringArguments& sArgs,
0104                                       const DDStringVectorArguments& vsArgs) {
0105   waferTypes_ = static_cast<int>(nArgs["WaferTypes"]);
0106   facingTypes_ = static_cast<int>(nArgs["FacingTypes"]);
0107   partialTypes_ = static_cast<int>(nArgs["PartialTypes"]);
0108   orientationTypes_ = static_cast<int>(nArgs["OrientationTypes"]);
0109 #ifdef EDM_ML_DEBUG
0110   edm::LogVerbatim("HGCalGeom") << "Number of types of wafers: " << waferTypes_ << " facings: " << facingTypes_
0111                                 << " partials: " << partialTypes_ << " Orientations: " << orientationTypes_;
0112 #endif
0113   firstLayer_ = static_cast<int>(nArgs["FirstLayer"]);
0114   absorbMode_ = static_cast<int>(nArgs["AbsorberMode"]);
0115   sensitiveMode_ = static_cast<int>(nArgs["SensitiveMode"]);
0116 #ifdef EDM_ML_DEBUG
0117   edm::LogVerbatim("HGCalGeom") << "First Layer " << firstLayer_ << " and "
0118                                 << "Absober:Sensitive mode " << absorbMode_ << ":" << sensitiveMode_;
0119 #endif
0120   zMinBlock_ = nArgs["zMinBlock"];
0121   waferSize_ = nArgs["waferSize"];
0122   waferSepar_ = nArgs["SensorSeparation"];
0123   sectors_ = (int)(nArgs["Sectors"]);
0124   alpha_ = (1._pi) / sectors_;
0125   cosAlpha_ = cos(alpha_);
0126   rotstr_ = sArgs["LayerRotation"];
0127 #ifdef EDM_ML_DEBUG
0128   edm::LogVerbatim("HGCalGeom") << "zStart " << zMinBlock_ << " wafer width " << waferSize_ << " separations "
0129                                 << waferSepar_ << " sectors " << sectors_ << ":" << convertRadToDeg(alpha_) << ":"
0130                                 << cosAlpha_ << " rotation matrix " << rotstr_;
0131 #endif
0132   waferFull_ = vsArgs["WaferNamesFull"];
0133   waferPart_ = vsArgs["WaferNamesPartial"];
0134 #ifdef EDM_ML_DEBUG
0135   edm::LogVerbatim("HGCalGeom") << "DDHGCalSiliconModule: " << waferFull_.size() << " full and " << waferPart_.size()
0136                                 << " partial modules\nDDHGCalSiliconModule:Full Modules:";
0137   unsigned int i1max = static_cast<unsigned int>(waferFull_.size());
0138   for (unsigned int i1 = 0; i1 < i1max; i1 += 2) {
0139     std::ostringstream st1;
0140     unsigned int i2 = std::min((i1 + 2), i1max);
0141     for (unsigned int i = i1; i < i2; ++i)
0142       st1 << " [" << i << "] " << waferFull_[i];
0143     edm::LogVerbatim("HGCalGeom") << st1.str() << std::endl;
0144   }
0145   edm::LogVerbatim("HGCalGeom") << "DDHGCalSiliconModule: Partial Modules:";
0146   i1max = static_cast<unsigned int>(waferPart_.size());
0147   for (unsigned int i1 = 0; i1 < i1max; i1 += 2) {
0148     std::ostringstream st1;
0149     unsigned int i2 = std::min((i1 + 2), i1max);
0150     for (unsigned int i = i1; i < i2; ++i)
0151       st1 << " [" << i << "] " << waferPart_[i];
0152     edm::LogVerbatim("HGCalGeom") << st1.str() << std::endl;
0153   }
0154 #endif
0155   materials_ = vsArgs["MaterialNames"];
0156   names_ = vsArgs["VolumeNames"];
0157   thick_ = vArgs["Thickness"];
0158   copyNumber_.resize(materials_.size(), 1);
0159 #ifdef EDM_ML_DEBUG
0160   edm::LogVerbatim("HGCalGeom") << "DDHGCalSiliconModule: " << materials_.size() << " types of volumes";
0161   for (unsigned int i = 0; i < names_.size(); ++i)
0162     edm::LogVerbatim("HGCalGeom") << "Volume [" << i << "] " << names_[i] << " of thickness " << thick_[i]
0163                                   << " filled with " << materials_[i] << " first copy number " << copyNumber_[i];
0164 #endif
0165   layers_ = dbl_to_int(vArgs["Layers"]);
0166   layerThick_ = vArgs["LayerThick"];
0167 #ifdef EDM_ML_DEBUG
0168   edm::LogVerbatim("HGCalGeom") << "There are " << layers_.size() << " blocks";
0169   for (unsigned int i = 0; i < layers_.size(); ++i)
0170     edm::LogVerbatim("HGCalGeom") << "Block [" << i << "] of thickness " << layerThick_[i] << " with " << layers_[i]
0171                                   << " layers";
0172 #endif
0173   layerType_ = dbl_to_int(vArgs["LayerType"]);
0174   layerSense_ = dbl_to_int(vArgs["LayerSense"]);
0175   layerOrient_ = dbl_to_int(vArgs["LayerTypes"]);
0176   for (unsigned int k = 0; k < layerOrient_.size(); ++k)
0177     layerOrient_[k] = HGCalTypes::layerType(layerOrient_[k]);
0178 #ifdef EDM_ML_DEBUG
0179   for (unsigned int i = 0; i < layerOrient_.size(); ++i)
0180     edm::LogVerbatim("HGCalGeom") << "LayerTypes [" << i << "] " << layerOrient_[i];
0181 #endif
0182   if (firstLayer_ > 0) {
0183     for (unsigned int i = 0; i < layerType_.size(); ++i) {
0184       if (layerSense_[i] > 0) {
0185         int ii = layerType_[i];
0186         copyNumber_[ii] = (layerSense_[i] == 1) ? firstLayer_ : (firstLayer_ + 1);
0187 #ifdef EDM_ML_DEBUG
0188         edm::LogVerbatim("HGCalGeom") << "First copy number for layer type " << i << ":" << ii << " with "
0189                                       << materials_[ii] << " changed to " << copyNumber_[ii];
0190 #endif
0191       }
0192     }
0193   } else {
0194     firstLayer_ = 1;
0195   }
0196 #ifdef EDM_ML_DEBUG
0197   edm::LogVerbatim("HGCalGeom") << "There are " << layerType_.size() << " layers";
0198   for (unsigned int i = 0; i < layerType_.size(); ++i)
0199     edm::LogVerbatim("HGCalGeom") << "Layer [" << i << "] with material type " << layerType_[i] << " sensitive class "
0200                                   << layerSense_[i];
0201 #endif
0202   slopeB_ = vArgs["SlopeBottom"];
0203   zFrontB_ = vArgs["ZFrontBottom"];
0204   rMinFront_ = vArgs["RMinFront"];
0205   slopeT_ = vArgs["SlopeTop"];
0206   zFrontT_ = vArgs["ZFrontTop"];
0207   rMaxFront_ = vArgs["RMaxFront"];
0208 #ifdef EDM_ML_DEBUG
0209   for (unsigned int i = 0; i < slopeB_.size(); ++i)
0210     edm::LogVerbatim("HGCalGeom") << "Bottom Block [" << i << "] Zmin " << zFrontB_[i] << " Rmin " << rMinFront_[i]
0211                                   << " Slope " << slopeB_[i];
0212   for (unsigned int i = 0; i < slopeT_.size(); ++i)
0213     edm::LogVerbatim("HGCalGeom") << "Top Block [" << i << "] Zmin " << zFrontT_[i] << " Rmax " << rMaxFront_[i]
0214                                   << " Slope " << slopeT_[i];
0215 #endif
0216   waferIndex_ = dbl_to_int(vArgs["WaferIndex"]);
0217   waferProperty_ = dbl_to_int(vArgs["WaferProperties"]);
0218   waferLayerStart_ = dbl_to_int(vArgs["WaferLayerStart"]);
0219 #ifdef EDM_ML_DEBUG
0220   edm::LogVerbatim("HGCalGeom") << "waferProperties with " << waferIndex_.size() << " entries in "
0221                                 << waferLayerStart_.size() << " layers";
0222   for (unsigned int k = 0; k < waferLayerStart_.size(); ++k)
0223     edm::LogVerbatim("HGCalGeom") << "LayerStart[" << k << "] " << waferLayerStart_[k];
0224   for (unsigned int k = 0; k < waferIndex_.size(); ++k)
0225     edm::LogVerbatim("HGCalGeom") << "Wafer[" << k << "] " << waferIndex_[k] << " ("
0226                                   << HGCalWaferIndex::waferLayer(waferIndex_[k]) << ", "
0227                                   << HGCalWaferIndex::waferU(waferIndex_[k]) << ", "
0228                                   << HGCalWaferIndex::waferV(waferIndex_[k]) << ") : ("
0229                                   << HGCalProperty::waferThick(waferProperty_[k]) << ":"
0230                                   << HGCalProperty::waferPartial(waferProperty_[k]) << ":"
0231                                   << HGCalProperty::waferOrient(waferProperty_[k]) << ")";
0232 #endif
0233   nameSpace_ = DDCurrentNamespace::ns();
0234 #ifdef EDM_ML_DEBUG
0235   edm::LogVerbatim("HGCalGeom") << "DDHGCalSiliconModule: NameSpace " << nameSpace_ << ":";
0236 #endif
0237 }
0238 
0239 ////////////////////////////////////////////////////////////////////
0240 // DDHGCalSiliconModule methods...
0241 ////////////////////////////////////////////////////////////////////
0242 
0243 void DDHGCalSiliconModule::execute(DDCompactView& cpv) {
0244 #ifdef EDM_ML_DEBUG
0245   edm::LogVerbatim("HGCalGeom") << "==>> Constructing DDHGCalSiliconModule...";
0246   copies_.clear();
0247 #endif
0248   constructLayers(parent(), cpv);
0249 #ifdef EDM_ML_DEBUG
0250   edm::LogVerbatim("HGCalGeom") << "DDHGCalSiliconModule: " << copies_.size() << " different wafer copy numbers";
0251   int k(0);
0252   for (std::unordered_set<int>::const_iterator itr = copies_.begin(); itr != copies_.end(); ++itr, ++k) {
0253     edm::LogVerbatim("HGCalGeom") << "Copy [" << k << "] : " << (*itr);
0254   }
0255   copies_.clear();
0256   edm::LogVerbatim("HGCalGeom") << "<<== End of DDHGCalSiliconModule construction...";
0257 #endif
0258 }
0259 
0260 void DDHGCalSiliconModule::constructLayers(const DDLogicalPart& module, DDCompactView& cpv) {
0261   double zi(zMinBlock_);
0262   int laymin(0);
0263   for (unsigned int i = 0; i < layers_.size(); i++) {
0264     double zo = zi + layerThick_[i];
0265     double routF = HGCalGeomTools::radius(zi, zFrontT_, rMaxFront_, slopeT_);
0266     int laymax = laymin + layers_[i];
0267     double zz = zi;
0268     double thickTot(0);
0269     for (int ly = laymin; ly < laymax; ++ly) {
0270       int ii = layerType_[ly];
0271       int copy = copyNumber_[ii];
0272       double hthick = 0.5 * thick_[ii];
0273       double rinB = HGCalGeomTools::radius(zo - tol1_, zFrontB_, rMinFront_, slopeB_);
0274       zz += hthick;
0275       thickTot += thick_[ii];
0276 
0277       std::string name = names_[ii] + std::to_string(copy);
0278 #ifdef EDM_ML_DEBUG
0279       edm::LogVerbatim("HGCalGeom") << "DDHGCalSiliconModule: Layer " << ly << ":" << ii << " Front " << zi << ", "
0280                                     << routF << " Back " << zo << ", " << rinB << " superlayer thickness "
0281                                     << layerThick_[i];
0282 #endif
0283       DDName matName(DDSplit(materials_[ii]).first, DDSplit(materials_[ii]).second);
0284       DDMaterial matter(matName);
0285       DDLogicalPart glog;
0286       if (layerSense_[ly] < 1) {
0287         std::vector<double> pgonZ, pgonRin, pgonRout;
0288         double rmax = routF * cosAlpha_ - tol1_;
0289         HGCalGeomTools::radius(zz - hthick,
0290                                zz + hthick,
0291                                zFrontB_,
0292                                rMinFront_,
0293                                slopeB_,
0294                                zFrontT_,
0295                                rMaxFront_,
0296                                slopeT_,
0297                                -layerSense_[ly],
0298                                pgonZ,
0299                                pgonRin,
0300                                pgonRout);
0301         for (unsigned int isec = 0; isec < pgonZ.size(); ++isec) {
0302           pgonZ[isec] -= zz;
0303           if (layerSense_[ly] == 0 || absorbMode_ == 0)
0304             pgonRout[isec] = rmax;
0305           else
0306             pgonRout[isec] = pgonRout[isec] * cosAlpha_ - tol1_;
0307         }
0308         DDSolid solid =
0309             DDSolidFactory::polyhedra(DDName(name, nameSpace_), sectors_, -alpha_, 2._pi, pgonZ, pgonRin, pgonRout);
0310         glog = DDLogicalPart(solid.ddname(), matter, solid);
0311 #ifdef EDM_ML_DEBUG
0312         edm::LogVerbatim("HGCalGeom") << "DDHGCalSiliconModule: " << solid.name() << " polyhedra of " << sectors_
0313                                       << " sectors covering " << convertRadToDeg(-alpha_) << ":"
0314                                       << convertRadToDeg(-alpha_ + 2._pi) << " with " << pgonZ.size()
0315                                       << " sections and filled with " << matName;
0316         for (unsigned int k = 0; k < pgonZ.size(); ++k)
0317           edm::LogVerbatim("HGCalGeom") << "[" << k << "] z " << pgonZ[k] << " R " << pgonRin[k] << ":" << pgonRout[k];
0318 #endif
0319       } else {
0320         double rins =
0321             (sensitiveMode_ < 1) ? rinB : HGCalGeomTools::radius(zz + hthick - tol1_, zFrontB_, rMinFront_, slopeB_);
0322         double routs =
0323             (sensitiveMode_ < 1) ? routF : HGCalGeomTools::radius(zz - hthick, zFrontT_, rMaxFront_, slopeT_);
0324         DDSolid solid = DDSolidFactory::tubs(DDName(name, nameSpace_), hthick, rins, routs, 0.0, 2._pi);
0325         glog = DDLogicalPart(solid.ddname(), matter, solid);
0326 #ifdef EDM_ML_DEBUG
0327         edm::LogVerbatim("HGCalGeom") << "DDHGCalSiliconModule: " << solid.name() << " Tubs made of " << matName
0328                                       << " of dimensions " << rinB << ":" << rins << ", " << routF << ":" << routs
0329                                       << ", " << hthick << ", 0.0, 360.0 and position " << glog.name() << " number "
0330                                       << copy << ":" << layerOrient_[copy - firstLayer_];
0331 #endif
0332         positionSensitive(glog, (copy - firstLayer_), cpv);
0333       }
0334       DDTranslation r1(0, 0, zz);
0335       DDRotation rot;
0336 #ifdef EDM_ML_DEBUG
0337       std::string rotName("Null");
0338 #endif
0339       if ((layerSense_[ly] > 0) && (layerOrient_[copy - firstLayer_] == HGCalTypes::WaferCenterR)) {
0340         rot = DDRotation(DDName(DDSplit(rotstr_).first, DDSplit(rotstr_).second));
0341 #ifdef EDM_ML_DEBUG
0342         rotName = rotstr_;
0343 #endif
0344       }
0345       cpv.position(glog, module, copy, r1, rot);
0346       int inc = ((layerSense_[ly] > 0) && (facingTypes_ > 1)) ? 2 : 1;
0347       copyNumber_[ii] = copy + inc;
0348 #ifdef EDM_ML_DEBUG
0349       edm::LogVerbatim("HGCalGeom") << "DDHGCalSiliconModule: " << glog.name() << " number " << copy
0350                                     << " positioned in " << module.name() << " at " << r1 << " with " << rotName
0351                                     << " rotation";
0352 #endif
0353       zz += hthick;
0354     }  // End of loop over layers in a block
0355     zi = zo;
0356     laymin = laymax;
0357     // Make consistency check of all the partitions of the block
0358     if (std::abs(thickTot - layerThick_[i]) >= tol2_) {
0359       if (thickTot > layerThick_[i]) {
0360         edm::LogError("HGCalGeom") << "Thickness of the partition " << layerThick_[i] << " is smaller than " << thickTot
0361                                    << ": thickness of all its components **** ERROR ****";
0362       } else {
0363         edm::LogWarning("HGCalGeom") << "Thickness of the partition " << layerThick_[i] << " does not match with "
0364                                      << thickTot << " of the components";
0365       }
0366     }
0367   }  // End of loop over blocks
0368 }
0369 
0370 void DDHGCalSiliconModule::positionSensitive(const DDLogicalPart& glog, int layer, DDCompactView& cpv) {
0371   static const double sqrt3 = std::sqrt(3.0);
0372   int layercenter = layerOrient_[layer];
0373   int layertype = (layerOrient_[layer] == HGCalTypes::WaferCenterB) ? 1 : 0;
0374   int firstWafer = waferLayerStart_[layer];
0375   int lastWafer = ((layer + 1 < static_cast<int>(waferLayerStart_.size())) ? waferLayerStart_[layer + 1]
0376                                                                            : static_cast<int>(waferIndex_.size()));
0377   double r = 0.5 * (waferSize_ + waferSepar_);
0378   double R = 2.0 * r / sqrt3;
0379   double dy = 0.75 * R;
0380   const auto& xyoff = geomTools_.shiftXY(layercenter, (waferSize_ + waferSepar_));
0381 #ifdef EDM_ML_DEBUG
0382   int ium(0), ivm(0), kount(0);
0383   std::vector<int> ntype(3, 0);
0384   edm::LogVerbatim("HGCalGeom") << "DDHGCalSiliconModule: " << glog.ddname() << "  r " << r << " R " << R << " dy "
0385                                 << dy << " Shift " << xyoff.first << ":" << xyoff.second << " WaferSize "
0386                                 << (waferSize_ + waferSepar_) << " index " << firstWafer << ":" << (lastWafer - 1);
0387 #endif
0388   for (int k = firstWafer; k < lastWafer; ++k) {
0389     int u = HGCalWaferIndex::waferU(waferIndex_[k]);
0390     int v = HGCalWaferIndex::waferV(waferIndex_[k]);
0391 #ifdef EDM_ML_DEBUG
0392     int iu = std::abs(u);
0393     int iv = std::abs(v);
0394 #endif
0395     int nr = 2 * v;
0396     int nc = -2 * u + v;
0397     double xpos = xyoff.first + nc * r;
0398     double ypos = xyoff.second + nr * dy;
0399     int type = HGCalProperty::waferThick(waferProperty_[k]);
0400     int part = HGCalProperty::waferPartial(waferProperty_[k]);
0401     int orien = HGCalProperty::waferOrient(waferProperty_[k]);
0402     std::string wafer;
0403     int i(999);
0404     if (part == HGCalTypes::WaferFull) {
0405       i = layertype * waferTypes_ + type;
0406       wafer = waferFull_[i];
0407     } else {
0408       i = (part - 1) * waferTypes_ * facingTypes_ * orientationTypes_ + layertype * waferTypes_ * orientationTypes_ +
0409           type * orientationTypes_ + orien;
0410 #ifdef EDM_ML_DEBUG
0411       edm::LogVerbatim("HGCalGeom") << " layertype:type:part:orien:ind " << layertype << ":" << type << ":" << part
0412                                     << ":" << orien << ":" << i << ":" << waferPart_.size();
0413 #endif
0414       wafer = waferPart_[i];
0415     }
0416     int copy = HGCalTypes::packTypeUV(type, u, v);
0417 #ifdef EDM_ML_DEBUG
0418     edm::LogVerbatim("HGCalGeom") << " DDHGCalSiliconModule: Layer " << HGCalWaferIndex::waferLayer(waferIndex_[k])
0419                                   << " Wafer " << wafer << " number " << copy << " type:part:orien:ind " << type << ":"
0420                                   << part << ":" << orien << ":" << i << " layer:u:v:indx " << (layer + firstLayer_)
0421                                   << ":" << u << ":" << v;
0422     if (iu > ium)
0423       ium = iu;
0424     if (iv > ivm)
0425       ivm = iv;
0426     kount++;
0427     if (copies_.count(copy) == 0)
0428       copies_.insert(copy);
0429 #endif
0430     DDTranslation tran(xpos, ypos, 0.0);
0431     DDRotation rotation;
0432     DDName name = DDName(DDSplit(wafer).first, DDSplit(wafer).second);
0433     cpv.position(name, glog.ddname(), copy, tran, rotation);
0434 #ifdef EDM_ML_DEBUG
0435     ++ntype[type];
0436     edm::LogVerbatim("HGCalGeom") << " DDHGCalSiliconModule: " << name << " number " << copy << " type " << layertype
0437                                   << ":" << type << " positioned in " << glog.ddname() << " at " << tran
0438                                   << " with no rotation";
0439 #endif
0440   }
0441 #ifdef EDM_ML_DEBUG
0442   edm::LogVerbatim("HGCalGeom") << "DDHGCalSiliconModule: Maximum # of u " << ium << " # of v " << ivm << " and "
0443                                 << kount << " wafers (" << ntype[0] << ":" << ntype[1] << ":" << ntype[2] << ") for "
0444                                 << glog.ddname();
0445 #endif
0446 }
0447 
0448 DEFINE_EDM_PLUGIN(DDAlgorithmFactory, DDHGCalSiliconModule, "hgcal:DDHGCalSiliconModule");