Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-06 08:10:43

0001 // This file is part of the Acts project.
0002 //
0003 // Copyright (C) 2016-2020 CERN for the benefit of the Acts project
0004 //
0005 // This Source Code Form is subject to the terms of the Mozilla Public
0006 // License, v. 2.0. If a copy of the MPL was not distributed with this
0007 // file, You can obtain one at http://mozilla.org/MPL/2.0/.
0008 
0009 #include "ActsExamples/TGeoDetector/TGeoITkModuleSplitter.hpp"
0010 
0011 #include "Acts/Definitions/Algebra.hpp"
0012 #include "Acts/Plugins/TGeo/TGeoDetectorElement.hpp"
0013 #include "Acts/Surfaces/AnnulusBounds.hpp"
0014 #include "Acts/Surfaces/RectangleBounds.hpp"
0015 #include "Acts/Surfaces/Surface.hpp"
0016 #include "Acts/Surfaces/SurfaceBounds.hpp"
0017 
0018 #include <algorithm>
0019 #include <array>
0020 #include <cstddef>
0021 #include <sstream>
0022 
0023 ActsExamples::TGeoITkModuleSplitter::TGeoITkModuleSplitter(
0024     const ActsExamples::TGeoITkModuleSplitter::Config& cfg,
0025     std::unique_ptr<const Acts::Logger> logger)
0026     : m_cfg(cfg), m_logger(std::move(logger)) {
0027   initSplitCategories();
0028 }
0029 
0030 void ActsExamples::TGeoITkModuleSplitter::initSplitCategories() {
0031   m_splitCategories.reserve(m_cfg.splitPatterns.size());
0032   for (const std::pair<const std::string, std::string>& pattern_split_category :
0033        m_cfg.splitPatterns) {
0034     // mark pattern for disc or barrel module splits:
0035     bool is_disk = false;
0036     if (m_cfg.discMap.find(pattern_split_category.second) !=
0037         m_cfg.discMap.end()) {
0038       is_disk = true;
0039     } else if (m_cfg.barrelMap.find(pattern_split_category.second) ==
0040                m_cfg.barrelMap.end()) {
0041       ACTS_ERROR(
0042           pattern_split_category.second +
0043           " is neither a category name for barrel or disk module splits.");
0044       continue;
0045     }
0046     m_splitCategories.push_back(
0047         std::make_tuple(std::regex(pattern_split_category.first),
0048                         pattern_split_category.second, is_disk));
0049   }
0050 }
0051 
0052 /// If applicable, returns a split detector element
0053 inline std::vector<std::shared_ptr<const Acts::TGeoDetectorElement>>
0054 ActsExamples::TGeoITkModuleSplitter::split(
0055     const Acts::GeometryContext& gctx,
0056     std::shared_ptr<const Acts::TGeoDetectorElement> detElement) const {
0057   // Is the current node covered by this splitter?
0058   const TGeoNode& node = detElement->tgeoNode();
0059   auto sensorName = std::string(node.GetName());
0060 
0061   static const char* category_names[2] = {"barrel", "disc"};
0062   for (const std::tuple<std::regex, std::string, bool>& split_category :
0063        m_splitCategories) {
0064     if (std::regex_match(sensorName, std::get<0>(split_category))) {
0065       ACTS_DEBUG("Splitting " +
0066                  std::string(category_names[std::get<2>(split_category)]) +
0067                  " node " + sensorName + " using split ranges of category " +
0068                  std::get<1>(split_category));
0069       if (!std::get<2>(split_category)) {
0070         return ActsExamples::TGeoITkModuleSplitter::splitBarrelModule(
0071             gctx, detElement, m_cfg.barrelMap.at(std::get<1>(split_category)));
0072       } else {
0073         return ActsExamples::TGeoITkModuleSplitter::splitDiscModule(
0074             gctx, detElement, m_cfg.discMap.at(std::get<1>(split_category)));
0075       }
0076     }
0077   }
0078   ACTS_DEBUG("No matching configuration found. Node " +
0079              std::string(detElement->tgeoNode().GetName()) +
0080              " will not be split.");
0081 
0082   return {detElement};
0083 }
0084 
0085 /// If applicable, returns a split detector element
0086 inline std::vector<std::shared_ptr<const Acts::TGeoDetectorElement>>
0087 ActsExamples::TGeoITkModuleSplitter::splitBarrelModule(
0088     const Acts::GeometryContext& gctx,
0089     const std::shared_ptr<const Acts::TGeoDetectorElement>& detElement,
0090     unsigned int nSegments) const {
0091   // Retrieve the surface
0092   auto identifier = detElement->identifier();
0093   const Acts::Surface& surface = detElement->surface();
0094   const Acts::SurfaceBounds& bounds = surface.bounds();
0095   if (bounds.type() != Acts::SurfaceBounds::eRectangle || nSegments <= 1u) {
0096     ACTS_WARNING("Invalid splitting config for barrel node: " +
0097                  std::string(detElement->tgeoNode().GetName()) +
0098                  "! Node will not be slpit.");
0099     return {detElement};
0100   }
0101 
0102   // Output container for the submodules
0103   std::vector<std::shared_ptr<const Acts::TGeoDetectorElement>> detElements =
0104       {};
0105   detElements.reserve(nSegments);
0106 
0107   // Get the geometric information
0108   double thickness = detElement->thickness();
0109   const Acts::Transform3& transform = surface.transform(gctx);
0110   // Determine the new bounds
0111   const std::vector<double> boundsValues = bounds.values();
0112   double lengthX = (boundsValues[Acts::RectangleBounds::eMaxX] -
0113                     boundsValues[Acts::RectangleBounds::eMinX]) /
0114                    nSegments;
0115   double lengthY = boundsValues[Acts::RectangleBounds::eMaxY] -
0116                    boundsValues[Acts::RectangleBounds::eMinY];
0117   auto rectBounds =
0118       std::make_shared<Acts::RectangleBounds>(0.5 * lengthX, 0.5 * lengthY);
0119   // Translation for every subelement
0120   auto localTranslation = Acts::Vector2(-0.5 * lengthX * (nSegments - 1), 0.);
0121   const auto step = Acts::Vector2(lengthX, 0.);
0122   ACTS_DEBUG("Rectangle bounds for new node (half length): " +
0123              std::to_string(rectBounds->halfLengthX()) + ", " +
0124              std::to_string(rectBounds->halfLengthY()));
0125 
0126   for (std::size_t i = 0; i < nSegments; i++) {
0127     Acts::Vector3 globalTranslation =
0128         surface.localToGlobal(gctx, localTranslation, {}) -
0129         transform.translation();
0130     auto elemTransform =
0131         Acts::Transform3(transform).pretranslate(globalTranslation);
0132     auto element = std::make_shared<const Acts::TGeoDetectorElement>(
0133         identifier, detElement->tgeoNode(), elemTransform, rectBounds,
0134         thickness);
0135     detElements.push_back(std::move(element));
0136 
0137     localTranslation += step;
0138   }
0139   return detElements;
0140 }
0141 
0142 /// If applicable, returns a split detector element
0143 inline std::vector<std::shared_ptr<const Acts::TGeoDetectorElement>>
0144 ActsExamples::TGeoITkModuleSplitter::splitDiscModule(
0145     const Acts::GeometryContext& gctx,
0146     const std::shared_ptr<const Acts::TGeoDetectorElement>& detElement,
0147     const std::vector<ActsExamples::TGeoITkModuleSplitter::SplitRange>&
0148         splitRanges) const {
0149   // Retrieve the surface
0150   auto identifier = detElement->identifier();
0151   const Acts::Surface& surface = detElement->surface();
0152   const Acts::SurfaceBounds& bounds = surface.bounds();
0153 
0154   // Check annulus bounds origin
0155   auto printOrigin = [&](const Acts::Surface& sf) {
0156     Acts::Vector3 discOrigin =
0157         sf.localToGlobal(gctx, Acts::Vector2(0., 0.), Acts::Vector3::Zero());
0158     std::string out =
0159         "Disc surface origin at: " + std::to_string(discOrigin[0]) + ", " +
0160         std::to_string(discOrigin[1]) + ", " + std::to_string(discOrigin[2]);
0161     return out;
0162   };
0163   ACTS_DEBUG(printOrigin(surface));
0164 
0165   if (bounds.type() != Acts::SurfaceBounds::eAnnulus || splitRanges.empty()) {
0166     ACTS_WARNING("Invalid splitting config for disk node: " +
0167                  std::string(detElement->tgeoNode().GetName()) +
0168                  "! Node will not be slpit.");
0169     return {detElement};
0170   }
0171 
0172   auto nSegments = splitRanges.size();
0173 
0174   // Output container for the submodules
0175   std::vector<std::shared_ptr<const Acts::TGeoDetectorElement>> detElements =
0176       {};
0177   detElements.reserve(nSegments);
0178 
0179   // Get the geometric information
0180   double thickness = detElement->thickness();
0181   const Acts::Transform3& transform = surface.transform(gctx);
0182   const std::vector<double> boundsValues = bounds.values();
0183   std::array<double, Acts::AnnulusBounds::eSize> values{};
0184   std::copy_n(boundsValues.begin(), Acts::AnnulusBounds::eSize, values.begin());
0185 
0186   for (std::size_t i = 0; i < nSegments; i++) {
0187     values[Acts::AnnulusBounds::eMinR] = splitRanges[i].first;
0188     values[Acts::AnnulusBounds::eMaxR] = splitRanges[i].second;
0189     auto annulusBounds = std::make_shared<Acts::AnnulusBounds>(values);
0190     ACTS_DEBUG(
0191         "New r bounds for node: " + std::to_string(annulusBounds->rMin()) +
0192         ", " + std::to_string(annulusBounds->rMax()));
0193 
0194     auto element = std::make_shared<const Acts::TGeoDetectorElement>(
0195         identifier, detElement->tgeoNode(), transform, annulusBounds,
0196         thickness);
0197     detElements.push_back(std::move(element));
0198   }
0199   return detElements;
0200 }