File indexing completed on 2025-10-15 08:04:58
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "ActsPlugins/Covfie/FieldConversion.hpp"
0010
0011 #include <cmath>
0012 #include <stdexcept>
0013 #include <type_traits>
0014
0015 using namespace Acts;
0016
0017 namespace ActsPlugins::Covfie {
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029 template <typename magnetic_field_t>
0030 auto newBuilder(const magnetic_field_t& magneticField,
0031 typename magnetic_field_t::Cache& cache,
0032 const std::array<std::size_t, 3>& nPoints, const Vector3& min,
0033 const Vector3& max) {
0034 using Field = covfie::field<BuilderBackend>;
0035
0036
0037
0038
0039
0040 std::vector<double> maxima = {
0041 std::nexttoward(max[0], -std::numeric_limits<double>::infinity()),
0042 std::nexttoward(max[1], -std::numeric_limits<double>::infinity()),
0043 std::nexttoward(max[1], -std::numeric_limits<double>::infinity()),
0044 };
0045
0046 Field field(covfie::make_parameter_pack(
0047 Field::backend_t::configuration_t{nPoints[0], nPoints[1], nPoints[2]}));
0048
0049 Field::view_t view(field);
0050
0051 std::array<double, 3> sampleSpacing = {
0052 (max.x() - min.x()) / (nPoints[0] - 1),
0053 (max.y() - min.y()) / (nPoints[1] - 1),
0054 (max.z() - min.z()) / (nPoints[2] - 1)};
0055
0056 for (std::size_t x = 0; x < nPoints[0]; x++) {
0057 for (std::size_t y = 0; y < nPoints[1]; y++) {
0058 for (std::size_t z = 0; z < nPoints[2]; z++) {
0059 Vector3 position{std::min(x * sampleSpacing[0] + min[0], maxima[0]),
0060 std::min(y * sampleSpacing[1] + min[1], maxima[1]),
0061 std::min(z * sampleSpacing[2] + min[2], maxima[2])};
0062
0063 Field::view_t::output_t& p = view.at(x, y, z);
0064 Result<Vector3> result = magneticField.getField(position, cache);
0065
0066 if (!result.ok()) {
0067 throw std::runtime_error("Field lookup failed!");
0068 }
0069
0070 Vector3 rv = *result;
0071 p[0] = static_cast<float>(rv[0]);
0072 p[1] = static_cast<float>(rv[1]);
0073 p[2] = static_cast<float>(rv[2]);
0074 }
0075 }
0076 }
0077
0078 return field;
0079 }
0080
0081
0082
0083
0084
0085
0086
0087
0088 template <typename backend_t>
0089 typename backend_t::configuration_t affineConfiguration(
0090 const std::array<std::size_t, 3>& nPoints, const Vector3& min,
0091 const Vector3& max) {
0092 auto scaling = covfie::algebra::affine<3>::scaling(
0093 static_cast<float>((nPoints[0] - 1) / (max[0] - min[0])),
0094 static_cast<float>((nPoints[1] - 1) / (max[1] - min[1])),
0095 static_cast<float>((nPoints[2] - 1) / (max[2] - min[2])));
0096
0097 auto translation = covfie::algebra::affine<3>::translation(
0098 static_cast<float>(-min[0]), static_cast<float>(-min[1]),
0099 static_cast<float>(-min[2]));
0100
0101 return {scaling * translation};
0102 }
0103
0104
0105
0106
0107
0108
0109 template <typename backend_t>
0110 typename backend_t::configuration_t clampConfigurationFloat(
0111 const Vector3& min, const Vector3& max) {
0112 return {{std::nextafter(static_cast<float>(min[0]),
0113 std::numeric_limits<float>::infinity()),
0114 std::nextafter(static_cast<float>(min[1]),
0115 std::numeric_limits<float>::infinity()),
0116 std::nextafter(static_cast<float>(min[2]),
0117 std::numeric_limits<float>::infinity())},
0118 {std::nextafter(static_cast<float>(max[0]),
0119 -std::numeric_limits<float>::infinity()),
0120 std::nextafter(static_cast<float>(max[1]),
0121 -std::numeric_limits<float>::infinity()),
0122 std::nextafter(static_cast<float>(max[2]),
0123 -std::numeric_limits<float>::infinity())}};
0124 }
0125
0126
0127
0128
0129
0130
0131
0132
0133 template <typename magnetic_field_t>
0134 InterpolatedField covfieFieldLinear(const magnetic_field_t& magneticField,
0135 typename magnetic_field_t::Cache& cache,
0136 const std::array<std::size_t, 3>& nPoints,
0137 const Vector3& min, const 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
0150
0151
0152
0153
0154
0155
0156 InterpolatedField covfieField(const MagneticFieldProvider& magneticField,
0157 MagneticFieldProvider::Cache& cache,
0158 const std::array<std::size_t, 3>& nPoints,
0159 const Vector3& min, const Vector3& max) {
0160 return covfieFieldLinear(magneticField, cache, nPoints, min, max);
0161 }
0162
0163
0164
0165
0166 InterpolatedField covfieField(const InterpolatedMagneticField& magneticField) {
0167 MagneticFieldContext ctx;
0168 auto cache = magneticField.makeCache(ctx);
0169 const std::vector<double>& old_min = magneticField.getMin();
0170 const std::vector<double>& old_max = magneticField.getMax();
0171 const std::vector<std::size_t>& old_nbins = magneticField.getNBins();
0172 Vector3 min{old_min.at(0), old_min.at(1), old_min.at(2)};
0173 Vector3 max{old_max.at(0), old_max.at(1), old_max.at(2)};
0174 std::array<std::size_t, 3> nPoints{old_nbins.at(0), old_nbins.at(1),
0175 old_nbins.at(2)};
0176 return covfieFieldLinear(magneticField, cache, nPoints, min, max);
0177 }
0178
0179
0180
0181
0182 ConstantField covfieField(const ConstantBField& magneticField) {
0183 auto B = magneticField.getField();
0184 ConstantField field(
0185 covfie::make_parameter_pack(ConstantField::backend_t::configuration_t{
0186 static_cast<float>(B[0]), static_cast<float>(B[1]),
0187 static_cast<float>(B[2])}));
0188 return field;
0189 }
0190
0191 }