Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /acts/Python/Plugins/src/TGeo.cpp was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

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/Root/TGeoDetectorElement.hpp"
0010 #include "ActsPlugins/Root/TGeoLayerBuilder.hpp"
0011 #include "ActsPlugins/Root/TGeoParser.hpp"
0012 #include "ActsPython/Utilities/Helpers.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 PYBIND11_MODULE(ActsPluginsPythonBindingsTGeo, tgeo) {
0026   using namespace Acts;
0027   using namespace ActsPlugins;
0028 
0029   // Basic bindings: detector element
0030   {
0031     py::class_<TGeoDetectorElement, std::shared_ptr<TGeoDetectorElement>>(
0032         tgeo, "TGeoDetectorElement")
0033         .def("surface", [](const TGeoDetectorElement& self) {
0034           return self.surface().getSharedPtr();
0035         });
0036   }
0037 
0038   {
0039     /// Helper function to test if the automatic geometry conversion works
0040     ///
0041     /// @param rootFileName is the name of the GDML file
0042     /// @param sensitiveMatches is a list of strings to match sensitive volumes
0043     /// @param localAxes is the TGeo->ACTS axis conversion convention
0044     /// @param scaleConversion is a unit scalor conversion factor
0045     tgeo.def(
0046         "convertToElements",
0047         [](const std::string& rootFileName,
0048            const std::vector<std::string>& sensitiveMatches,
0049            const std::string& localAxes, double scaleConversion) {
0050           // Return vector
0051           std::vector<std::shared_ptr<const TGeoDetectorElement>> tgElements;
0052           // TGeo import
0053           TGeoManager::Import(rootFileName.c_str());
0054           if (gGeoManager != nullptr) {
0055             auto tVolume = gGeoManager->GetTopVolume();
0056             if (tVolume != nullptr) {
0057               TGeoHMatrix gmatrix = TGeoIdentity(tVolume->GetName());
0058 
0059               TGeoParser::Options tgpOptions;
0060               tgpOptions.volumeNames = {tVolume->GetName()};
0061               tgpOptions.targetNames = sensitiveMatches;
0062               tgpOptions.unit = scaleConversion;
0063               TGeoParser::State tgpState;
0064               tgpState.volume = tVolume;
0065               tgpState.onBranch = true;
0066 
0067               TGeoParser::select(tgpState, tgpOptions, gmatrix);
0068               tgElements.reserve(tgpState.selectedNodes.size());
0069 
0070               for (const auto& snode : tgpState.selectedNodes) {
0071                 auto identifier = TGeoDetectorElement::Identifier();
0072                 auto tgElement = TGeoLayerBuilder::defaultElementFactory(
0073                     identifier, *snode.node, *snode.transform, localAxes,
0074                     scaleConversion, nullptr);
0075                 tgElements.emplace_back(tgElement);
0076               }
0077             }
0078           }
0079           // Return the elements
0080           return tgElements;
0081         });
0082   }
0083 }