File indexing completed on 2025-10-25 07:56:55
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include <boost/test/unit_test.hpp>
0010
0011 #include "Acts/Definitions/Algebra.hpp"
0012 #include "Acts/Definitions/Units.hpp"
0013 #include "Acts/MagneticField/MagneticFieldContext.hpp"
0014 #include "Acts/MagneticField/SolenoidBField.hpp"
0015 #include "Acts/Utilities/Result.hpp"
0016 #include "ActsTests/CommonHelpers/FloatComparisons.hpp"
0017
0018 #include <cstddef>
0019
0020 using namespace Acts;
0021 using namespace Acts::UnitLiterals;
0022
0023 namespace ActsTests {
0024
0025 BOOST_AUTO_TEST_SUITE(MagneticFieldSuite)
0026
0027 BOOST_AUTO_TEST_CASE(TestSolenoidBField) {
0028
0029 MagneticFieldContext mfContext = MagneticFieldContext();
0030
0031 SolenoidBField::Config cfg{};
0032 cfg.length = 5.8_m;
0033 cfg.radius = (2.56 + 2.46) * 0.5 * 0.5_m;
0034 cfg.nCoils = 1154;
0035 cfg.bMagCenter = 2_T;
0036 SolenoidBField bField(cfg);
0037
0038 auto cache = bField.makeCache(mfContext);
0039 CHECK_CLOSE_ABS(bField.getField({0, 0, 0}, cache).value(),
0040 Vector3(0, 0, 2.0_T), 1e-6_T);
0041
0042
0043
0044
0045 double tol = 1e-6;
0046 double tol_B = 1e-6_T;
0047 std::size_t steps = 20;
0048 for (std::size_t i = 0; i < steps; i++) {
0049 double r = 1.5 * cfg.radius / steps * i;
0050 BOOST_TEST_CONTEXT("r=" << r) {
0051 Vector3 B1 = bField.getField({r, 0, 0}, cache).value();
0052 Vector3 B2 = bField.getField({-r, 0, 0}, cache).value();
0053 CHECK_SMALL(B1.x(), tol);
0054 CHECK_SMALL(B1.y(), tol);
0055 BOOST_CHECK_GT(std::abs(B1.z()), tol_B);
0056
0057 CHECK_CLOSE_ABS(B1, B2, tol_B);
0058
0059
0060 for (std::size_t j = 0; j <= steps; j++) {
0061
0062 double z = (1.5 * cfg.length / 2.) / steps * j;
0063 BOOST_TEST_CONTEXT("z=" << z) {
0064 Vector3 B_zp_rp = bField.getField({r, 0, z}, cache).value();
0065 Vector3 B_zn_rp = bField.getField({r, 0, -z}, cache).value();
0066 Vector3 B_zp_rn = bField.getField({-r, 0, z}, cache).value();
0067 Vector3 B_zn_rn = bField.getField({-r, 0, -z}, cache).value();
0068
0069
0070
0071
0072
0073
0074
0075
0076
0077
0078
0079
0080
0081
0082
0083
0084
0085 BOOST_CHECK_GT(std::abs(B_zp_rp.z()), tol_B);
0086 BOOST_CHECK_GT(std::abs(B_zn_rp.z()), tol_B);
0087 BOOST_CHECK_GT(std::abs(B_zn_rn.z()), tol_B);
0088 BOOST_CHECK_GT(std::abs(B_zp_rn.z()), tol_B);
0089 if (i > 0) {
0090
0091 CHECK_CLOSE_ABS(B_zp_rp.z(), B_zp_rn.z(), tol_B);
0092 CHECK_CLOSE_ABS(B_zn_rp.z(), B_zn_rn.z(), tol_B);
0093
0094 CHECK_CLOSE_ABS(B_zp_rp.x(), -B_zp_rn.x(), tol_B);
0095 CHECK_CLOSE_ABS(B_zn_rp.x(), -B_zn_rn.x(), tol_B);
0096 }
0097 if (j > 0) {
0098
0099 CHECK_CLOSE_ABS(B_zp_rp.z(), B_zn_rp.z(), tol_B);
0100 CHECK_CLOSE_ABS(B_zp_rn.z(), B_zn_rn.z(), tol_B);
0101
0102 CHECK_CLOSE_ABS(B_zp_rp.x(), -B_zn_rp.x(), tol_B);
0103 CHECK_CLOSE_ABS(B_zp_rn.x(), -B_zn_rn.x(), tol_B);
0104 }
0105 }
0106 }
0107 }
0108 }
0109
0110 }
0111
0112 BOOST_AUTO_TEST_SUITE_END()
0113
0114 }