Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:12:04

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/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 "Acts/Utilities/Helpers.hpp"
0024 #include "ActsExamples/EventData/GeometryContainers.hpp"
0025 #include "ActsExamples/EventData/SimHit.hpp"
0026 #include "ActsExamples/EventData/SimSpacePoint.hpp"
0027 #include "ActsExamples/Io/Svg/SvgPointWriter.hpp"
0028 #include "ActsExamples/Io/Svg/SvgTrackingGeometryWriter.hpp"
0029 
0030 #include <memory>
0031 #include <string>
0032 #include <tuple>
0033 #include <vector>
0034 
0035 #include <pybind11/pybind11.h>
0036 #include <pybind11/stl.h>
0037 
0038 namespace py = pybind11;
0039 using namespace pybind11::literals;
0040 
0041 using namespace Acts;
0042 using namespace ActsExamples;
0043 
0044 namespace {
0045 
0046 // A cache object
0047 using PortalCache = std::list<std::string>;
0048 
0049 /// @brief helper tuple to define which views and which range to be used
0050 ///
0051 /// @param view is the view type, e.g. 'xy', 'zr', ...
0052 /// @param selection is the selection of the volume, e.g. 'all', 'sensitives', 'portals', 'materials'
0053 using ViewAndRange =
0054     std::tuple<std::string, std::vector<std::string>, Acts::Extent>;
0055 
0056 /// Helper function to be picked in different access patterns
0057 ///
0058 /// @param pVolume is the proto volume to be drawn
0059 /// @param identification is the identification of the volume
0060 /// @param viewAndRange is the view, selection and range to be drawn
0061 /// @param portalCache is a portal cache to avoid multiple drawings of the same portal
0062 ///
0063 /// Returns an svg object in the right view
0064 actsvg::svg::object drawDetectorVolume(const Svg::ProtoVolume& pVolume,
0065                                        const std::string& identification,
0066                                        const ViewAndRange& viewAndRange,
0067                                        PortalCache& portalCache) {
0068   actsvg::svg::object svgDet;
0069   svgDet._id = identification;
0070   svgDet._tag = "g";
0071 
0072   const auto& [view, selection, viewRange] = viewAndRange;
0073 
0074   // Translate selection into booleans
0075   const bool all = rangeContainsValue(selection, "all");
0076   const bool sensitives = rangeContainsValue(selection, "sensitives");
0077   const bool portals = rangeContainsValue(selection, "portals");
0078   const bool materials = rangeContainsValue(selection, "materials");
0079 
0080   // Helper lambda for material selection
0081   auto materialSel = [&](const Svg::ProtoSurface& s) -> bool {
0082     return (materials && s._decorations.contains("material"));
0083   };
0084 
0085   // Helper lambda for view range selection
0086   auto viewRangeSel = [](const Svg::ProtoSurface& s,
0087                          const Extent& vRange) -> bool {
0088     bool accept = false;
0089     for (const auto& v : s._vertices) {
0090       if (vRange.contains(v)) {
0091         accept = true;
0092         break;
0093       }
0094     }
0095     return accept;
0096   };
0097 
0098   // -------------------- surface section
0099   // The surfaces to be drawn
0100   std::vector<Svg::ProtoVolume::surface_type> sSurfaces;
0101   sSurfaces.reserve(pVolume._v_surfaces.size());
0102   for (const auto& s : pVolume._v_surfaces) {
0103     if ((all || sensitives || materialSel(s)) && viewRangeSel(s, viewRange)) {
0104       sSurfaces.push_back(s);
0105     }
0106   }
0107 
0108   // Now draw all the surfaces
0109   for (const auto& vs : sSurfaces) {
0110     if (view == "xy") {
0111       svgDet.add_object(Svg::View::xy(vs, identification));
0112     } else if (view == "zr") {
0113       svgDet.add_object(Svg::View::zr(vs, identification));
0114     } else {
0115       throw std::invalid_argument("Unknown view type");
0116     }
0117   }
0118 
0119   // -------------------- portal section
0120   std::vector<Svg::ProtoPortal> sPortals;
0121   sPortals.reserve(pVolume._portals.size());
0122   for (const auto& vp : pVolume._portals) {
0123     if ((all || portals || materialSel(vp._surface)) &&
0124         viewRangeSel(vp._surface, viewRange)) {
0125       sPortals.push_back(vp);
0126     }
0127   }
0128 
0129   // Now draw all the portals - if not already in the cache
0130   for (const auto& vp : sPortals) {
0131     auto pgID = vp._surface._decorations.find("geo_id");
0132     std::string gpIDs = "";
0133     if (pgID != vp._surface._decorations.end()) {
0134       gpIDs = pgID->second._id;
0135     }
0136 
0137     if (rangeContainsValue(portalCache, gpIDs)) {
0138       continue;
0139     }
0140 
0141     // Register this portal to the cache
0142     portalCache.insert(portalCache.begin(), gpIDs);
0143 
0144     if (view == "xy") {
0145       svgDet.add_object(Svg::View::xy(vp, identification));
0146     } else if (view == "zr") {
0147       svgDet.add_object(Svg::View::zr(vp, identification));
0148     } else {
0149       throw std::invalid_argument("Unknown view type");
0150     }
0151   }
0152   return svgDet;
0153 }
0154 
0155 // Helper function to be picked in different access patterns
0156 std::vector<actsvg::svg::object> drawDetector(
0157     const Acts::GeometryContext& gctx,
0158     const Acts::Experimental::Detector& detector,
0159     const std::string& identification,
0160     const std::vector<std::tuple<int, Svg::DetectorVolumeConverter::Options>>&
0161         volumeIdxOpts,
0162     const std::vector<ViewAndRange>& viewAndRanges) {
0163   PortalCache portalCache;
0164 
0165   // The svg object to be returned
0166   std::vector<actsvg::svg::object> svgDetViews;
0167   svgDetViews.reserve(viewAndRanges.size());
0168   for (unsigned int i = 0; i < viewAndRanges.size(); ++i) {
0169     actsvg::svg::object svgDet;
0170     svgDet._id = identification;
0171     svgDet._tag = "g";
0172     svgDetViews.push_back(svgDet);
0173   }
0174 
0175   for (const auto& [vidx, vopts] : volumeIdxOpts) {
0176     // Get the volume and convert it
0177     const auto& v = detector.volumes()[vidx];
0178     auto [pVolume, pGrid] =
0179         Svg::DetectorVolumeConverter::convert(gctx, *v, vopts);
0180 
0181     for (auto [iv, var] : Acts::enumerate(viewAndRanges)) {
0182       auto [view, selection, range] = var;
0183       // Get the view and the range
0184       auto svgVolView = drawDetectorVolume(
0185           pVolume, identification + "_vol" + std::to_string(vidx) + "_" + view,
0186           var, portalCache);
0187       svgDetViews[iv].add_object(svgVolView);
0188     }
0189   }
0190   return svgDetViews;
0191 }
0192 
0193 }  // namespace
0194 
0195 namespace Acts::Python {
0196 void addSvg(Context& ctx) {
0197   auto [m, mex] = ctx.get("main", "examples");
0198 
0199   auto svg = m.def_submodule("svg");
0200 
0201   // Some basics
0202   py::class_<actsvg::svg::object>(svg, "object");
0203 
0204   py::class_<actsvg::svg::file>(svg, "file")
0205       .def(py::init<>())
0206       .def("addObject", &actsvg::svg::file::add_object)
0207       .def("addObjects", &actsvg::svg::file::add_objects)
0208       .def("clip",
0209            [](actsvg::svg::file& self, std::array<actsvg::scalar, 4> box) {
0210              self.set_view_box(box);
0211            })
0212       .def("write", [](actsvg::svg::file& self, const std::string& filename) {
0213         std::ofstream file(filename);
0214         file << self;
0215         file.close();
0216       });
0217 
0218   // Core components, added as an acts.svg submodule
0219   {
0220     auto c = py::class_<Svg::Style>(svg, "Style").def(py::init<>());
0221     ACTS_PYTHON_STRUCT_BEGIN(c, Svg::Style);
0222     ACTS_PYTHON_MEMBER(fillColor);
0223     ACTS_PYTHON_MEMBER(fillOpacity);
0224     ACTS_PYTHON_MEMBER(highlightColor);
0225     ACTS_PYTHON_MEMBER(highlights);
0226     ACTS_PYTHON_MEMBER(strokeWidth);
0227     ACTS_PYTHON_MEMBER(strokeColor);
0228     ACTS_PYTHON_MEMBER(highlightStrokeWidth);
0229     ACTS_PYTHON_MEMBER(highlightStrokeColor);
0230     ACTS_PYTHON_MEMBER(fontSize);
0231     ACTS_PYTHON_MEMBER(fontColor);
0232     ACTS_PYTHON_MEMBER(quarterSegments);
0233     ACTS_PYTHON_STRUCT_END();
0234   }
0235 
0236   // How surfaces should be drawn: Svg Surface options & drawning
0237   {
0238     auto c = py::class_<Svg::SurfaceConverter::Options>(svg, "SurfaceOptions")
0239                  .def(py::init<>());
0240     ACTS_PYTHON_STRUCT_BEGIN(c, Svg::SurfaceConverter::Options);
0241     ACTS_PYTHON_MEMBER(style);
0242     ACTS_PYTHON_MEMBER(templateSurface);
0243     ACTS_PYTHON_STRUCT_END();
0244 
0245     // Define the proto surface
0246     py::class_<Svg::ProtoSurface>(svg, "ProtoSurface");
0247     // Convert an Acts::Surface object into an acts::svg::proto::surface
0248     svg.def("convertSurface", &Svg::SurfaceConverter::convert);
0249 
0250     // Define the view functions
0251     svg.def("viewSurface", [](const Svg::ProtoSurface& pSurface,
0252                               const std::string& identification,
0253                               const std::string& view = "xy") {
0254       if (view == "xy") {
0255         return Svg::View::xy(pSurface, identification);
0256       } else if (view == "zr") {
0257         return Svg::View::zr(pSurface, identification);
0258       } else if (view == "zphi") {
0259         return Svg::View::zphi(pSurface, identification);
0260       } else if (view == "zrphi") {
0261         return Svg::View::zrphi(pSurface, identification);
0262       } else {
0263         throw std::invalid_argument("Unknown view type");
0264       }
0265     });
0266   }
0267 
0268   // How portals should be drawn: Svg Portal options & drawning
0269   {
0270     auto c = py::class_<Svg::PortalConverter::Options>(svg, "PortalOptions")
0271                  .def(py::init<>());
0272 
0273     ACTS_PYTHON_STRUCT_BEGIN(c, Svg::PortalConverter::Options);
0274     ACTS_PYTHON_MEMBER(surfaceOptions);
0275     ACTS_PYTHON_MEMBER(linkLength);
0276     ACTS_PYTHON_MEMBER(volumeIndices);
0277     ACTS_PYTHON_STRUCT_END();
0278 
0279     // Define the proto portal
0280     py::class_<Svg::ProtoPortal>(svg, "ProtoPortal");
0281     // Convert an Acts::Experimental::Portal object into an
0282     // acts::svg::proto::portal
0283     svg.def("convertPortal", &Svg::PortalConverter::convert);
0284 
0285     // Define the view functions
0286     svg.def("viewPortal", [](const Svg::ProtoPortal& pPortal,
0287                              const std::string& identification,
0288                              const std::string& view = "xy") {
0289       if (view == "xy") {
0290         return Svg::View::xy(pPortal, identification);
0291       } else if (view == "zr") {
0292         return Svg::View::zr(pPortal, identification);
0293       } else {
0294         throw std::invalid_argument("Unknown view type");
0295       }
0296     });
0297   }
0298 
0299   // Draw primitives
0300   {
0301     svg.def("drawArrow", &actsvg::draw::arrow);
0302 
0303     svg.def("drawText", &actsvg::draw::text);
0304 
0305     svg.def("drawInfoBox", &Svg::infoBox);
0306   }
0307 
0308   // Draw Eta Lines
0309   {
0310     svg.def(
0311         "drawEtaLines",
0312         [](const std::string& id, actsvg ::scalar z, actsvg::scalar r,
0313            const std::vector<actsvg::scalar>& etaMain,
0314            actsvg::scalar strokeWidthMain, unsigned int sizeMain,
0315            bool labelMain, const std::vector<actsvg::scalar>& etaSub,
0316            actsvg::scalar strokeWidthSub, const std::vector<int> strokeDashSub,
0317            unsigned int sizeSub, bool labelSub) {
0318           // The main eta lines
0319           actsvg::style::stroke strokeMain;
0320           strokeMain._width = strokeWidthMain;
0321           actsvg::style::font fontMain;
0322           fontMain._size = sizeMain;
0323 
0324           actsvg::style::stroke strokeSub;
0325           strokeSub._width = strokeWidthSub;
0326           strokeSub._dasharray = strokeDashSub;
0327           actsvg::style::font fontSub;
0328           fontSub._size = sizeSub;
0329 
0330           return actsvg::display::eta_lines(
0331               id, z, r,
0332               {std::tie(etaMain, strokeMain, labelMain, fontMain),
0333                std::tie(etaSub, strokeSub, labelSub, fontSub)});
0334         });
0335   }
0336 
0337   {
0338     auto gco = py::class_<Svg::GridConverter::Options>(svg, "GridOptions")
0339                    .def(py::init<>());
0340     ACTS_PYTHON_STRUCT_BEGIN(gco, Svg::GridConverter::Options);
0341     ACTS_PYTHON_MEMBER(style);
0342     ACTS_PYTHON_STRUCT_END();
0343 
0344     auto isco = py::class_<Svg::IndexedSurfacesConverter::Options>(
0345                     svg, "IndexedSurfacesOptions")
0346                     .def(py::init<>());
0347     ACTS_PYTHON_STRUCT_BEGIN(isco, Svg::IndexedSurfacesConverter::Options);
0348     ACTS_PYTHON_MEMBER(gridOptions);
0349     ACTS_PYTHON_STRUCT_END();
0350   }
0351 
0352   // How detector volumes are drawn: Svg DetectorVolume options & drawning
0353   {
0354     auto c = py::class_<Svg::DetectorVolumeConverter::Options>(
0355                  svg, "DetectorVolumeOptions")
0356                  .def(py::init<>());
0357 
0358     ACTS_PYTHON_STRUCT_BEGIN(c, Svg::DetectorVolumeConverter::Options);
0359     ACTS_PYTHON_MEMBER(portalIndices);
0360     ACTS_PYTHON_MEMBER(portalOptions);
0361     ACTS_PYTHON_MEMBER(surfaceOptions);
0362     ACTS_PYTHON_MEMBER(indexedSurfacesOptions);
0363     ACTS_PYTHON_STRUCT_END();
0364 
0365     // Define the proto volume & indexed surface grid
0366     py::class_<Svg::ProtoVolume>(svg, "ProtoVolume");
0367     py::class_<Svg::ProtoIndexedSurfaceGrid>(svg, "ProtoIndexedSurfaceGrid");
0368 
0369     // Define the proto grid
0370     py::class_<Svg::ProtoGrid>(svg, "ProtoGrid");
0371 
0372     // Convert an Acts::Experimental::DetectorVolume object into an
0373     // acts::svg::proto::volume
0374     svg.def("convertDetectorVolume", &Svg::DetectorVolumeConverter::convert);
0375 
0376     // Define the view functions
0377     svg.def("drawDetectorVolume", &drawDetectorVolume);
0378   }
0379 
0380   // Draw the ProtoIndexedSurfaceGrid
0381   {
0382     svg.def("drawIndexedSurfaces",
0383             [](const Svg::ProtoIndexedSurfaceGrid& pIndexedSurfaceGrid,
0384                const std::string& identification) {
0385               return Svg::View::xy(pIndexedSurfaceGrid, identification);
0386             });
0387   }
0388 
0389   // How a detector is drawn: Svg Detector options & drawning
0390   { svg.def("drawDetector", &drawDetector); }
0391 
0392   // Legacy geometry drawing
0393   {
0394     using DefinedStyle = std::pair<GeometryIdentifier, Svg::Style>;
0395     using DefinedStyleSet = std::vector<DefinedStyle>;
0396 
0397     auto sm = py::class_<GeometryHierarchyMap<Svg::Style>>(svg, "StyleMap")
0398                   .def(py::init<DefinedStyleSet>(), py::arg("elements"));
0399 
0400     auto c = py::class_<Svg::LayerConverter::Options>(svg, "LayerOptions")
0401                  .def(py::init<>());
0402     ACTS_PYTHON_STRUCT_BEGIN(c, Svg::LayerConverter::Options);
0403     ACTS_PYTHON_MEMBER(name);
0404     ACTS_PYTHON_MEMBER(surfaceStyles);
0405     ACTS_PYTHON_MEMBER(zRange);
0406     ACTS_PYTHON_MEMBER(phiRange);
0407     ACTS_PYTHON_MEMBER(gridInfo);
0408     ACTS_PYTHON_MEMBER(moduleInfo);
0409     ACTS_PYTHON_MEMBER(projectionInfo);
0410     ACTS_PYTHON_MEMBER(labelProjection);
0411     ACTS_PYTHON_MEMBER(labelGauge);
0412     ACTS_PYTHON_STRUCT_END();
0413   }
0414 
0415   {
0416     using DefinedLayerOptions =
0417         std::pair<GeometryIdentifier, Svg::LayerConverter::Options>;
0418     using DefinedLayerOptionsSet = std::vector<DefinedLayerOptions>;
0419 
0420     auto lom =
0421         py::class_<GeometryHierarchyMap<Svg::LayerConverter::Options>>(
0422             svg, "LayerOptionMap")
0423             .def(py::init<DefinedLayerOptionsSet>(), py::arg("elements"));
0424 
0425     auto c = py::class_<Svg::TrackingGeometryConverter::Options>(
0426                  svg, "TrackingGeometryOptions")
0427                  .def(py::init<>());
0428     ACTS_PYTHON_STRUCT_BEGIN(c, Svg::TrackingGeometryConverter::Options);
0429     ACTS_PYTHON_MEMBER(prefix);
0430     ACTS_PYTHON_MEMBER(layerOptions);
0431     ACTS_PYTHON_STRUCT_END();
0432   }
0433 
0434   // Components from the ActsExamples - part of acts.examples
0435 
0436   {
0437     using Writer = ActsExamples::SvgTrackingGeometryWriter;
0438     auto w = py::class_<Writer, std::shared_ptr<Writer>>(
0439                  mex, "SvgTrackingGeometryWriter")
0440                  .def(py::init<const Writer::Config&, Acts::Logging::Level>(),
0441                       py::arg("config"), py::arg("level"))
0442                  .def("write", py::overload_cast<const AlgorithmContext&,
0443                                                  const Acts::TrackingGeometry&>(
0444                                    &Writer::write));
0445 
0446     auto c = py::class_<Writer::Config>(w, "Config").def(py::init<>());
0447     ACTS_PYTHON_STRUCT_BEGIN(c, Writer::Config);
0448     ACTS_PYTHON_MEMBER(outputDir);
0449     ACTS_PYTHON_MEMBER(converterOptions);
0450     ACTS_PYTHON_STRUCT_END();
0451   }
0452   {
0453     using Writer = ActsExamples::SvgPointWriter<ActsExamples::SimSpacePoint>;
0454     auto w =
0455         py::class_<Writer, ActsExamples::IWriter, std::shared_ptr<Writer>>(
0456             mex, "SvgSimSpacePointWriter")
0457             .def(py::init<const Writer::Config&, Acts::Logging::Level>(),
0458                  py::arg("config"), py::arg("level"))
0459             .def("write",
0460                  py::overload_cast<const AlgorithmContext&>(&Writer::write));
0461 
0462     auto c = py::class_<Writer::Config>(w, "Config").def(py::init<>());
0463     ACTS_PYTHON_STRUCT_BEGIN(c, Writer::Config);
0464     ACTS_PYTHON_MEMBER(writerName);
0465     ACTS_PYTHON_MEMBER(trackingGeometry);
0466     ACTS_PYTHON_MEMBER(inputCollection);
0467     ACTS_PYTHON_MEMBER(infoBoxTitle);
0468     ACTS_PYTHON_MEMBER(outputDir);
0469     ACTS_PYTHON_STRUCT_END();
0470   }
0471   {
0472     using Writer =
0473         ActsExamples::SvgPointWriter<ActsExamples::SimHit,
0474                                      ActsExamples::AccessorPositionXYZ>;
0475     auto w =
0476         py::class_<Writer, ActsExamples::IWriter, std::shared_ptr<Writer>>(
0477             mex, "SvgSimHitWriter")
0478             .def(py::init<const Writer::Config&, Acts::Logging::Level>(),
0479                  py::arg("config"), py::arg("level"))
0480             .def("write",
0481                  py::overload_cast<const AlgorithmContext&>(&Writer::write));
0482 
0483     auto c = py::class_<Writer::Config>(w, "Config").def(py::init<>());
0484     ACTS_PYTHON_STRUCT_BEGIN(c, Writer::Config);
0485     ACTS_PYTHON_MEMBER(writerName);
0486     ACTS_PYTHON_MEMBER(trackingGeometry);
0487     ACTS_PYTHON_MEMBER(inputCollection);
0488     ACTS_PYTHON_MEMBER(infoBoxTitle);
0489     ACTS_PYTHON_MEMBER(outputDir);
0490     ACTS_PYTHON_STRUCT_END();
0491   }
0492 }
0493 }  // namespace Acts::Python