HGCalMixRotatedLayer

ddcms_det_element_DDCMS_hgcal_DDHGCalMixRotatedLayer

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612
///////////////////////////////////////////////////////////////////////////////
// File: DDHGCalMixRotatedLayer.cc
// Description: Geometry factory class for HGCal (Mix) adopted for DD4hep
///////////////////////////////////////////////////////////////////////////////

#include <cmath>
#include <memory>
#include <string>
#include <unordered_set>
#include <vector>

#include "Geometry/HGCalCommonData/interface/HGCalCell.h"
#include "Geometry/HGCalCommonData/interface/HGCalCassette.h"
#include "Geometry/HGCalCommonData/interface/HGCalGeomTools.h"
#include "Geometry/HGCalCommonData/interface/HGCalParameters.h"
#include "Geometry/HGCalCommonData/interface/HGCalProperty.h"
#include "Geometry/HGCalCommonData/interface/HGCalTileIndex.h"
#include "Geometry/HGCalCommonData/interface/HGCalTypes.h"
#include "Geometry/HGCalCommonData/interface/HGCalWaferIndex.h"
#include "Geometry/HGCalCommonData/interface/HGCalWaferType.h"
#include "DD4hep/DetFactoryHelper.h"
#include "DataFormats/Math/interface/angle_units.h"
#include "DetectorDescription/DDCMS/interface/DDPlugins.h"
#include "DetectorDescription/DDCMS/interface/DDutils.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"

//#define EDM_ML_DEBUG
using namespace angle_units::operators;

