Analysis Software
Documentation for sPHENIX simulation software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
FieldMapTextIo.cpp
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file FieldMapTextIo.cpp
1 // This file is part of the Acts project.
2 //
3 // Copyright (C) 2017-2020 CERN for the benefit of the Acts project
4 //
5 // This Source Code Form is subject to the terms of the Mozilla Public
6 // License, v. 2.0. If a copy of the MPL was not distributed with this
7 // file, You can obtain one at http://mozilla.org/MPL/2.0/.
8 
10 
12 
13 #include <fstream>
14 #include <vector>
15 
16 namespace {
17 constexpr size_t kDefaultSize = 1 << 15;
18 }
19 
22  const std::function<size_t(std::array<size_t, 2> binsRZ,
23  std::array<size_t, 2> nBinsRZ)>&
24  localToGlobalBin,
25  const std::string& fieldMapFile, Acts::ActsScalar lengthUnit,
26  Acts::ActsScalar BFieldUnit, bool firstQuadrant) {
28  // Grid position points in r and z
29  std::vector<double> rPos;
30  std::vector<double> zPos;
31  // components of magnetic field on grid points
32  std::vector<Acts::Vector2> bField;
33  // reserve estimated size
34  rPos.reserve(kDefaultSize);
35  zPos.reserve(kDefaultSize);
36  bField.reserve(kDefaultSize);
37  // [1] Read in file and fill values
38  std::ifstream map_file(fieldMapFile.c_str(), std::ios::in);
40  double r = 0., z = 0.;
41  double br = 0., bz = 0.;
42  while (std::getline(map_file, line)) {
43  if (line.empty() || line[0] == '%' || line[0] == '#' ||
44  line.find_first_not_of(' ') == std::string::npos) {
45  continue;
46  }
47 
48  std::istringstream tmp(line);
49  tmp >> r >> z >> br >> bz;
50  rPos.push_back(r);
51  zPos.push_back(z);
52  bField.push_back(Acts::Vector2(br, bz));
53  }
54  map_file.close();
55  rPos.shrink_to_fit();
56  zPos.shrink_to_fit();
57  bField.shrink_to_fit();
59  return Acts::fieldMapRZ(localToGlobalBin, rPos, zPos, bField, lengthUnit,
60  BFieldUnit, firstQuadrant);
61 }
62 
65  const std::function<size_t(std::array<size_t, 3> binsXYZ,
66  std::array<size_t, 3> nBinsXYZ)>&
67  localToGlobalBin,
68  const std::string& fieldMapFile, Acts::ActsScalar lengthUnit,
69  Acts::ActsScalar BFieldUnit, bool firstOctant) {
71  // Grid position points in x, y and z
72  std::vector<double> xPos;
73  std::vector<double> yPos;
74  std::vector<double> zPos;
75  // components of magnetic field on grid points
76  std::vector<Acts::Vector3> bField;
77  // reserve estimated size
78  xPos.reserve(kDefaultSize);
79  yPos.reserve(kDefaultSize);
80  zPos.reserve(kDefaultSize);
81  bField.reserve(kDefaultSize);
82  // [1] Read in file and fill values
83  std::ifstream map_file(fieldMapFile.c_str(), std::ios::in);
85  double x = 0., y = 0., z = 0.;
86  double bx = 0., by = 0., bz = 0.;
87  while (std::getline(map_file, line)) {
88  if (line.empty() || line[0] == '%' || line[0] == '#' ||
89  line.find_first_not_of(' ') == std::string::npos) {
90  continue;
91  }
92 
93  std::istringstream tmp(line);
94  tmp >> x >> y >> z >> bx >> by >> bz;
95  xPos.push_back(x);
96  yPos.push_back(y);
97  zPos.push_back(z);
98  bField.push_back(Acts::Vector3(bx, by, bz));
99  }
100  map_file.close();
101  xPos.shrink_to_fit();
102  yPos.shrink_to_fit();
103  zPos.shrink_to_fit();
104  bField.shrink_to_fit();
105  return Acts::fieldMapXYZ(localToGlobalBin, xPos, yPos, zPos, bField,
106  lengthUnit, BFieldUnit, firstOctant);
107 }