|
||||
File indexing completed on 2025-01-18 09:27:39
0001 // This file is part of the Acts project. 0002 // 0003 // Copyright (C) 2017-2018 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 http://mozilla.org/MPL/2.0/. 0008 0009 #pragma once 0010 0011 #include "Acts/Definitions/Algebra.hpp" 0012 #include "Acts/MagneticField/MagneticFieldContext.hpp" 0013 #include "Acts/MagneticField/MagneticFieldProvider.hpp" 0014 #include "Acts/Utilities/Result.hpp" 0015 0016 #include <cstddef> 0017 #include <functional> 0018 0019 namespace Acts { 0020 0021 /// @ingroup MagneticField 0022 // 0023 /// @class SolenoidBField 0024 /// Implements a multi-coil solenoid magnetic field. On every call, the field 0025 /// is evaluated at that exact position. The field has radially symmetry, the 0026 /// field vectors point in +z direction. 0027 /// The config exposes a target field value in the center. This value is used 0028 /// to empirically determine a scale factor which reproduces this field value 0029 /// in the center. 0030 /// 0031 /// E_1(k^2) = complete elliptic integral of the 1st kind 0032 /// E_2(k^2) = complete elliptic integral of the 2nd kind 0033 /// 0034 /// E_1(k^2) and E_2(k^2) are usually indicated as K(k^2) and E(k^2) in 0035 /// literature, 0036 /// respectively 0037 /// _ 0038 /// 2 / pi / 2 2 2 - 1 / 2 0039 /// E (k ) = | ( 1 - k sin {theta} ) dtheta 0040 /// 1 _/ 0 0041 /// 0042 /// _ ____________________ 0043 /// 2 / pi / 2| / 2 2 0044 /// E (k ) = | |/ 1 - k sin {theta} dtheta 0045 /// 2 _/ 0 0046 /// 0047 /// k^2 = is a function of the point (r, z) and of the radius of the coil R 0048 /// 0049 /// 2 4Rr 0050 /// k = --------------- 0051 /// 2 2 0052 /// (R + r) + z 0053 /// Using these, you can evaluate the two components B_r and B_z of the 0054 /// magnetic field: 0055 /// _ _ 0056 /// mu I | / 2 \ | 0057 /// 0 kz | |2 - k | 2 2 | 0058 /// B (r, z) = ----- ------ | |-------|E (k ) - E (k ) | 0059 /// r 4pi ___ | | 2| 2 1 | 0060 /// | / 3 |_ \2 - 2k / _| 0061 /// |/ Rr 0062 /// 0063 /// _ _ 0064 /// mu I | / 2 \ | 0065 /// 0 k | | (R + r)k - 2r | 2 2 | 0066 /// B (r,z) = ----- ---- | | -------------- | E (k ) + E (k ) | 0067 /// z 4pi __ | | 2 | 2 1 | 0068 /// |/Rr |_ \ 2r(1 - k ) / _| 0069 /// 0070 class SolenoidBField final : public MagneticFieldProvider { 0071 public: 0072 struct Cache { 0073 /// @brief Constructor with magnetic field context 0074 Cache(const MagneticFieldContext& /*mctx*/) {} 0075 }; 0076 0077 /// Config struct for the SolenoidBfield. 0078 struct Config { 0079 /// Radius at which the coils are located. 0080 double radius; 0081 /// Extent of the solenoid in z. It goes from 0082 /// -length/2 to +length/2 by convention 0083 double length; 0084 /// The number of coils that make up the solenoid 0085 std::size_t nCoils; 0086 /// The target magnetic field strength at the center. 0087 /// This will be used to scale coefficients 0088 double bMagCenter; 0089 }; 0090 0091 /// @brief the constructor with a shared pointer 0092 /// @note since it is a shared field, we enforce it to be const 0093 /// @tparam bField is the shared BField to be stored 0094 SolenoidBField(Config config); 0095 0096 /// @brief Retrieve magnetic field value in local (r,z) coordinates 0097 /// 0098 /// @param [in] position local 2D position 0099 Vector2 getField(const Vector2& position) const; 0100 0101 /// @copydoc MagneticFieldProvider::makeCache(const MagneticFieldContext&) const 0102 MagneticFieldProvider::Cache makeCache( 0103 const MagneticFieldContext& mctx) const override; 0104 0105 /// @brief Get the B field at a position 0106 /// 0107 /// @param position The position to query at 0108 Vector3 getField(const Vector3& position) const; 0109 0110 /// @copydoc MagneticFieldProvider::getField(const Vector3&,MagneticFieldProvider::Cache&) const 0111 Result<Vector3> getField(const Vector3& position, 0112 MagneticFieldProvider::Cache& cache) const override; 0113 0114 /// @copydoc MagneticFieldProvider::getFieldGradient(const Vector3&,ActsMatrix<3,3>&,MagneticFieldProvider::Cache&) const 0115 /// 0116 /// @note currently the derivative is not calculated 0117 /// @todo return derivative 0118 Result<Vector3> getFieldGradient( 0119 const Vector3& position, ActsMatrix<3, 3>& derivative, 0120 MagneticFieldProvider::Cache& cache) const override; 0121 0122 private: 0123 Config m_cfg; 0124 double m_scale; 0125 double m_dz; 0126 double m_R2; 0127 0128 Vector2 multiCoilField(const Vector2& pos, double scale) const; 0129 0130 Vector2 singleCoilField(const Vector2& pos, double scale) const; 0131 0132 double B_r(const Vector2& pos, double scale) const; 0133 0134 double B_z(const Vector2& pos, double scale) const; 0135 0136 double k2(double r, double z) const; 0137 }; 0138 0139 } // namespace Acts
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |