Back to home page

sPhenix code displayed by LXR

 
 

    


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

0001 // This file is part of the Acts project.
0002 //
0003 // Copyright (C) 2022-2024 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 "Acts/Detector/Detector.hpp"
0010 #include "Acts/Detector/DetectorVolume.hpp"
0011 #include "Acts/Detector/Portal.hpp"
0012 #include "Acts/Geometry/GeometryContext.hpp"
0013 #include "Acts/Geometry/TrackingGeometry.hpp"
0014 #include "Acts/Plugins/ActSVG/DetectorVolumeSvgConverter.hpp"
0015 #include "Acts/Plugins/ActSVG/IndexedSurfacesSvgConverter.hpp"
0016 #include "Acts/Plugins/ActSVG/LayerSvgConverter.hpp"
0017 #include "Acts/Plugins/ActSVG/PortalSvgConverter.hpp"
0018 #include "Acts/Plugins/ActSVG/SurfaceSvgConverter.hpp"
0019 #include "Acts/Plugins/ActSVG/SvgUtils.hpp"
0020 #include "Acts/Plugins/ActSVG/TrackingGeometrySvgConverter.hpp"
0021 #include "Acts/Plugins/Python/Utilities.hpp"
0022 #include "Acts/Utilities/Enumerate.hpp"
0023 #include "ActsExamples/EventData/GeometryContainers.hpp"
0024 #include "ActsExamples/EventData/SimHit.hpp"
0025 #include "ActsExamples/EventData/SimSpacePoint.hpp"
0026 #include "ActsExamples/Io/Svg/SvgPointWriter.hpp"
0027 #include "ActsExamples/Io/Svg/SvgTrackingGeometryWriter.hpp"
0028 
0029 #include <memory>
0030 #include <string>
0031 #include <tuple>
0032 #include <vector>
0033 
0034 #include <pybind11/pybind11.h>
0035 #include <pybind11/stl.h>
0036 
0037 namespace py = pybind11;
0038 using namespace pybind11::literals;
0039 
0040 using namespace Acts;
0041 using namespace ActsExamples;
0042 
0043 namespace {
0044 
0045 // A cache object
0046 using PortalCache = std::list<std::string>;
0047 
0048 /// @brief helper tuple to define which views and which range to be used
0049 ///
0050 /// @param view is the view type, e.g. 'xy', 'zr', ...
0051 /// @param selection is the selection of the volume, e.g. 'all', 'sensitives', 'portals', 'materials'
0052 using ViewAndRange =
0053     std::tuple<std::string, std::vector<std::string>, Acts::Extent>;
0054 
0055 /// Helper function to be picked in different access patterns
0056 ///
0057 /// @param pVolume is the proto volume to be drawn
0058 /// @param identification is the identification of the volume
0059 /// @param viewAndRange is the view, selection and range to be drawn
0060 /// @param portalCache is a portal cache to avoid multiple drawings of the same portal
0061 ///
0062 /// Returns an svg object in the right view
0063 actsvg::svg::object viewDetectorVolume(const Svg::ProtoVolume& pVolume,
0064                                        const std::string& identification,
0065                                        const ViewAndRange& viewAndRange,
0066                                        PortalCache& portalCache) {
0067   actsvg::svg::object svgDet;
0068   svgDet._id = identification;
0069   svgDet._tag = "g";
0070 
0071   auto [view, selection, viewRange] = viewAndRange;
0072 
0073   // Translate selection into booleans
0074   bool all =
0075       std::find(selection.begin(), selection.end(), "all") != selection.end();
0076   bool sensitives = std::find(selection.begin(), selection.end(),
0077                               "sensitives") != selection.end();
0078   bool portals = std::find(selection.begin(), selection.end(), "portals") !=
0079                  selection.end();
0080   bool materials = std::find(selection.begin(), selection.end(), "materials") !=
0081                    selection.end();
0082 
0083   // Helper lambda for material selection
0084   auto materialSel = [&](const Svg::ProtoSurface& s) -> bool {
0085     return (materials &&
0086             s._decorations.find("material") != s._decorations.end());
0087   };
0088 
0089   // Helper lambda for view range selection
0090   auto viewRangeSel = [](const Svg::ProtoSurface& s,
0091                          const Extent& vRange) -> bool {
0092     bool accept = false;
0093     for (const auto& v : s._vertices) {
0094       if (vRange.contains(v)) {
0095         accept = true;
0096         break;
0097       }
0098     }
0099     return accept;
0100   };
0101 
0102   // -------------------- surface section
0103   // The surfaces to be drawn
0104   std::vector<Svg::ProtoVolume::surface_type> sSurfaces;
0105   sSurfaces.reserve(pVolume._v_surfaces.size());
0106   for (const auto& s : pVolume._v_surfaces) {
0107     if ((all || sensitives || materialSel(s)) && viewRangeSel(s, viewRange)) {
0108       sSurfaces.push_back(s);
0109     }
0110   }
0111 
0112   // Now draw all the surfaces
0113   for (const auto& vs : sSurfaces) {
0114     if (view == "xy") {
0115       svgDet.add_object(Svg::View::xy(vs, identification));
0116     } else if (view == "zr") {
0117       svgDet.add_object(Svg::View::zr(vs, identification));
0118     } else {
0119       throw std::invalid_argument("Unknown view type");
0120     }
0121   }
0122 
0123   // -------------------- portal section
0124   std::vector<Svg::ProtoPortal> sPortals;
0125   sPortals.reserve(pVolume._portals.size());
0126   for (const auto& vp : pVolume._portals) {
0127     if ((all || portals || materialSel(vp._surface)) &&
0128         viewRangeSel(vp._surface, viewRange)) {
0129       sPortals.push_back(vp);
0130     }
0131   }
0132 
0133   // Now draw all the portals - if not already in the cache
0134   for (const auto& vp : sPortals) {
0135     auto pgID = vp._surface._decorations.find("geo_id");
0136     std::string gpIDs = "";
0137     if (pgID != vp._surface._decorations.end()) {
0138       gpIDs = pgID->second._id;
0139     }
0140 
0141     if (std::find(portalCache.begin(), portalCache.end(), gpIDs) !=
0142         portalCache.end()) {
0143       continue;
0144     }
0145 
0146     // Register this portal to the cache
0147     portalCache.insert(portalCache.begin(), gpIDs);
0148 
0149     if (view == "xy") {
0150       svgDet.add_object(Svg::View::xy(vp, identification));
0151     } else if (view == "zr") {
0152       svgDet.add_object(Svg::View::zr(vp, identification));
0153     } else {
0154       throw std::invalid_argument("Unknown view type");
0155     }
0156   }
0157   return svgDet;
0158 }
0159 
0160 // Helper function to be picked in different access patterns
0161 void viewDetector(
0162     const Acts::GeometryContext& gctx,
0163     const Acts::Experimental::Detector& detector,
0164     const std::string& identification,
0165     const std::vector<std::tuple<int, Svg::DetectorVolumeConverter::Options>>&
0166         volumeIdxOpts,
0167     const std::vector<ViewAndRange>& viewAndRanges, const std::string& saveAs) {
0168   PortalCache portalCache;
0169 
0170   // The svg object to be returned
0171   std::vector<actsvg::svg::object> svgDetViews;
0172   svgDetViews.reserve(viewAndRanges.size());
0173   for (unsigned int i = 0; i < viewAndRanges.size(); ++i) {
0174     actsvg::svg::object svgDet;
0175     svgDet._id = identification;
0176     svgDet._tag = "g";
0177     svgDetViews.push_back(svgDet);
0178   }
0179 
0180   for (const auto& [vidx, vopts] : volumeIdxOpts) {
0181     // Get the volume and convert it
0182     const auto& v = detector.volumes()[vidx];
0183     auto [pVolume, pGrid] =
0184         Svg::DetectorVolumeConverter::convert(gctx, *v, vopts);
0185 
0186     for (auto [iv, var] : Acts::enumerate(viewAndRanges)) {
0187       auto [view, selection, range] = var;
0188       // Get the view and the range
0189       auto svgVolView = viewDetectorVolume(
0190           pVolume, identification + "_vol" + std::to_string(vidx) + "_" + view,
0191           var, portalCache);
0192       svgDetViews[iv].add_object(svgVolView);
0193     }
0194   }
0195 
0196   for (auto [iv, var] : Acts::enumerate(viewAndRanges)) {
0197     auto [view, selection, range] = var;
0198     Svg::toFile({svgDetViews[iv]}, saveAs + "_" + view + ".svg");
0199   }
0200 }
0201 
0202 }  // namespace
0203 
0204 namespace Acts::Python {
0205 void addSvg(Context& ctx) {
0206   auto [m, mex] = ctx.get("main", "examples");
0207 
0208   auto svg = m.def_submodule("svg");
0209 
0210   // Some basics
0211   py::class_<actsvg::svg::object>(svg, "object");
0212 
0213   // Core components, added as an acts.svg submodule
0214   {
0215     auto c = py::class_<Svg::Style>(svg, "Style").def(py::init<>());
0216     ACTS_PYTHON_STRUCT_BEGIN(c, Svg::Style);
0217     ACTS_PYTHON_MEMBER(fillColor);
0218     ACTS_PYTHON_MEMBER(fillOpacity);
0219     ACTS_PYTHON_MEMBER(highlightColor);
0220     ACTS_PYTHON_MEMBER(highlights);
0221     ACTS_PYTHON_MEMBER(strokeWidth);
0222     ACTS_PYTHON_MEMBER(strokeColor);
0223     ACTS_PYTHON_MEMBER(nSegments);
0224     ACTS_PYTHON_STRUCT_END();
0225   }
0226 
0227   // How surfaces should be drawn: Svg Surface options & drawning
0228   {
0229     auto c = py::class_<Svg::SurfaceConverter::Options>(svg, "SurfaceOptions")
0230                  .def(py::init<>());
0231     ACTS_PYTHON_STRUCT_BEGIN(c, Svg::SurfaceConverter::Options);
0232     ACTS_PYTHON_MEMBER(style);
0233     ACTS_PYTHON_MEMBER(templateSurface);
0234     ACTS_PYTHON_STRUCT_END();
0235 
0236     // Define the proto surface
0237     py::class_<Svg::ProtoSurface>(svg, "ProtoSurface");
0238     // Convert an Acts::Surface object into an acts::svg::proto::surface
0239     svg.def("convertSurface", &Svg::SurfaceConverter::convert);
0240 
0241     // Define the view functions
0242     svg.def("viewSurface", [](const Svg::ProtoSurface& pSurface,
0243                               const std::string& identification,
0244                               const std::string& view = "xy") {
0245       if (view == "xy") {
0246         return Svg::View::xy(pSurface, identification);
0247       } else if (view == "zr") {
0248         return Svg::View::zr(pSurface, identification);
0249       } else if (view == "zphi") {
0250         return Svg::View::zphi(pSurface, identification);
0251       } else if (view == "zrphi") {
0252         return Svg::View::zrphi(pSurface, identification);
0253       } else {
0254         throw std::invalid_argument("Unknown view type");
0255       }
0256     });
0257   }
0258 
0259   // How portals should be drawn: Svg Portal options & drawning
0260   {
0261     auto c = py::class_<Svg::PortalConverter::Options>(svg, "PortalOptions")
0262                  .def(py::init<>());
0263 
0264     ACTS_PYTHON_STRUCT_BEGIN(c, Svg::PortalConverter::Options);
0265     ACTS_PYTHON_MEMBER(surfaceOptions);
0266     ACTS_PYTHON_MEMBER(linkLength);
0267     ACTS_PYTHON_MEMBER(volumeIndices);
0268     ACTS_PYTHON_STRUCT_END();
0269 
0270     // Define the proto portal
0271     py::class_<Svg::ProtoPortal>(svg, "ProtoPortal");
0272     // Convert an Acts::Experimental::Portal object into an
0273     // acts::svg::proto::portal
0274     svg.def("convertPortal", &Svg::PortalConverter::convert);
0275 
0276     // Define the view functions
0277     svg.def("viewPortal", [](const Svg::ProtoPortal& pPortal,
0278                              const std::string& identification,
0279                              const std::string& view = "xy") {
0280       if (view == "xy") {
0281         return Svg::View::xy(pPortal, identification);
0282       } else if (view == "zr") {
0283         return Svg::View::zr(pPortal, identification);
0284       } else {
0285         throw std::invalid_argument("Unknown view type");
0286       }
0287     });
0288   }
0289 
0290   // How detector volumes are drawn: Svg DetectorVolume options & drawning
0291   {
0292     auto c = py::class_<Svg::DetectorVolumeConverter::Options>(
0293                  svg, "DetectorVolumeOptions")
0294                  .def(py::init<>());
0295 
0296     ACTS_PYTHON_STRUCT_BEGIN(c, Svg::DetectorVolumeConverter::Options);
0297     ACTS_PYTHON_MEMBER(portalIndices);
0298     ACTS_PYTHON_MEMBER(portalOptions);
0299     ACTS_PYTHON_MEMBER(surfaceOptions);
0300     ACTS_PYTHON_STRUCT_END();
0301 
0302     // Define the proto volume & indexed surface grid
0303     py::class_<Svg::ProtoVolume>(svg, "ProtoVolume");
0304     py::class_<Svg::ProtoIndexedSurfaceGrid>(svg, "ProtoIndexedSurfaceGrid");
0305 
0306     // Convert an Acts::Experimental::DetectorVolume object into an
0307     // acts::svg::proto::volume
0308     svg.def("convertDetectorVolume", &Svg::DetectorVolumeConverter::convert);
0309 
0310     // Define the view functions
0311     svg.def("viewDetectorVolume", &viewDetectorVolume);
0312   }
0313 
0314   // How a detector is drawn: Svg Detector options & drawning
0315   { svg.def("viewDetector", &viewDetector); }
0316 
0317   // Legacy geometry drawing
0318   {
0319     using DefinedStyle = std::pair<GeometryIdentifier, Svg::Style>;
0320     using DefinedStyleSet = std::vector<DefinedStyle>;
0321 
0322     auto sm = py::class_<GeometryHierarchyMap<Svg::Style>>(svg, "StyleMap")
0323                   .def(py::init<DefinedStyleSet>(), py::arg("elements"));
0324 
0325     auto c = py::class_<Svg::LayerConverter::Options>(svg, "LayerOptions")
0326                  .def(py::init<>());
0327     ACTS_PYTHON_STRUCT_BEGIN(c, Svg::LayerConverter::Options);
0328     ACTS_PYTHON_MEMBER(name);
0329     ACTS_PYTHON_MEMBER(surfaceStyles);
0330     ACTS_PYTHON_MEMBER(zRange);
0331     ACTS_PYTHON_MEMBER(phiRange);
0332     ACTS_PYTHON_MEMBER(gridInfo);
0333     ACTS_PYTHON_MEMBER(moduleInfo);
0334     ACTS_PYTHON_MEMBER(projectionInfo);
0335     ACTS_PYTHON_MEMBER(labelProjection);
0336     ACTS_PYTHON_MEMBER(labelGauge);
0337     ACTS_PYTHON_STRUCT_END();
0338   }
0339 
0340   {
0341     using DefinedLayerOptions =
0342         std::pair<GeometryIdentifier, Svg::LayerConverter::Options>;
0343     using DefinedLayerOptionsSet = std::vector<DefinedLayerOptions>;
0344 
0345     auto lom =
0346         py::class_<GeometryHierarchyMap<Svg::LayerConverter::Options>>(
0347             svg, "LayerOptionMap")
0348             .def(py::init<DefinedLayerOptionsSet>(), py::arg("elements"));
0349 
0350     auto c = py::class_<Svg::TrackingGeometryConverter::Options>(
0351                  svg, "TrackingGeometryOptions")
0352                  .def(py::init<>());
0353     ACTS_PYTHON_STRUCT_BEGIN(c, Svg::TrackingGeometryConverter::Options);
0354     ACTS_PYTHON_MEMBER(prefix);
0355     ACTS_PYTHON_MEMBER(layerOptions);
0356     ACTS_PYTHON_STRUCT_END();
0357   }
0358 
0359   // Components from the ActsExamples - part of acts.examples
0360 
0361   {
0362     using Writer = ActsExamples::SvgTrackingGeometryWriter;
0363     auto w = py::class_<Writer, std::shared_ptr<Writer>>(
0364                  mex, "SvgTrackingGeometryWriter")
0365                  .def(py::init<const Writer::Config&, Acts::Logging::Level>(),
0366                       py::arg("config"), py::arg("level"))
0367                  .def("write", py::overload_cast<const AlgorithmContext&,
0368                                                  const Acts::TrackingGeometry&>(
0369                                    &Writer::write));
0370 
0371     auto c = py::class_<Writer::Config>(w, "Config").def(py::init<>());
0372     ACTS_PYTHON_STRUCT_BEGIN(c, Writer::Config);
0373     ACTS_PYTHON_MEMBER(outputDir);
0374     ACTS_PYTHON_MEMBER(converterOptions);
0375     ACTS_PYTHON_STRUCT_END();
0376   }
0377   {
0378     using Writer = ActsExamples::SvgPointWriter<ActsExamples::SimSpacePoint>;
0379     auto w =
0380         py::class_<Writer, ActsExamples::IWriter, std::shared_ptr<Writer>>(
0381             mex, "SvgSimSpacePointWriter")
0382             .def(py::init<const Writer::Config&, Acts::Logging::Level>(),
0383                  py::arg("config"), py::arg("level"))
0384             .def("write",
0385                  py::overload_cast<const AlgorithmContext&>(&Writer::write));
0386 
0387     auto c = py::class_<Writer::Config>(w, "Config").def(py::init<>());
0388     ACTS_PYTHON_STRUCT_BEGIN(c, Writer::Config);
0389     ACTS_PYTHON_MEMBER(writerName);
0390     ACTS_PYTHON_MEMBER(trackingGeometry);
0391     ACTS_PYTHON_MEMBER(inputCollection);
0392     ACTS_PYTHON_MEMBER(infoBoxTitle);
0393     ACTS_PYTHON_MEMBER(outputDir);
0394     ACTS_PYTHON_STRUCT_END();
0395   }
0396   {
0397     using Writer =
0398         ActsExamples::SvgPointWriter<ActsExamples::SimHit,
0399                                      ActsExamples::AccessorPositionXYZ>;
0400     auto w =
0401         py::class_<Writer, ActsExamples::IWriter, std::shared_ptr<Writer>>(
0402             mex, "SvgSimHitWriter")
0403             .def(py::init<const Writer::Config&, Acts::Logging::Level>(),
0404                  py::arg("config"), py::arg("level"))
0405             .def("write",
0406                  py::overload_cast<const AlgorithmContext&>(&Writer::write));
0407 
0408     auto c = py::class_<Writer::Config>(w, "Config").def(py::init<>());
0409     ACTS_PYTHON_STRUCT_BEGIN(c, Writer::Config);
0410     ACTS_PYTHON_MEMBER(writerName);
0411     ACTS_PYTHON_MEMBER(trackingGeometry);
0412     ACTS_PYTHON_MEMBER(inputCollection);
0413     ACTS_PYTHON_MEMBER(infoBoxTitle);
0414     ACTS_PYTHON_MEMBER(outputDir);
0415     ACTS_PYTHON_STRUCT_END();
0416   }
0417 }
0418 }  // namespace Acts::Python