Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /include/ActsExamples/ITkModuleSplitting/ITkModuleSplitting.hpp was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

0001 // This file is part of the ACTS project.
0002 //
0003 // Copyright (C) 2016 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 https://mozilla.org/MPL/2.0/.
0008 
0009 #pragma once
0010 
0011 #include "Acts/Definitions/Algebra.hpp"
0012 #include "Acts/Surfaces/AnnulusBounds.hpp"
0013 #include "Acts/Surfaces/RectangleBounds.hpp"
0014 #include "Acts/Surfaces/Surface.hpp"
0015 #include "Acts/Surfaces/SurfaceBounds.hpp"
0016 
0017 #include <algorithm>
0018 #include <array>
0019 #include <cstddef>
0020 
0021 namespace ActsExamples::ITk {
0022 
0023 template <typename detector_element_t, typename element_factory_t>
0024 inline std::vector<std::shared_ptr<const detector_element_t>> splitBarrelModule(
0025     const Acts::GeometryContext& gctx,
0026     const std::shared_ptr<const detector_element_t>& detElement,
0027     unsigned int nSegments, const element_factory_t& factory,
0028     const std::string& name,
0029     const Acts::Logger& logger = Acts::getDummyLogger()) {
0030   // Retrieve the surface
0031   const Acts::Surface& surface = detElement->surface();
0032   const Acts::SurfaceBounds& bounds = surface.bounds();
0033   if (bounds.type() != Acts::SurfaceBounds::eRectangle || nSegments <= 1u) {
0034     ACTS_WARNING("Invalid splitting config for barrel node: "
0035                  << name << "! Node will not be slpit.");
0036     return {detElement};
0037   }
0038 
0039   // Output container for the submodules
0040   std::vector<std::shared_ptr<const detector_element_t>> detElements = {};
0041   detElements.reserve(nSegments);
0042 
0043   // Get the geometric information
0044   const Acts::Transform3& transform = surface.transform(gctx);
0045   // Determine the new bounds
0046   const std::vector<double> boundsValues = bounds.values();
0047 
0048   double lengthX = (boundsValues[Acts::RectangleBounds::eMaxX] -
0049                     boundsValues[Acts::RectangleBounds::eMinX]) /
0050                    nSegments;
0051   double lengthY = boundsValues[Acts::RectangleBounds::eMaxY] -
0052                    boundsValues[Acts::RectangleBounds::eMinY];
0053   auto rectBounds =
0054       std::make_shared<Acts::RectangleBounds>(0.5 * lengthX, 0.5 * lengthY);
0055   // Translation for every subelement
0056   auto localTranslation = Acts::Vector2(-0.5 * lengthX * (nSegments - 1), 0.);
0057   const auto step = Acts::Vector2(lengthX, 0.);
0058   ACTS_DEBUG("Rectangle bounds for new node (half length): " +
0059              std::to_string(rectBounds->halfLengthX()) + ", " +
0060              std::to_string(rectBounds->halfLengthY()));
0061 
0062   for (std::size_t i = 0; i < nSegments; i++) {
0063     Acts::Vector3 globalTranslation =
0064         surface.localToGlobal(gctx, localTranslation, {}) -
0065         transform.translation();
0066     auto elemTransform =
0067         Acts::Transform3(transform).pretranslate(globalTranslation);
0068     detElements.emplace_back(factory(elemTransform, rectBounds));
0069 
0070     localTranslation += step;
0071   }
0072   return detElements;
0073 }
0074 
0075 template <typename detector_element_t, typename element_factory_t>
0076 inline std::vector<std::shared_ptr<detector_element_t>> splitDiscModule(
0077     const Acts::GeometryContext& gctx,
0078     const std::shared_ptr<detector_element_t>& detElement,
0079     const std::vector<std::pair<double, double>>& splitRanges,
0080     const element_factory_t& factory, const std::string& name,
0081     const Acts::Logger& logger = Acts::getDummyLogger()) {
0082   // Retrieve the surface
0083   const Acts::Surface& surface = detElement->surface();
0084   const Acts::SurfaceBounds& bounds = surface.bounds();
0085 
0086   // Check annulus bounds origin
0087   auto printOrigin = [&](const Acts::Surface& sf) {
0088     Acts::Vector3 discOrigin =
0089         sf.localToGlobal(gctx, Acts::Vector2(0., 0.), Acts::Vector3::Zero());
0090     std::string out =
0091         "Disc surface origin at: " + std::to_string(discOrigin[0]) + ", " +
0092         std::to_string(discOrigin[1]) + ", " + std::to_string(discOrigin[2]);
0093     return out;
0094   };
0095   ACTS_DEBUG(printOrigin(surface));
0096 
0097   if (bounds.type() != Acts::SurfaceBounds::eAnnulus || splitRanges.empty()) {
0098     ACTS_WARNING("Invalid splitting config for disk node: "
0099                  << name << "! Node will not be slpit.");
0100     return {detElement};
0101   }
0102 
0103   auto nSegments = splitRanges.size();
0104 
0105   // Output container for the submodules
0106   std::vector<std::shared_ptr<const detector_element_t>> detElements = {};
0107   detElements.reserve(nSegments);
0108 
0109   // Get the geometric information
0110   const Acts::Transform3& transform = surface.transform(gctx);
0111   const std::vector<double> boundsValues = bounds.values();
0112   std::array<double, Acts::AnnulusBounds::eSize> values{};
0113 
0114   std::copy_n(boundsValues.begin(), Acts::AnnulusBounds::eSize, values.begin());
0115 
0116   for (std::size_t i = 0; i < nSegments; i++) {
0117     if (boundsValues[Acts::AnnulusBounds::eMinR] > splitRanges[i].first ||
0118         boundsValues[Acts::AnnulusBounds::eMaxR] < splitRanges[i].second) {
0119       ACTS_WARNING(
0120           "Radius pattern not within the original bounds, node will not be "
0121           "split!");
0122       return {detElement};
0123     }
0124 
0125     values[Acts::AnnulusBounds::eMinR] = splitRanges[i].first;
0126     values[Acts::AnnulusBounds::eMaxR] = splitRanges[i].second;
0127     auto annulusBounds = std::make_shared<Acts::AnnulusBounds>(values);
0128     ACTS_DEBUG(
0129         "New r bounds for node: " + std::to_string(annulusBounds->rMin()) +
0130         ", " + std::to_string(annulusBounds->rMax()));
0131 
0132     auto element = factory(transform, annulusBounds);
0133     detElements.push_back(std::move(element));
0134   }
0135   return detElements;
0136 }
0137 
0138 }  // namespace ActsExamples::ITk