struct HGCalMixRotatedLayer {
  HGCalMixRotatedLayer() { throw cms::Exception("HGCalGeom") << "Wrong initialization to HGCalMixRotatedLayer"; }
  HGCalMixRotatedLayer(cms::DDParsingContext& ctxt, xml_h e) {
    cms::DDNamespace ns(ctxt, e, true);
    cms::DDAlgoArguments args(ctxt, e);

#ifdef EDM_ML_DEBUG
    edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: Creating an instance";
#endif

    static constexpr double tol1 = 0.01 * dd4hep::mm;
    dd4hep::Volume mother = ns.volume(args.parentName());

    waferTypes_ = args.value<int>("WaferTypes");
    facingTypes_ = args.value<int>("FacingTypes");
    orientationTypes_ = args.value<int>("OrientationTypes");
    placeOffset_ = args.value<int>("PlaceOffset");
    phiBinsScint_ = args.value<int>("NPhiBinScint");
    forFireworks_ = args.value<int>("ForFireWorks");
#ifdef EDM_ML_DEBUG
    edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer::Number of types of wafers: " << waferTypes_
                                  << " facings: " << facingTypes_ << " Orientations: " << orientationTypes_
                                  << " PlaceOffset: " << placeOffset_ << "; number of cells along phi " << phiBinsScint_
                                  << " forFireworks_: " << forFireworks_;
#endif
    firstLayer_ = args.value<int>("FirstLayer");
    absorbMode_ = args.value<int>("AbsorberMode");
    sensitiveMode_ = args.value<int>("SensitiveMode");
#ifdef EDM_ML_DEBUG
    edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer::First Layer " << firstLayer_ << " and "
                                  << "Absober:Sensitive mode " << absorbMode_ << ":" << sensitiveMode_;
#endif
    zMinBlock_ = args.value<double>("zMinBlock");
    waferSize_ = args.value<double>("waferSize");
    waferSepar_ = args.value<double>("SensorSeparation");
    sectors_ = args.value<int>("Sectors");
    cassettes_ = args.value<int>("Cassettes");
    alpha_ = (1._pi) / sectors_;
    cosAlpha_ = cos(alpha_);
#ifdef EDM_ML_DEBUG
    edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: zStart " << cms::convert2mm(zMinBlock_) << " wafer width "
                                  << cms::convert2mm(waferSize_) << " separations " << cms::convert2mm(waferSepar_)
                                  << " sectors " << sectors_ << ":" << convertRadToDeg(alpha_) << ":" << cosAlpha_
                                  << " with " << cassettes_ << " cassettes";
#endif
    slopeB_ = args.value<std::vector<double>>("SlopeBottom");
    zFrontB_ = args.value<std::vector<double>>("ZFrontBottom");
    rMinFront_ = args.value<std::vector<double>>("RMinFront");
    slopeT_ = args.value<std::vector<double>>("SlopeTop");
    zFrontT_ = args.value<std::vector<double>>("ZFrontTop");
    rMaxFront_ = args.value<std::vector<double>>("RMaxFront");
#ifdef EDM_ML_DEBUG
    for (unsigned int i = 0; i < slopeB_.size(); ++i)
      edm::LogVerbatim("HGCalGeom") << "Bottom Block [" << i << "] Zmin " << cms::convert2mm(zFrontB_[i]) << " Rmin "
                                    << cms::convert2mm(rMinFront_[i]) << " Slope " << slopeB_[i];
    for (unsigned int i = 0; i < slopeT_.size(); ++i)
      edm::LogVerbatim("HGCalGeom") << "Top Block [" << i << "] Zmin " << cms::convert2mm(zFrontT_[i]) << " Rmax "
                                    << cms::convert2mm(rMaxFront_[i]) << " Slope " << slopeT_[i];
#endif

    waferFull_ = args.value<std::vector<std::string>>("WaferNamesFull");
    waferPart_ = args.value<std::vector<std::string>>("WaferNamesPartial");
#ifdef EDM_ML_DEBUG
    edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << waferFull_.size() << " full and "
                                  << waferPart_.size() << " partial modules\nDDHGCalMixRotatedLayer:Full Modules:";
    unsigned int i1max = static_cast<unsigned int>(waferFull_.size());
    for (unsigned int i1 = 0; i1 < i1max; i1 += 2) {
      std::ostringstream st1;
      unsigned int i2 = std::min((i1 + 2), i1max);
      for (unsigned int i = i1; i < i2; ++i)
        st1 << " [" << i << "] " << waferFull_[i];
      edm::LogVerbatim("HGCalGeom") << st1.str();
    }
    edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: Partial Modules:";
    i1max = static_cast<unsigned int>(waferPart_.size());
    for (unsigned int i1 = 0; i1 < i1max; i1 += 2) {
      std::ostringstream st1;
      unsigned int i2 = std::min((i1 + 2), i1max);
      for (unsigned int i = i1; i < i2; ++i)
        st1 << " [" << i << "] " << waferPart_[i];
      edm::LogVerbatim("HGCalGeom") << st1.str();
    }
#endif

    materials_ = args.value<std::vector<std::string>>("MaterialNames");
    names_ = args.value<std::vector<std::string>>("VolumeNames");
    thick_ = args.value<std::vector<double>>("Thickness");
    copyNumber_.resize(materials_.size(), 1);
#ifdef EDM_ML_DEBUG
    edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << materials_.size() << " types of volumes";
    for (unsigned int i = 0; i < names_.size(); ++i)
      edm::LogVerbatim("HGCalGeom") << "Volume [" << i << "] " << names_[i] << " of thickness "
                                    << cms::convert2mm(thick_[i]) << " filled with " << materials_[i]
                                    << " first copy number " << copyNumber_[i];
#endif
    layers_ = args.value<std::vector<int>>("Layers");
    layerThick_ = args.value<std::vector<double>>("LayerThick");
#ifdef EDM_ML_DEBUG
    edm::LogVerbatim("HGCalGeom") << "There are " << layers_.size() << " blocks";
    for (unsigned int i = 0; i < layers_.size(); ++i)
      edm::LogVerbatim("HGCalGeom") << "Block [" << i << "] of thickness " << cms::convert2mm(layerThick_[i])
                                    << " with " << layers_[i] << " layers";
#endif
    layerType_ = args.value<std::vector<int>>("LayerType");
    layerSense_ = args.value<std::vector<int>>("LayerSense");
    layerOrient_ = args.value<std::vector<int>>("LayerTypes");
    for (unsigned int k = 0; k < layerOrient_.size(); ++k)
      layerOrient_[k] = HGCalTypes::layerType(layerOrient_[k]);
#ifdef EDM_ML_DEBUG
    for (unsigned int i = 0; i < layerOrient_.size(); ++i)
      edm::LogVerbatim("HGCalGeom") << "LayerTypes [" << i << "] " << layerOrient_[i];
#endif
    if (firstLayer_ > 0) {
      for (unsigned int i = 0; i < layerType_.size(); ++i) {
        if (layerSense_[i] > 0) {
          int ii = layerType_[i];
          copyNumber_[ii] = firstLayer_;
#ifdef EDM_ML_DEBUG
          edm::LogVerbatim("HGCalGeom") << "First copy number for layer type " << i << ":" << ii << " with "
                                        << materials_[ii] << " changed to " << copyNumber_[ii];
#endif
          break;
        }
      }
    } else {
      firstLayer_ = 1;
    }
#ifdef EDM_ML_DEBUG
    edm::LogVerbatim("HGCalGeom") << "There are " << layerType_.size() << " layers";
    for (unsigned int i = 0; i < layerType_.size(); ++i)
      edm::LogVerbatim("HGCalGeom") << "Layer [" << i << "] with material type " << layerType_[i] << " sensitive class "
                                    << layerSense_[i];
#endif
    materialTop_ = args.value<std::vector<std::string>>("TopMaterialNames");
    namesTop_ = args.value<std::vector<std::string>>("TopVolumeNames");
    layerThickTop_ = args.value<std::vector<double>>("TopLayerThickness");
    layerTypeTop_ = args.value<std::vector<int>>("TopLayerType");
    copyNumberTop_.resize(materialTop_.size(), firstLayer_);
#ifdef EDM_ML_DEBUG
    edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << materialTop_.size()
                                  << " types of volumes in the top part";
    for (unsigned int i = 0; i < materialTop_.size(); ++i)
      edm::LogVerbatim("HGCalGeom") << "Volume [" << i << "] " << namesTop_[i] << " of thickness "
                                    << cms::convert2mm(layerThickTop_[i]) << " filled with " << materialTop_[i]
                                    << " first copy number " << copyNumberTop_[i];
    edm::LogVerbatim("HGCalGeom") << "There are " << layerTypeTop_.size() << " layers in the top part";
    for (unsigned int i = 0; i < layerTypeTop_.size(); ++i)
      edm::LogVerbatim("HGCalGeom") << "Layer [" << i << "] with material type " << layerTypeTop_[i];
#endif
    waferIndex_ = args.value<std::vector<int>>("WaferIndex");
    waferProperty_ = args.value<std::vector<int>>("WaferProperties");
    waferLayerStart_ = args.value<std::vector<int>>("WaferLayerStart");
    cassetteShift_ = args.value<std::vector<double>>("CassetteShift");
#ifdef EDM_ML_DEBUG
    edm::LogVerbatim("HGCalGeom") << "waferProperties with " << waferIndex_.size() << " entries in "
                                  << waferLayerStart_.size() << " layers";
    for (unsigned int k = 0; k < waferLayerStart_.size(); ++k)
      edm::LogVerbatim("HGCalGeom") << "LayerStart[" << k << "] " << waferLayerStart_[k];
    for (unsigned int k = 0; k < waferIndex_.size(); ++k)
      edm::LogVerbatim("HGCalGeom") << "[" << k << "] " << waferIndex_[k] << " ("
                                    << HGCalWaferIndex::waferLayer(waferIndex_[k]) << ", "
                                    << HGCalWaferIndex::waferU(waferIndex_[k]) << ", "
                                    << HGCalWaferIndex::waferV(waferIndex_[k]) << ") : ("
                                    << HGCalProperty::waferThick(waferProperty_[k]) << ":"
                                    << HGCalProperty::waferPartial(waferProperty_[k]) << ":"
                                    << HGCalProperty::waferOrient(waferProperty_[k]) << ")";
    edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << cassetteShift_.size()
                                  << " elements for cassette shifts";
    unsigned int j1max = cassetteShift_.size();
    for (unsigned int j1 = 0; j1 < j1max; j1 += 6) {
      std::ostringstream st1;
      unsigned int j2 = std::min((j1 + 6), j1max);
      for (unsigned int j = j1; j < j2; ++j)
        st1 << " [" << j << "] " << std::setw(9) << cms::convert2mm(cassetteShift_[j]);
      edm::LogVerbatim("HGCalGeom") << st1.str();
    }
#endif
    tileRMin_ = args.value<std::vector<double>>("TileRMin");
    tileRMax_ = args.value<std::vector<double>>("TileRMax");
    tileIndex_ = args.value<std::vector<int>>("TileLayerRings");
    tilePhis_ = args.value<std::vector<int>>("TilePhiRange");
    tileLayerStart_ = args.value<std::vector<int>>("TileLayerStart");
#ifdef EDM_ML_DEBUG
    edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer:: with " << tileRMin_.size() << " rings";
    for (unsigned int k = 0; k < tileRMin_.size(); ++k)
      edm::LogVerbatim("HGCalGeom") << "Ring[" << k << "] " << cms::convert2mm(tileRMin_[k]) << " : "
                                    << cms::convert2mm(tileRMax_[k]);
    edm::LogVerbatim("HGCalGeom") << "TileProperties with " << tileIndex_.size() << " entries in "
                                  << tileLayerStart_.size() << " layers";
    for (unsigned int k = 0; k < tileLayerStart_.size(); ++k)
      edm::LogVerbatim("HGCalGeom") << "LayerStart[" << k << "] " << tileLayerStart_[k];
    for (unsigned int k = 0; k < tileIndex_.size(); ++k)
      edm::LogVerbatim("HGCalGeom") << "[" << k << "] " << tileIndex_[k] << " ("
                                    << "Layer " << std::get<0>(HGCalTileIndex::tileUnpack(tileIndex_[k])) << " Ring "
                                    << std::get<1>(HGCalTileIndex::tileUnpack(tileIndex_[k])) << ":"
                                    << std::get<2>(HGCalTileIndex::tileUnpack(tileIndex_[k])) << ") Phi "
                                    << std::get<1>(HGCalTileIndex::tileUnpack(tilePhis_[k])) << ":"
                                    << std::get<2>(HGCalTileIndex::tileUnpack(tilePhis_[k]));
#endif
    cassette_.setParameter(cassettes_, cassetteShift_);

#ifdef EDM_ML_DEBUG
    edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: NameSpace " << ns.name();

    edm::LogVerbatim("HGCalGeom") << "==>> Constructing DDHGCalMixRotatedLayer...";
    copies_.clear();
#endif

    double zi(zMinBlock_);
    int laymin(0);
    for (unsigned int i = 0; i < layers_.size(); i++) {
      double zo = zi + layerThick_[i];
      double routF = HGCalGeomTools::radius(zi, zFrontT_, rMaxFront_, slopeT_);
      int laymax = laymin + layers_[i];
      double zz = zi;
      double thickTot(0);
      for (int ly = laymin; ly < laymax; ++ly) {
        int ii = layerType_[ly];
        int copy = copyNumber_[ii];
        double hthick = 0.5 * thick_[ii];
        double rinB = HGCalGeomTools::radius(zo, zFrontB_, rMinFront_, slopeB_);
        zz += hthick;
        thickTot += thick_[ii];

        std::string name = names_[ii] + std::to_string(copy);
#ifdef EDM_ML_DEBUG
        edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: Layer " << ly << ":" << ii << " Front "
                                      << cms::convert2mm(zi) << ", " << cms::convert2mm(routF) << " Back "
                                      << cms::convert2mm(zo) << ", " << cms::convert2mm(rinB)
                                      << " superlayer thickness " << cms::convert2mm(layerThick_[i]);
#endif

        dd4hep::Material matter = ns.material(materials_[ii]);
        dd4hep::Volume glog;

        if (layerSense_[ly] < 1) {
          std::vector<double> pgonZ, pgonRin, pgonRout;
          double rmax =
              (std::min(routF, HGCalGeomTools::radius(zz + hthick, zFrontT_, rMaxFront_, slopeT_)) * cosAlpha_) - tol1;
          HGCalGeomTools::radius(zz - hthick,
                                 zz + hthick,
                                 zFrontB_,
                                 rMinFront_,
                                 slopeB_,
                                 zFrontT_,
                                 rMaxFront_,
                                 slopeT_,
                                 -layerSense_[ly],
                                 pgonZ,
                                 pgonRin,
                                 pgonRout);
          for (unsigned int isec = 0; isec < pgonZ.size(); ++isec) {
            pgonZ[isec] -= zz;
            if (layerSense_[ly] == 0 || absorbMode_ == 0)
              pgonRout[isec] = rmax;
            else
              pgonRout[isec] = pgonRout[isec] * cosAlpha_ - tol1;
          }

          dd4hep::Solid solid = dd4hep::Polyhedra(sectors_, -alpha_, 2._pi, pgonZ, pgonRin, pgonRout);
          ns.addSolidNS(ns.prepend(name), solid);
          glog = dd4hep::Volume(solid.name(), solid, matter);
          ns.addVolumeNS(glog);
#ifdef EDM_ML_DEBUG
          edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << solid.name() << " polyhedra of " << sectors_
                                        << " sectors covering " << convertRadToDeg(-alpha_) << ":"
                                        << convertRadToDeg(-alpha_ + 2._pi) << " with " << pgonZ.size() << " sections";
          for (unsigned int k = 0; k < pgonZ.size(); ++k)
            edm::LogVerbatim("HGCalGeom") << "[" << k << "] z " << cms::convert2mm(pgonZ[k]) << " R "
                                          << cms::convert2mm(pgonRin[k]) << ":" << cms::convert2mm(pgonRout[k]);
#endif
        } else {
          double rins =
              (sensitiveMode_ < 1) ? rinB : HGCalGeomTools::radius(zz + hthick, zFrontB_, rMinFront_, slopeB_);
          double routs =
              (sensitiveMode_ < 1) ? routF : HGCalGeomTools::radius(zz - hthick, zFrontT_, rMaxFront_, slopeT_);
          dd4hep::Solid solid = dd4hep::Tube(rins, routs, hthick, 0.0, 2._pi);
          ns.addSolidNS(ns.prepend(name), solid);
          glog = dd4hep::Volume(solid.name(), solid, matter);
          ns.addVolumeNS(glog);

#ifdef EDM_ML_DEBUG
          edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << solid.name() << " Tubs made of "
                                        << matter.name() << " of dimensions " << cms::convert2mm(rinB) << ":"
                                        << cms::convert2mm(rins) << ", " << cms::convert2mm(routF) << ":"
                                        << cms::convert2mm(routs) << ", " << cms::convert2mm(hthick)
                                        << ", 0.0, 360.0 and positioned in: " << glog.name() << " number " << copy;
#endif
          positionMix(ctxt, e, glog, name, copy, thick_[ii], matter);
        }

        dd4hep::Position r1(0, 0, zz);
        mother.placeVolume(glog, copy, r1);
        ++copyNumber_[ii];
#ifdef EDM_ML_DEBUG
        edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << glog.name() << " number " << copy
                                      << " positioned in " << mother.name() << " at (0,0," << cms::convert2mm(zz)
                                      << ") with no rotation";
#endif
        zz += hthick;
      }  // End of loop over layers in a block
      zi = zo;
      laymin = laymax;
      if (std::abs(thickTot - layerThick_[i]) > tol2_) {
        if (thickTot > layerThick_[i]) {
          edm::LogError("HGCalGeom") << "Thickness of the partition " << cms::convert2mm(layerThick_[i])
                                     << " is smaller than " << cms::convert2mm(thickTot)
                                     << ": thickness of all its components **** ERROR ****";
        } else {
          edm::LogWarning("HGCalGeom") << "Thickness of the partition " << cms::convert2mm(layerThick_[i])
                                       << " does not match with " << cms::convert2mm(thickTot) << " of the components";
        }
      }
    }  // End of loop over blocks

#ifdef EDM_ML_DEBUG
    edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << copies_.size() << " different wafer copy numbers";
    int k(0);
    for (std::unordered_set<int>::const_iterator itr = copies_.begin(); itr != copies_.end(); ++itr, ++k) {
      edm::LogVerbatim("HGCalGeom") << "Copy [" << k << "] : " << (*itr);
    }
    copies_.clear();
    edm::LogVerbatim("HGCalGeom") << "<<== End of DDHGCalMixRotatedLayer construction...";
#endif
  }

