Skip to content

Commit

Permalink
Adding G4Trap converter in Geant4Converters and its test-case
Browse files Browse the repository at this point in the history
  • Loading branch information
cms committed Oct 25, 2024
1 parent aff6611 commit f16ac93
Show file tree
Hide file tree
Showing 2 changed files with 82 additions and 12 deletions.
24 changes: 12 additions & 12 deletions Plugins/Geant4/src/Geant4Converters.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -37,8 +37,8 @@
#include "G4RotationMatrix.hh"
#include "G4ThreeVector.hh"
#include "G4Transform3D.hh"
#include "G4Trd.hh"
#include "G4Trap.hh"
#include "G4Trd.hh"
#include "G4Tubs.hh"
#include "G4VPhysicalVolume.hh"
#include "G4VSolid.hh"
Expand Down Expand Up @@ -249,10 +249,10 @@ Acts::Geant4ShapeConverter::trapezoidBounds(const G4Trap& g4Trap) {
ActsScalar theta = static_cast<ActsScalar>(g4Trap.GetTheta());
ActsScalar z = static_cast<ActsScalar>(g4Trap.GetZHalfLength());

ActsScalar hlX0 = (x1 + x2)*0.5;
ActsScalar hlX1 = 2*z*std::tan(theta)*std::cos(phi) + (x3+x4)*0.5;
ActsScalar hlX0 = (x1 + x2) * 0.5;
ActsScalar hlX1 = 2 * z * std::tan(theta) * std::cos(phi) + (x3 + x4) * 0.5;
ActsScalar hlY0 = y1;
ActsScalar hlY1 = y2 + 2*z*std::tan(theta)*std::sin(phi);
ActsScalar hlY1 = y2 + 2 * z * std::tan(theta) * std::sin(phi);
ActsScalar hlZ = z;

std::vector<ActsScalar> dXYZ = {(hlX0 + hlX1) * 0.5, (hlY0 + hlY1) * 0.5,
Expand All @@ -269,26 +269,26 @@ Acts::Geant4ShapeConverter::trapezoidBounds(const G4Trap& g4Trap) {
std::array<int, 2u> rAxes = {};
switch (minPos) {
case 0: {
halfLengthXminY = std::min(hlY0,hlY1);
halfLengthXmaxY = std::max(hlY0,hlY1);
halfLengthXminY = std::min(hlY0, hlY1);
halfLengthXmaxY = std::max(hlY0, hlY1);
halfLengthY = hlZ;
rAxes = {1, 2};
} break;
case 1: {
halfLengthXminY = std::min(hlX0,hlX1);
halfLengthXmaxY = std::max(hlX0,hlX1);
halfLengthXminY = std::min(hlX0, hlX1);
halfLengthXmaxY = std::max(hlX0, hlX1);
halfLengthY = hlZ;
rAxes = {0, -2};
} break;
case 2: {
if (std::abs(hlY0 - hlY1) < std::abs(hlX0 - hlX1)) {
halfLengthXminY = std::min(hlX0,hlX1);
halfLengthXmaxY = std::max(hlX0,hlX1);
halfLengthXminY = std::min(hlX0, hlX1);
halfLengthXmaxY = std::max(hlX0, hlX1);
halfLengthY = (hlY0 + hlY1) * 0.5;
rAxes = {0, 1};
} else {
halfLengthXminY = std::min(hlY0,hlY1);
halfLengthXmaxY = std::max(hlY0,hlY1);
halfLengthXminY = std::min(hlY0, hlY1);
halfLengthXmaxY = std::max(hlY0, hlY1);
halfLengthY = (hlX0 + hlX1) * 0.5;
rAxes = {-1, 0};
}
Expand Down
70 changes: 70 additions & 0 deletions Tests/UnitTests/Plugins/Geant4/Geant4ConvertersTests.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -210,6 +210,76 @@ BOOST_AUTO_TEST_CASE(Geant4TrapzoidConversion) {
CHECK_CLOSE_ABS(thicknessY, 2., 10e-10);
}

BOOST_AUTO_TEST_CASE(Geant4TrapzoidConversion) {
// Standard Trap: XY are already well defined
G4Trap trapXY("trapXY", 2, 0.523599, 0.785398, 125, 200, 125, 0.174533, 50, 125, 50, 0.174533);
auto [boundsXY, axesXY, thicknessZ] =
Acts::Geant4ShapeConverter{}.TrapezoidBounds(trapXY);
CHECK_CLOSE_ABS(
boundsXY->get(Acts::TrapezoidBounds::BoundValues::eHalfLengthXnegY), 51,
10e-10);
CHECK_CLOSE_ABS(
boundsXY->get(Acts::TrapezoidBounds::BoundValues::eHalfLengthXposY), 125,
10e-10);
CHECK_CLOSE_ABS(
boundsXY->get(Acts::TrapezoidBounds::BoundValues::eHalfLengthY), 125,
10e-10);
auto refXY = std::array<int, 2u>{0, 1};
BOOST_CHECK(axesXY == refXY);
CHECK_CLOSE_ABS(thicknessZ, 4., 10e-10);

// Flipped, yX are the coordinates
G4Trap trapyX("trapyX", 2, 0.523599, 0.785398, 50, 125, 50, 0.174533, 125, 200, 125, 0.174533);
auto [boundsyX, axesyX, thicknessZ2] =
Acts::Geant4ShapeConverter{}.trapezoidBounds(trapyX);
CHECK_CLOSE_ABS(
boundsyX->get(Acts::TrapezoidBounds::BoundValues::eHalfLengthXnegY), 87,
10e-10);
CHECK_CLOSE_ABS(
boundsyX->get(Acts::TrapezoidBounds::BoundValues::eHalfLengthXposY), 164,
10e-10);
CHECK_CLOSE_ABS(
boundsyX->get(Acts::TrapezoidBounds::BoundValues::eHalfLengthY), 88,
10e-10);
auto refyX = std::array<int, 2u>{-1, 0};
BOOST_CHECK(axesyX == refyX);
CHECK_CLOSE_ABS(thicknessZ2, 4., 10e-10);

// YZ span the trapezoid
G4Trap trapYZ("trapYZ", 200, 0.523599, 0.785398, 140, 2, 2, 0.174533, 120, 2, 2, 0.174533);
auto [boundsYZ, axesYZ, thicknessX] =
Acts::Geant4ShapeConverter{}.trapezoidBounds(trapYZ);
CHECK_CLOSE_ABS(
boundsYZ->get(Acts::TrapezoidBounds::BoundValues::eHalfLengthXnegY), 140.,
10e-10);
CHECK_CLOSE_ABS(
boundsYZ->get(Acts::TrapezoidBounds::BoundValues::eHalfLengthXposY), 283.,
10e-10);
CHECK_CLOSE_ABS(
boundsYZ->get(Acts::TrapezoidBounds::BoundValues::eHalfLengthY), 200.,
10e-10);
auto refYZ = std::array<int, 2u>{1, 2};
BOOST_CHECK(axesYZ == refYZ);
CHECK_CLOSE_ABS(thicknessX, 166, 10e-10);

// Xz span the trapezoid
G4Trap trapXz("trapXz", 200, 0.523599, 0.785398, 2, 150, 100, 0.174533, 2, 150, 100, 0.174533);
auto [boundsXz, axesXz, thicknessY] =
Acts::Geant4ShapeConverter{}.trapezoidBounds(trapXz);
CHECK_CLOSE_ABS(
boundsXz->get(Acts::TrapezoidBounds::BoundValues::eHalfLengthXnegY), 125.,
10e-10);
CHECK_CLOSE_ABS(
boundsXz->get(Acts::TrapezoidBounds::BoundValues::eHalfLengthXposY), 288.,
10e-10);
CHECK_CLOSE_ABS(
boundsXz->get(Acts::TrapezoidBounds::BoundValues::eHalfLengthY), 200.,
10e-10);
auto refXz = std::array<int, 2u>{0, -2};
BOOST_CHECK(axesXz == refXz);
CHECK_CLOSE_ABS(thicknessY, 166, 10e-10);
}

BOOST_AUTO_TEST_CASE(Geant4PlanarConversion) {
G4Box boxXY("boxXY", 23., 34., 1.);
auto pBoundsBox =
Expand Down

0 comments on commit f16ac93

Please sign in to comment.