Skip to content

Commit

Permalink
refactor!: remove ActsScalar (#3873)
Browse files Browse the repository at this point in the history
We concluded, that there is not really a use for `ActsScalar` as it is
now:
- CPUs should have the same performance for `float` and `double`
- If we wanted to check performance, we would want separate scalar types
for all components
- It is too complicated to adapt everything correctly.

blocked by:
- #3894
  • Loading branch information
AJPfleger authored Nov 26, 2024
1 parent a0aed8c commit 622902b
Show file tree
Hide file tree
Showing 298 changed files with 2,131 additions and 2,333 deletions.
7 changes: 0 additions & 7 deletions CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -550,13 +550,6 @@ if(ACTS_BUILD_DOCS)
find_package(Sphinx REQUIRED)
endif()

if(ACTS_CUSTOM_SCALARTYPE)
message(
STATUS
"Building Acts with custom scalar type: ${ACTS_CUSTOM_SCALARTYPE}"
)
endif()

# core library, core plugins, and other components
add_component(Core Core)
add_subdirectory(Plugins)
Expand Down
7 changes: 0 additions & 7 deletions Core/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -35,13 +35,6 @@ if(ACTS_SOURCELINK_SBO_SIZE)
)
endif()

if(ACTS_CUSTOM_SCALARTYPE)
target_compile_definitions(
ActsCore
PUBLIC -DACTS_CUSTOM_SCALARTYPE=${ACTS_CUSTOM_SCALARTYPE}
)
endif()