  void positionMix(cms::DDParsingContext& ctxt,
                   xml_h e,
                   const dd4hep::Volume& glog,
                   const std::string& nameM,
                   int copyM,
                   double thick,
                   const dd4hep::Material& matter) {
    cms::DDNamespace ns(ctxt, e, true);

    // Make the top part first
    for (unsigned int ly = 0; ly < layerTypeTop_.size(); ++ly) {
      int ii = layerTypeTop_[ly];
      copyNumberTop_[ii] = copyM;
    }
    double hthick = 0.5 * thick;
    double dphi = (2._pi) / phiBinsScint_;
    double thickTot(0), zpos(-hthick);
    for (unsigned int ly = 0; ly < layerTypeTop_.size(); ++ly) {
      int ii = layerTypeTop_[ly];
      int copy = copyNumberTop_[ii];
      int layer = copy - firstLayer_;
      double hthickl = 0.5 * layerThickTop_[ii];
      thickTot += layerThickTop_[ii];
      zpos += hthickl;
      dd4hep::Material matter1 = ns.material(materialTop_[ii]);
      unsigned int k = 0;
      int firstTile = tileLayerStart_[layer];
      int lastTile = ((layer + 1 < static_cast<int>(tileLayerStart_.size())) ? tileLayerStart_[layer + 1]
                                                                             : static_cast<int>(tileIndex_.size()));
#ifdef EDM_ML_DEBUG
      edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: Layer " << ly << ":" << ii << " Copy " << copy
                                    << " Tiles " << firstTile << ":" << lastTile;
#endif
      for (int ti = firstTile; ti < lastTile; ++ti) {
        double r1 = tileRMin_[std::get<1>(HGCalTileIndex::tileUnpack(tileIndex_[ti])) - 1];
        double r2 = tileRMax_[std::get<2>(HGCalTileIndex::tileUnpack(tileIndex_[ti])) - 1];
        int cassette = std::get<0>(HGCalTileIndex::tileUnpack(tilePhis_[ti]));
        int fimin = std::get<1>(HGCalTileIndex::tileUnpack(tilePhis_[ti]));
        int fimax = std::get<2>(HGCalTileIndex::tileUnpack(tilePhis_[ti]));
        double phi1 = dphi * (fimin - 1);
        double phi2 = (forFireworks_ == 1) ? (dphi * (fimax - fimin + 1)) : (dphi * fimax);
        auto cshift = cassette_.getShift(layer + 1, 1, cassette);
#ifdef EDM_ML_DEBUG
        int cassette0 = HGCalCassette::cassetteType(2, 1, cassette);  //
        edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: Layer " << copy << ":" << (layer + 1) << " iR "
                                      << std::get<1>(HGCalTileIndex::tileUnpack(tileIndex_[ti])) << ":"
                                      << std::get<2>(HGCalTileIndex::tileUnpack(tileIndex_[ti])) << " R "
                                      << cms::convert2mm(r1) << ":" << cms::convert2mm(r2) << " Thick "
                                      << cms::convert2mm((2.0 * hthickl)) << " phi " << fimin << ":" << fimax << ":"
                                      << convertRadToDeg(phi1) << ":" << convertRadToDeg(phi2) << " cassette "
                                      << cassette << ":" << cassette0 << " Shift " << cms::convert2mm(cshift.first)
                                      << ":" << cms::convert2mm(cshift.second);
#endif
        std::string name = namesTop_[ii] + "L" + std::to_string(copy) + "F" + std::to_string(k);
        ++k;
        dd4hep::Solid solid = dd4hep::Tube(r1, r2, hthickl, phi1, phi2);
        ns.addSolidNS(ns.prepend(name), solid);
        dd4hep::Volume glog1 = dd4hep::Volume(solid.name(), solid, matter1);
        ns.addVolumeNS(glog1);
#ifdef EDM_ML_DEBUG
        edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << glog1.name() << " Tubs made of "
                                      << materialTop_[ii] << " of dimensions " << cms::convert2mm(r1) << ", "
                                      << cms::convert2mm(r2) << ", " << cms::convert2mm(hthickl) << ", "
                                      << convertRadToDeg(phi1) << ", " << convertRadToDeg(phi2);
#endif
        dd4hep::Position tran(-cshift.first, cshift.second, zpos);
        glog.placeVolume(glog1, copy, tran);
#ifdef EDM_ML_DEBUG
        edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: Position " << glog1.name() << " number " << copy
                                      << " in " << glog.name() << " at (" << cms::convert2mm(cshift.first) << ", "
                                      << cms::convert2mm(cshift.second) << ", " << cms::convert2mm(zpos)
                                      << ") with no rotation";
#endif
      }
      ++copyNumberTop_[ii];
      zpos += hthickl;
    }
    if (std::abs(thickTot - thick) > tol2_) {
      if (thickTot > thick) {
        edm::LogError("HGCalGeom") << "DDHGCalMixRotatedLayer: Thickness of the partition " << cms::convert2mm(thick)
                                   << " is smaller than " << cms::convert2mm(thickTot)
                                   << ": thickness of all its components in the top part **** ERROR ****";
      } else {
        edm::LogWarning("HGCalGeom") << "DDHGCalMixRotatedLayer: Thickness of the partition " << cms::convert2mm(thick)
                                     << " does not match with " << cms::convert2mm(thickTot)
                                     << " of the components in top part";
      }
    }

    // Make the bottom part next
    int layer = (copyM - firstLayer_);
    static const double sqrt3 = std::sqrt(3.0);
    int layercenter = layerOrient_[layer];
    int layertype = HGCalTypes::layerFrontBack(layerOrient_[layer]);
    int firstWafer = waferLayerStart_[layer];
    int lastWafer = ((layer + 1 < static_cast<int>(waferLayerStart_.size())) ? waferLayerStart_[layer + 1]
                                                                             : static_cast<int>(waferIndex_.size()));
    double delx = 0.5 * (waferSize_ + waferSepar_);
    double dely = 2.0 * delx / sqrt3;
    double dy = 0.75 * dely;
    const auto& xyoff = geomTools_.shiftXY(layercenter, (waferSize_ + waferSepar_));
#ifdef EDM_ML_DEBUG
    int ium(0), ivm(0), kount(0);
    std::vector<int> ntype(3, 0);
    edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: " << glog.name() << "  r " << cms::convert2mm(delx)
                                  << " R " << cms::convert2mm(dely) << " dy " << cms::convert2mm(dy) << " Shift "
                                  << cms::convert2mm(xyoff.first) << ":" << cms::convert2mm(xyoff.second)
                                  << " WaferSize " << cms::convert2mm((waferSize_ + waferSepar_)) << " index "
                                  << firstWafer << ":" << (lastWafer - 1);
#endif
    for (int k = firstWafer; k < lastWafer; ++k) {
      int u = HGCalWaferIndex::waferU(waferIndex_[k]);
      int v = HGCalWaferIndex::waferV(waferIndex_[k]);
#ifdef EDM_ML_DEBUG
      int iu = std::abs(u);
      int iv = std::abs(v);
#endif
      int nr = 2 * v;
      int nc = -2 * u + v;
      int type = HGCalProperty::waferThick(waferProperty_[k]);
      int part = HGCalProperty::waferPartial(waferProperty_[k]);
      int orien = HGCalProperty::waferOrient(waferProperty_[k]);
      int cassette = HGCalProperty::waferCassette(waferProperty_[k]);
      int place = HGCalCell::cellPlacementIndex(1, layertype, orien);
#ifdef EDM_ML_DEBUG
      edm::LogVerbatim("HGCalGeom")
          << "DDHGCalMixRotatedLayer::index:Property:layertype:type:part:orien:cassette:place:offsets:ind " << k << ":"
          << waferProperty_[k] << ":" << layertype << ":" << type << ":" << part << ":" << orien << ":" << cassette
          << ":" << place;
#endif
      auto cshift = cassette_.getShift(layer + 1, -1, cassette);
      double xpos = xyoff.first - cshift.first + nc * delx;
      double ypos = xyoff.second + cshift.second + nr * dy;
#ifdef EDM_ML_DEBUG
      double xorig = xyoff.first + nc * delx;
      double yorig = xyoff.second + nr * dy;
      double angle = std::atan2(yorig, xorig);
      edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer::Wafer: layer " << layer + 1 << " cassette " << cassette
                                    << " Shift " << cms::convert2mm(cshift.first) << ":"
                                    << cms::convert2mm(cshift.second) << " Original " << cms::convert2mm(xorig) << ":"
                                    << cms::convert2mm(yorig) << ":" << convertRadToDeg(angle) << " Final "
                                    << cms::convert2mm(xpos) << ":" << cms::convert2mm(ypos);
#endif
      std::string wafer;
      int i(999);
      if (part == HGCalTypes::WaferFull) {
        i = type * facingTypes_ * orientationTypes_ + place - placeOffset_;
#ifdef EDM_ML_DEBUG
        edm::LogVerbatim("HGCalGeom") << " FullWafer type:place:ind " << type << ":" << place << ":" << i << ":"
                                      << waferFull_.size();
#endif
        wafer = waferFull_[i];
      } else {
        int partoffset =
            (part >= HGCalTypes::WaferHDTop) ? HGCalTypes::WaferPartHDOffset : HGCalTypes::WaferPartLDOffset;
        i = (part - partoffset) * facingTypes_ * orientationTypes_ +
            HGCalTypes::WaferTypeOffset[type] * facingTypes_ * orientationTypes_ + place - placeOffset_;
#ifdef EDM_ML_DEBUG
        edm::LogVerbatim("HGCalGeom") << " layertype:type:part:orien:cassette:place:offsets:ind " << layertype << ":"
                                      << type << ":" << part << ":" << orien << ":" << cassette << ":" << place << ":"
                                      << partoffset << ":" << HGCalTypes::WaferTypeOffset[type] << ":" << i << ":"
                                      << waferPart_.size();
#endif
        wafer = waferPart_[i];
      }
      int copy = HGCalTypes::packTypeUV(type, u, v);
#ifdef EDM_ML_DEBUG
      edm::LogVerbatim("HGCalGeom") << " DDHGCalMixRotatedLayer: Layer " << HGCalWaferIndex::waferLayer(waferIndex_[k])
                                    << " Wafer " << wafer << " number " << copy << " type :part:orien:ind " << type
                                    << ":" << part << ":" << orien << ":" << i << " layer:u:v " << (layer + firstLayer_)
                                    << ":" << u << ":" << v;
      if (iu > ium)
        ium = iu;
      if (iv > ivm)
        ivm = iv;
      kount++;
      if (copies_.count(copy) == 0)
        copies_.insert(copy);
#endif
      dd4hep::Position tran(xpos, ypos, 0.0);
      glog.placeVolume(ns.volume(wafer), copy, tran);
#ifdef EDM_ML_DEBUG
      ++ntype[type];
      edm::LogVerbatim("HGCalGeom") << " DDHGCalMixRotatedLayer: " << wafer << " number " << copy << " type "
                                    << layertype << ":" << type << " positioned in " << glog.name() << " at ("
                                    << cms::convert2mm(xpos) << "," << cms::convert2mm(ypos) << ",0) with no rotation";
#endif
    }

#ifdef EDM_ML_DEBUG
    edm::LogVerbatim("HGCalGeom") << "DDHGCalMixRotatedLayer: Maximum # of u " << ium << " # of v " << ivm << " and "
                                  << kount << " wafers (" << ntype[0] << ":" << ntype[1] << ":" << ntype[2] << ") for "
                                  << glog.name();
#endif
  }

