Class Reference for E1039 Core & Analysis Software
SQField3DCartesian.cc
Go to the documentation of this file.
1 #include "SQField3DCartesian.h"
2 
3 #include <TFile.h>
4 #include <TTree.h>
5 
6 #include <iostream>
7 #include <fstream>
8 #include <sstream>
9 
10 SQField3DCartesian::SQField3DCartesian(const std::string& fname, const float magfield_rescale)
11  : filename(fname)
12  , fieldstr(magfield_rescale)
13  , xvals()
14  , yvals()
15  , zvals()
16  , xsteps(0)
17  , ysteps(0)
18  , zsteps(0)
19  , xmin(1000000)
20  , xmax(-1000000)
21  , ymin(1000000)
22  , ymax(-1000000)
23  , zmin(1000000)
24  , zmax(-1000000)
25  , xstepsize(-1.)
26  , ystepsize(-1.)
27  , zstepsize(-1.)
28 {
29  std::cout << "\n================ Begin Construct Mag Field =====================" << std::endl;
30  std::cout << "\n-----------------------------------------------------------"
31  << "\n Magnetic field Module - Verbosity:"
32  << "\n-----------------------------------------------------------\n";
33 
34  //load from input file first
35  fpoints.clear();
36  if(filename.find(".root") != std::string::npos)
37  {
38  TFile* inputFile = TFile::Open(filename.c_str());
39  if(!inputFile)
40  {
41  std::cout << "SQField3DCartesian: cannot find the input ROOT field map, abort " << filename << std::endl;
42  exit(EXIT_FAILURE);
43  }
44  std::cout << "SQField3DCartesian: loading field map from ROOT file " << filename << std::endl;
45 
46  float x, y, z, Bx, By, Bz;
47  TTree* inputTree = (TTree*)inputFile->Get("fieldmap");
48  inputTree->SetBranchAddress("x", &x);
49  inputTree->SetBranchAddress("y", &y);
50  inputTree->SetBranchAddress("z", &z);
51  inputTree->SetBranchAddress("Bx", &Bx);
52  inputTree->SetBranchAddress("By", &By);
53  inputTree->SetBranchAddress("Bz", &Bz);
54 
55  unsigned int nRecords = inputTree->GetEntries();
56  fpoints.reserve(nRecords);
57  for(unsigned int i = 0; i < nRecords; ++i)
58  {
59  inputTree->GetEntry(i);
60 
61  FieldPoint p;
62  p.x = x*cm; p.y = y*cm; p.z = z*cm;
63  p.B.SetXYZ(fieldstr*Bx*tesla, fieldstr*By*tesla, fieldstr*Bz*tesla);
64  fpoints.push_back(p);
65 
66  xvals.insert(p.x);
67  yvals.insert(p.y);
68  zvals.insert(p.z);
69  }
70 
71  inputFile->Close();
72  }
73  else //load from ascii input
74  {
75  std::ifstream fin(filename.c_str());
76  if(!fin)
77  {
78  std::cout << "SQField3DCartesian: cannot find the input ASCII field map, abort " << filename << std::endl;
79  exit(EXIT_FAILURE);
80  }
81  std::cout << "SQField3DCartesian: loading field map from ASCII file " << filename << std::endl;
82 
83  std::string line;
84  getline(fin, line); //header line
85  while(getline(fin, line))
86  {
87  float x, y, z, Bx, By, Bz;
88  std::stringstream ss(line);
89  ss >> x >> y >> z >> Bx >> By >> Bz;
90 
91  FieldPoint p;
92  p.x = x*cm; p.y = y*cm; p.z = z*cm;
93  p.B.SetXYZ(fieldstr*Bx*tesla, fieldstr*By*tesla, fieldstr*Bz*tesla);
94  fpoints.push_back(p);
95 
96  xvals.insert(p.x);
97  yvals.insert(p.y);
98  zvals.insert(p.z);
99  }
100 
101  fin.close();
102  }
103 
104  //Extract information on each dimension
105  xmin = *(xvals.begin());
106  xmax = *(xvals.rbegin());
107  ymin = *(yvals.begin());
108  ymax = *(yvals.rbegin());
109  zmin = *(zvals.begin());
110  zmax = *(zvals.rbegin());
111 
112  xsteps = xvals.size();
113  ysteps = yvals.size();
114  zsteps = zvals.size();
115 
116  xstepsize = (xmax - xmin)/(xsteps - 1);
117  ystepsize = (ymax - ymin)/(ysteps - 1);
118  zstepsize = (zmax - zmin)/(zsteps - 1);
119 
120  //NOTE here we assume the input files (text or root) are pre-sorted. Otherwise we need to sort the fpoints vector
121  //using something like this:
122  //std::sort(fpoints.begin(), fpoints.end(), [&](FieldPoint a, FieldPoint b)
123  //{ return (int((a.z-zmin)/zstepsize) + zsteps*(int((a.y-ymin)/ystepsize)) + ysteps*(int((a.x-xmin)/xstepsize))) <
124  // (int((b.z-zmin)/zstepsize) + zsteps*(int((b.y-ymin)/ystepsize)) + ysteps*(int((b.x-xmin)/xstepsize))); })
125 
126  std::cout << " its demensions and ranges of the field: " << fpoints.size() << std::endl;
127  std::cout << " X: " << xsteps << " bins, " << xmin/cm << " cm -- " << xmax/cm << " cm, step size = " << xstepsize/cm << " cm." << std::endl;
128  std::cout << " Y: " << ysteps << " bins, " << ymin/cm << " cm -- " << ymax/cm << " cm, step size = " << ystepsize/cm << " cm." << std::endl;
129  std::cout << " Z: " << zsteps << " bins, " << zmin/cm << " cm -- " << zmax/cm << " cm, step size = " << zstepsize/cm << " cm." << std::endl;
130 
131  std::cout << "\n================= End Construct Mag Field ======================\n" << std::endl;
132 }
133 
135 {}
136 
137 int SQField3DCartesian::GetGlobalIndex(int xIdx, int yIdx, int zIdx) const
138 {
139  return zIdx + zsteps*(yIdx + ysteps*xIdx);
140 }
141 
142 void SQField3DCartesian::GetFieldValue(const double point[4], double* Bfield) const
143 {
144  //This 3D intepolation algorithm is based on the wiki link http://en.wikipedia.org/wiki/Trilinear_interpolation
145 
146  for(int i = 0; i < 3; ++i) Bfield[i] = 0.;
147  double x = point[0];
148  double y = point[1];
149  double z = point[2];
150 
151  int ubx = int((xsteps-1)*(x - xmin)/(xmax - xmin));
152  int uby = int((ysteps-1)*(y - ymin)/(ymax - ymin));
153  int ubz = int((zsteps-1)*(z - zmin)/(zmax - zmin));
154  int obx = ubx + 1;
155  int oby = uby + 1;
156  int obz = ubz + 1;
157 
158  if(ubx < 0 || uby < 0 || ubz < 0 ||
159  obx >= xsteps || oby >= ysteps || obz >= zsteps)
160  {
161  /*
162  static int ifirst = 0;
163  if(ifirst < 10)
164  {
165  ++ifirst;
166  std::cout << "SQField3DCartesian::GetFieldValue: WARNING " << filename << " "
167  << "Out-of-boundary coordinates: "
168  << "x: " << x/cm
169  << ", y: " << y/cm
170  << ", z: " << z/cm
171  << " bailing out returning zero bfield" << std::endl;
172  }
173  */
174  return;
175  }
176 
177  int idx000 = GetGlobalIndex(ubx, uby, ubz);
178  int idx001 = GetGlobalIndex(ubx, uby, obz);
179  int idx010 = GetGlobalIndex(ubx, oby, ubz);
180  int idx011 = GetGlobalIndex(ubx, oby, obz);
181  int idx100 = GetGlobalIndex(obx, uby, ubz);
182  int idx101 = GetGlobalIndex(obx, uby, obz);
183  int idx110 = GetGlobalIndex(obx, oby, ubz);
184  int idx111 = GetGlobalIndex(obx, oby, obz);
185 
186  double xp = (x - fpoints[idx000].x)/xstepsize;
187  double yp = (y - fpoints[idx000].y)/ystepsize;
188  double zp = (z - fpoints[idx000].z)/zstepsize;
189 
190  TVector3 i1 = fpoints[idx000].B*(1. - zp) + fpoints[idx001].B*zp;
191  TVector3 i2 = fpoints[idx010].B*(1. - zp) + fpoints[idx011].B*zp;
192  TVector3 j1 = fpoints[idx100].B*(1. - zp) + fpoints[idx101].B*zp;
193  TVector3 j2 = fpoints[idx110].B*(1. - zp) + fpoints[idx111].B*zp;
194 
195  TVector3 w1 = i1*(1. - yp) + i2*yp;
196  TVector3 w2 = j1*(1. - yp) + j2*yp;
197 
198  TVector3 res = w1*(1. - xp) + w2*xp;
199  Bfield[0] = res.X();
200  Bfield[1] = res.Y();
201  Bfield[2] = res.Z();
202 }
std::set< double > zvals
std::vector< FieldPoint > fpoints
void GetFieldValue(const double Point[4], double *Bfield) const
int GetGlobalIndex(int xIdx, int yIdx, int zIdx) const
return the index of the global index based on the local index
std::set< double > xvals
std::set< double > yvals
SQField3DCartesian(const std::string &fname, const float magfield_rescale=1.0)