File indexing completed on 2025-10-16 08:04:39
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010 #include "Acts/Definitions/Units.hpp"
0011 #include "Acts/MagneticField/ConstantBField.hpp"
0012 #include "Acts/MagneticField/MagneticFieldContext.hpp"
0013 #include "Acts/MagneticField/MagneticFieldProvider.hpp"
0014 #include "Acts/MagneticField/SolenoidBField.hpp"
0015
0016
0017 #include "ActsPlugins/Covfie/FieldConversion.hpp"
0018
0019
0020 #include <array>
0021 #include <cmath>
0022 #include <sstream>
0023 #include <vector>
0024
0025
0026 #include <boost/test/unit_test.hpp>
0027
0028 using namespace Acts;
0029 using namespace ActsPlugins;
0030 using namespace Acts::UnitLiterals;
0031
0032 template <typename view_t, typename iterator_t>
0033 void checkMagneticFieldEqual(const MagneticFieldProvider& fieldProvider,
0034 MagneticFieldProvider::Cache& cache, view_t view,
0035 iterator_t points, float error_margin_half_width) {
0036 for (auto point : points) {
0037 auto x = point[0], y = point[1], z = point[2];
0038
0039 auto lookupResult = fieldProvider.getField(Vector3{x, y, z}, cache);
0040 if (!lookupResult.ok()) {
0041 throw std::runtime_error{"Field lookup failure"};
0042 }
0043 auto actsValueX = (*lookupResult)[0], actsValueY = (*lookupResult)[1],
0044 actsValueZ = (*lookupResult)[2];
0045
0046 auto covfieValues = view.at(x, y, z);
0047 auto covfieValueX = covfieValues[0], covfieValueY = covfieValues[1],
0048 covfieValueZ = covfieValues[2];
0049
0050 auto isEqual =
0051 std::abs(covfieValueX - actsValueX) <= error_margin_half_width &&
0052 std::abs(covfieValueY - actsValueY) <= error_margin_half_width &&
0053 std::abs(covfieValueZ - actsValueZ) <= error_margin_half_width;
0054
0055 std::stringstream ss;
0056 ss << "Fields are not equal at position (" << x << ", " << y << ", " << z
0057 << "). Acts: (" << actsValueX << ", " << actsValueY << ", " << actsValueZ
0058 << "), Covfie: (" << covfieValueX << ", " << covfieValueY << ", "
0059 << covfieValueZ << ")" << std::endl;
0060
0061 BOOST_CHECK_MESSAGE(isEqual, ss.str());
0062 }
0063 }
0064
0065 namespace ActsTests {
0066
0067 BOOST_AUTO_TEST_SUITE(CovfieSuite)
0068
0069 BOOST_AUTO_TEST_CASE(InterpolatedMagneticField1) {
0070 auto localToGlobalBin_xyz = [](std::array<std::size_t, 3> binsXYZ,
0071 std::array<std::size_t, 3> nBinsXYZ) {
0072 return (binsXYZ.at(0) * (nBinsXYZ.at(1) * nBinsXYZ.at(2)) +
0073 binsXYZ.at(1) * nBinsXYZ.at(2) + binsXYZ.at(2));
0074 };
0075
0076 std::vector<double> xPos = {0., 1., 2., 3.};
0077 std::vector<double> yPos = {0., 1., 2., 3.};
0078 std::vector<double> zPos = {0., 1., 2., 3.};
0079
0080 std::vector<Vector3> bField_xyz;
0081 for (int i = 0; i < 64; i++) {
0082 bField_xyz.push_back(Vector3(i, i, i));
0083 }
0084
0085 MagneticFieldContext fieldContext;
0086 auto actsField = fieldMapXYZ(localToGlobalBin_xyz, xPos, yPos, zPos,
0087 bField_xyz, 1, 1, false);
0088 MagneticFieldProvider::Cache cache = actsField.makeCache(fieldContext);
0089
0090 Covfie::InterpolatedField field = Covfie::covfieField(actsField);
0091 typename Covfie::InterpolatedField::view_t view(field);
0092
0093 std::array<std::array<float, 3>, 14> points = {{
0094 {0.f, 0.f, 0.f},
0095 {1.f, 1.f, 1.f},
0096 {2.f, 2.f, 2.f},
0097 {2.9f, 2.9f, 2.9f},
0098 {1.2f, 2.5f, 0.8f},
0099 {0.7f, 1.9f, 2.3f},
0100 {2.1f, 0.3f, 1.5f},
0101 {0.4f, 2.8f, 2.9f},
0102 {1.6f, 1.2f, 0.5f},
0103 {2.3f, 0.6f, 2.2f},
0104 {1.1f, 2.7f, 1.3f},
0105 {0.9f, 1.4f, 2.7f},
0106 {2.4f, 1.8f, 0.9f},
0107 {0.6f, 2.2f, 2.1f},
0108 }};
0109
0110 checkMagneticFieldEqual(actsField, cache, view, points, 0.0001);
0111 }
0112
0113 BOOST_AUTO_TEST_CASE(InterpolatedMagneticField2) {
0114 auto localToGlobalBin_xyz = [](std::array<std::size_t, 3> binsXYZ,
0115 std::array<std::size_t, 3> nBinsXYZ) {
0116 return (binsXYZ.at(0) * (nBinsXYZ.at(1) * nBinsXYZ.at(2)) +
0117 binsXYZ.at(1) * nBinsXYZ.at(2) + binsXYZ.at(2));
0118 };
0119
0120 std::vector<double> xPos = {8., 12., 16., 20.};
0121 std::vector<double> yPos = {8., 12., 16., 20.};
0122 std::vector<double> zPos = {8., 12., 16., 20.};
0123
0124 std::vector<Vector3> bField_xyz;
0125 for (int i = 0; i < 64; i++) {
0126 bField_xyz.push_back(Vector3(i, i * i * 0.01, i));
0127 }
0128
0129 MagneticFieldContext fieldContext;
0130 auto actsField = fieldMapXYZ(localToGlobalBin_xyz, xPos, yPos, zPos,
0131 bField_xyz, 1, 1, false);
0132 MagneticFieldProvider::Cache cache = actsField.makeCache(fieldContext);
0133
0134 Covfie::InterpolatedField field = Covfie::covfieField(actsField);
0135 typename Covfie::InterpolatedField::view_t view(field);
0136
0137 std::array<std::array<float, 3>, 14> points = {{
0138 {8.f, 8.f, 8.f},
0139 {12.f, 12.f, 12.f},
0140 {16.f, 16.f, 16.f},
0141 {19.9, 19.9, 19.9},
0142 {8.1f, 10.2f, 12.3f},
0143 {9.4f, 11.5f, 13.6f},
0144 {10.7f, 12.8f, 14.9f},
0145 {11.0f, 13.1f, 15.2f},
0146 {12.3f, 14.4f, 16.5f},
0147 {13.6f, 15.7f, 17.8f},
0148 {14.9f, 16.0f, 18.1f},
0149 {16.2f, 17.3f, 19.4f},
0150 {17.5f, 18.6f, 19.7f},
0151 {18.8f, 19.9f, 14.0f},
0152 }};
0153
0154 checkMagneticFieldEqual(actsField, cache, view, points, 0.0001f);
0155 }
0156
0157 BOOST_AUTO_TEST_CASE(ConstantMagneticField1) {
0158 ConstantBField actsField(Vector3{1.3f, 2.5f, 2.f});
0159 MagneticFieldContext ctx;
0160 MagneticFieldProvider::Cache cache = actsField.makeCache(ctx);
0161
0162 Covfie::ConstantField field = Covfie::covfieField(actsField);
0163 typename Covfie::ConstantField::view_t view(field);
0164
0165 std::array<std::array<float, 3>, 13> points = {{
0166 {8.f, 8.f, 8.f},
0167 {12.f, 12.f, 12.f},
0168 {16.f, 16.f, 16.f},
0169 {8.1f, 10.2f, 12.3f},
0170 {9.4f, 11.5f, 13.6f},
0171 {10.7f, 12.8f, 14.9f},
0172 {11.0f, 13.1f, 15.2f},
0173 {12.3f, 14.4f, 16.5f},
0174 {13.6f, 15.7f, 17.8f},
0175 {14.9f, 16.0f, 18.1f},
0176 {16.2f, 17.3f, 19.4f},
0177 {17.5f, 18.6f, 19.7f},
0178 {18.8f, 19.9f, 14.0f},
0179 }};
0180
0181 checkMagneticFieldEqual(actsField, cache, view, points, 0.0001f);
0182 }
0183
0184 BOOST_AUTO_TEST_CASE(SolenoidBField1) {
0185 SolenoidBField::Config cfg{};
0186 cfg.length = 5.8_m;
0187 cfg.radius = (2.56 + 2.46) * 0.5 * 0.5_m;
0188 cfg.nCoils = 1154;
0189 cfg.bMagCenter = 2_T;
0190 SolenoidBField actsField(cfg);
0191 MagneticFieldContext ctx;
0192 MagneticFieldProvider::Cache cache = actsField.makeCache(ctx);
0193
0194 Covfie::InterpolatedField field = Covfie::covfieField(
0195 actsField, cache, {21UL, 21UL, 21UL}, {0., 0., 0.}, {20., 20., 20.});
0196 typename Covfie::InterpolatedField::view_t view(field);
0197
0198 std::array<std::array<float, 3>, 13> points = {{
0199 {8.f, 8.f, 8.f},
0200 {12.f, 12.f, 12.f},
0201 {16.f, 16.f, 16.f},
0202 {8.1f, 10.2f, 12.3f},
0203 {9.4f, 11.5f, 13.6f},
0204 {10.7f, 12.8f, 14.9f},
0205 {11.0f, 13.1f, 15.2f},
0206 {12.3f, 14.4f, 16.5f},
0207 {13.6f, 15.7f, 17.8f},
0208 {14.9f, 16.0f, 18.1f},
0209 {16.2f, 17.3f, 19.4f},
0210 {17.5f, 18.6f, 19.7f},
0211 {18.8f, 19.9f, 14.0f},
0212 }};
0213
0214 checkMagneticFieldEqual(actsField, cache, view, points, 0.0001);
0215 }
0216
0217 BOOST_AUTO_TEST_SUITE_END()
0218
0219 }