  //Required data members to cache the values from XML file
  HGCalGeomTools geomTools_;
  HGCalCassette cassette_;

  static constexpr double tol2_ = 0.00001 * dd4hep::mm;

  int waferTypes_;                        // Number of wafer types
  int facingTypes_;                       // Types of facings of modules toward IP
  int orientationTypes_;                  // Number of partial wafer orienations
  int placeOffset_;                       // Offset for placement
  int phiBinsScint_;                      // Maximum number of cells along phi
  int forFireworks_;                      // Needed for Fireworks(1)/Geant4(0)
  int firstLayer_;                        // Copy # of the first sensitive layer
  int absorbMode_;                        // Absorber mode
  int sensitiveMode_;                     // Sensitive mode
  double zMinBlock_;                      // Starting z-value of the block
  double waferSize_;                      // Width of the wafer
  double waferSepar_;                     // Sensor separation
  int sectors_;                           // Sectors
  int cassettes_;                         // Cassettes
  std::vector<double> slopeB_;            // Slope at the lower R
  std::vector<double> zFrontB_;           // Starting Z values for the slopes
  std::vector<double> rMinFront_;         // Corresponding rMin's
  std::vector<double> slopeT_;            // Slopes at the larger R
  std::vector<double> zFrontT_;           // Starting Z values for the slopes
  std::vector<double> rMaxFront_;         // Corresponding rMax's
  std::vector<std::string> waferFull_;    // Names of full wafer modules
  std::vector<std::string> waferPart_;    // Names of partial wafer modules
  std::vector<std::string> materials_;    // Materials
  std::vector<std::string> names_;        // Names
  std::vector<double> thick_;             // Thickness of the material
  std::vector<int> copyNumber_;           // Initial copy numbers
  std::vector<int> layers_;               // Number of layers in a section
  std::vector<double> layerThick_;        // Thickness of each section
  std::vector<int> layerType_;            // Type of the layer
  std::vector<int> layerSense_;           // Content of a layer (sensitive?)
  std::vector<std::string> materialTop_;  // Materials of top layers
  std::vector<std::string> namesTop_;     // Names of top layers
  std::vector<double> layerThickTop_;     // Thickness of the top sections
  std::vector<int> layerTypeTop_;         // Type of the Top layer
  std::vector<int> copyNumberTop_;        // Initial copy numbers (top section)
  std::vector<int> layerOrient_;          // Layer orientation for the silicon component
  std::vector<int> waferIndex_;           // Wafer index for the types
  std::vector<int> waferProperty_;        // Wafer property
  std::vector<int> waferLayerStart_;      // Start index of wafers in each layer
  std::vector<double> cassetteShift_;     // Shifts of the cassetes
  std::vector<double> tileRMin_;          // Minimum radius of each ring
  std::vector<double> tileRMax_;          // Maximum radius of each ring
  std::vector<int> tileIndex_;            // Index of tile (layer/start|end ring)
  std::vector<int> tilePhis_;             // Tile phi range for each index
  std::vector<int> tileLayerStart_;       // Start index of tiles in each layer
  std::unordered_set<int> copies_;        // List of copy #'s
  double alpha_, cosAlpha_;
};

static long algorithm(dd4hep::Detector& /* description */, cms::DDParsingContext& ctxt, xml_h e) {
  HGCalMixRotatedLayer mixRotatedLayerAlgo(ctxt, e);
  return cms::s_executed;
}

DECLARE_DDCMS_DETELEMENT(DDCMS_hgcal_DDHGCalMixRotatedLayer, algorithm)