File indexing completed on 2025-01-18 09:12:04
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "Acts/Plugins/Python/Utilities.hpp"
0010 #include "Acts/Plugins/TGeo/TGeoDetectorElement.hpp"
0011 #include "Acts/Plugins/TGeo/TGeoLayerBuilder.hpp"
0012 #include "Acts/Plugins/TGeo/TGeoParser.hpp"
0013
0014 #include <vector>
0015
0016 #include <TGeoManager.h>
0017 #include <TGeoVolume.h>
0018 #include <pybind11/pybind11.h>
0019 #include <pybind11/stl.h>
0020 #include <pybind11/stl/filesystem.h>
0021
0022 namespace py = pybind11;
0023 using namespace pybind11::literals;
0024
0025 namespace Acts::Python {
0026 void addTGeo(Context& ctx) {
0027 auto [m, mex] = ctx.get("main", "examples");
0028
0029 auto tgeo = mex.def_submodule("tgeo");
0030
0031 {
0032 py::class_<Acts::TGeoDetectorElement,
0033 std::shared_ptr<Acts::TGeoDetectorElement>>(
0034 tgeo, "TGeoDetectorElement")
0035 .def("surface", [](const Acts::TGeoDetectorElement& self) {
0036 return self.surface().getSharedPtr();
0037 });
0038 }
0039
0040 {
0041
0042
0043
0044
0045
0046
0047 tgeo.def("_convertToElements",
0048 [](const std::string& rootFileName,
0049 const std::vector<std::string>& sensitiveMatches,
0050 const std::string& localAxes, double scaleConversion) {
0051
0052 std::vector<std::shared_ptr<const Acts::TGeoDetectorElement>>
0053 tgElements;
0054
0055 TGeoManager::Import(rootFileName.c_str());
0056 if (gGeoManager != nullptr) {
0057 auto tVolume = gGeoManager->GetTopVolume();
0058 if (tVolume != nullptr) {
0059 TGeoHMatrix gmatrix = TGeoIdentity(tVolume->GetName());
0060
0061 TGeoParser::Options tgpOptions;
0062 tgpOptions.volumeNames = {tVolume->GetName()};
0063 tgpOptions.targetNames = sensitiveMatches;
0064 tgpOptions.unit = scaleConversion;
0065 TGeoParser::State tgpState;
0066 tgpState.volume = tVolume;
0067 tgpState.onBranch = true;
0068
0069 TGeoParser::select(tgpState, tgpOptions, gmatrix);
0070 tgElements.reserve(tgpState.selectedNodes.size());
0071
0072 for (const auto& snode : tgpState.selectedNodes) {
0073 auto identifier = Acts::TGeoDetectorElement::Identifier();
0074 auto tgElement = TGeoLayerBuilder::defaultElementFactory(
0075 identifier, *snode.node, *snode.transform, localAxes,
0076 scaleConversion, nullptr);
0077 tgElements.emplace_back(tgElement);
0078 }
0079 }
0080 }
0081
0082 return tgElements;
0083 });
0084 }
0085 }
0086
0087 }