Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-12-06 08:57:28

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 #include "ActsPlugins/GeoModel/GeoModelDetectorObjectFactory.hpp"
0010 
0011 #include "Acts/Geometry/CuboidVolumeBounds.hpp"
0012 #include "Acts/Geometry/CutoutCylinderVolumeBounds.hpp"
0013 #include "Acts/Geometry/CylinderVolumeBounds.hpp"
0014 #include "Acts/Geometry/GeometryIdentifier.hpp"
0015 #include "Acts/Geometry/TrapezoidVolumeBounds.hpp"
0016 #include "ActsPlugins/GeoModel/GeoModelConverters.hpp"
0017 #include "ActsPlugins/GeoModel/IGeoShapeConverter.hpp"
0018 
0019 #include <algorithm>
0020 #include <iostream>
0021 #include <typeinfo>
0022 
0023 #include <GeoModelHelpers/GeoShapeUtils.h>
0024 #include <GeoModelKernel/GeoBox.h>
0025 #include <GeoModelKernel/GeoPcon.h>
0026 #include <GeoModelKernel/GeoShapeShift.h>
0027 #include <GeoModelKernel/GeoShapeSubtraction.h>
0028 #include <GeoModelKernel/GeoShapeUnion.h>
0029 #include <GeoModelKernel/GeoSimplePolygonBrep.h>
0030 #include <GeoModelKernel/GeoTrd.h>
0031 #include <GeoModelKernel/GeoTube.h>
0032 #include <GeoModelKernel/GeoTubs.h>
0033 
0034 using namespace Acts;
0035 
0036 namespace ActsPlugins {
0037 
0038 GeoModelDetectorObjectFactory::GeoModelDetectorObjectFactory(
0039     const Config &cfg, std::unique_ptr<const Logger> mlogger)
0040     : m_logger(std::move(mlogger)), m_cfg(cfg) {}
0041 
0042 void GeoModelDetectorObjectFactory::construct(Cache &cache,
0043                                               const GeometryContext &gctx,
0044                                               const GeoModelTree &geoModelTree,
0045                                               const Options &options) {
0046   for (const std::string &q : options.queries) {
0047     ACTS_VERBOSE("Constructing detector elements for query " << q);
0048     // load data from database according to querie (Muon)
0049     auto qFPV = geoModelTree.publisher->getPublishedVol(q);
0050 
0051     /** Full physical volumes represent  logical detector units.*/
0052     for (const auto &[name, physVol] : qFPV) {
0053       ACTS_DEBUG("Convert volume " << name);
0054       convertFpv(name, physVol, cache, gctx);
0055     }
0056   }
0057 }
0058 
0059 void GeoModelDetectorObjectFactory::convertSensitive(
0060     const PVConstLink &geoPV, const Transform3 &transform,
0061     SurfaceBoundFactory &boundFactory,
0062     std::vector<GeoModelSensitiveSurface> &sensitives) {
0063   const GeoLogVol *logVol = geoPV->getLogVol();
0064   const GeoShape *shape = logVol->getShape();
0065   int shapeId = shape->typeID();
0066   const std::string &name = logVol->getName();
0067   std::shared_ptr<const IGeoShapeConverter> converter =
0068       geoShapesConverters(shapeId);
0069   if (converter == nullptr) {
0070     throw std::runtime_error("The converter for " + printGeoShape(shape) +
0071                              " is a nullptr");
0072   }
0073 
0074   auto converted =
0075       converter->toSensitiveSurface(geoPV, transform, boundFactory);
0076   if (converted.ok()) {
0077     const auto &[el, sf] = converted.value();
0078 
0079     if (!el || !sf) {
0080       throw std::runtime_error(
0081           "The Detector Element or the Surface is a nullptr");
0082     }
0083     sensitives.push_back(converted.value());
0084     ACTS_VERBOSE("(successfully converted: "
0085                  << name << " / " << printGeoShape(shape) << " / "
0086                  << logVol->getMaterial()->getName() << ")");
0087     return;
0088   }
0089   ACTS_ERROR(name << " / " << printGeoShape(shape)
0090                   << " could not be converted by any converter");
0091 }
0092 
0093 std::vector<GeoChildNodeWithTrf>
0094 GeoModelDetectorObjectFactory::findAllSubVolumes(const PVConstLink &vol) const {
0095   /// Fetch the direct children of the volume
0096   std::vector<GeoChildNodeWithTrf> subvolumes = getChildrenWithRef(vol, false);
0097   std::vector<GeoChildNodeWithTrf> sensitives;
0098   sensitives.reserve(subvolumes.size());
0099   for (auto &subvolume : subvolumes) {
0100     /// Check whether material or GeoNameTag satisfy the user defined patterns
0101     if (matches(subvolume.nodeName, subvolume.volume)) {
0102       sensitives.push_back(subvolume);
0103     }
0104     /// If the volume has no children nothing can be done further
0105     if (subvolume.volume->getNChildVols() == 0) {
0106       continue;
0107     }
0108     /// Enter the next recursion level to check whether there're sensitive
0109     /// children
0110     std::vector<GeoChildNodeWithTrf> senssubsubvolumes =
0111         findAllSubVolumes(subvolume.volume);
0112     /* Append the found volumes to the output, but  update the transforms
0113      * of the nodes before. They're expressed with respect to their parent,
0114      * which is the grand child of the volume passed to the function call
0115      * -> apply on each grand child also the transform of the child. */
0116     std::transform(std::make_move_iterator(senssubsubvolumes.begin()),
0117                    std::make_move_iterator(senssubsubvolumes.end()),
0118                    std::back_inserter(sensitives),
0119                    [&subvolume](GeoChildNodeWithTrf &&volume) {
0120                      volume.transform = subvolume.transform * volume.transform;
0121                      return volume;
0122                    });
0123   }
0124   return sensitives;
0125 }
0126 
0127 bool GeoModelDetectorObjectFactory::convertBox(const std::string &name) const {
0128   auto convB = std::ranges::any_of(m_cfg.convertBox, [&](const auto &n) {
0129     return name.find(n) != std::string::npos;
0130   });
0131   return convB;
0132 }
0133 
0134 void GeoModelDetectorObjectFactory::convertFpv(
0135     const std::string &name, const FpvConstLink &fpv, Cache &cache,
0136     const GeometryContext & /*gctx*/) {
0137   const std::size_t prevSize = cache.sensitiveSurfaces.size();
0138   {
0139     /** Search all subvolumes that may be converted to sensitive surfaces */
0140     std::vector<GeoChildNodeWithTrf> subVolToTrf = findAllSubVolumes(fpv);
0141 
0142     std::vector<GeoModelSensitiveSurface> sensitives;
0143     sensitives.reserve(subVolToTrf.size());
0144 
0145     for (const auto &convertMe : subVolToTrf) {
0146       /** Align the surface with the global position of the detector */
0147       const Transform3 transform =
0148           fpv->getAbsoluteTransform() * convertMe.transform;
0149       convertSensitive(convertMe.volume, transform, *cache.surfBoundFactory,
0150                        sensitives);
0151     }
0152 
0153     if (sensitives.empty() && matches(name, fpv)) {
0154       convertSensitive(fpv, fpv->getAbsoluteTransform(),
0155                        *cache.surfBoundFactory, cache.sensitiveSurfaces);
0156     }
0157     cache.sensitiveSurfaces.insert(cache.sensitiveSurfaces.end(),
0158                                    std::make_move_iterator(sensitives.begin()),
0159                                    std::make_move_iterator(sensitives.end()));
0160     // Set the corresponding database entry name to all sensitive surfaces
0161     for (auto i = prevSize; i < cache.sensitiveSurfaces.size(); ++i) {
0162       const auto &detEl = std::get<0>(cache.sensitiveSurfaces[i]);
0163       detEl->setDatabaseEntryName(name);
0164       ACTS_VERBOSE("Set database name of the DetectorElement to "
0165                    << detEl->databaseEntryName());
0166     }
0167   }
0168   // Extract the bounding box surrounding the surface
0169   if (convertBox(name)) {
0170     ConvertedGeoVol &convEnvelope = cache.volumeBoxFPVs.emplace_back();
0171     convEnvelope.name = name;
0172     convEnvelope.fullPhysVol = fpv;
0173     convEnvelope.volume = GeoModel::convertVolume(fpv->getAbsoluteTransform(),
0174                                                   fpv->getLogVol()->getShape(),
0175                                                   *cache.volumeBoundFactory);
0176     std::transform(cache.sensitiveSurfaces.begin() + prevSize,
0177                    cache.sensitiveSurfaces.end(),
0178                    std::back_inserter(convEnvelope.surfaces),
0179                    [](const GeoModelSensitiveSurface &sensitive) {
0180                      return std::get<1>(sensitive);
0181                    });
0182   }
0183 }
0184 // function to determine if object fits query
0185 bool GeoModelDetectorObjectFactory::matches(const std::string &name,
0186                                             const PVConstLink &physvol) const {
0187   if (m_cfg.nameList.empty() && m_cfg.materialList.empty()) {
0188     return true;
0189   }
0190 
0191   auto matchName = std::ranges::any_of(m_cfg.nameList, [&](const auto &n) {
0192     return name.find(n) != std::string::npos;
0193   });
0194 
0195   std::string matStr = physvol->getLogVol()->getMaterial()->getName();
0196 
0197   auto matchMaterial = std::ranges::any_of(
0198       m_cfg.materialList,
0199       [&](const auto &m) { return matStr.find(m) != std::string::npos; });
0200 
0201   bool match = matchMaterial && matchName;
0202 
0203   // for the fullphysvol we only check the name
0204   if (m_cfg.nameList.empty()) {
0205     return matchMaterial;
0206   }
0207 
0208   // if no material specified or we're looking at fpv judge by name only
0209   if (m_cfg.materialList.empty() ||
0210       dynamic_pointer_cast<const GeoVFullPhysVol>(physvol)) {
0211     return matchName;
0212   }
0213   return match;
0214 }
0215 }  // namespace ActsPlugins