if(ACTS_LOG_FAILURE_THRESHOLD)
message(
STATUS
Expand Down
8 changes: 4 additions & 4 deletions Core/include/Acts/Clusterization/TimedClusterization.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -17,16 +17,16 @@ namespace Acts::Ccl {

template <typename Cell>
concept HasRetrievableTimeInfo = requires(Cell cell) {
{ getCellTime(cell) } -> std::same_as<Acts::ActsScalar>;
{ getCellTime(cell) } -> std::same_as<double>;
};

template <Acts::Ccl::HasRetrievableTimeInfo Cell, std::size_t N>
struct TimedConnect : public Acts::Ccl::DefaultConnect<Cell, N> {
Acts::ActsScalar timeTolerance{std::numeric_limits<Acts::ActsScalar>::max()};
double timeTolerance{std::numeric_limits<double>::max()};

TimedConnect() = default;
TimedConnect(Acts::ActsScalar time);
TimedConnect(Acts::ActsScalar time, bool commonCorner)
TimedConnect(double time);
TimedConnect(double time, bool commonCorner)
requires(N == 2);
~TimedConnect() override = default;

Expand Down
5 changes: 2 additions & 3 deletions Core/include/Acts/Clusterization/TimedClusterization.ipp
Original file line number Diff line number Diff line change
Expand Up @@ -9,11 +9,10 @@
namespace Acts::Ccl {

template <Acts::Ccl::HasRetrievableTimeInfo Cell, std::size_t N>
TimedConnect<Cell, N>::TimedConnect(Acts::ActsScalar time)
: timeTolerance(time) {}
TimedConnect<Cell, N>::TimedConnect(double time) : timeTolerance(time) {}

template <Acts::Ccl::HasRetrievableTimeInfo Cell, std::size_t N>
TimedConnect<Cell, N>::TimedConnect(Acts::ActsScalar time, bool commonCorner)
TimedConnect<Cell, N>::TimedConnect(double time, bool commonCorner)
requires(N == 2)
: Acts::Ccl::DefaultConnect<Cell, N>(commonCorner), timeTolerance(time) {}

Expand Down
33 changes: 11 additions & 22 deletions Core/include/Acts/Definitions/Algebra.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -41,30 +41,19 @@ namespace Acts {
/// any conditions e.g. square size, for the dynamic-sized case. Consequently,
/// no dynamic-sized symmetric matrix type is defined. Use the
/// `ActsDynamicMatrix` instead.
///

/// Common scalar (floating point type used for the default algebra types.
///
/// Defaults to `double` but can be customized by the user.
#ifdef ACTS_CUSTOM_SCALARTYPE
using ActsScalar = ACTS_CUSTOM_SCALARTYPE;
#else
using ActsScalar = double;
#endif

template <unsigned int kSize>
using ActsVector = Eigen::Matrix<ActsScalar, kSize, 1>;
using ActsVector = Eigen::Matrix<double, kSize, 1>;

template <unsigned int kRows, unsigned int kCols>
using ActsMatrix = Eigen::Matrix<ActsScalar, kRows, kCols>;
using ActsMatrix = Eigen::Matrix<double, kRows, kCols>;

template <unsigned int kSize>
using ActsSquareMatrix = Eigen::Matrix<ActsScalar, kSize, kSize>;
using ActsSquareMatrix = Eigen::Matrix<double, kSize, kSize>;

using ActsDynamicVector = Eigen::Matrix<ActsScalar, Eigen::Dynamic, 1>;
using ActsDynamicVector = Eigen::Matrix<double, Eigen::Dynamic, 1>;

using ActsDynamicMatrix =
Eigen::Matrix<ActsScalar, Eigen::Dynamic, Eigen::Dynamic>;
using ActsDynamicMatrix = Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>;

/// @defgroup coordinates-types Fixed-size vector/matrix types for coordinates
///
Expand All @@ -84,22 +73,22 @@ using SquareMatrix3 = ActsSquareMatrix<3>;
using SquareMatrix4 = ActsSquareMatrix<4>;

// pure translation transformations
using Translation2 = Eigen::Translation<ActsScalar, 2>;
using Translation3 = Eigen::Translation<ActsScalar, 3>;
using Translation2 = Eigen::Translation<double, 2>;
using Translation3 = Eigen::Translation<double, 3>;

// linear (rotation) matrices
using RotationMatrix2 = ActsMatrix<2, 2>;
using RotationMatrix3 = ActsMatrix<3, 3>;

// pure rotation defined by a rotation angle around a rotation axis
using AngleAxis3 = Eigen::AngleAxis<ActsScalar>;
using AngleAxis3 = Eigen::AngleAxis<double>;

// combined affine transformations. types are chosen for better data alignment:
// - 2d affine compact stored as 2x3 matrix
// - 3d affine stored as 4x4 matrix
using Transform2 = Eigen::Transform<ActsScalar, 2, Eigen::AffineCompact>;
using Transform3 = Eigen::Transform<ActsScalar, 3, Eigen::Affine>;
using Transform2 = Eigen::Transform<double, 2, Eigen::AffineCompact>;
using Transform3 = Eigen::Transform<double, 3, Eigen::Affine>;

constexpr ActsScalar s_transformEquivalentTolerance = 1e-9;
constexpr double s_transformEquivalentTolerance = 1e-9;

} // namespace Acts
4 changes: 2 additions & 2 deletions Core/include/Acts/Definitions/Direction.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -41,7 +41,7 @@ class Direction final {
/// @param scalar is the signed value
///
/// @return a direction enum
static constexpr Direction fromScalar(ActsScalar scalar) {
static constexpr Direction fromScalar(double scalar) {
assert(scalar != 0);
return scalar >= 0 ? Value::Positive : Value::Negative;
}
Expand All @@ -53,7 +53,7 @@ class Direction final {
/// @param scalar is the signed value
///
/// @return a direction enum
static constexpr Direction fromScalarZeroAsPositive(ActsScalar scalar) {
static constexpr Direction fromScalarZeroAsPositive(double scalar) {
return scalar >= 0 ? Value::Positive : Value::Negative;
}

Expand Down
7 changes: 3 additions & 4 deletions Core/include/Acts/Definitions/Tolerance.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -15,19 +15,18 @@
namespace Acts {

/// Tolerance for being numerical equal for geometry building
static constexpr ActsScalar s_epsilon =
3 * std::numeric_limits<ActsScalar>::epsilon();
static constexpr double s_epsilon = 3 * std::numeric_limits<double>::epsilon();

/// Tolerance for being on Surface
///
/// @note This is intentionally given w/o an explicit unit to avoid having
/// to include the units header unnecessarily. With the native length
/// unit of mm this corresponds to 0.1um.
static constexpr ActsScalar s_onSurfaceTolerance = 1e-4;
static constexpr double s_onSurfaceTolerance = 1e-4;

/// Tolerance for not being within curvilinear projection
/// this allows using the same curvilinear frame to eta = 6,
/// validity tested with IntegrationTests/PropagationTest
static constexpr ActsScalar s_curvilinearProjTolerance = 0.999995;
static constexpr double s_curvilinearProjTolerance = 0.999995;

} // namespace Acts
6 changes: 3 additions & 3 deletions Core/include/Acts/Detector/Blueprint.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -51,7 +51,7 @@ struct Node final {
/// @param cs the children of the node
/// @param e the estimated extent of the node (optional)
Node(const std::string& n, const Transform3& t, VolumeBounds::BoundsType bt,
const std::vector<ActsScalar>& bv, const std::vector<BinningValue>& bss,
const std::vector<double>& bv, const std::vector<BinningValue>& bss,
std::vector<std::unique_ptr<Node>> cs = {}, const Extent& e = Extent())
: name(n),
transform(t),
Expand All @@ -73,7 +73,7 @@ struct Node final {
/// @param isb the internal structure builder (optional)
/// @param e the estimated extent of the node (optional)
Node(const std::string& n, const Transform3& t, VolumeBounds::BoundsType bt,
const std::vector<ActsScalar>& bv,
const std::vector<double>& bv,
std::shared_ptr<const IInternalStructureBuilder> isb = nullptr,
const Extent& e = Extent())
: name(n),
Expand All @@ -90,7 +90,7 @@ struct Node final {
/// The boundary type
VolumeBounds::BoundsType boundsType = VolumeBounds::eOther;
/// The associated values
std::vector<ActsScalar> boundaryValues = {};
std::vector<double> boundaryValues = {};
/// Parent node - nullptr for root only
const Node* parent = nullptr;
/// Branch definitions: children
Expand Down
4 changes: 2 additions & 2 deletions Core/include/Acts/Detector/DetectorVolume.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -61,8 +61,8 @@ class Detector;
class DetectorVolume : public std::enable_shared_from_this<DetectorVolume> {
public:
using BoundingBox =
Acts::AxisAlignedBoundingBox<Acts::Experimental::DetectorVolume,
Acts::ActsScalar, 3>;
Acts::AxisAlignedBoundingBox<Acts::Experimental::DetectorVolume, double,
3>;

friend class DetectorVolumeFactory;

Expand Down
10 changes: 5 additions & 5 deletions Core/include/Acts/Detector/KdtSurfacesProvider.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -34,10 +34,10 @@ class KdtSurfaces {
public:
/// Broadcast the surface KDT type
using KDTS =
KDTree<kDIM, std::shared_ptr<Surface>, ActsScalar, std::array, bSize>;
KDTree<kDIM, std::shared_ptr<Surface>, double, std::array, bSize>;

/// Broadcast the query definition
using Query = std::array<ActsScalar, kDIM>;
using Query = std::array<double, kDIM>;

/// Broadcast the entry
using Entry = std::pair<Query, std::shared_ptr<Surface>>;
Expand Down Expand Up @@ -86,7 +86,7 @@ class KdtSurfaces {
///
/// @return the matching surfaces from the KDT structure
std::vector<std::shared_ptr<Surface>> surfaces(
const RangeXD<kDIM, ActsScalar>& range) const {
const RangeXD<kDIM, double>& range) const {
// Strip the surfaces
std::vector<std::shared_ptr<Surface>> surfacePtrs;
auto surfaceQuery = m_kdt->rangeSearchWithKey(range);
Expand All @@ -101,7 +101,7 @@ class KdtSurfaces {
///
/// @return the matching surfaces fpulled from the KDT structure
std::vector<std::shared_ptr<Surface>> surfaces(const Extent& extent) const {
RangeXD<kDIM, ActsScalar> qRange;
RangeXD<kDIM, double> qRange;
for (auto [ibv, v] : enumerate(m_casts)) {
qRange[ibv] = extent.range(v);
}
Expand Down Expand Up @@ -148,7 +148,7 @@ class KdtSurfaces {
float weight = 1. / cQueries.size();
for (auto& q : cQueries) {
std::transform(c.begin(), c.end(), q.begin(), c.begin(),
std::plus<ActsScalar>());
std::plus<double>());
}
std::for_each(c.begin(), c.end(), [&](auto& v) { v *= weight; });
return c;
Expand Down
2 changes: 1 addition & 1 deletion Core/include/Acts/Detector/MultiWireStructureBuilder.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -38,7 +38,7 @@ class MultiWireStructureBuilder {
Transform3 transform = Transform3::Identity();

/// The bounds of the multi-wire volume
std::vector<ActsScalar> mlBounds = {};
std::vector<double> mlBounds = {};

// The binning of the multi wire structure
std::vector<ProtoBinning> mlBinning = {};
Expand Down
17 changes: 8 additions & 9 deletions Core/include/Acts/Detector/ProtoBinning.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -34,7 +34,7 @@ struct ProtoBinning {
/// The axis boundary type: Open, Bound or Closed
Acts::AxisBoundaryType boundaryType = Acts::AxisBoundaryType::Bound;
/// The binning edges
std::vector<ActsScalar> edges = {};
std::vector<double> edges = {};
/// An expansion for the filling (in bins)
std::size_t expansion = 0u;
/// Indication if this is an auto-range binning
Expand All @@ -47,7 +47,7 @@ struct ProtoBinning {
/// @param e the bin edges (variable binning)
/// @param exp the expansion (in bins)
ProtoBinning(BinningValue bValue, Acts::AxisBoundaryType bType,
const std::vector<ActsScalar>& e, std::size_t exp = 0u)
const std::vector<double>& e, std::size_t exp = 0u)
: binValue(bValue),
axisType(Acts::AxisType::Variable),
boundaryType(bType),
Expand All @@ -67,9 +67,8 @@ struct ProtoBinning {
/// @param maxE the highest edge of the binning
/// @param nbins the number of bins
/// @param exp the expansion (in bins)
ProtoBinning(BinningValue bValue, Acts::AxisBoundaryType bType,
ActsScalar minE, ActsScalar maxE, std::size_t nbins,
std::size_t exp = 0u)
ProtoBinning(BinningValue bValue, Acts::AxisBoundaryType bType, double minE,
double maxE, std::size_t nbins, std::size_t exp = 0u)
: binValue(bValue), boundaryType(bType), expansion(exp) {
if (minE >= maxE) {
std::string msg = "ProtoBinning: Invalid binning for value '";
Expand All @@ -84,7 +83,7 @@ struct ProtoBinning {
"ProtoBinning: Invalid binning, at least one bin is needed.");
}

ActsScalar stepE = (maxE - minE) / nbins;
double stepE = (maxE - minE) / nbins;
edges.reserve(nbins + 1);
for (std::size_t i = 0; i <= nbins; i++) {
edges.push_back(minE + i * stepE);
Expand Down Expand Up @@ -140,15 +139,15 @@ struct BinningDescription {
Acts::AxisBoundaryType boundaryType =
bData.option == open ? Acts::AxisBoundaryType::Bound
: Acts::AxisBoundaryType::Closed;
std::vector<ActsScalar> edges;
std::vector<double> edges;
if (bData.type == equidistant) {
bDesc.binning.push_back(ProtoBinning(bData.binvalue, boundaryType,
bData.min, bData.max, bData.bins(),
0u));

} else {
std::for_each(bData.boundaries().begin(), bData.boundaries().end(),
[&](ActsScalar edge) { edges.push_back(edge); });
[&](double edge) { edges.push_back(edge); });
bDesc.binning.push_back(
ProtoBinning(bData.binvalue, boundaryType, edges, 0u));
}
Expand All @@ -170,7 +169,7 @@ struct BinningDescription {
} else {
std::vector<float> edges;
std::for_each(b.edges.begin(), b.edges.end(),
[&](ActsScalar edge) { edges.push_back(edge); });
[&](double edge) { edges.push_back(edge); });
binUtility += BinUtility(edges, bOption, b.binValue);
}
}
Expand Down
2 changes: 1 addition & 1 deletion Core/include/Acts/Detector/ProtoSupport.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -38,7 +38,7 @@ struct ProtoSupport {
Surface::SurfaceType type = Surface::SurfaceType::Other;

/// The offset of the support to an estimated position (e.g. from an extent)
ActsScalar offset = 0.;
double offset = 0.;

/// A given extent from the volume, this allows to set support surfaces
/// to fit into given volume extensions (flagged by the binning value
Expand Down
2 changes: 1 addition & 1 deletion Core/include/Acts/Detector/VolumeStructureBuilder.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -42,7 +42,7 @@ class VolumeStructureBuilder : public IExternalStructureBuilder {
/// The starting transform
Transform3 transform = Transform3::Identity();
/// The values (if already defined)
std::vector<ActsScalar> boundValues = {};
std::vector<double> boundValues = {};
/// The optional extent to feed into the values
std::optional<Extent> extent = std::nullopt;
/// Some auxiliary information
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -73,7 +73,7 @@ DetectorComponent::PortalContainer connect(
/// @param logLevel is the screen logging level
///
/// @return extracted boundary values
std::array<std::vector<ActsScalar>, 3u> xyzBoundaries(
std::array<std::vector<double>, 3u> xyzBoundaries(
const GeometryContext& gctx,
const std::vector<const DetectorVolume*>& volumes,
Acts::Logging::Level logLevel = Acts::Logging::INFO);
Expand Down
Loading

0 comments on commit 622902b

Please sign in to comment.