Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /////////////////////////////////////////////////////////////////////////////
0002 // File: DDHCalBarrelAlgo.cc
0003 //   adapted from CCal(G4)HcalBarrel.cc
0004 // Description: Geometry factory class for Hcal Barrel
0005 ///////////////////////////////////////////////////////////////////////////////
0006 
0007 #include <cmath>
0008 #include <algorithm>
0009 #include <map>
0010 #include <string>
0011 #include <vector>
0012 
0013 #include "DataFormats/Math/interface/angle_units.h"
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015 #include "DetectorDescription/Core/interface/DDSplit.h"
0016 #include "DetectorDescription/DDCMS/interface/DDPlugins.h"
0017 #include "DetectorDescription/DDCMS/interface/DDutils.h"
0018 #include "DD4hep/DetFactoryHelper.h"
0019 
0020 //#define EDM_ML_DEBUG
0021 using namespace angle_units::operators;
0022 
0023 struct HcalBarrelAlgo {
0024   //General Volume
0025   //      <----- Zmax ------>
0026   //Router************************-------
0027   //      *                      *Rstep2|        Theta angle w.r.t. vertical
0028   //      *                      *---------------
0029   //      *                     *               |
0030   //      *                    *Theta[i]        Rmax[i]
0031   //      *                   *---------------  |
0032   //                        *Theta[0] Rmax[0]|  |
0033   //Rinner*****************----------------------
0034 
0035   std::string genMaterial;     //General material
0036   int nsectors;                //Number of potenital straight edges
0037   int nsectortot;              //Number of straight edges (actual)
0038   int nhalf;                   //Number of half modules
0039   double rinner, router;       //See picture
0040   int rzones;                  //  ....
0041   std::vector<double> theta;   //  .... (in degrees)
0042   std::vector<double> rmax;    //  ....
0043   std::vector<double> zoff;    //  ....
0044   std::vector<double> ttheta;  //tan(theta)
0045   std::string rotHalf;         //Rotation matrix of the second half
0046   std::string rotns;           //Name space for Rotation matrices
0047 
0048   //Upper layers inside general volume
0049   //     <---- Zout ---->
0050   //  |  ****************     |
0051   //  |  *              *     Wstep
0052   //  W  *              ***** |
0053   //  |  *                  *
0054   //  |  ********************
0055   //     <------ Zin ------->
0056   // Middle layers inside general volume
0057   //     <------ Zout ------>         Zout = Full sector Z at position
0058   //  |  ********************         Zin  = Full sector Z at position
0059   //  |  *                 *
0060   //  W  *                * Angle = Theta sector
0061   //  |  *               *  )
0062   //  |  ****************--------
0063   //     <------ Zin ------->
0064 
0065   // Lower layers
0066   //     <------ Zout ------>         Zin(i)=Zout(i-1)
0067   //  |  ********************         Zout(i)=Zin(i)+W(i)/tan(Theta(i))
0068   //  |  *                 *
0069   //  W  *                *  Theta
0070   //  |  *               *
0071   //  |  ****************--------
0072   //     <--- Zin ------>
0073 
0074   int nLayers;                          //Number of layers
0075   std::vector<int> layerId;             //Number identification
0076   std::vector<std::string> layerLabel;  //String identification
0077   std::vector<std::string> layerMat;    //Material
0078   std::vector<double> layerWidth;       //W in picture
0079   std::vector<double> layerD1;          //d1 in front picture
0080   std::vector<double> layerD2;          //d2 in front picture
0081   std::vector<double> layerAlpha;       //Angular width of the middle tiles
0082   std::vector<double> layerT1;          //t in front picture (side)
0083   std::vector<double> layerT2;          //t in front picture (front)
0084   std::vector<int> layerAbsorb;         //Absorber flag
0085   std::vector<double> layerGap;         //Gap at the edge
0086 
0087   int nAbsorber;                        //Number of absorber layers in middle
0088   std::vector<std::string> absorbName;  //Absorber name
0089   std::vector<std::string> absorbMat;   //Absorber material
0090   std::vector<double> absorbD;          //Distance from the bottom surface
0091   std::vector<double> absorbT;          //Thickness
0092   std::string middleMat;                //Material of the detector layer
0093   double middleD;                       //Distance from the bottom surface
0094   double middleW;                       //Half width
0095   int nMidAbs;                          //Number of absorbers in front part
0096   std::vector<std::string> midName;     //Absorber names in the front part
0097   std::vector<std::string> midMat;      //Absorber material
0098   std::vector<double> midW;             //Half width
0099   std::vector<double> midT;             //Thickness
0100 
0101   std::vector<std::string> sideMat;      //Material for special side layers
0102   std::vector<double> sideD;             //Depth from bottom surface
0103   std::vector<double> sideT;             //Thickness
0104   int nSideAbs;                          //Number of absorbers in special side
0105   std::vector<std::string> sideAbsName;  //Absorber name
0106   std::vector<std::string> sideAbsMat;   //Absorber material
0107   std::vector<double> sideAbsW;          //Half width
0108 
0109   // Detectors. Each volume inside the layer has the shape:
0110   //
0111   // ******************************* |
0112   // *\\\\\\\Plastic\\\\\\\\\\\\\\\* T2
0113   // ******************************* |
0114   // *////Scintillator/////////////* Tsc
0115   // ******************************* |
0116   // *\\\\\\\Plastic\\\\\\\\\\\\\\\* T1
0117   // ******************************* |   |
0118   // *         Air                 *     dP1
0119   // *******************************     |
0120   //
0121   std::string detMat;    //fill material
0122   std::string detRot;    //Rotation matrix for the 2nd
0123   std::string detMatPl;  //Plastic material
0124   std::string detMatSc;  //Scintillator material
0125   std::vector<int> detType;
0126   std::vector<double> detdP1;     //Air gap (side)
0127   std::vector<double> detdP2;     //Air gap (centre)
0128   std::vector<double> detT11;     //Back plastic thickness (side)
0129   std::vector<double> detT12;     //Back plastic thickness (centre)
0130   std::vector<double> detTsc;     //Scintillator
0131   std::vector<double> detT21;     //Front plastic thickness (side)
0132   std::vector<double> detT22;     //Front plastic thickness (centre)
0133   std::vector<double> detWidth1;  //Width of phi(1,4) megatiles
0134   std::vector<double> detWidth2;  //Width of phi(2,3) megatiles
0135   std::vector<int> detPosY;       //Positioning of phi(1,4) tiles - 0 centre
0136 
0137   std::string idName;       //Name of the "parent" volume.
0138   std::string idNameSpace;  //Namespace of this and ALL sub-parts
0139   int idOffset;             // Geant4 ID's...    = 3000;
0140 
0141   HcalBarrelAlgo() = delete;
0142 
0143   HcalBarrelAlgo(cms::DDParsingContext& ctxt, xml_h& e) {
0144     cms::DDNamespace ns(ctxt, e, true);
0145     cms::DDAlgoArguments args(ctxt, e);
0146 #ifdef EDM_ML_DEBUG
0147     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: Creating an instance";
0148 #endif
0149 
0150     genMaterial = args.value<std::string>("MaterialName");
0151     nsectors = args.value<int>("NSector");
0152     nsectortot = args.value<int>("NSectorTot");
0153     nhalf = args.value<int>("NHalf");
0154     rinner = args.value<double>("RIn");
0155     router = args.value<double>("ROut");
0156     rzones = args.value<int>("RZones");
0157     rotHalf = args.value<std::string>("RotHalf");
0158     rotns = args.value<std::string>("RotNameSpace");
0159 
0160     theta = args.value<std::vector<double> >("Theta");
0161     rmax = args.value<std::vector<double> >("RMax");
0162     zoff = args.value<std::vector<double> >("ZOff");
0163     for (int i = 0; i < rzones; i++) {
0164       ttheta.emplace_back(tan(theta[i]));  //*deg already done in XML
0165     }
0166     if (rzones > 3)
0167       rmax[2] = (zoff[3] - zoff[2]) / ttheta[2];
0168 
0169 #ifdef EDM_ML_DEBUG
0170     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: General material " << genMaterial << "\tSectors " << nsectors
0171                                  << ", " << nsectortot << "\tHalves " << nhalf << "\tRotation matrix " << rotns << ":"
0172                                  << rotHalf << "\n\t\t" << cms::convert2mm(rinner) << "\t" << cms::convert2mm(router)
0173                                  << "\t" << rzones;
0174     for (int i = 0; i < rzones; i++)
0175       edm::LogVerbatim("HCalGeom") << "\tTheta[" << i << "] = " << theta[i] << "\trmax[" << i
0176                                    << "] = " << cms::convert2mm(rmax[i]) << "\tzoff[" << i
0177                                    << "] = " << cms::convert2mm(zoff[i]);
0178 #endif
0179     ///////////////////////////////////////////////////////////////
0180     //Layers
0181     nLayers = args.value<int>("NLayers");
0182 #ifdef EDM_ML_DEBUG
0183     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: Layer\t" << nLayers;
0184 #endif
0185     layerId = args.value<std::vector<int> >("Id");
0186     layerLabel = args.value<std::vector<std::string> >("LayerLabel");
0187     layerMat = args.value<std::vector<std::string> >("LayerMat");
0188     layerWidth = args.value<std::vector<double> >("LayerWidth");
0189     layerD1 = args.value<std::vector<double> >("D1");
0190     layerD2 = args.value<std::vector<double> >("D2");
0191     layerAlpha = args.value<std::vector<double> >("Alpha2");
0192     layerT1 = args.value<std::vector<double> >("T1");
0193     layerT2 = args.value<std::vector<double> >("T2");
0194     layerAbsorb = args.value<std::vector<int> >("AbsL");
0195     layerGap = args.value<std::vector<double> >("Gap");
0196 #ifdef EDM_ML_DEBUG
0197     for (int i = 0; i < nLayers; i++)
0198       edm::LogVerbatim("HCalGeom") << layerLabel[i] << "\t" << layerId[i] << "\t" << layerMat[i] << "\t"
0199                                    << cms::convert2mm(layerWidth[i]) << "\t" << cms::convert2mm(layerD1[i]) << "\t"
0200                                    << cms::convert2mm(layerD2[i]) << "\t" << layerAlpha[i] << "\t"
0201                                    << cms::convert2mm(layerT1[i]) << "\t" << cms::convert2mm(layerT2[i]) << "\t"
0202                                    << layerAbsorb[i] << "\t" << cms::convert2mm(layerGap[i]);
0203 #endif
0204 
0205     ///////////////////////////////////////////////////////////////
0206     //Absorber Layers and middle part
0207     absorbName = args.value<std::vector<std::string> >("AbsorbName");
0208     absorbMat = args.value<std::vector<std::string> >("AbsorbMat");
0209     absorbD = args.value<std::vector<double> >("AbsorbD");
0210     absorbT = args.value<std::vector<double> >("AbsorbT");
0211     nAbsorber = absorbName.size();
0212 #ifdef EDM_ML_DEBUG
0213     for (int i = 0; i < nAbsorber; i++)
0214       edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << absorbName[i] << " Material " << absorbMat[i] << " d "
0215                                    << cms::convert2mm(absorbD[i]) << " t " << cms::convert2mm(absorbT[i]);
0216 #endif
0217     middleMat = args.value<std::string>("MiddleMat");
0218     middleD = args.value<double>("MiddleD");
0219     middleW = args.value<double>("MiddleW");
0220 #ifdef EDM_ML_DEBUG
0221     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: Middle material " << middleMat << " d "
0222                                  << cms::convert2mm(middleD) << " w " << cms::convert2mm(middleW);
0223 #endif
0224     midName = args.value<std::vector<std::string> >("MidAbsName");
0225     midMat = args.value<std::vector<std::string> >("MidAbsMat");
0226     midW = args.value<std::vector<double> >("MidAbsW");
0227     midT = args.value<std::vector<double> >("MidAbsT");
0228     nMidAbs = midName.size();
0229 #ifdef EDM_ML_DEBUG
0230     for (int i = 0; i < nMidAbs; i++)
0231       edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << midName[i] << " Material " << midMat[i] << " W "
0232                                    << cms::convert2mm(midW[i]) << " T " << cms::convert2mm(midT[i]);
0233 #endif
0234 
0235     //Absorber layers in the side part
0236     sideMat = args.value<std::vector<std::string> >("SideMat");
0237     sideD = args.value<std::vector<double> >("SideD");
0238     sideT = args.value<std::vector<double> >("SideT");
0239 #ifdef EDM_ML_DEBUG
0240     for (unsigned int i = 0; i < sideMat.size(); i++)
0241       edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: Side material " << sideMat[i] << " d "
0242                                    << cms::convert2mm(sideD[i]) << " t " << cms::convert2mm(sideT[i]);
0243 #endif
0244     sideAbsName = args.value<std::vector<std::string> >("SideAbsName");
0245     sideAbsMat = args.value<std::vector<std::string> >("SideAbsMat");
0246     sideAbsW = args.value<std::vector<double> >("SideAbsW");
0247     nSideAbs = sideAbsName.size();
0248 #ifdef EDM_ML_DEBUG
0249     for (int i = 0; i < nSideAbs; i++)
0250       edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << sideAbsName[i] << " Material " << sideAbsMat[i] << " W "
0251                                    << cms::convert2mm(sideAbsW[i]);
0252 #endif
0253 
0254     ///////////////////////////////////////////////////////////////
0255     // Detectors
0256 
0257     detMat = args.value<std::string>("DetMat");
0258     detRot = args.value<std::string>("DetRot");
0259     detMatPl = args.value<std::string>("DetMatPl");
0260     detMatSc = args.value<std::string>("DetMatSc");
0261 #ifdef EDM_ML_DEBUG
0262     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: Detector (" << nLayers << ") Rotation matrix " << rotns << ":"
0263                                  << detRot << "\n\t\t" << detMat << "\t" << detMatPl << "\t" << detMatSc;
0264 #endif
0265     detType = args.value<std::vector<int> >("DetType");
0266     detdP1 = args.value<std::vector<double> >("DetdP1");
0267     detdP2 = args.value<std::vector<double> >("DetdP2");
0268     detT11 = args.value<std::vector<double> >("DetT11");
0269     detT12 = args.value<std::vector<double> >("DetT12");
0270     detTsc = args.value<std::vector<double> >("DetTsc");
0271     detT21 = args.value<std::vector<double> >("DetT21");
0272     detT22 = args.value<std::vector<double> >("DetT22");
0273     detWidth1 = args.value<std::vector<double> >("DetWidth1");
0274     detWidth2 = args.value<std::vector<double> >("DetWidth2");
0275     detPosY = args.value<std::vector<int> >("DetPosY");
0276 #ifdef EDM_ML_DEBUG
0277     for (int i = 0; i < nLayers; i++)
0278       edm::LogVerbatim("HCalGeom") << i + 1 << "\t" << detType[i] << "\t" << cms::convert2mm(detdP1[i]) << ", "
0279                                    << cms::convert2mm(detdP2[i]) << "\t" << cms::convert2mm(detT11[i]) << ", "
0280                                    << cms::convert2mm(detT12[i]) << "\t" << cms::convert2mm(detTsc[i]) << "\t"
0281                                    << cms::convert2mm(detT21[i]) << ", " << cms::convert2mm(detT22[i]) << "\t"
0282                                    << cms::convert2mm(detWidth1[i]) << "\t" << cms::convert2mm(detWidth2[i]) << "\t"
0283                                    << detPosY[i];
0284 #endif
0285 
0286     //  idName = parentName.name();
0287     idName = args.value<std::string>("MotherName");
0288     idNameSpace = ns.name();
0289     idOffset = args.value<int>("IdOffset");
0290 #ifdef EDM_ML_DEBUG
0291     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: Parent " << args.parentName() << " idName " << idName
0292                                  << " NameSpace " << idNameSpace << " Offset " << idOffset;
0293     edm::LogVerbatim("HCalGeom") << "==>> Constructing DDHCalBarrelAlgo...\n"
0294                                  << "DDHCalBarrelAlgo: General volume...";
0295 #endif
0296 
0297     double alpha = (1._pi) / nsectors;
0298     double dphi = nsectortot * (2._pi) / nsectors;
0299     int nsec, ntot(15);
0300     if (nhalf == 1)
0301       nsec = 8;
0302     else
0303       nsec = 15;
0304     int nf = ntot - nsec;
0305 
0306     //Calculate zmin... see HCalBarrel.hh picture. For polyhedra
0307     //Rmin and Rmax are distances to vertex
0308     double zmax = zoff[3];
0309     double zstep5 = zoff[4];
0310     double zstep4 = (zoff[1] + rmax[1] * ttheta[1]);
0311     if ((zoff[2] + rmax[1] * ttheta[2]) > zstep4)
0312       zstep4 = (zoff[2] + rmax[1] * ttheta[2]);
0313     double zstep3 = (zoff[1] + rmax[0] * ttheta[1]);
0314     double zstep2 = (zoff[0] + rmax[0] * ttheta[0]);
0315     double zstep1 = (zoff[0] + rinner * ttheta[0]);
0316     double rout = router;
0317     double rout1 = rmax[3];
0318     double rin = rinner;
0319     double rmid1 = rmax[0];
0320     double rmid2 = rmax[1];
0321     double rmid3 = (zoff[4] - zoff[2]) / ttheta[2];
0322     double rmid4 = rmax[2];
0323 
0324     std::vector<double> pgonZ = {-zmax,
0325                                  -zstep5,
0326                                  -zstep5,
0327                                  -zstep4,
0328                                  -zstep3,
0329                                  -zstep2,
0330                                  -zstep1,
0331                                  0,
0332                                  zstep1,
0333                                  zstep2,
0334                                  zstep3,
0335                                  zstep4,
0336                                  zstep5,
0337                                  zstep5,
0338                                  zmax};
0339 
0340     std::vector<double> pgonRmin = {
0341         rmid4, rmid3, rmid3, rmid2, rmid1, rmid1, rin, rin, rin, rmid1, rmid1, rmid2, rmid3, rmid3, rmid4};
0342 
0343     std::vector<double> pgonRmax = {
0344         rout1, rout1, rout, rout, rout, rout, rout, rout, rout, rout, rout, rout, rout, rout1, rout1};
0345 
0346     std::vector<double> pgonZHalf = {0, zstep1, zstep2, zstep3, zstep4, zstep5, zstep5, zmax};
0347 
0348     std::vector<double> pgonRminHalf = {rin, rin, rmid1, rmid1, rmid2, rmid3, rmid3, rmid4};
0349 
0350     std::vector<double> pgonRmaxHalf = {rout, rout, rout, rout, rout, rout, rout1, rout1};
0351 
0352     std::string name("Null");
0353     dd4hep::Solid solid;
0354     if (nf == 0) {
0355       solid = dd4hep::Polyhedra(ns.prepend(idName), nsectortot, -alpha, dphi, pgonZ, pgonRmin, pgonRmax);
0356 #ifdef EDM_ML_DEBUG
0357       edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Polyhedra made of " << genMaterial
0358                                    << " with " << nsectortot << " sectors from " << convertRadToDeg(-alpha) << " to "
0359                                    << convertRadToDeg(-alpha + dphi) << " and with " << nsec << " sections ";
0360       for (unsigned int i = 0; i < pgonZ.size(); i++)
0361         edm::LogVerbatim("HCalGeom") << "\t"
0362                                      << "\tZ = " << cms::convert2mm(pgonZ[i])
0363                                      << "\tRmin = " << cms::convert2mm(pgonRmin[i])
0364                                      << "\tRmax = " << cms::convert2mm(pgonRmax[i]);
0365 #endif
0366     } else {
0367       solid = dd4hep::Polyhedra(ns.prepend(idName), nsectortot, -alpha, dphi, pgonZHalf, pgonRminHalf, pgonRmaxHalf);
0368 #ifdef EDM_ML_DEBUG
0369       edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Polyhedra made of " << genMaterial
0370                                    << " with " << nsectortot << " sectors from " << convertRadToDeg(-alpha) << " to "
0371                                    << convertRadToDeg(-alpha + dphi) << " and with " << nsec << " sections ";
0372       for (unsigned int i = 0; i < pgonZHalf.size(); i++)
0373         edm::LogVerbatim("HCalGeom") << "\t\tZ = " << cms::convert2mm(pgonZHalf[i])
0374                                      << "\tRmin = " << cms::convert2mm(pgonRminHalf[i])
0375                                      << "\tRmax = " << cms::convert2mm(pgonRmaxHalf[i]);
0376 #endif
0377     }
0378 
0379     dd4hep::Material matter = ns.material(genMaterial);
0380     dd4hep::Volume genlogic(solid.name(), solid, matter);
0381     dd4hep::Volume parentName = ns.volume(args.parentName());
0382     parentName.placeVolume(genlogic, 1);
0383 #ifdef EDM_ML_DEBUG
0384     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << genlogic.name() << " number 1 positioned in "
0385                                  << parentName.name() << " at (0, 0, 0) with no rotation";
0386 #endif
0387     //Forward and backwards halfs
0388     name = idName + "Half";
0389     nf = (ntot + 1) / 2;
0390     solid = dd4hep::Polyhedra(ns.prepend(name), nsectortot, -alpha, dphi, pgonZHalf, pgonRminHalf, pgonRmaxHalf);
0391 #ifdef EDM_ML_DEBUG
0392     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Polyhedra made of " << genMaterial
0393                                  << " with " << nsectortot << " sectors from " << convertRadToDeg(-alpha) << " to "
0394                                  << convertRadToDeg(-alpha + dphi) << " and with " << nf << " sections ";
0395     for (unsigned int i = 0; i < pgonZHalf.size(); i++)
0396       edm::LogVerbatim("HCalGeom") << "\t\tZ = " << cms::convert2mm(pgonZHalf[i])
0397                                    << "\tRmin = " << cms::convert2mm(pgonRminHalf[i])
0398                                    << "\tRmax = " << cms::convert2mm(pgonRmaxHalf[i]);
0399 #endif
0400     dd4hep::Volume genlogich(solid.name(), solid, matter);
0401     genlogic.placeVolume(genlogich, 1);
0402 #ifdef EDM_ML_DEBUG
0403     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << genlogich.name() << " number 1 positioned in "
0404                                  << genlogic.name() << " at (0, 0, 0) with no rotation";
0405 #endif
0406     if (nhalf != 1) {
0407       dd4hep::Rotation3D rot = getRotation(rotHalf, rotns, ns);
0408       genlogic.placeVolume(genlogich, 2, rot);
0409 #ifdef EDM_ML_DEBUG
0410       edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << genlogich.name() << " number 2 positioned in "
0411                                    << genlogic.name() << " at (0, 0, 0) with " << rot;
0412 #endif
0413     }  //end if (getNhalf...
0414 
0415     //Construct sector (from -alpha to +alpha)
0416     name = idName + "Module";
0417     solid = dd4hep::Polyhedra(ns.prepend(name), 1, -alpha, 2 * alpha, pgonZHalf, pgonRminHalf, pgonRmaxHalf);
0418     dd4hep::Volume seclogic(solid.name(), solid, matter);
0419 #ifdef EDM_ML_DEBUG
0420     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Polyhedra made of " << matter.name()
0421                                  << " with 1 sector from " << convertRadToDeg(-alpha) << " to "
0422                                  << convertRadToDeg(alpha) << " and with " << nf << " sections";
0423     for (unsigned int i = 0; i < pgonZHalf.size(); i++)
0424       edm::LogVerbatim("HCalGeom") << "\t\tZ = " << cms::convert2mm(pgonZHalf[i])
0425                                    << "\tRmin = " << cms::convert2mm(pgonRminHalf[i])
0426                                    << "\tRmax = " << cms::convert2mm(pgonRmaxHalf[i]);
0427 #endif
0428 
0429     for (int ii = 0; ii < nsectortot; ii++) {
0430       double phi = ii * 2 * alpha;
0431       dd4hep::Rotation3D rotation;
0432       if (phi != 0) {
0433         rotation = dd4hep::RotationZ(phi);
0434 #ifdef EDM_ML_DEBUG
0435         edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: Creating a new "
0436                                      << "rotation around Z " << convertRadToDeg(phi);
0437 #endif
0438       }  //if phideg!=0
0439       genlogich.placeVolume(seclogic, ii + 1, rotation);
0440 #ifdef EDM_ML_DEBUG
0441       edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << seclogic.name() << " number " << ii + 1
0442                                    << " positioned in " << genlogich.name() << " at (0, 0, 0) with " << rotation;
0443 #endif
0444     }
0445 
0446     //Construct the things inside the sector
0447 #ifdef EDM_ML_DEBUG
0448     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: Layers (" << nLayers << ") ...";
0449 #endif
0450     rin = rinner;
0451     for (int i = 0; i < nLayers; i++) {
0452       std::string name = idName + layerLabel[i];
0453       dd4hep::Material matter = ns.material(layerMat[i]);
0454       double width = layerWidth[i];
0455       double rout = rin + width;
0456 
0457       int in = 0, out = 0;
0458       for (int j = 0; j < rzones - 1; j++) {
0459         if (rin >= rmax[j])
0460           in = j + 1;
0461         if (rout > rmax[j])
0462           out = j + 1;
0463       }
0464       double zout = zoff[in] + rin * ttheta[in];
0465 
0466       //!!!!!!!!!!!!!!!!!Should be zero. And removed as soon as
0467       //vertical walls are allowed in SolidPolyhedra
0468 #ifdef EDM_ML_DEBUG
0469       int nsec = 2;
0470 #endif
0471       std::vector<double> pgonZ, pgonRmin, pgonRmax;
0472       // index 0
0473       pgonZ.emplace_back(0);
0474       pgonRmin.emplace_back(rin);
0475       pgonRmax.emplace_back(rout);
0476       // index 1
0477       pgonZ.emplace_back(zout);
0478       pgonRmin.emplace_back(rin);
0479       pgonRmax.emplace_back(rout);
0480       if (in == out) {
0481         if (in <= 3) {
0482           //index 2
0483           pgonZ.emplace_back(zoff[in] + rout * ttheta[in]);
0484           pgonRmin.emplace_back(pgonRmax[1]);
0485           pgonRmax.emplace_back(pgonRmax[1]);
0486 #ifdef EDM_ML_DEBUG
0487           nsec++;
0488 #endif
0489         }
0490       } else {
0491         if (in == 3) {
0492           //redo index 1, add index 2
0493           pgonZ[1] = (zoff[out] + rmax[out] * ttheta[out]);
0494           pgonZ.emplace_back(pgonZ[1]);
0495           pgonRmin.emplace_back(pgonRmin[1]);
0496           pgonRmax.emplace_back(rmax[in]);
0497           //index 3
0498           pgonZ.emplace_back(zoff[in] + rmax[in] * ttheta[in]);
0499           pgonRmin.emplace_back(pgonRmin[2]);
0500           pgonRmax.emplace_back(pgonRmax[2]);
0501 #ifdef EDM_ML_DEBUG
0502           nsec += 2;
0503 #endif
0504         } else {
0505           //index 2
0506           pgonZ.emplace_back(zoff[in] + rmax[in] * ttheta[in]);
0507           pgonRmin.emplace_back(rmax[in]);
0508           pgonRmax.emplace_back(pgonRmax[1]);
0509 #ifdef EDM_ML_DEBUG
0510           nsec++;
0511 #endif
0512           if (in == 0) {
0513             pgonZ.emplace_back(zoff[out] + rmax[in] * ttheta[out]);
0514             pgonRmin.emplace_back(pgonRmin[2]);
0515             pgonRmax.emplace_back(pgonRmax[2]);
0516 #ifdef EDM_ML_DEBUG
0517             nsec++;
0518 #endif
0519           }
0520           if (in <= 1) {
0521             pgonZ.emplace_back(zoff[out] + rout * ttheta[out]);
0522             pgonRmin.emplace_back(rout);
0523             pgonRmax.emplace_back(rout);
0524 #ifdef EDM_ML_DEBUG
0525             nsec++;
0526 #endif
0527           }
0528         }
0529       }
0530       //Solid & volume
0531       dd4hep::Solid solid;
0532       double alpha1 = alpha;
0533       if (layerGap[i] > 1.e-6) {
0534         double rmid = 0.5 * (rin + rout);
0535         double width = rmid * tan(alpha) - layerGap[i];
0536         alpha1 = atan(width / rmid);
0537 #ifdef EDM_ML_DEBUG
0538         edm::LogVerbatim("HCalGeom") << "\tAlpha_1 modified from " << convertRadToDeg(alpha) << " to "
0539                                      << convertRadToDeg(alpha1) << " Rmid " << cms::convert2mm(rmid)
0540                                      << " Reduced width " << cms::convert2mm(width);
0541 #endif
0542       }
0543       solid = dd4hep::Polyhedra(ns.prepend(name), 1, -alpha1, 2 * alpha1, pgonZ, pgonRmin, pgonRmax);
0544       dd4hep::Volume glog(solid.name(), solid, matter);
0545 #ifdef EDM_ML_DEBUG
0546       edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " (Layer " << i << ") Polyhedra made of "
0547                                    << matter.name() << " with 1 sector from " << convertRadToDeg(-alpha1) << " to "
0548                                    << convertRadToDeg(alpha1) << " and with " << nsec << " sections";
0549       for (unsigned int k = 0; k < pgonZ.size(); k++)
0550         edm::LogVerbatim("HCalGeom") << "\t\t" << cms::convert2mm(pgonZ[k]) << "\t" << cms::convert2mm(pgonRmin[k])
0551                                      << "\t" << cms::convert2mm(pgonRmax[k]);
0552 #endif
0553 
0554       seclogic.placeVolume(glog, layerId[i]);
0555 #ifdef EDM_ML_DEBUG
0556       edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << glog.name() << " number " << layerId[i]
0557                                    << " positioned in " << seclogic.name() << " at (0,0,0) with no rotation";
0558 #endif
0559       constructInsideLayers(glog,
0560                             layerLabel[i],
0561                             layerId[i],
0562                             layerAbsorb[i],
0563                             rin,
0564                             layerD1[i],
0565                             alpha1,
0566                             layerD2[i],
0567                             layerAlpha[i],
0568                             layerT1[i],
0569                             layerT2[i],
0570                             ns);
0571       rin = rout;
0572     }
0573   }
0574 
0575   void constructInsideLayers(dd4hep::Volume& laylog,
0576                              const std::string& nm,
0577                              int id,
0578                              int nAbs,
0579                              double rin,
0580                              double d1,
0581                              double alpha1,
0582                              double d2,
0583                              double alpha2,
0584                              double t1,
0585                              double t2,
0586                              cms::DDNamespace& ns) {
0587 #ifdef EDM_ML_DEBUG
0588     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: \t\tInside layer " << id << "...";
0589 #endif
0590 
0591     ///////////////////////////////////////////////////////////////
0592     //Pointers to the Rotation Matrices and to the Materials
0593     dd4hep::Rotation3D rot = getRotation(detRot, rotns, ns);
0594 
0595     std::string nam0 = nm + "In";
0596     std::string name = idName + nam0;
0597     dd4hep::Material matter = ns.material(detMat);
0598 
0599     dd4hep::Solid solid;
0600     dd4hep::Volume glog, mother;
0601     double rsi, dx, dy, dz, x, y;
0602     int i, in;
0603     //Two lower volumes
0604     if (alpha1 > 0) {
0605       rsi = rin + d1;
0606       in = 0;
0607       for (i = 0; i < rzones - 1; i++) {
0608         if (rsi >= rmax[i])
0609           in = i + 1;
0610       }
0611       dx = 0.5 * t1;
0612       dy = 0.5 * rsi * (tan(alpha1) - tan(alpha2));
0613       dz = 0.5 * (zoff[in] + rsi * ttheta[in]);
0614       x = rsi + dx;
0615       y = 0.5 * rsi * (tan(alpha1) + tan(alpha2));
0616       dd4hep::Position r11(x, y, dz);
0617       dd4hep::Position r12(x, -y, dz);
0618 
0619       solid = dd4hep::Box(ns.prepend(name + "1"), dx, dy, dz);
0620 #ifdef EDM_ML_DEBUG
0621       edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Box made of " << matter.name()
0622                                    << " of dimensions " << cms::convert2mm(dx) << ", " << cms::convert2mm(dy) << ", "
0623                                    << cms::convert2mm(dz);
0624 #endif
0625       glog = dd4hep::Volume(solid.name(), solid, matter);
0626 
0627       if (nAbs != 0) {
0628         mother = constructSideLayer(laylog, name, nAbs, rin, alpha1, ns);
0629       } else {
0630         mother = laylog;
0631       }
0632       mother.placeVolume(glog, idOffset + 1, r11);
0633       mother.placeVolume(glog, idOffset + 2, dd4hep::Transform3D(rot, r12));
0634 #ifdef EDM_ML_DEBUG
0635       edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << glog.name() << " Number " << (idOffset + 1)
0636                                    << " positioned in " << mother.name() << " at (" << cms::convert2mm(x) << ","
0637                                    << cms::convert2mm(y) << "," << cms::convert2mm(dz) << ") with no rotation\n"
0638                                    << "DDHCalBarrelAlgo: " << glog.name() << " Number " << (idOffset + 2)
0639                                    << " positioned in " << mother.name() << " at (" << cms::convert2mm(x) << ","
0640                                    << -cms::convert2mm(y) << "," << cms::convert2mm(dz) << ") with " << rot;
0641 #endif
0642       //Constructin the plastics and scintillators inside
0643       constructInsideDetectors(glog, nam0 + "1", id, dx, dy, dz, 1, ns);
0644     }
0645 
0646     //Upper volume
0647     rsi = rin + d2;
0648     in = 0;
0649     for (i = 0; i < rzones - 1; i++) {
0650       if (rsi >= rmax[i])
0651         in = i + 1;
0652     }
0653     dx = 0.5 * t2;
0654     dy = 0.5 * rsi * tan(alpha2);
0655     dz = 0.5 * (zoff[in] + rsi * ttheta[in]);
0656     x = rsi + dx;
0657     dd4hep::Position r21(x, dy, dz);
0658     dd4hep::Position r22(x, -dy, dz);
0659 
0660     solid = dd4hep::Box(ns.prepend(name + "2"), dx, dy, dz);
0661 #ifdef EDM_ML_DEBUG
0662     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Box made of " << matter.name()
0663                                  << " of dimensions " << cms::convert2mm(dx) << ", " << cms::convert2mm(dy) << ", "
0664                                  << cms::convert2mm(dz);
0665 #endif
0666     glog = dd4hep::Volume(solid.name(), solid, matter);
0667 
0668     if (nAbs < 0) {
0669       mother = constructMidLayer(laylog, name, rin, alpha1, ns);
0670     } else {
0671       mother = laylog;
0672     }
0673     mother.placeVolume(glog, idOffset + 3, r21);
0674     mother.placeVolume(glog, idOffset + 4, dd4hep::Transform3D(rot, r22));
0675 #ifdef EDM_ML_DEBUG
0676     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << glog.name() << " Number " << (idOffset + 3)
0677                                  << " positioned in " << mother.name() << " at (" << cms::convert2mm(x) << ","
0678                                  << cms::convert2mm(dy) << "," << cms::convert2mm(dz)
0679                                  << ") with no rotation\nDDHCalBarrelAlgo: " << glog.name() << " Number "
0680                                  << (idOffset + 4) << " positioned in " << mother.name() << " at ("
0681                                  << cms::convert2mm(x) << "," << -cms::convert2mm(dy) << "," << cms::convert2mm(dz)
0682                                  << ") with " << rot;
0683 #endif
0684     //Constructin the plastics and scintillators inside
0685     constructInsideDetectors(glog, nam0 + "2", id, dx, dy, dz, 2, ns);
0686   }
0687 
0688   dd4hep::Volume constructSideLayer(
0689       dd4hep::Volume& laylog, const std::string& nm, int nAbs, double rin, double alpha, cms::DDNamespace& ns) {
0690     //Extra absorber layer
0691     int k = abs(nAbs) - 1;
0692     std::string namek = nm + "Side";
0693     double rsi = rin + sideD[k];
0694     int in = 0;
0695     for (int i = 0; i < rzones - 1; i++) {
0696       if (rsi >= rmax[i])
0697         in = i + 1;
0698     }
0699     std::vector<double> pgonZ, pgonRmin, pgonRmax;
0700     // index 0
0701     pgonZ.emplace_back(0.0);
0702     pgonRmin.emplace_back(rsi);
0703     pgonRmax.emplace_back(rsi + sideT[k]);
0704     // index 1
0705     pgonZ.emplace_back(zoff[in] + rsi * ttheta[in]);
0706     pgonRmin.emplace_back(rsi);
0707     pgonRmax.emplace_back(pgonRmax[0]);
0708     // index 2
0709     pgonZ.emplace_back(zoff[in] + pgonRmax[0] * ttheta[in]);
0710     pgonRmin.emplace_back(pgonRmax[1]);
0711     pgonRmax.emplace_back(pgonRmax[1]);
0712     dd4hep::Solid solid = dd4hep::Polyhedra(ns.prepend(namek), 1, -alpha, 2 * alpha, pgonZ, pgonRmin, pgonRmax);
0713     dd4hep::Material matter = ns.material(sideMat[k]);
0714     dd4hep::Volume glog(solid.name(), solid, matter);
0715 #ifdef EDM_ML_DEBUG
0716     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Polyhedra made of " << sideMat[k]
0717                                  << " with 1 sector from " << convertRadToDeg(-alpha) << " to "
0718                                  << convertRadToDeg(alpha) << " and with " << pgonZ.size() << " sections";
0719     for (unsigned int ii = 0; ii < pgonZ.size(); ii++)
0720       edm::LogVerbatim("HCalGeom") << "\t\tZ = " << cms::convert2mm(pgonZ[ii])
0721                                    << "\tRmin = " << cms::convert2mm(pgonRmin[ii])
0722                                    << "\tRmax = " << cms::convert2mm(pgonRmax[ii]);
0723 #endif
0724 
0725     laylog.placeVolume(glog, 1);
0726 #ifdef EDM_ML_DEBUG
0727     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << glog.name() << " Number 1 positioned in " << laylog.name()
0728                                  << " at (0,0,0) with no rotation";
0729 #endif
0730     if (nAbs < 0) {
0731       dd4hep::Volume mother = glog;
0732       double rmid = pgonRmax[0];
0733       for (int i = 0; i < nSideAbs; i++) {
0734         double alpha1 = atan(sideAbsW[i] / rmid);
0735         if (alpha1 > 0) {
0736           std::string name = namek + sideAbsName[i];
0737           solid = dd4hep::Polyhedra(ns.prepend(name), 1, -alpha1, 2 * alpha1, pgonZ, pgonRmin, pgonRmax);
0738           dd4hep::Material matter = ns.material(sideAbsMat[i]);
0739           dd4hep::Volume log(solid.name(), solid, matter);
0740 #ifdef EDM_ML_DEBUG
0741           edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Polyhedra made of " << sideAbsMat[i]
0742                                        << " with 1 sector from " << convertRadToDeg(-alpha1) << " to "
0743                                        << convertRadToDeg(alpha1) << " and with " << pgonZ.size() << " sections";
0744           for (unsigned int ii = 0; ii < pgonZ.size(); ii++)
0745             edm::LogVerbatim("HCalGeom") << "\t\tZ = " << cms::convert2mm(pgonZ[ii])
0746                                          << "\tRmin = " << cms::convert2mm(pgonRmin[ii])
0747                                          << "\tRmax = " << cms::convert2mm(pgonRmax[ii]);
0748 #endif
0749 
0750           mother.placeVolume(log, 1);
0751 #ifdef EDM_ML_DEBUG
0752           edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << log.name() << " Number 1 positioned in "
0753                                        << mother.name() << " at (0,0,0) with no rotation";
0754 #endif
0755           mother = log;
0756         }
0757       }
0758     }
0759     return glog;
0760   }
0761 
0762   dd4hep::Volume constructMidLayer(
0763       dd4hep::Volume laylog, const std::string& nm, double rin, double alpha, cms::DDNamespace& ns) {
0764     dd4hep::Solid solid;
0765     dd4hep::Volume log, glog;
0766     std::string name = nm + "Mid";
0767     for (int k = 0; k < nAbsorber; k++) {
0768       std::string namek = name + absorbName[k];
0769       double rsi = rin + absorbD[k];
0770       int in = 0;
0771       for (int i = 0; i < rzones - 1; i++) {
0772         if (rsi >= rmax[i])
0773           in = i + 1;
0774       }
0775       std::vector<double> pgonZ, pgonRmin, pgonRmax;
0776       // index 0
0777       pgonZ.emplace_back(0.0);
0778       pgonRmin.emplace_back(rsi);
0779       pgonRmax.emplace_back(rsi + absorbT[k]);
0780       // index 1
0781       pgonZ.emplace_back(zoff[in] + rsi * ttheta[in]);
0782       pgonRmin.emplace_back(rsi);
0783       pgonRmax.emplace_back(pgonRmax[0]);
0784       // index 2
0785       pgonZ.emplace_back(zoff[in] + pgonRmax[0] * ttheta[in]);
0786       pgonRmin.emplace_back(pgonRmax[1]);
0787       pgonRmax.emplace_back(pgonRmax[1]);
0788       solid = dd4hep::Polyhedra(ns.prepend(namek), 1, -alpha, 2 * alpha, pgonZ, pgonRmin, pgonRmax);
0789       dd4hep::Material matter = ns.material(absorbMat[k]);
0790       log = dd4hep::Volume(solid.name(), solid, matter);
0791 #ifdef EDM_ML_DEBUG
0792       edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Polyhedra made of " << matter.name()
0793                                    << " with 1 sector from " << convertRadToDeg(-alpha) << " to "
0794                                    << convertRadToDeg(alpha) << " and with " << pgonZ.size() << " sections";
0795       for (unsigned int ii = 0; ii < pgonZ.size(); ii++)
0796         edm::LogVerbatim("HCalGeom") << "\t\tZ = " << cms::convert2mm(pgonZ[ii])
0797                                      << "\tRmin = " << cms::convert2mm(pgonRmin[ii])
0798                                      << "\tRmax = " << cms::convert2mm(pgonRmax[ii]);
0799 #endif
0800 
0801       laylog.placeVolume(log, 1);
0802 #ifdef EDM_ML_DEBUG
0803       edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << log.name() << " Number 1 positioned in " << laylog.name()
0804                                    << " at (0,0,0) with no rotation";
0805 #endif
0806       if (k == 0) {
0807         double rmin = pgonRmin[0];
0808         double rmax = pgonRmax[0];
0809         dd4hep::Volume mother = log;
0810         for (int i = 0; i < 1; i++) {
0811           double alpha1 = atan(midW[i] / rmin);
0812           std::string namek = name + midName[i];
0813           solid = dd4hep::Polyhedra(ns.prepend(namek), 1, -alpha1, 2 * alpha1, pgonZ, pgonRmin, pgonRmax);
0814           dd4hep::Material matter1 = ns.material(midMat[i]);
0815           log = dd4hep::Volume(solid.name(), solid, matter1);
0816 #ifdef EDM_ML_DEBUG
0817           edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Polyhedra made of "
0818                                        << matter1.name() << " with 1 sector from " << convertRadToDeg(-alpha1) << " to "
0819                                        << convertRadToDeg(alpha1) << " and with " << pgonZ.size() << " sections";
0820           for (unsigned int ii = 0; ii < pgonZ.size(); ii++)
0821             edm::LogVerbatim("HCalGeom") << "\t\tZ = " << cms::convert2mm(pgonZ[ii])
0822                                          << "\tRmin = " << cms::convert2mm(pgonRmin[ii])
0823                                          << "\tRmax = " << cms::convert2mm(pgonRmax[ii]);
0824 #endif
0825 
0826           mother.placeVolume(log, 1);
0827 #ifdef EDM_ML_DEBUG
0828           edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << log.name() << " Number 1 positioned in "
0829                                        << mother.name() << " at (0,0,0) with no rotation";
0830 #endif
0831           mother = log;
0832         }
0833 
0834         // Now the layer with detectors
0835         double rmid = rmin + middleD;
0836         pgonRmin[0] = rmid;
0837         pgonRmax[0] = rmax;
0838         pgonRmin[1] = rmid;
0839         pgonRmax[1] = rmax;
0840         pgonZ[1] = zoff[in] + rmid * ttheta[in];
0841         pgonRmin[2] = rmax;
0842         pgonRmax[2] = rmax;
0843         pgonZ[2] = zoff[in] + rmax * ttheta[in];
0844         double alpha1 = atan(middleW / rmin);
0845         solid = dd4hep::Polyhedra(ns.prepend(name), 1, -alpha1, 2 * alpha1, pgonZ, pgonRmin, pgonRmax);
0846         dd4hep::Material matter1 = ns.material(middleMat);
0847         glog = dd4hep::Volume(solid.name(), solid, matter1);
0848 #ifdef EDM_ML_DEBUG
0849         edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Polyhedra made of " << matter1.name()
0850                                      << " with 1 sector from " << convertRadToDeg(-alpha1) << " to "
0851                                      << convertRadToDeg(alpha1) << " and with " << pgonZ.size() << " sections";
0852         for (unsigned int ii = 0; ii < pgonZ.size(); ii++)
0853           edm::LogVerbatim("HCalGeom") << "\t\tZ = " << cms::convert2mm(pgonZ[ii])
0854                                        << "\tRmin = " << cms::convert2mm(pgonRmin[ii])
0855                                        << "\tRmax = " << cms::convert2mm(pgonRmax[ii]);
0856 #endif
0857 
0858         mother.placeVolume(glog, 1);
0859 #ifdef EDM_ML_DEBUG
0860         edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << glog.name() << " Number 1 positioned in "
0861                                      << mother.name() << " at (0,0,0) with no rotation";
0862 #endif
0863         // Now the remaining absorber layers
0864         for (int i = 1; i < nMidAbs; i++) {
0865           namek = name + midName[i];
0866           rmid = rmin + midT[i];
0867           pgonRmin[0] = rmin;
0868           pgonRmax[0] = rmid;
0869           pgonRmin[1] = rmin;
0870           pgonRmax[1] = rmid;
0871           pgonZ[1] = zoff[in] + rmin * ttheta[in];
0872           pgonRmin[2] = rmid;
0873           pgonRmax[2] = rmid;
0874           pgonZ[2] = zoff[in] + rmid * ttheta[in];
0875           alpha1 = atan(midW[i] / rmin);
0876           solid = dd4hep::Polyhedra(ns.prepend(namek), 1, -alpha1, 2 * alpha1, pgonZ, pgonRmin, pgonRmax);
0877           dd4hep::Material matter2 = ns.material(midMat[i]);
0878           log = dd4hep::Volume(solid.name(), solid, matter2);
0879 #ifdef EDM_ML_DEBUG
0880           edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Polyhedra made of "
0881                                        << matter2.name() << " with 1 sector from " << convertRadToDeg(-alpha1) << " to "
0882                                        << convertRadToDeg(alpha1) << " and with " << pgonZ.size() << " sections";
0883           for (unsigned int ii = 0; ii < pgonZ.size(); ii++)
0884             edm::LogVerbatim("HCalGeom") << "\t\tZ = " << cms::convert2mm(pgonZ[ii])
0885                                          << "\tRmin = " << cms::convert2mm(pgonRmin[ii])
0886                                          << "\tRmax = " << cms::convert2mm(pgonRmax[ii]);
0887 #endif
0888 
0889           mother.placeVolume(log, i);
0890 #ifdef EDM_ML_DEBUG
0891           edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << log.name() << " Number " << i << " positioned in "
0892                                        << mother.name() << " at (0,0,0) with no rotation";
0893 #endif
0894           mother = log;
0895         }
0896       }
0897     }
0898     return glog;
0899   }
0900 
0901   void constructInsideDetectors(dd4hep::Volume& detector,
0902                                 const std::string& name,
0903                                 int id,
0904                                 double dx,
0905                                 double dy,
0906                                 double dz,
0907                                 int type,
0908                                 cms::DDNamespace& ns) {
0909 #ifdef EDM_ML_DEBUG
0910     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: \t\tInside detector " << id << "...";
0911 #endif
0912 
0913     dd4hep::Material plmatter = ns.material(detMatPl);
0914     dd4hep::Material scmatter = ns.material(detMatSc);
0915     std::string plname = DDSplit(detector.name()).first + "Plastic_";
0916     std::string scname = idName + "Scintillator" + name;
0917 
0918     id--;
0919     dd4hep::Solid solid;
0920     dd4hep::Volume glog;
0921     double wid, y = 0;
0922     double dx1, dx2, shiftX;
0923 
0924     if (type == 1) {
0925       wid = 0.5 * detWidth1[id];
0926       dx1 = 0.5 * detT11[id];
0927       dx2 = 0.5 * detT21[id];
0928       shiftX = detdP1[id];
0929       if (detPosY[id] > 0)
0930         y = -dy + wid;
0931     } else {
0932       wid = 0.5 * detWidth2[id];
0933       dx1 = 0.5 * detT12[id];
0934       dx2 = 0.5 * detT22[id];
0935       shiftX = detdP2[id];
0936     }
0937 
0938     solid = dd4hep::Box(ns.prepend(plname + "1"), dx1, wid, dz);
0939     glog = dd4hep::Volume(solid.name(), solid, plmatter);
0940 #ifdef EDM_ML_DEBUG
0941     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Box made of " << plmatter.name()
0942                                  << " of dimensions " << cms::convert2mm(dx1) << ", " << cms::convert2mm(wid) << ", "
0943                                  << cms::convert2mm(dz);
0944 #endif
0945 
0946     double x = shiftX + dx1 - dx;
0947     detector.placeVolume(glog, 1, dd4hep::Position(x, y, 0));
0948 #ifdef EDM_ML_DEBUG
0949     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << glog.name() << " Number 1 positioned in " << detector.name()
0950                                  << " at (" << cms::convert2mm(x) << "," << cms::convert2mm(y)
0951                                  << ",0) with no rotation";
0952 #endif
0953     solid = dd4hep::Box(ns.prepend(scname), 0.5 * detTsc[id], wid, dz);
0954     glog = dd4hep::Volume(solid.name(), solid, scmatter);
0955 #ifdef EDM_ML_DEBUG
0956     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Box made of " << scmatter.name()
0957                                  << " of dimensions " << cms::convert2mm(0.5 * detTsc[id]) << ", "
0958                                  << cms::convert2mm(wid) << ", " << cms::convert2mm(dz);
0959 #endif
0960 
0961     x += dx1 + 0.5 * detTsc[id];
0962     int copyNo = id * 10 + detType[id];
0963     detector.placeVolume(glog, copyNo, dd4hep::Position(x, y, 0));
0964 #ifdef EDM_ML_DEBUG
0965     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << glog.name() << " Number " << copyNo << " positioned in "
0966                                  << detector.name() << " at (" << cms::convert2mm(x) << "," << cms::convert2mm(y)
0967                                  << ",0) with no rotation";
0968 #endif
0969     solid = dd4hep::Box(ns.prepend(plname + "2"), dx2, wid, dz);
0970     glog = dd4hep::Volume(solid.name(), solid, plmatter);
0971 #ifdef EDM_ML_DEBUG
0972     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Box made of " << plmatter.name()
0973                                  << " of dimensions " << cms::convert2mm(dx2) << ", " << cms::convert2mm(wid) << ", "
0974                                  << cms::convert2mm(dz);
0975 #endif
0976 
0977     x += 0.5 * detTsc[id] + dx2;
0978     detector.placeVolume(glog, 1, dd4hep::Position(x, y, 0));
0979 #ifdef EDM_ML_DEBUG
0980     edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << glog.name() << " Number 1 positioned in " << detector.name()
0981                                  << " at (" << cms::convert2mm(x) << "," << cms::convert2mm(y)
0982                                  << ",0) with no rotation";
0983 #endif
0984   }
0985 
0986   dd4hep::Rotation3D getRotation(std::string& rotation, std::string& rotns, cms::DDNamespace& ns) {
0987     std::string rot = (strchr(rotation.c_str(), NAMESPACE_SEP) == nullptr) ? (rotns + ":" + rotation) : rotation;
0988 #ifdef EDM_ML_DEBUG
0989     edm::LogVerbatim("HCalGeom") << "getRotation: " << rotation << ":" << rot << ":" << ns.rotation(rot);
0990 #endif
0991     return ns.rotation(rot);
0992   }
0993 };
0994 
0995 static long algorithm(dd4hep::Detector& /* description */, cms::DDParsingContext& ctxt, xml_h e) {
0996   HcalBarrelAlgo hcalbarrelalgo(ctxt, e);
0997 #ifdef EDM_ML_DEBUG
0998   edm::LogVerbatim("HCalGeom") << "<<== End of DDHCalBarrelAlgo construction";
0999 #endif
1000   return cms::s_executed;
1001 }
1002 
1003 // first argument is the type from the xml file
1004 DECLARE_DDCMS_DETELEMENT(DDCMS_hcal_DDHCalBarrelAlgo, algorithm)