Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-02-22 09:35:01

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 "Acts/Plugins/GeoModel/GeoModelDetectorObjectFactory.hpp"
0010 
0011 #include "Acts/Detector/GeometryIdGenerator.hpp"
0012 #include "Acts/Detector/PortalGenerators.hpp"
0013 #include "Acts/Geometry/CuboidVolumeBounds.hpp"
0014 #include "Acts/Geometry/CutoutCylinderVolumeBounds.hpp"
0015 #include "Acts/Geometry/CylinderVolumeBounds.hpp"
0016 #include "Acts/Geometry/GeometryIdentifier.hpp"
0017 #include "Acts/Geometry/TrapezoidVolumeBounds.hpp"
0018 #include "Acts/Navigation/InternalNavigation.hpp"
0019 #include "Acts/Plugins/GeoModel/GeoModelConverters.hpp"
0020 #include "Acts/Plugins/GeoModel/IGeoShapeConverter.hpp"
0021 
0022 #include <algorithm>
0023 #include <iostream>
0024 #include <typeinfo>
0025 
0026 #include <GeoModelKernel/GeoBox.h>
0027 #include <GeoModelKernel/GeoPcon.h>
0028 #include <GeoModelKernel/GeoShapeShift.h>
0029 #include <GeoModelKernel/GeoShapeSubtraction.h>
0030 #include <GeoModelKernel/GeoShapeUnion.h>
0031 #include <GeoModelKernel/GeoSimplePolygonBrep.h>
0032 #include <GeoModelKernel/GeoTrd.h>
0033 #include <GeoModelKernel/GeoTube.h>
0034 #include <GeoModelKernel/GeoTubs.h>
0035 
0036 namespace {
0037 std::string recType(const GeoShapeShift &gshift);
0038 std::string recType(const GeoShapeUnion &gunion);
0039 std::string recType(const GeoShape &gshape);
0040 
0041 std::string recType(const GeoShapeShift &gshift) {
0042   return "Shift[" + recType(*gshift.getOp()) + "]";
0043 }
0044 std::string recType(const GeoShapeUnion &gunion) {
0045   return "Union[" + recType(*gunion.getOpA()) + ", " +
0046          recType(*gunion.getOpB()) + "]";
0047 }
0048 std::string recType(const GeoShape &gshape) {
0049   if (auto ptr = dynamic_cast<const GeoShapeUnion *>(&gshape); ptr != nullptr) {
0050     return recType(*ptr);
0051   }
0052   if (auto ptr = dynamic_cast<const GeoShapeShift *>(&gshape); ptr != nullptr) {
0053     return recType(*ptr);
0054   }
0055   return gshape.type();
0056 }
0057 
0058 }  // namespace
0059 
0060 namespace Acts {
0061 Acts::GeoModelDetectorObjectFactory::GeoModelDetectorObjectFactory(
0062     const Config &cfg, std::unique_ptr<const Logger> mlogger)
0063     : m_logger(std::move(mlogger)), m_cfg(cfg) {}
0064 
0065 void Acts::GeoModelDetectorObjectFactory::construct(
0066     Cache &cache, const GeometryContext &gctx, const GeoModelTree &geoModelTree,
0067     const Options &options) {
0068   if (geoModelTree.geoReader == nullptr) {
0069     throw std::invalid_argument("GeoModelTree has no GeoModelReader");
0070   }
0071   for (const auto &q : options.queries) {
0072     ACTS_VERBOSE("Constructing detector elements for query " << q);
0073     // load data from database according to querie (Muon)
0074     auto qFPV = geoModelTree.geoReader
0075                     ->getPublishedNodes<std::string, GeoFullPhysVol *>(q);
0076 
0077     // go through each fpv
0078     for (const auto &[name, fpv] : qFPV) {
0079       PVConstLink physVol{fpv};
0080       // if the match lambda returns false skip the rest of the loop
0081       if (!matches(name, physVol)) {
0082         continue;
0083       }
0084       convertFpv(name, fpv, cache, gctx);
0085     }
0086   }
0087 }
0088 
0089 void Acts::GeoModelDetectorObjectFactory::convertSensitive(
0090     const PVConstLink &geoPV, const Acts::Transform3 &transform,
0091     std::vector<GeoModelSensitiveSurface> &sensitives) {
0092   const GeoLogVol *logVol = geoPV->getLogVol();
0093   const GeoShape *shape = logVol->getShape();
0094   int shapeId = shape->typeID();
0095   std::string name = logVol->getName();
0096   std::shared_ptr<const Acts::IGeoShapeConverter> converter =
0097       Acts::geoShapesConverters(shapeId);
0098   if (converter == nullptr) {
0099     throw std::runtime_error("The converter for " + recType(*shape) +
0100                              " is nullptr");
0101   }
0102 
0103   auto converted = converter->toSensitiveSurface(geoPV, transform);
0104   if (converted.ok()) {
0105     sensitives.push_back(converted.value());
0106     const auto &[el, sf] = converted.value();
0107     auto ss = std::get<1>(converted.value());
0108 
0109     ACTS_VERBOSE("(successfully converted: "
0110                  << name << " / " << recType(*shape) << " / "
0111                  << logVol->getMaterial()->getName() << ")");
0112 
0113     if (!el || !sf) {
0114       throw std::runtime_error(
0115           "The Detector Element or the Surface is nullptr");
0116     }
0117     return;
0118   }
0119   ACTS_ERROR(name << " / " << recType(*shape)
0120                   << ") could not be converted by any converter");
0121 }
0122 
0123 std::vector<GeoChildNodeWithTrf>
0124 Acts::GeoModelDetectorObjectFactory::findAllSubVolumes(const PVConstLink &vol) {
0125   std::vector<GeoChildNodeWithTrf> subvolumes = getChildrenWithRef(vol, false);
0126   std::vector<GeoChildNodeWithTrf> sensitives;
0127   for (auto subvolume : subvolumes) {
0128     if (matches(subvolume.nodeName, subvolume.volume)) {
0129       sensitives.push_back(subvolume);
0130     }
0131     std::vector<GeoChildNodeWithTrf> senssubsubvolumes =
0132         findAllSubVolumes(subvolume.volume);
0133     std::transform(std::make_move_iterator(senssubsubvolumes.begin()),
0134                    std::make_move_iterator(senssubsubvolumes.end()),
0135                    std::back_inserter(sensitives),
0136                    [&subvolume](GeoChildNodeWithTrf &&volume) {
0137                      volume.transform = subvolume.transform * volume.transform;
0138                      return volume;
0139                    });
0140   }
0141   return sensitives;
0142 }
0143 
0144 bool Acts::GeoModelDetectorObjectFactory::convertBox(std::string name) {
0145   auto convB = std::ranges::any_of(m_cfg.convertBox, [&](const auto &n) {
0146     return name.find(n) != std::string::npos;
0147   });
0148   return convB;
0149 }
0150 
0151 void Acts::GeoModelDetectorObjectFactory::convertFpv(
0152     const std::string &name, GeoFullPhysVol *fpv, Cache &cache,
0153     const GeometryContext &gctx) {
0154   const auto prevSize = cache.sensitiveSurfaces.size();
0155   PVConstLink physVol{fpv};
0156 
0157   // get children
0158   std::vector<GeoChildNodeWithTrf> subvolumes =
0159       getChildrenWithRef(physVol, false);
0160   // vector containing all subvolumes to be converted to surfaces
0161   std::vector<GeoChildNodeWithTrf> surfaces = findAllSubVolumes(physVol);
0162   std::vector<GeoModelSensitiveSurface> sensitives;
0163 
0164   for (const auto &surface : surfaces) {
0165     const Transform3 &transform =
0166         fpv->getAbsoluteTransform() * surface.transform;
0167     convertSensitive(surface.volume, transform, sensitives);
0168   }
0169   cache.sensitiveSurfaces.insert(cache.sensitiveSurfaces.end(),
0170                                  sensitives.begin(), sensitives.end());
0171   if (convertBox(name)) {
0172     const GeoLogVol *logVol =
0173         physVol->getLogVol();  // get logVol for the shape of the volume
0174     const GeoShape *shape = logVol->getShape();  // get shape
0175     const Acts::Transform3 &fpvtransform = fpv->getAbsoluteTransform(nullptr);
0176 
0177     // convert bounding boxes with surfaces inside
0178     std::shared_ptr<Experimental::DetectorVolume> box =
0179         Acts::GeoModel::convertDetectorVolume(gctx, *shape, name, fpvtransform,
0180                                               sensitives);
0181     cache.boundingBoxes.push_back(box);
0182   }
0183   // If fpv has no subs and should not be converted to volume convert to surface
0184   else if (subvolumes.empty()) {
0185     // convert fpvs to surfaces
0186     const Transform3 &transform = fpv->getAbsoluteTransform();
0187     convertSensitive(fpv, transform, cache.sensitiveSurfaces);
0188   }
0189 
0190   // Set the corresponding database entry name to all sensitive surfaces
0191   for (auto i = prevSize; i < cache.sensitiveSurfaces.size(); ++i) {
0192     auto &[detEl, _] = cache.sensitiveSurfaces[i];
0193     detEl->setDatabaseEntryName(name);
0194   }
0195 }
0196 // function to determine if object fits query
0197 bool Acts::GeoModelDetectorObjectFactory::matches(const std::string &name,
0198                                                   const PVConstLink &physvol) {
0199   if (m_cfg.nameList.empty() && m_cfg.materialList.empty()) {
0200     return true;
0201   }
0202 
0203   auto matchName = std::ranges::any_of(m_cfg.nameList, [&](const auto &n) {
0204     return name.find(n) != std::string::npos;
0205   });
0206 
0207   std::string matStr = physvol->getLogVol()->getMaterial()->getName();
0208 
0209   auto matchMaterial = std::ranges::any_of(
0210       m_cfg.materialList,
0211       [&](const auto &m) { return matStr.find(m) != std::string::npos; });
0212 
0213   bool match = matchMaterial && matchName;
0214   GeoIntrusivePtr<const GeoVFullPhysVol> fullVol =
0215       dynamic_pointer_cast<const GeoVFullPhysVol>(physvol);
0216 
0217   // for the fullphysvol we only check the name
0218   if (m_cfg.nameList.empty()) {
0219     return matchMaterial;
0220   }
0221 
0222   // if no material specified or we're looking at fpv judge by name only
0223   if (m_cfg.materialList.empty() || !(fullVol == nullptr)) {
0224     return matchName;
0225   }
0226   return match;
0227 }
0228 }  // namespace Acts