Class Reference for E1039 Core & Analysis Software
PHField3DCartesian.cc
Go to the documentation of this file.
1 #include "PHField3DCartesian.h"
2 
3 #include <TFile.h>
4 #include <TNtuple.h>
5 
6 #include <boost/tuple/tuple.hpp>
7 #include <boost/tuple/tuple_comparison.hpp>
8 
9 #include <iostream>
10 #include <cassert>
11 #include <fstream>
12 
13 using namespace std;
14 using namespace CLHEP; // units
15 
16 #define UNIT_LENGTH cm
17 #define UNIT_FIELD tesla
18 
19 //#define _DEBUG_ON
20 
21 namespace {
22  typedef boost::tuple<double, double, double> trio;
23 }
24 
25 PHField3DCartesian::PHField3DCartesian(const string &fname, const float magfield_rescale)
26  : filename(fname)
27  , xvals()
28  , yvals()
29  , zvals()
30  , xmin(1000000)
31  , xmax(-1000000)
32  , ymin(1000000)
33  , ymax(-1000000)
34  , zmin(1000000)
35  , zmax(-1000000)
36  , xstepsize(NAN)
37  , ystepsize(NAN)
38  , zstepsize(NAN)
39  , xkey_save(NAN)
40  , ykey_save(NAN)
41  , zkey_save(NAN)
42  , cache_hits(0)
43  , cache_misses(0)
44 {
45  for (int i = 0; i < 2; i++)
46  {
47  for (int j = 0; j < 2; j++)
48  {
49  for (int k = 0; k < 2; k++)
50  {
51  for (int l = 0; l < 3; l++)
52  {
53  xyz[i][j][k][l] = NAN;
54  bf[i][j][k][l] = NAN;
55  }
56  }
57  }
58  }
59  cout << "\n================ Begin Construct Mag Field =====================" << endl;
60  cout << "\n-----------------------------------------------------------"
61  << "\n Magnetic field Module - Verbosity:"
62  << "\n-----------------------------------------------------------\n";
63 
64  bool root_input = false;
65  if(root_input) {
66  // open file
67  TFile *rootinput = TFile::Open(filename.c_str());
68  if (!rootinput)
69  {
70  cout << "\n could not open " << filename << " exiting now" << endl;
71  exit(1);
72  }
73  cout << "\n ---> "
74  "Reading the field grid from "
75  << filename << " ... " << endl;
76  rootinput->cd();
77 
78  // get root NTuple objects
79  TNtuple *field_map = (TNtuple *) gDirectory->Get("ntuple");
80  Float_t ROOT_X, ROOT_Y, ROOT_Z;
81  Float_t ROOT_BX, ROOT_BY, ROOT_BZ;
82  field_map->SetBranchAddress("x", &ROOT_X);
83  field_map->SetBranchAddress("y", &ROOT_Y);
84  field_map->SetBranchAddress("z", &ROOT_Z);
85  field_map->SetBranchAddress("Bx", &ROOT_BX);
86  field_map->SetBranchAddress("By", &ROOT_BY);
87  field_map->SetBranchAddress("Bz", &ROOT_BZ);
88 
89  for (int i = 0; i < field_map->GetEntries(); i++)
90  {
91  field_map->GetEntry(i);
92  trio coord_key(ROOT_X * UNIT_LENGTH, ROOT_Y * UNIT_LENGTH, ROOT_Z * UNIT_LENGTH);
93  trio field_val(ROOT_BX * UNIT_FIELD, ROOT_BY * UNIT_FIELD, ROOT_BZ * UNIT_FIELD);
94  xvals.insert(ROOT_X * UNIT_LENGTH);
95  yvals.insert(ROOT_Y * UNIT_LENGTH);
96  zvals.insert(ROOT_Z * UNIT_LENGTH);
97  fieldmap[coord_key] = field_val;
98 
99  if (rootinput) rootinput->Close();
100  }
101  } else {
102  ifstream file( filename.c_str() );
103  if (!file.good())
104  cout << "Field map input file error." << endl;
105 
106  char buffer[256];
107  file.getline(buffer,256);
108 
109  float ROOT_X, ROOT_Y, ROOT_Z;
110  float ROOT_BX, ROOT_BY, ROOT_BZ;
111  while (file >> ROOT_X >> ROOT_Y >> ROOT_Z >> ROOT_BX >> ROOT_BY >> ROOT_BZ) {
112  trio coord_key(ROOT_X * UNIT_LENGTH, ROOT_Y * UNIT_LENGTH, ROOT_Z * UNIT_LENGTH);
113  trio field_val(ROOT_BX * magfield_rescale * UNIT_FIELD, ROOT_BY * magfield_rescale * UNIT_FIELD, ROOT_BZ * magfield_rescale * UNIT_FIELD);
114  //cout << ROOT_X << " " << ROOT_Y << " " << ROOT_Z << endl;
115  //cout << ROOT_BX << " " << ROOT_BY << " " << ROOT_BZ << endl;
116  xvals.insert(ROOT_X * UNIT_LENGTH);
117  yvals.insert(ROOT_Y * UNIT_LENGTH);
118  zvals.insert(ROOT_Z * UNIT_LENGTH);
119  fieldmap[coord_key] = field_val;
120  }
121  file.close();
122  }
123 
124  xmin = *(xvals.begin());
125  xmax = *(xvals.rbegin());
126 
127  ymin = *(yvals.begin());
128  ymax = *(yvals.rbegin());
129  //TODO comment out for now, figure out why
130 // if (ymin != xmin || ymax != xmax)
131 // {
132 // cout << "PHField3DCartesian: Compiler bug!!!!!!!! Do not use inlining!!!!!!" << endl;
133 // cout << "exiting now - recompile with -fno-inline" << endl;
134 // exit(1);
135 // }
136 
137 // yuhw @ 10/10/18
138 // if (magfield_rescale != 1.0)
139 // {
140 // cout << "PHField3DCartesian: Rescale not implemented" << endl;
141 // exit(1);
142 // }
143 
144  zmin = *(zvals.begin());
145  zmax = *(zvals.rbegin());
146 
147  xstepsize = (xmax - xmin) / (xvals.size() - 1);
148  ystepsize = (ymax - ymin) / (yvals.size() - 1);
149  zstepsize = (zmax - zmin) / (zvals.size() - 1);
150 
151  std::cout << " its demensions and ranges of the field: " << fieldmap.size() << std::endl;
152  std::cout << " X: " << xvals.size() << " bins, " << xmin/cm << " cm -- " << xmax/cm << " cm." << std::endl;
153  std::cout << " Y: " << xvals.size() << " bins, " << ymin/cm << " cm -- " << ymax/cm << " cm." << std::endl;
154  std::cout << " Z: " << xvals.size() << " bins, " << zmin/cm << " cm -- " << zmax/cm << " cm." << std::endl;
155 
156 #ifdef _DEBUG_ON
157  auto key = boost::make_tuple(0, 0, 0);
158  auto magval = fieldmap.find(key);
159  if(magval!=fieldmap.end()) {
160  cout
161  << "DEBUG: "
162  << (magval->second).get<0>() << ", "
163  << (magval->second).get<1>() << ", "
164  << (magval->second).get<2>() << ", "
165  << endl;
166  }
167 #endif
168 
169 
170  cout << "\n================= End Construct Mag Field ======================\n"
171  << endl;
172 }
173 
175 {
176  // cout << "PHField3DCartesian: cache hits: " << cache_hits
177  // << " cache misses: " << cache_misses
178  // << endl;
179 }
180 
181 void PHField3DCartesian::GetFieldValue(const double point[4], double *Bfield) const
182 {
183  //TODO Test!
184 // cout << "untested code - I don't know if this is being used, drop me a line (with the field) and I test this"
185 // << endl
186 // << " Chris P." << endl;
187 // assert(0);
188 
189  static double xsav = -1000000.;
190  static double ysav = -1000000.;
191  static double zsav = -1000000.;
192 
193  double x = point[0];
194  double y = point[1];
195  double z = point[2];
196 
197  Bfield[0] = 0.0;
198  Bfield[1] = 0.0;
199  Bfield[2] = 0.0;
200  if (!isfinite(x) || !isfinite(y) || !isfinite(z))
201  {
202  static int ifirst = 0;
203  if (ifirst < 10)
204  {
205  cout << "PHField3DCartesian::GetFieldValue: "
206  << "Invalid coordinates: "
207  << "x: " << x
208  << ", y: " << y
209  << ", z: " << z
210  << " bailing out returning zero bfield"
211  << endl;
212  cout << "previous point: "
213  << "x: " << xsav
214  << ", y: " << ysav
215  << ", z: " << zsav
216  << endl;
217  ifirst++;
218  }
219  return;
220  }
221  xsav = x;
222  ysav = y;
223  zsav = z;
224 
225  if (point[0] < xmin || point[0] > xmax ||
226  point[1] < ymin || point[1] > ymax ||
227  point[2] < zmin || point[2] > zmax)
228  {
229  //std::cout << " sdsadsadsadasdasdasdas" << std::endl;
230  return;
231  }
232  set<double>::const_iterator it = xvals.lower_bound(x);
233  if (it == xvals.begin())
234  {
235  cout << "x too small - outside range: " << x << endl;
236  return;
237  }
238  double xkey[2];
239  xkey[0] = *it;
240  --it;
241  xkey[1] = *it;
242 
243  it = yvals.lower_bound(y);
244  if (it == yvals.begin())
245  {
246  cout << "y too small - outside range: " << y << endl;
247  return;
248  }
249  double ykey[2];
250  ykey[0] = *it;
251  --it;
252  ykey[1] = *it;
253 
254  it = zvals.lower_bound(z);
255  if (it == zvals.begin())
256  {
257  cout << "z too small - outside range: " << z << endl;
258  return;
259  }
260  double zkey[2];
261  zkey[0] = *it;
262  --it;
263  zkey[1] = *it;
264 
265  if (xkey_save != xkey[0] ||
266  ykey_save != ykey[0] ||
267  zkey_save != zkey[0])
268  {
269  cache_misses++;
270  xkey_save = xkey[0];
271  ykey_save = ykey[0];
272  zkey_save = zkey[0];
273 
274  map<boost::tuple<double, double, double>, boost::tuple<double, double, double> >::const_iterator magval;
275  trio key;
276  for (int i = 0; i < 2; i++)
277  {
278  for (int j = 0; j < 2; j++)
279  {
280  for (int k = 0; k < 2; k++)
281  {
282  key = boost::make_tuple(xkey[i], ykey[j], zkey[k]);
283  magval = fieldmap.find(key);
284  if (magval == fieldmap.end())
285  {
286  cout << "could not locate key, x: " << xkey[i] / UNIT_LENGTH
287  << ", y: " << ykey[j] / UNIT_LENGTH
288  << ", z: " << zkey[k] / UNIT_LENGTH << endl;
289  return;
290  }
291  xyz[i][j][k][0] = (magval->first).get<0>();
292  xyz[i][j][k][1] = (magval->first).get<1>();
293  xyz[i][j][k][2] = (magval->first).get<2>();
294  bf[i][j][k][0] = (magval->second).get<0>();
295  bf[i][j][k][1] = (magval->second).get<1>();
296  bf[i][j][k][2] = (magval->second).get<2>();
297 #ifdef _DEBUG_ON
298  cout << "read x/y/z: "
299  << xyz[i][j][k][0]/cm << "/"
300  << xyz[i][j][k][1]/cm << "/"
301  << xyz[i][j][k][2]/cm << " cm"
302  << " bx/by/bz: "
303  << bf[i][j][k][0]/tesla << "/"
304  << bf[i][j][k][1]/tesla << "/"
305  << bf[i][j][k][2]/tesla << " tesla"
306  << endl;
307 #endif
308  }
309  }
310  }
311  }
312  else
313  {
314  cache_hits++;
315  }
316  // linear extrapolation in cube:
317  //Vxyz = V000 (1 - x) (1 - y) (1 - z) +
318  //V100 x (1 - y) (1 - z) +
319  //V010 (1 - x) y (1 - z) +
320  //V001 (1 - x) (1 - y) z +
321  //V101 x (1 - y) z +
322  //V011 (1 - x) y z +
323  //V110 x y (1 - z) +
324  //V111 x y z
325 
326  double xinblock = point[0] - xkey[0];
327  double yinblock = point[1] - ykey[0];
328  double zinblock = point[2] - zkey[0];
329  // cout << "x/y/z stepsize: " << xstepsize << "/" << ystepsize << "/" << zstepsize << endl;
330  // cout << "x/y/z inblock: " << xinblock << "/" << yinblock << "/" << zinblock << endl;
331  for (int i = 0; i < 3; i++)
332  {
333  Bfield[i] = bf[0][0][0][i] * (xstepsize - xinblock) * (ystepsize - yinblock) * (zstepsize - zinblock) +
334  bf[1][0][0][i] * xinblock * (ystepsize - yinblock) * (zstepsize - zinblock) +
335  bf[0][1][0][i] * (xstepsize - xinblock) * yinblock * (zstepsize - zinblock) +
336  bf[0][0][1][i] * (xstepsize - xinblock) * (ystepsize - yinblock) * zinblock +
337  bf[1][0][1][i] * xinblock * (ystepsize - yinblock) * zinblock +
338  bf[0][1][1][i] * (xstepsize - xinblock) * yinblock * zinblock +
339  bf[1][1][0][i] * xinblock * yinblock * (zstepsize - zinblock) +
340  bf[1][1][1][i] * xinblock * yinblock * zinblock;
341  Bfield[i] /= (xstepsize*ystepsize*zstepsize);
342  }
343 
344  return;
345 }
#define UNIT_LENGTH
#define UNIT_FIELD
std::map< boost::tuple< double, double, double >, boost::tuple< double, double, double > > fieldmap
std::set< double > zvals
double xyz[2][2][2][3]
std::set< double > yvals
std::set< double > xvals
double bf[2][2][2][3]
PHField3DCartesian(const std::string &fname, const float magfield_rescale=1.0)
void GetFieldValue(const double Point[4], double *Bfield) const