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