Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:12:15

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 "Acts/Plugins/Covfie/FieldConversion.hpp"
0010 
0011 #include <cmath>
0012 #include <stdexcept>
0013 #include <type_traits>
0014 
0015 namespace Acts::CovfiePlugin {
0016 
0017 /// @brief Creates a strided covfie field that stores the values of the
0018 /// magnetic field in the volume given by min and max using a fixed sample
0019 /// spacing (determined by nPoints).
0020 ///
0021 /// @param magneticField The acts magnetic field.
0022 /// @param cache The acts cache.
0023 /// @param nPoints 3D array of containing the number of bins for each axis.
0024 /// @param min (min_x, min_y, min_z)
0025 /// @param max (max_x, max_y, max_z)
0026 /// @return A strided covfie field.
0027 template <typename magnetic_field_t>
0028 auto newBuilder(const magnetic_field_t& magneticField,
0029                 typename magnetic_field_t::Cache& cache,
0030                 const std::array<std::size_t, 3>& nPoints,
0031                 const Acts::Vector3& min, const Acts::Vector3& max) {
0032   using Field = covfie::field<BuilderBackend>;
0033 
0034   // Hack to avoid the fact that the domain of ACTS magnetic fields is defined
0035   // as a half-open interval. Has the potential to introduce very minor
0036   // floating point errors, but no easy way to fix this right now.
0037   // TODO: Fix the aforementioned problem.
0038   std::vector<double> maxima = {
0039       std::nexttoward(max[0], -std::numeric_limits<double>::infinity()),
0040       std::nexttoward(max[1], -std::numeric_limits<double>::infinity()),
0041       std::nexttoward(max[1], -std::numeric_limits<double>::infinity()),
0042   };
0043 
0044   Field field(covfie::make_parameter_pack(
0045       Field::backend_t::configuration_t{nPoints[0], nPoints[1], nPoints[2]}));
0046 
0047   Field::view_t view(field);
0048 
0049   std::array<double, 3> sampleSpacing = {
0050       (max.x() - min.x()) / (nPoints[0] - 1),
0051       (max.y() - min.y()) / (nPoints[1] - 1),
0052       (max.z() - min.z()) / (nPoints[2] - 1)};
0053 
0054   for (std::size_t x = 0; x < nPoints[0]; x++) {
0055     for (std::size_t y = 0; y < nPoints[1]; y++) {
0056       for (std::size_t z = 0; z < nPoints[2]; z++) {
0057         Acts::Vector3 position{
0058             std::min(x * sampleSpacing[0] + min[0], maxima[0]),
0059             std::min(y * sampleSpacing[1] + min[1], maxima[1]),
0060             std::min(z * sampleSpacing[2] + min[2], maxima[2])};
0061 
0062         Field::view_t::output_t& p = view.at(x, y, z);
0063         Result<Vector3> result = magneticField.getField(position, cache);
0064 
0065         if (!result.ok()) {
0066           throw std::runtime_error("Field lookup failed!");
0067         }
0068 
0069         Acts::Vector3 rv = *result;
0070         p[0] = static_cast<float>(rv[0]);
0071         p[1] = static_cast<float>(rv[1]);
0072         p[2] = static_cast<float>(rv[2]);
0073       }
0074     }
0075   }
0076 
0077   return field;
0078 }
0079 
0080 /// @brief Generate the affine covfie configuration (scaling and rotation)
0081 /// given the size of the field (min and max)
0082 ///
0083 /// @param nPoints 3D array of containing the number of bins for each axis.
0084 /// @param min (min_x, min_y, min_z)
0085 /// @param max (max_x, max_y, max_z)
0086 /// @return The affine field configuration.
0087 template <typename backend_t>
0088 typename backend_t::configuration_t affineConfiguration(
0089     const std::array<std::size_t, 3>& nPoints, const Acts::Vector3& min,
0090     const Acts::Vector3& max) {
0091   auto scaling = covfie::algebra::affine<3>::scaling(
0092       static_cast<float>((nPoints[0] - 1) / (max[0] - min[0])),
0093       static_cast<float>((nPoints[1] - 1) / (max[1] - min[1])),
0094       static_cast<float>((nPoints[2] - 1) / (max[2] - min[2])));
0095 
0096   auto translation = covfie::algebra::affine<3>::translation(
0097       static_cast<float>(-min[0]), static_cast<float>(-min[1]),
0098       static_cast<float>(-min[2]));
0099 
0100   return {scaling * translation};
0101 }
0102 
0103 /// @brief Uses std::nextafter to generates a clamp backend
0104 /// configuration where arguments min and max hold floating point values.
0105 /// @param min (min_x, min_y, min_z)
0106 /// @param max (max_x, max_y, max_z)
0107 /// @return The clamp field configuration.
0108 template <typename backend_t>
0109 typename backend_t::configuration_t clampConfigurationFloat(
0110     const Acts::Vector3& min, const Acts::Vector3& max) {
0111   return {{std::nextafter(static_cast<float>(min[0]),
0112                           std::numeric_limits<float>::infinity()),
0113            std::nextafter(static_cast<float>(min[1]),
0114                           std::numeric_limits<float>::infinity()),
0115            std::nextafter(static_cast<float>(min[2]),
0116                           std::numeric_limits<float>::infinity())},
0117           {std::nextafter(static_cast<float>(max[0]),
0118                           -std::numeric_limits<float>::infinity()),
0119            std::nextafter(static_cast<float>(max[1]),
0120                           -std::numeric_limits<float>::infinity()),
0121            std::nextafter(static_cast<float>(max[2]),
0122                           -std::numeric_limits<float>::infinity())}};
0123 }
0124 
0125 /// @brief Creates a covfie field from a generic magnetic field.
0126 /// @param magneticField The generic magnetic field.
0127 /// @param cache The cache.
0128 /// @param nPoints 3D array of containing the number of bins for each axis.
0129 /// @param min (min_x, min_y, min_z)
0130 /// @param max (max_x, max_y, max_z)
0131 /// @return A clamp affine linear strided covfie field.
0132 template <typename magnetic_field_t>
0133 InterpolatedField covfieFieldLinear(const magnetic_field_t& magneticField,
0134                                     typename magnetic_field_t::Cache& cache,
0135                                     const std::array<std::size_t, 3>& nPoints,
0136                                     const Acts::Vector3& min,
0137                                     const Acts::Vector3& max) {
0138   auto builder = newBuilder(magneticField, cache, nPoints, min, max);
0139   InterpolatedField field(covfie::make_parameter_pack(
0140       clampConfigurationFloat<InterpolatedField::backend_t>(min, max),
0141       affineConfiguration<InterpolatedField::backend_t::backend_t>(nPoints, min,
0142                                                                    max),
0143       InterpolatedField::backend_t::backend_t::backend_t::configuration_t{},
0144       builder.backend()));
0145 
0146   return field;
0147 }
0148 
0149 /// @brief Creates a covfie field from a magnetic field provider by sampling it.
0150 /// @param magneticField The acts magnetic field provider.
0151 /// @param cache The acts cache.
0152 /// @param nPoints 3D array of containing the number of bins for each axis.
0153 /// @param min (min_x, min_y, min_z)
0154 /// @param max (max_x, max_y, max_z)
0155 /// @return A clamp affine linear strided covfie field.
0156 InterpolatedField covfieField(const Acts::MagneticFieldProvider& magneticField,
0157                               Acts::MagneticFieldProvider::Cache& cache,
0158                               const std::array<std::size_t, 3>& nPoints,
0159                               const Acts::Vector3& min,
0160                               const Acts::Vector3& max) {
0161   return covfieFieldLinear(magneticField, cache, nPoints, min, max);
0162 }
0163 
0164 /// @brief Creates a covfie field from an interpolated magnetic field.
0165 /// @param magneticField The acts interpolated magnetic field.
0166 /// @return A clamp affine linear strided covfie field.
0167 InterpolatedField covfieField(
0168     const Acts::InterpolatedMagneticField& magneticField) {
0169   Acts::MagneticFieldContext ctx;
0170   auto cache = magneticField.makeCache(ctx);
0171   const std::vector<double>& old_min = magneticField.getMin();
0172   const std::vector<double>& old_max = magneticField.getMax();
0173   const std::vector<std::size_t>& old_nbins = magneticField.getNBins();
0174   Acts::Vector3 min{old_min.at(0), old_min.at(1), old_min.at(2)};
0175   Acts::Vector3 max{old_max.at(0), old_max.at(1), old_max.at(2)};
0176   std::array<std::size_t, 3> nPoints{old_nbins.at(0), old_nbins.at(1),
0177                                      old_nbins.at(2)};
0178   return covfieFieldLinear(magneticField, cache, nPoints, min, max);
0179 }
0180 
0181 /// @brief Creates a covfie field from a constant B field.
0182 /// @param magneticField The acts constant magnetic field.
0183 /// @return A constant covfie field.
0184 ConstantField covfieField(const Acts::ConstantBField& magneticField) {
0185   auto B = magneticField.getField();
0186   ConstantField field(
0187       covfie::make_parameter_pack(ConstantField::backend_t::configuration_t{
0188           static_cast<float>(B[0]), static_cast<float>(B[1]),
0189           static_cast<float>(B[2])}));
0190   return field;
0191 }
0192 
0193 }  // namespace Acts::CovfiePlugin