#include <ActsPlugins/Geant4/Geant4Converters.hpp>
|
| std::tuple< std::shared_ptr< Acts::CylinderBounds >, double > | cylinderBounds (const G4Tubs &g4Tubs) |
| | Convert to cylinder bounds.
|
| std::shared_ptr< Acts::LineBounds > | lineBounds (const G4Tubs &g4Tubs) |
| | Convert to line/straw bounds.
|
| std::tuple< std::shared_ptr< Acts::PlanarBounds >, std::array< int, 2u >, double > | planarBounds (const G4VSolid &g4Solid) |
| | Convert to general solid into a planar shape.
|
| std::tuple< std::shared_ptr< Acts::RadialBounds >, double > | radialBounds (const G4Tubs &g4Tubs) |
| | Convert to radial bounds.
|
| std::tuple< std::shared_ptr< Acts::RectangleBounds >, std::array< int, 2u >, double > | rectangleBounds (const G4Box &g4Box) |
| | Convert to rectangle bounds.
|
| std::tuple< std::shared_ptr< Acts::TrapezoidBounds >, std::array< int, 2u >, double > | trapezoidBounds (const G4Trap &g4Trap) |
| | Convert to trapezoid bounds - from Trap.
|
| std::tuple< std::shared_ptr< Acts::TrapezoidBounds >, std::array< int, 2u >, double > | trapezoidBounds (const G4Trd &g4Trd) |
| | Convert to trapezoid bounds - from Trd.
|
|
| bool | keepAxisOrder = false |
| | A description to keep the axis order, if it is set to false cyclic re-ordering will happen, otherwise axis flip if needed in order to keep the system right-handed.
|
| double | scale = 1. |
| | A scale between Geant4 and ACTS.
|
◆ cylinderBounds()
| std::tuple< std::shared_ptr< Acts::CylinderBounds >, double > ActsPlugins::Geant4ShapeConverter::cylinderBounds |
( |
const G4Tubs & | g4Tubs | ) |
|
Convert to cylinder bounds.
- Parameters
-
| g4Tubs | a Geant4 tube shape |
- Returns
- an Acts Cylinder bounds object, and thickness
◆ lineBounds()
| std::shared_ptr< Acts::LineBounds > ActsPlugins::Geant4ShapeConverter::lineBounds |
( |
const G4Tubs & | g4Tubs | ) |
|
Convert to line/straw bounds.
- Parameters
-
| g4Tubs | a Geant4 tube shape |
- Returns
- an Acts line bounds object and thickness
◆ planarBounds()
| std::tuple< std::shared_ptr< Acts::PlanarBounds >, std::array< int, 2u >, double > ActsPlugins::Geant4ShapeConverter::planarBounds |
( |
const G4VSolid & | g4Solid | ) |
|
Convert to general solid into a planar shape.
- Parameters
-
| g4Solid | a Geant4 solid shape |
- Returns
- an ACTS Planar bounds object, the axes, and the thickness of the compressed dimension
◆ radialBounds()
| std::tuple< std::shared_ptr< Acts::RadialBounds >, double > ActsPlugins::Geant4ShapeConverter::radialBounds |
( |
const G4Tubs & | g4Tubs | ) |
|
Convert to radial bounds.
- Parameters
-
| g4Tubs | a Geant4 tube shape |
- Returns
- an Acts Radial bounds object and thickness
◆ rectangleBounds()
| std::tuple< std::shared_ptr< Acts::RectangleBounds >, std::array< int, 2u >, double > ActsPlugins::Geant4ShapeConverter::rectangleBounds |
( |
const G4Box & | g4Box | ) |
|
Convert to rectangle bounds.
- Parameters
-
- Returns
- an ACTS Rectangle bounds shape, axis orientation, and thickness
◆ trapezoidBounds() [1/2]
| std::tuple< std::shared_ptr< Acts::TrapezoidBounds >, std::array< int, 2u >, double > ActsPlugins::Geant4ShapeConverter::trapezoidBounds |
( |
const G4Trap & | g4Trap | ) |
|
Convert to trapezoid bounds - from Trap.
- Parameters
-
| g4Trap | a Geant4 trapezoid shape |
- Returns
- an ACTS Trapezoid bounds object, axis orientation, and thickness
◆ trapezoidBounds() [2/2]
| std::tuple< std::shared_ptr< Acts::TrapezoidBounds >, std::array< int, 2u >, double > ActsPlugins::Geant4ShapeConverter::trapezoidBounds |
( |
const G4Trd & | g4Trd | ) |
|
Convert to trapezoid bounds - from Trd.
- Parameters
-
| g4Trd | a Geant4 trapezoid shape |
- Returns
- an ACTS Trapezoid bounds object, axis orientation, and thickness
◆ keepAxisOrder
| bool ActsPlugins::Geant4ShapeConverter::keepAxisOrder = false |
A description to keep the axis order, if it is set to false cyclic re-ordering will happen, otherwise axis flip if needed in order to keep the system right-handed.
◆ scale
| double ActsPlugins::Geant4ShapeConverter::scale = 1. |
A scale between Geant4 and ACTS.