ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BFieldUtils.cpp
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file BFieldUtils.cpp
1 // This file is part of the Acts project.
2 //
3 // Copyright (C) 2017 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 
11 #include <fstream>
12 
16 #include "TFile.h"
17 #include "TROOT.h"
18 #include "TTree.h"
19 
24  std::function<size_t(std::array<size_t, 2> binsRZ,
25  std::array<size_t, 2> nBinsRZ)>
26  localToGlobalBin,
27  std::string fieldMapFile, double lengthUnit, double BFieldUnit,
28  size_t nPoints, bool firstQuadrant) {
30  // Grid position points in r and z
31  std::vector<double> rPos;
32  std::vector<double> zPos;
33  // components of magnetic field on grid points
34  std::vector<Acts::Vector2D> bField;
35  // reserve estimated size
36  rPos.reserve(nPoints);
37  zPos.reserve(nPoints);
38  bField.reserve(nPoints);
39  // [1] Read in file and fill values
40  std::ifstream map_file(fieldMapFile.c_str(), std::ios::in);
41  std::string line;
42  double r = 0., z = 0.;
43  double br = 0., bz = 0.;
44  while (std::getline(map_file, line)) {
45  if (line.empty() || line[0] == '%' || line[0] == '#' ||
46  line.find_first_not_of(' ') == std::string::npos)
47  continue;
48 
49  std::istringstream tmp(line);
50  tmp >> r >> z >> br >> bz;
51  rPos.push_back(r);
52  zPos.push_back(z);
53  bField.push_back(Acts::Vector2D(br, bz));
54  }
55  map_file.close();
57  return Acts::fieldMapperRZ(localToGlobalBin, rPos, zPos, bField, lengthUnit,
58  BFieldUnit, firstQuadrant);
59 }
60 
63  Acts::detail::EquidistantAxis, Acts::detail::EquidistantAxis>>
65  std::function<size_t(std::array<size_t, 3> binsXYZ,
66  std::array<size_t, 3> nBinsXYZ)>
67  localToGlobalBin,
68  std::string fieldMapFile, double lengthUnit, double BFieldUnit,
69  size_t nPoints, 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::Vector3D> bField;
77  // reserve estimated size
78  xPos.reserve(nPoints);
79  yPos.reserve(nPoints);
80  zPos.reserve(nPoints);
81  bField.reserve(nPoints);
82  // [1] Read in file and fill values
83  std::ifstream map_file(fieldMapFile.c_str(), std::ios::in);
84  std::string line;
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  std::istringstream tmp(line);
93  tmp >> x >> y >> z >> bx >> by >> bz;
94  xPos.push_back(x);
95  yPos.push_back(y);
96  zPos.push_back(z);
97  bField.push_back(Acts::Vector3D(bx, by, bz));
98  }
99  map_file.close();
100 
101  return Acts::fieldMapperXYZ(localToGlobalBin, xPos, yPos, zPos, bField,
102  lengthUnit, BFieldUnit, firstOctant);
103 }
104 
107  Acts::detail::EquidistantAxis>>
109  std::function<size_t(std::array<size_t, 2> binsRZ,
110  std::array<size_t, 2> nBinsRZ)>
111  localToGlobalBin,
112  std::string fieldMapFile, std::string treeName, double lengthUnit,
113  double BFieldUnit, bool firstQuadrant) {
115  // Grid position points in r and z
116  std::vector<double> rPos;
117  std::vector<double> zPos;
118  // components of magnetic field on grid points
119  std::vector<Acts::Vector2D> bField;
120  // [1] Read in file and fill values
121  TFile* inputFile = TFile::Open(fieldMapFile.c_str());
122  TTree* tree = (TTree*)inputFile->Get(treeName.c_str());
123  Int_t entries = tree->GetEntries();
124 
125  double r, z;
126  double Br, Bz;
127 
128  tree->SetBranchAddress("r", &r);
129  tree->SetBranchAddress("z", &z);
130 
131  tree->SetBranchAddress("Br", &Br);
132  tree->SetBranchAddress("Bz", &Bz);
133 
134  // reserve size
135  rPos.reserve(entries);
136  zPos.reserve(entries);
137  bField.reserve(entries);
138 
139  for (int i = 0; i < entries; i++) {
140  tree->GetEvent(i);
141  rPos.push_back(r);
142  zPos.push_back(z);
143  bField.push_back(Acts::Vector2D(Br, Bz));
144  }
145  inputFile->Close();
147  return Acts::fieldMapperRZ(localToGlobalBin, rPos, zPos, bField, lengthUnit,
148  BFieldUnit, firstQuadrant);
149 }
150 
153  Acts::detail::EquidistantAxis, Acts::detail::EquidistantAxis>>
155  std::function<size_t(std::array<size_t, 3> binsXYZ,
156  std::array<size_t, 3> nBinsXYZ)>
157  localToGlobalBin,
158  std::string fieldMapFile, std::string treeName, double lengthUnit,
159  double BFieldUnit, bool firstOctant) {
161  // Grid position points in x, y and z
162  std::vector<double> xPos;
163  std::vector<double> yPos;
164  std::vector<double> zPos;
165  // components of magnetic field on grid points
166  std::vector<Acts::Vector3D> bField;
167  // [1] Read in file and fill values
168  TFile* inputFile = TFile::Open(fieldMapFile.c_str());
169  TTree* tree = (TTree*)inputFile->Get(treeName.c_str());
170  Int_t entries = tree->GetEntries();
171 
172  double x, y, z;
173  double Bx, By, Bz;
174 
175  tree->SetBranchAddress("x", &x);
176  tree->SetBranchAddress("y", &y);
177  tree->SetBranchAddress("z", &z);
178 
179  tree->SetBranchAddress("Bx", &Bx);
180  tree->SetBranchAddress("By", &By);
181  tree->SetBranchAddress("Bz", &Bz);
182 
183  // reserve size
184  xPos.reserve(entries);
185  yPos.reserve(entries);
186  zPos.reserve(entries);
187  bField.reserve(entries);
188 
189  for (int i = 0; i < entries; i++) {
190  tree->GetEvent(i);
191  xPos.push_back(x);
192  yPos.push_back(y);
193  zPos.push_back(z);
194  bField.push_back(Acts::Vector3D(Bx, By, Bz));
195  }
196  inputFile->Close();
197 
198  return Acts::fieldMapperXYZ(localToGlobalBin, xPos, yPos, zPos, bField,
199  lengthUnit, BFieldUnit, firstOctant);
200 }