Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:13:07

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 // Acts include(s)
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 // Covfie Plugin include(s)
0017 #include "Acts/Plugins/Covfie/FieldConversion.hpp"
0018 
0019 // System include(s)
0020 #include <array>
0021 #include <cmath>
0022 #include <sstream>
0023 #include <vector>
0024 
0025 // Boost include(s)
0026 #include <boost/test/unit_test.hpp>
0027 
0028 using namespace Acts::UnitLiterals;
0029 
0030 template <typename view_t, typename iterator_t>
0031 void checkMagneticFieldEqual(const Acts::MagneticFieldProvider& fieldProvider,
0032                              Acts::MagneticFieldProvider::Cache& cache,
0033                              view_t view, iterator_t points,
0034                              float error_margin_half_width) {
0035   for (auto point : points) {
0036     auto x = point[0], y = point[1], z = point[2];
0037 
0038     auto lookupResult = fieldProvider.getField(Acts::Vector3{x, y, z}, cache);
0039     if (!lookupResult.ok()) {
0040       throw std::runtime_error{"Field lookup failure"};
0041     }
0042     auto actsValueX = (*lookupResult)[0], actsValueY = (*lookupResult)[1],
0043          actsValueZ = (*lookupResult)[2];
0044 
0045     auto covfieValues = view.at(x, y, z);
0046     auto covfieValueX = covfieValues[0], covfieValueY = covfieValues[1],
0047          covfieValueZ = covfieValues[2];
0048 
0049     auto isEqual =
0050         std::abs(covfieValueX - actsValueX) <= error_margin_half_width &&
0051         std::abs(covfieValueY - actsValueY) <= error_margin_half_width &&
0052         std::abs(covfieValueZ - actsValueZ) <= error_margin_half_width;
0053 
0054     std::stringstream ss;
0055     ss << "Fields are not equal at position (" << x << ", " << y << ", " << z
0056        << "). Acts: (" << actsValueX << ", " << actsValueY << ", " << actsValueZ
0057        << "), Covfie: (" << covfieValueX << ", " << covfieValueY << ", "
0058        << covfieValueZ << ")" << std::endl;
0059 
0060     BOOST_CHECK_MESSAGE(isEqual, ss.str());
0061   }
0062 }
0063 
0064 BOOST_AUTO_TEST_SUITE(CovfiePlugin)
0065 
0066 BOOST_AUTO_TEST_CASE(InterpolatedMagneticField1) {
0067   auto localToGlobalBin_xyz = [](std::array<std::size_t, 3> binsXYZ,
0068                                  std::array<std::size_t, 3> nBinsXYZ) {
0069     return (binsXYZ.at(0) * (nBinsXYZ.at(1) * nBinsXYZ.at(2)) +
0070             binsXYZ.at(1) * nBinsXYZ.at(2) + binsXYZ.at(2));
0071   };
0072 
0073   std::vector<double> xPos = {0., 1., 2., 3.};
0074   std::vector<double> yPos = {0., 1., 2., 3.};
0075   std::vector<double> zPos = {0., 1., 2., 3.};
0076 
0077   std::vector<Acts::Vector3> bField_xyz;
0078   for (int i = 0; i < 64; i++) {
0079     bField_xyz.push_back(Acts::Vector3(i, i, i));
0080   }
0081 
0082   Acts::MagneticFieldContext fieldContext;
0083   auto actsField = Acts::fieldMapXYZ(localToGlobalBin_xyz, xPos, yPos, zPos,
0084                                      bField_xyz, 1, 1, false);
0085   Acts::MagneticFieldProvider::Cache cache = actsField.makeCache(fieldContext);
0086 
0087   Acts::CovfiePlugin::InterpolatedField field =
0088       Acts::CovfiePlugin::covfieField(actsField);
0089   typename Acts::CovfiePlugin::InterpolatedField::view_t view(field);
0090 
0091   std::array<std::array<float, 3>, 14> points = {{
0092       {0.f, 0.f, 0.f},
0093       {1.f, 1.f, 1.f},
0094       {2.f, 2.f, 2.f},
0095       {2.9f, 2.9f, 2.9f},
0096       {1.2f, 2.5f, 0.8f},
0097       {0.7f, 1.9f, 2.3f},
0098       {2.1f, 0.3f, 1.5f},
0099       {0.4f, 2.8f, 2.9f},
0100       {1.6f, 1.2f, 0.5f},
0101       {2.3f, 0.6f, 2.2f},
0102       {1.1f, 2.7f, 1.3f},
0103       {0.9f, 1.4f, 2.7f},
0104       {2.4f, 1.8f, 0.9f},
0105       {0.6f, 2.2f, 2.1f},
0106   }};
0107 
0108   checkMagneticFieldEqual(actsField, cache, view, points, 0.0001);
0109 }
0110 
0111 BOOST_AUTO_TEST_CASE(InterpolatedMagneticField2) {
0112   auto localToGlobalBin_xyz = [](std::array<std::size_t, 3> binsXYZ,
0113                                  std::array<std::size_t, 3> nBinsXYZ) {
0114     return (binsXYZ.at(0) * (nBinsXYZ.at(1) * nBinsXYZ.at(2)) +
0115             binsXYZ.at(1) * nBinsXYZ.at(2) + binsXYZ.at(2));
0116   };
0117 
0118   std::vector<double> xPos = {8., 12., 16., 20.};
0119   std::vector<double> yPos = {8., 12., 16., 20.};
0120   std::vector<double> zPos = {8., 12., 16., 20.};
0121 
0122   std::vector<Acts::Vector3> bField_xyz;
0123   for (int i = 0; i < 64; i++) {
0124     bField_xyz.push_back(Acts::Vector3(i, i * i * 0.01, i));
0125   }
0126 
0127   Acts::MagneticFieldContext fieldContext;
0128   auto actsField = Acts::fieldMapXYZ(localToGlobalBin_xyz, xPos, yPos, zPos,
0129                                      bField_xyz, 1, 1, false);
0130   Acts::MagneticFieldProvider::Cache cache = actsField.makeCache(fieldContext);
0131 
0132   Acts::CovfiePlugin::InterpolatedField field =
0133       Acts::CovfiePlugin::covfieField(actsField);
0134   typename Acts::CovfiePlugin::InterpolatedField::view_t view(field);
0135 
0136   std::array<std::array<float, 3>, 14> points = {{
0137       {8.f, 8.f, 8.f},
0138       {12.f, 12.f, 12.f},
0139       {16.f, 16.f, 16.f},
0140       {19.9, 19.9, 19.9},
0141       {8.1f, 10.2f, 12.3f},
0142       {9.4f, 11.5f, 13.6f},
0143       {10.7f, 12.8f, 14.9f},
0144       {11.0f, 13.1f, 15.2f},
0145       {12.3f, 14.4f, 16.5f},
0146       {13.6f, 15.7f, 17.8f},
0147       {14.9f, 16.0f, 18.1f},
0148       {16.2f, 17.3f, 19.4f},
0149       {17.5f, 18.6f, 19.7f},
0150       {18.8f, 19.9f, 14.0f},
0151   }};
0152 
0153   checkMagneticFieldEqual(actsField, cache, view, points, 0.0001f);
0154 }
0155 
0156 BOOST_AUTO_TEST_CASE(ConstantMagneticField1) {
0157   Acts::ConstantBField actsField(Acts::Vector3{1.3f, 2.5f, 2.f});
0158   Acts::MagneticFieldContext ctx;
0159   Acts::MagneticFieldProvider::Cache cache = actsField.makeCache(ctx);
0160 
0161   Acts::CovfiePlugin::ConstantField field =
0162       Acts::CovfiePlugin::covfieField(actsField);
0163   typename Acts::CovfiePlugin::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   Acts::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   Acts::SolenoidBField actsField(cfg);
0191   Acts::MagneticFieldContext ctx;
0192   Acts::MagneticFieldProvider::Cache cache = actsField.makeCache(ctx);
0193 
0194   Acts::CovfiePlugin::InterpolatedField field = Acts::CovfiePlugin::covfieField(
0195       actsField, cache, {21UL, 21UL, 21UL}, {0., 0., 0.}, {20., 20., 20.});
0196   typename Acts::CovfiePlugin::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()