Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /*
0002 //\class ME0GeometryBuilder
0003 
0004  Description: ME0 Geometry builder from DD & DD4hep
0005               DD4hep part added to the original old file (DD version) made by M. Maggi (INFN Bari)
0006               Sergio Lo Meo (sergio.lo.meo@cern.ch) following what Ianna Osborne made for DTs (DD4hep migration)
0007               Updated by Sunanda Banerjee (Fermilab) to make it work for DDD/DD4hep
0008             Updated:  7 August 2020
0009 */
0010 #include "Geometry/GEMGeometryBuilder/src/ME0GeometryBuilder.h"
0011 #include "Geometry/GEMGeometry/interface/ME0Geometry.h"
0012 #include "Geometry/GEMGeometry/interface/ME0EtaPartitionSpecs.h"
0013 
0014 #include "DetectorDescription/Core/interface/DDFilter.h"
0015 #include "DetectorDescription/Core/interface/DDFilteredView.h"
0016 #include "DetectorDescription/Core/interface/DDSolid.h"
0017 #include "DetectorDescription/DDCMS/interface/DDFilteredView.h"
0018 #include "DetectorDescription/DDCMS/interface/DDCompactView.h"
0019 #include "DetectorDescription/DDCMS/interface/DDSpecParRegistry.h"
0020 
0021 #include "Geometry/MuonNumbering/interface/MuonGeometryConstants.h"
0022 #include "Geometry/MuonNumbering/interface/MuonGeometryNumbering.h"
0023 #include "Geometry/MuonNumbering/interface/MuonBaseNumber.h"
0024 #include "Geometry/MuonNumbering/interface/ME0NumberingScheme.h"
0025 
0026 #include "DataFormats/GeometrySurface/interface/TrapezoidalPlaneBounds.h"
0027 #include "DataFormats/GeometryVector/interface/Basic3DVector.h"
0028 #include "DataFormats/Math/interface/GeantUnits.h"
0029 
0030 #include <algorithm>
0031 #include <iostream>
0032 #include <string>
0033 
0034 using namespace geant_units::operators;
0035 //#define EDM_ML_DEBUG
0036 
0037 ME0GeometryBuilder::ME0GeometryBuilder() {}
0038 
0039 ME0GeometryBuilder::~ME0GeometryBuilder() {}
0040 
0041 ME0Geometry* ME0GeometryBuilder::build(const DDCompactView* cview, const MuonGeometryConstants& muonConstants) {
0042   std::string attribute = "MuStructure";
0043   std::string value = "MuonEndCapME0";
0044   DDSpecificsMatchesValueFilter filter{DDValue(attribute, value, 0.0)};
0045   DDFilteredView fview(*cview, filter);
0046   return this->buildGeometry(fview, muonConstants);
0047 }
0048 
0049 // for DD4hep
0050 ME0Geometry* ME0GeometryBuilder::build(const cms::DDCompactView* cview, const MuonGeometryConstants& muonConstants) {
0051   std::string attribute = "MuStructure";
0052   std::string value = "MuonEndCapME0";
0053   const cms::DDFilter filter(attribute, value);
0054   cms::DDFilteredView fview(*cview, filter);
0055   return this->buildGeometry(fview, muonConstants);
0056 }
0057 
0058 ME0Geometry* ME0GeometryBuilder::buildGeometry(DDFilteredView& fv, const MuonGeometryConstants& muonConstants) {
0059   ME0Geometry* geometry = new ME0Geometry();
0060   MuonGeometryNumbering mdddnum(muonConstants);
0061   ME0NumberingScheme me0Num(muonConstants);
0062 
0063   LogTrace("ME0Geometry") << "Building the geometry service";
0064   LogTrace("ME0Geometry") << "About to run through the ME0 structure\n"
0065                           << "Top level logical part: " << fv.logicalPart().name().name();
0066 
0067 // ==========================================
0068 // ===  Test to understand the structure  ===
0069 // ==========================================
0070 #ifdef EDM_ML_DEBUG
0071   bool testChambers = fv.firstChild();
0072   LogTrace("ME0Geometry") << "doChamber = fv.firstChild() = " << testChambers;
0073 
0074   while (testChambers) {
0075     // to etapartitions
0076     LogTrace("ME0Geometry") << "to layer " << fv.firstChild();
0077     LogTrace("ME0Geometry") << "to etapt " << fv.firstChild();
0078     int rawId = me0Num.baseNumberToUnitNumber(mdddnum.geoHistoryToBaseNumber(fv.geoHistory()));
0079     ME0DetId detId = ME0DetId(rawId);
0080     ME0DetId detIdCh = detId.chamberId();
0081     // back to chambers
0082     LogTrace("ME0Geometry") << "back to layer " << fv.parent();
0083     LogTrace("ME0Geometry") << "back to chamb " << fv.parent();
0084 
0085     LogTrace("ME0Geometry") << "In DoChambers Loop :: ME0DetId " << detId << " = " << detId.rawId()
0086                             << " (which belongs to ME0Chamber " << detIdCh << " = " << detIdCh.rawId() << ")";
0087     LogTrace("ME0Geometry") << "Second level logical part: " << fv.logicalPart().name().name();
0088     DDBooleanSolid solid2 = (DDBooleanSolid)(fv.logicalPart().solid());
0089     std::vector<double> dpar2 = solid2.parameters();
0090     std::stringstream parameters2;
0091     for (unsigned int i = 0; i < dpar2.size(); ++i) {
0092       parameters2 << " dpar[" << i << "]=" << convertMmToCm(dpar2[i]) << "cm ";
0093     }
0094     LogTrace("ME0Geometry") << "Second level parameters: vector with size = " << dpar2.size() << " and elements "
0095                             << parameters2.str();
0096 
0097     bool doLayers = fv.firstChild();
0098 
0099     LogTrace("ME0Geometry") << "doLayer = fv.firstChild() = " << doLayers;
0100     while (doLayers) {
0101       // to etapartitions
0102       LogTrace("ME0Geometry") << "to etapt " << fv.firstChild();
0103       int rawId = me0Num.baseNumberToUnitNumber(mdddnum.geoHistoryToBaseNumber(fv.geoHistory()));
0104       ME0DetId detId = ME0DetId(rawId);
0105       ME0DetId detIdLa = detId.layerId();
0106       // back to layers
0107       LogTrace("ME0Geometry") << "back to layer " << fv.parent();
0108       LogTrace("ME0Geometry") << "In DoLayers Loop :: ME0DetId " << detId << " = " << detId.rawId()
0109                               << " (which belongs to ME0Layer " << detIdLa << " = " << detIdLa.rawId() << ")";
0110       LogTrace("ME0Geometry") << "Third level logical part: " << fv.logicalPart().name().name();
0111       DDBooleanSolid solid3 = (DDBooleanSolid)(fv.logicalPart().solid());
0112       std::vector<double> dpar3 = solid3.parameters();
0113       std::stringstream parameters3;
0114       for (unsigned int i = 0; i < dpar3.size(); ++i) {
0115         parameters3 << " dpar[" << i << "]=" << convertMmToCm(dpar3[i]) << "cm ";
0116       }
0117       LogTrace("ME0Geometry") << "Third level parameters: vector with size = " << dpar3.size() << " and elements "
0118                               << parameters3.str();
0119       bool doEtaParts = fv.firstChild();
0120 
0121       LogTrace("ME0Geometry") << "doEtaPart = fv.firstChild() = " << doEtaParts;
0122       while (doEtaParts) {
0123         LogTrace("ME0Geometry") << "In DoEtaParts Loop :: ME0DetId " << detId << " = " << detId.rawId();
0124         LogTrace("ME0Geometry") << "Fourth level logical part: " << fv.logicalPart().name().name();
0125         DDBooleanSolid solid4 = (DDBooleanSolid)(fv.logicalPart().solid());
0126         std::vector<double> dpar4 = solid4.parameters();
0127         std::stringstream parameters4;
0128         for (unsigned int i = 0; i < dpar4.size(); ++i) {
0129           parameters4 << " dpar[" << i << "]=" << convertMmToCm(dpar4[i]) << "cm ";
0130         }
0131         LogTrace("ME0Geometry") << "Fourth level parameters: vector with size = " << dpar4.size() << " and elements "
0132                                 << parameters4.str();
0133 
0134         doEtaParts = fv.nextSibling();
0135         LogTrace("ME0Geometry") << "doEtaPart = fv.nextSibling() = " << doEtaParts;
0136       }
0137       fv.parent();
0138       LogTrace("ME0Geometry") << "went back to parent :: name = " << fv.logicalPart().name().name()
0139                               << " will now ask for nextSibling";
0140       doLayers = fv.nextSibling();
0141       LogTrace("ME0Geometry") << "doLayer = fv.nextSibling() = " << doLayers;
0142     }
0143     fv.parent();
0144     LogTrace("ME0Geometry") << "went back to parent :: name = " << fv.logicalPart().name().name()
0145                             << " will now ask for nextSibling";
0146     testChambers = fv.nextSibling();
0147     LogTrace("ME0Geometry") << "doChamber = fv.nextSibling() = " << testChambers;
0148   }
0149   fv.parent();
0150 #endif
0151 
0152   // ==========================================
0153   // === Here the Real ME0 Geometry Builder ===
0154   // ==========================================
0155   bool doChambers = fv.firstChild();
0156 
0157   while (doChambers) {
0158     // to etapartitions and back again to pick up DetId
0159     fv.firstChild();
0160     fv.firstChild();
0161 
0162 #ifdef EDM_ML_DEBUG
0163     edm::LogVerbatim("ME0Geometry") << "MuonGeometry 1 " << fv.geoHistory() << " Levels "
0164                                     << mdddnum.geoHistoryToBaseNumber(fv.geoHistory()).getLevels();
0165 #endif
0166     int rawId = me0Num.baseNumberToUnitNumber(mdddnum.geoHistoryToBaseNumber(fv.geoHistory()));
0167     ME0DetId detId = ME0DetId(rawId);
0168     ME0DetId detIdCh = detId.chamberId();
0169 
0170     fv.parent();
0171     fv.parent();
0172 #ifdef EDM_ML_DEBUG
0173     edm::LogVerbatim("ME0Geometry") << "MuonGeometry 2 " << fv.geoHistory() << " Levels "
0174                                     << mdddnum.geoHistoryToBaseNumber(fv.geoHistory()).getLevels();
0175 #endif
0176     // build chamber
0177     ME0Chamber* me0Chamber = buildChamber(fv, detIdCh);
0178     geometry->add(me0Chamber);
0179 
0180     // loop over layers of the chamber
0181     bool doLayers = fv.firstChild();
0182 
0183     while (doLayers) {
0184       // to etapartitions and back again to pick up DetId
0185       fv.firstChild();
0186       int rawId = me0Num.baseNumberToUnitNumber(mdddnum.geoHistoryToBaseNumber(fv.geoHistory()));
0187       ME0DetId detId = ME0DetId(rawId);
0188       ME0DetId detIdLa = detId.layerId();
0189       fv.parent();
0190 #ifdef EDM_ML_DEBUG
0191       edm::LogVerbatim("ME0Geometry") << "MuonGeometry 3 " << fv.geoHistory() << " Levels "
0192                                       << mdddnum.geoHistoryToBaseNumber(fv.geoHistory()).getLevels();
0193 #endif
0194       // build layer
0195       ME0Layer* me0Layer = buildLayer(fv, detIdLa);
0196       me0Chamber->add(me0Layer);
0197       geometry->add(me0Layer);
0198 
0199       // loop over etapartitions of the layer
0200       bool doEtaParts = fv.firstChild();
0201 
0202       while (doEtaParts) {
0203         // pick up DetId
0204 #ifdef EDM_ML_DEBUG
0205         edm::LogVerbatim("ME0Geometry") << "MuonGeometry 4 " << fv.geoHistory() << " Levels "
0206                                         << mdddnum.geoHistoryToBaseNumber(fv.geoHistory()).getLevels();
0207 #endif
0208         int rawId = me0Num.baseNumberToUnitNumber(mdddnum.geoHistoryToBaseNumber(fv.geoHistory()));
0209         ME0DetId detId = ME0DetId(rawId);
0210 
0211         // build etapartition
0212         ME0EtaPartition* etaPart = buildEtaPartition(fv, detId);
0213         me0Layer->add(etaPart);
0214         geometry->add(etaPart);
0215 
0216         doEtaParts = fv.nextSibling();
0217       }
0218       fv.parent();
0219 
0220       doLayers = fv.nextSibling();
0221     }
0222     fv.parent();
0223 
0224     doChambers = fv.nextSibling();
0225   }
0226   return geometry;
0227 }
0228 
0229 ME0Chamber* ME0GeometryBuilder::buildChamber(DDFilteredView& fv, ME0DetId detId) const {
0230   LogTrace("ME0Geometry") << "buildChamber " << fv.logicalPart().name().name() << " " << detId << std::endl;
0231   DDBooleanSolid solid = (DDBooleanSolid)(fv.logicalPart().solid());
0232 
0233   std::vector<double> dpar = solid.parameters();
0234 
0235   double L = convertMmToCm(dpar[0]);  // length is along local Y
0236   double T = convertMmToCm(dpar[3]);  // thickness is long local Z
0237   double b = convertMmToCm(dpar[4]);  // bottom width is along local X
0238   double B = convertMmToCm(dpar[8]);  // top width is along local X
0239 
0240 #ifdef EDM_ML_DEBUG
0241   LogTrace("ME0Geometry") << " name of logical part = " << fv.logicalPart().name().name() << std::endl;
0242   LogTrace("ME0Geometry") << " dpar is vector with size = " << dpar.size() << std::endl;
0243   for (unsigned int i = 0; i < dpar.size(); ++i) {
0244     LogTrace("ME0Geometry") << " dpar [" << i << "] = " << convertMmToCm(dpar[i]) << " cm " << std::endl;
0245   }
0246   LogTrace("ME0Geometry") << "size  b: " << b << "cm, B: " << B << "cm,  L: " << L << "cm, T: " << T << "cm "
0247                           << std::endl;
0248 #endif
0249 
0250   bool isOdd = false;  // detId.chamber()%2;
0251   ME0BoundPlane surf(boundPlane(fv, new TrapezoidalPlaneBounds(b, B, L, T), isOdd));
0252   ME0Chamber* chamber = new ME0Chamber(detId.chamberId(), surf);
0253   return chamber;
0254 }
0255 
0256 ME0Layer* ME0GeometryBuilder::buildLayer(DDFilteredView& fv, ME0DetId detId) const {
0257   LogTrace("ME0Geometry") << "buildLayer " << fv.logicalPart().name().name() << " " << detId << std::endl;
0258 
0259   DDBooleanSolid solid = (DDBooleanSolid)(fv.logicalPart().solid());
0260 
0261   std::vector<double> dpar = solid.parameters();
0262   double L = convertMmToCm(dpar[0]);  // length is along local Y
0263   double t = convertMmToCm(dpar[3]);  // thickness is long local Z
0264   double b = convertMmToCm(dpar[4]);  // bottom width is along local X
0265   double B = convertMmToCm(dpar[8]);  // top width is along local X
0266 
0267 #ifdef EDM_ML_DEBUG
0268   LogTrace("ME0Geometry") << " name of logical part = " << fv.logicalPart().name().name() << std::endl;
0269   LogTrace("ME0Geometry") << " dpar is vector with size = " << dpar.size() << std::endl;
0270   for (unsigned int i = 0; i < dpar.size(); ++i) {
0271     LogTrace("ME0Geometry") << " dpar [" << i << "] = " << convertMmToCm(dpar[i]) << " cm " << std::endl;
0272   }
0273   LogTrace("ME0Geometry") << "size  b: " << b << "cm, B: " << B << "cm,  L: " << L << "cm, t: " << t << "cm "
0274                           << std::endl;
0275 #endif
0276 
0277   bool isOdd = false;  // detId.chamber()%2;
0278   ME0BoundPlane surf(boundPlane(fv, new TrapezoidalPlaneBounds(b, B, L, t), isOdd));
0279   ME0Layer* layer = new ME0Layer(detId.layerId(), surf);
0280   return layer;
0281 }
0282 
0283 ME0EtaPartition* ME0GeometryBuilder::buildEtaPartition(DDFilteredView& fv, ME0DetId detId) const {
0284   LogTrace("ME0Geometry") << "buildEtaPartition " << fv.logicalPart().name().name() << " " << detId << std::endl;
0285 
0286   // EtaPartition specific parameter (nstrips and npads)
0287   DDValue numbOfStrips("nStrips");
0288   DDValue numbOfPads("nPads");
0289   std::vector<const DDsvalues_type*> specs(fv.specifics());
0290   double nStrips = 0., nPads = 0.;
0291   for (const auto& is : specs) {
0292     if (DDfetch(is, numbOfStrips))
0293       nStrips = numbOfStrips.doubles()[0];
0294     if (DDfetch(is, numbOfPads))
0295       nPads = numbOfPads.doubles()[0];
0296   }
0297 
0298   LogTrace("ME0Geometry") << ((nStrips == 0.) ? ("No nStrips found!!")
0299                                               : ("Number of strips: " + std::to_string(nStrips)));
0300   LogTrace("ME0Geometry") << ((nPads == 0.) ? ("No nPads found!!") : ("Number of pads: " + std::to_string(nPads)));
0301 
0302   // EtaPartition specific parameter (size)
0303   std::vector<double> dpar = fv.logicalPart().solid().parameters();
0304   double b = convertMmToCm(dpar[4]);  // half bottom edge
0305   double B = convertMmToCm(dpar[8]);  // half top edge
0306   double L = convertMmToCm(dpar[0]);  // half apothem
0307   double t = convertMmToCm(dpar[3]);  // half thickness
0308 
0309 #ifdef EDM_ML_DEBUG
0310   LogTrace("ME0Geometry") << " name of logical part = " << fv.logicalPart().name().name() << std::endl;
0311   LogTrace("ME0Geometry") << " dpar is vector with size = " << dpar.size() << std::endl;
0312   for (unsigned int i = 0; i < dpar.size(); ++i) {
0313     LogTrace("ME0Geometry") << " dpar [" << i << "] = " << convertMmToCm(dpar[i]) << " cm " << std::endl;
0314   }
0315   LogTrace("ME0Geometry") << "size  b: " << b << "cm, B: " << B << "cm,  L: " << L << "cm, t: " << t << "cm "
0316                           << std::endl;
0317 #endif
0318 
0319   std::vector<float> pars;
0320   pars.emplace_back(b);
0321   pars.emplace_back(B);
0322   pars.emplace_back(L);
0323   pars.emplace_back(nStrips);
0324   pars.emplace_back(nPads);
0325 
0326   bool isOdd = false;  // detId.chamber()%2;
0327   ME0BoundPlane surf(boundPlane(fv, new TrapezoidalPlaneBounds(b, B, L, t), isOdd));
0328   std::string name = fv.logicalPart().name().name();
0329   ME0EtaPartitionSpecs* e_p_specs = new ME0EtaPartitionSpecs(GeomDetEnumerators::ME0, name, pars);
0330 
0331   ME0EtaPartition* etaPartition = new ME0EtaPartition(detId, surf, e_p_specs);
0332   return etaPartition;
0333 }
0334 
0335 ME0GeometryBuilder::ME0BoundPlane ME0GeometryBuilder::boundPlane(const DDFilteredView& fv,
0336                                                                  Bounds* bounds,
0337                                                                  bool isOddChamber) const {
0338   // extract the position
0339   const DDTranslation& trans(fv.translation());
0340   const Surface::PositionType posResult(
0341       float(convertMmToCm(trans.x())), float(convertMmToCm(trans.y())), float(convertMmToCm(trans.z())));
0342 
0343   const DDRotationMatrix& rotation = fv.rotation();
0344   DD3Vector x, y, z;
0345   rotation.GetComponents(x, y, z);
0346 
0347   Surface::RotationType rotResult(float(x.X()),
0348                                   float(x.Y()),
0349                                   float(x.Z()),
0350                                   float(y.X()),
0351                                   float(y.Y()),
0352                                   float(y.Z()),
0353                                   float(z.X()),
0354                                   float(z.Y()),
0355                                   float(z.Z()));
0356 
0357   //Change of axes for the forward
0358   Basic3DVector<float> newX(1., 0., 0.);
0359   Basic3DVector<float> newY(0., 0., 1.);
0360   Basic3DVector<float> newZ(0., 1., 0.);
0361   newY *= -1;
0362 
0363   rotResult.rotateAxes(newX, newY, newZ);
0364 
0365   return ME0BoundPlane(new BoundPlane(posResult, rotResult, bounds));
0366 }
0367 
0368 // dd4hep
0369 
0370 ME0Geometry* ME0GeometryBuilder::buildGeometry(cms::DDFilteredView& fv, const MuonGeometryConstants& muonConstants) {
0371   ME0Geometry* geometry = new ME0Geometry();
0372   MuonGeometryNumbering mdddnum(muonConstants);
0373   ME0NumberingScheme me0Num(muonConstants);
0374 
0375   static constexpr uint32_t levelChamber = 7;
0376   static constexpr uint32_t levelLayer = 8;
0377   uint32_t theLevelPart = muonConstants.getValue("level");
0378   uint32_t theSectorLevel = muonConstants.getValue("m0_sector") / theLevelPart;
0379   std::vector<ME0Chamber*> chambers;
0380   std::vector<ME0Layer*> layers;
0381 
0382   while (fv.firstChild()) {
0383     const auto& history = fv.history();
0384     MuonBaseNumber num(mdddnum.geoHistoryToBaseNumber(history));
0385     ME0DetId detId(me0Num.baseNumberToUnitNumber(num));
0386 #ifdef EDM_ML_DEBUG
0387     edm::LogVerbatim("ME0Geometry") << fv.name() << " with " << history.tags.size() << " Levels and ID " << detId
0388                                     << " Mask " << std::hex << ME0DetId::chamberIdMask_ << ":" << ME0DetId::layerIdMask_
0389                                     << std::dec << " and " << ME0DetId(((detId.rawId()) & ME0DetId::chamberIdMask_))
0390                                     << ":" << ME0DetId(((detId.rawId()) & ME0DetId::layerIdMask_)) << " Sector Level "
0391                                     << theSectorLevel << ":" << history.tags.size() << ":" << fv.level();
0392     for (unsigned int k = 0; k < history.tags.size(); ++k)
0393       edm::LogVerbatim("ME0Geometry") << "[" << k << "] Tag " << history.tags[k] << " Offset " << history.offsets[k]
0394                                       << " copy " << history.copyNos[k];
0395 #endif
0396 
0397     if (fv.level() == levelChamber) {
0398       // build chamber
0399       ME0Chamber* me0Chamber = buildChamber(fv, detId);
0400       chambers.emplace_back(me0Chamber);
0401     } else if (fv.level() == levelLayer) {
0402       // build layer
0403       ME0Layer* me0Layer = buildLayer(fv, detId);
0404       layers.emplace_back(me0Layer);
0405     } else if (history.tags.size() > theSectorLevel) {
0406       // build first eta partition
0407       ME0EtaPartition* etaPart = buildEtaPartition(fv, detId);
0408       geometry->add(etaPart);
0409     }
0410   }
0411 
0412   auto const& partitions = geometry->etaPartitions();
0413   for (auto& layer : layers) {
0414     uint32_t id0 = ((layer->id().rawId()) & ME0DetId::layerIdMask_);
0415     for (auto& etaPart : partitions) {
0416       if (((etaPart->id().rawId()) & ME0DetId::layerIdMask_) == id0) {
0417         layer->add(etaPart);
0418       }
0419     }
0420     geometry->add(layer);
0421   }
0422   for (auto& chamber : chambers) {
0423     uint32_t id0 = ((chamber->id().rawId()) & ME0DetId::chamberIdMask_);
0424     for (auto& layer : layers) {
0425       if (((layer->id().rawId()) & ME0DetId::chamberIdMask_) == id0) {
0426         chamber->add(layer);
0427       }
0428     }
0429     geometry->add(chamber);
0430   }
0431   return geometry;
0432 }
0433 
0434 ME0Chamber* ME0GeometryBuilder::buildChamber(cms::DDFilteredView& fv, ME0DetId detId) const {
0435   std::vector<double> dpar = fv.parameters();
0436 
0437   double L = k_ScaleFromDD4hep * dpar[3];
0438   double T = k_ScaleFromDD4hep * dpar[2];
0439   double b = k_ScaleFromDD4hep * dpar[0];
0440   double B = k_ScaleFromDD4hep * dpar[1];
0441   bool isOdd = false;
0442   ME0BoundPlane surf(boundPlane(fv, new TrapezoidalPlaneBounds(b, B, L, T), isOdd));
0443   ME0Chamber* chamber = new ME0Chamber(detId.chamberId(), surf);
0444 
0445   return chamber;
0446 }
0447 
0448 ME0Layer* ME0GeometryBuilder::buildLayer(cms::DDFilteredView& fv, ME0DetId detId) const {
0449   std::vector<double> dpar = fv.parameters();
0450 
0451   double L = k_ScaleFromDD4hep * dpar[3];
0452   double t = k_ScaleFromDD4hep * dpar[2];
0453   double b = k_ScaleFromDD4hep * dpar[0];
0454   double B = k_ScaleFromDD4hep * dpar[1];
0455   bool isOdd = false;
0456   ME0BoundPlane surf(boundPlane(fv, new TrapezoidalPlaneBounds(b, B, L, t), isOdd));
0457   ME0Layer* layer = new ME0Layer(detId.layerId(), surf);
0458 
0459   return layer;
0460 }
0461 
0462 ME0EtaPartition* ME0GeometryBuilder::buildEtaPartition(cms::DDFilteredView& fv, ME0DetId detId) const {
0463   //      auto nStrips = fv.get<double>("nStrips"); //it doesn't work
0464   //    auto nPads = fv.get<double>("nPads"); //it doesn't work
0465 
0466   auto nStrips = 384;  // from GEMSpecs
0467   auto nPads = 192;    // from GEMSpecs
0468 
0469   std::vector<double> dpar = fv.parameters();
0470 
0471   double b = k_ScaleFromDD4hep * dpar[0];
0472   double B = k_ScaleFromDD4hep * dpar[1];
0473   double L = k_ScaleFromDD4hep * dpar[3];
0474   double t = k_ScaleFromDD4hep * dpar[2];
0475 
0476   const std::vector<float> pars{float(dpar[0]), float(dpar[1]), float(dpar[3]), float(nStrips), float(nPads)};
0477 
0478   bool isOdd = false;
0479   ME0BoundPlane surf(boundPlane(fv, new TrapezoidalPlaneBounds(b, B, L, t), isOdd));
0480 
0481   std::string_view name = fv.name();
0482 
0483   ME0EtaPartitionSpecs* e_p_specs = new ME0EtaPartitionSpecs(GeomDetEnumerators::ME0, std::string(name), pars);
0484 
0485   ME0EtaPartition* etaPartition = new ME0EtaPartition(detId, surf, e_p_specs);
0486 
0487   return etaPartition;
0488 }
0489 
0490 ME0GeometryBuilder::ME0BoundPlane ME0GeometryBuilder::boundPlane(const cms::DDFilteredView& fv,
0491                                                                  Bounds* bounds,
0492                                                                  bool isOddChamber) const {
0493   // extract the position
0494   const Double_t* trans = fv.trans();
0495   Surface::PositionType posResult(
0496       k_ScaleFromDD4hep * trans[0], k_ScaleFromDD4hep * trans[1], k_ScaleFromDD4hep * trans[2]);
0497 
0498   // now the rotation
0499   DDRotationMatrix rotation;
0500   fv.rot(rotation);
0501   DD3Vector x, y, z;
0502   rotation.GetComponents(x, y, z);
0503   Surface::RotationType rotResult(float(x.X()),
0504                                   float(x.Y()),
0505                                   float(x.Z()),
0506                                   float(y.X()),
0507                                   float(y.Y()),
0508                                   float(y.Z()),
0509                                   float(z.X()),
0510                                   float(z.Y()),
0511                                   float(z.Z()));
0512 
0513   //Change of axes for the forward
0514   Basic3DVector<float> newX(1., 0., 0.);
0515   Basic3DVector<float> newY(0., 0., 1.);
0516   Basic3DVector<float> newZ(0., 1., 0.);
0517   newY *= -1;
0518 
0519   rotResult.rotateAxes(newX, newY, newZ);
0520 
0521   return ME0BoundPlane(new BoundPlane(posResult, rotResult, bounds));
0522 }