Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:11:45

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 #include "ActsExamples/MagneticField/FieldMapTextIo.hpp"
0010 
0011 #include "Acts/MagneticField/BFieldMapUtils.hpp"
0012 
0013 #include <fstream>
0014 #include <vector>
0015 
0016 namespace {
0017 constexpr std::size_t kDefaultSize = 1 << 15;
0018 }
0019 
0020 ActsExamples::detail::InterpolatedMagneticField2
0021 ActsExamples::makeMagneticFieldMapRzFromText(
0022     const std::function<std::size_t(std::array<std::size_t, 2> binsRZ,
0023                                     std::array<std::size_t, 2> nBinsRZ)>&
0024         localToGlobalBin,
0025     const std::string& fieldMapFile, double lengthUnit, double BFieldUnit,
0026     bool firstQuadrant) {
0027   /// [1] Read in field map file
0028   // Grid position points in r and z
0029   std::vector<double> rPos;
0030   std::vector<double> zPos;
0031   // components of magnetic field on grid points
0032   std::vector<Acts::Vector2> bField;
0033   // reserve estimated size
0034   rPos.reserve(kDefaultSize);
0035   zPos.reserve(kDefaultSize);
0036   bField.reserve(kDefaultSize);
0037   // [1] Read in file and fill values
0038   std::ifstream map_file(fieldMapFile.c_str(), std::ios::in);
0039   std::string line;
0040   double r = 0., z = 0.;
0041   double br = 0., bz = 0.;
0042   while (std::getline(map_file, line)) {
0043     if (line.empty() || line[0] == '%' || line[0] == '#' ||
0044         line.find_first_not_of(' ') == std::string::npos) {
0045       continue;
0046     }
0047 
0048     std::istringstream tmp(line);
0049     tmp >> r >> z >> br >> bz;
0050     rPos.push_back(r);
0051     zPos.push_back(z);
0052     bField.push_back(Acts::Vector2(br, bz));
0053   }
0054   map_file.close();
0055   rPos.shrink_to_fit();
0056   zPos.shrink_to_fit();
0057   bField.shrink_to_fit();
0058   /// [2] use helper function in core
0059   return Acts::fieldMapRZ(localToGlobalBin, rPos, zPos, bField, lengthUnit,
0060                           BFieldUnit, firstQuadrant);
0061 }
0062 
0063 ActsExamples::detail::InterpolatedMagneticField3
0064 ActsExamples::makeMagneticFieldMapXyzFromText(
0065     const std::function<std::size_t(std::array<std::size_t, 3> binsXYZ,
0066                                     std::array<std::size_t, 3> nBinsXYZ)>&
0067         localToGlobalBin,
0068     const std::string& fieldMapFile, double lengthUnit, double BFieldUnit,
0069     bool firstOctant) {
0070   /// [1] Read in field map file
0071   // Grid position points in x, y and z
0072   std::vector<double> xPos;
0073   std::vector<double> yPos;
0074   std::vector<double> zPos;
0075   // components of magnetic field on grid points
0076   std::vector<Acts::Vector3> bField;
0077   // reserve estimated size
0078   xPos.reserve(kDefaultSize);
0079   yPos.reserve(kDefaultSize);
0080   zPos.reserve(kDefaultSize);
0081   bField.reserve(kDefaultSize);
0082   // [1] Read in file and fill values
0083   std::ifstream map_file(fieldMapFile.c_str(), std::ios::in);
0084   std::string line;
0085   double x = 0., y = 0., z = 0.;
0086   double bx = 0., by = 0., bz = 0.;
0087   while (std::getline(map_file, line)) {
0088     if (line.empty() || line[0] == '%' || line[0] == '#' ||
0089         line.find_first_not_of(' ') == std::string::npos) {
0090       continue;
0091     }
0092 
0093     std::istringstream tmp(line);
0094     tmp >> x >> y >> z >> bx >> by >> bz;
0095     xPos.push_back(x);
0096     yPos.push_back(y);
0097     zPos.push_back(z);
0098     bField.push_back(Acts::Vector3(bx, by, bz));
0099   }
0100   map_file.close();
0101   xPos.shrink_to_fit();
0102   yPos.shrink_to_fit();
0103   zPos.shrink_to_fit();
0104   bField.shrink_to_fit();
0105   return Acts::fieldMapXYZ(localToGlobalBin, xPos, yPos, zPos, bField,
0106                            lengthUnit, BFieldUnit, firstOctant);
0107 }