Class Reference for E1039 Core & Analysis Software
PHField3DCylindrical.cc
Go to the documentation of this file.
1 
2 #include "PHField3DCylindrical.h"
3 
4 #include <TFile.h>
5 #include <TNtuple.h>
6 
7 #include <set>
8 #include <algorithm>
9 
10 using namespace std;
11 using namespace CLHEP; // units
12 
13 PHField3DCylindrical::PHField3DCylindrical(const string &filename, const int verb, const float magfield_rescale)
14  : PHField(verb)
15 {
16  cout << "\n================ Begin Construct Mag Field =====================" << endl;
17  cout << "\n-----------------------------------------------------------"
18  << "\n Magnetic field Module - Verbosity:" << verb_
19  << "\n-----------------------------------------------------------";
20 
21  // open file
22  TFile *rootinput = TFile::Open(filename.c_str());
23  if (!rootinput)
24  {
25  cout << "\n could not open " << filename << " exiting now" << endl;
26  exit(1);
27  }
28  cout << "\n ---> "
29  "Reading the field grid from "
30  << filename << " ... " << endl;
31  rootinput->cd();
32 
33  // get root NTuple objects
34  TNtuple *field_map = (TNtuple *) gDirectory->Get("map");
35  Float_t ROOT_Z, ROOT_R, ROOT_PHI;
36  Float_t ROOT_BZ, ROOT_BR, ROOT_BPHI;
37  field_map->SetBranchAddress("z", &ROOT_Z);
38  field_map->SetBranchAddress("r", &ROOT_R);
39  field_map->SetBranchAddress("phi", &ROOT_PHI);
40  field_map->SetBranchAddress("bz", &ROOT_BZ);
41  field_map->SetBranchAddress("br", &ROOT_BR);
42  field_map->SetBranchAddress("bphi", &ROOT_BPHI);
43 
44  // get the number of entries in the tree
45  int nz, nr, nphi;
46  nz = field_map->GetEntries("z>-1e6");
47  nr = field_map->GetEntries("r>-1e6");
48  nphi = field_map->GetEntries("phi>-1e6");
49  static const int NENTRIES = field_map->GetEntries();
50 
51  // run checks on entries
52  cout << " ---> The field grid contained " << NENTRIES << " entries" << endl;
53  if (verb_ > 0)
54  {
55  cout << "\n NENTRIES should be the same as the following values:"
56  << "\n [ Number of values r,phi,z: "
57  << nr << " " << nphi << " " << nz << " ]! " << endl;
58  }
59 
60  if (nz != nr || nz != nphi || nr != nphi)
61  {
62  cout << "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!"
63  << "\n The file you entered is not a \"table\" of values"
64  << "\n Something very likely went oh so wrong"
65  << "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!" << endl;
66  }
67 
68  // Keep track of the unique z, r, phi values in the grid using sets
69  std::set<float> z_set, r_set, phi_set;
70 
71  // Sort the entries to get rid of any stupid ordering problems...
72 
73  // We copy the TNtuple into a std::map (which is always sorted)
74  // using a 3-tuple of (z, r, phi) so it is sorted in z, then r, then
75  // phi.
76  if (verb_ > 0)
77  {
78  cout << " --> Sorting Entries..." << endl;
79  }
80  std::map<trio, trio> sorted_map;
81  for (int i = 0; i < field_map->GetEntries(); i++)
82  {
83  field_map->GetEntry(i);
84  trio coord_key(ROOT_Z, ROOT_R, ROOT_PHI);
85  trio field_val(ROOT_BZ, ROOT_BR, ROOT_BPHI);
86  sorted_map[coord_key] = field_val;
87 
88  z_set.insert(ROOT_Z * cm);
89  r_set.insert(ROOT_R * cm);
90  phi_set.insert(ROOT_PHI * deg);
91  }
92 
93  // couts for assurance
94  if (verb_ > 4)
95  {
96  map<trio, trio>::iterator it = sorted_map.begin();
97  print_map(it);
98  float last_z = it->first.get<0>();
99  for (it = sorted_map.begin(); it != sorted_map.end(); ++it)
100  {
101  if (it->first.get<0>() != last_z)
102  {
103  last_z = it->first.get<0>();
104  print_map(it);
105  }
106  }
107  }
108 
109  if (verb_ > 0)
110  {
111  cout << " --> Putting entries into containers... " << endl;
112  }
113 
114  // grab the minimum and maximum z values
115  minz_ = *(z_set.begin());
116  std::set<float>::iterator ziter = z_set.end();
117  --ziter;
118  maxz_ = *ziter;
119 
120  // initialize maps
121  nz = z_set.size();
122  nr = r_set.size();
123  nphi = phi_set.size();
124  r_map_.resize(nr, 0);
125  z_map_.resize(nz, 0);
126  phi_map_.resize(nphi, 0);
127 
128  std::copy(z_set.begin(), z_set.end(), z_map_.begin());
129  std::copy(phi_set.begin(), phi_set.end(), phi_map_.begin());
130  std::copy(r_set.begin(), r_set.end(), r_map_.begin());
131 
132  // initialize the field map vectors to the correct sizes
133  BFieldR_.resize(nz, vector<vector<float> >(nr, vector<float>(nphi, 0)));
134  BFieldPHI_.resize(nz, vector<vector<float> >(nr, vector<float>(nphi, 0)));
135  BFieldZ_.resize(nz, vector<vector<float> >(nr, vector<float>(nphi, 0)));
136 
137  // all of this assumes that z_prev < z , i.e. the table is ordered (as of right now)
138  unsigned int ir = 0, iphi = 0, iz = 0; // useful indexes to keep track of
139  map<trio, trio>::iterator iter = sorted_map.begin();
140  for (; iter != sorted_map.end(); ++iter)
141  {
142  // equivalent to ->GetEntry(iter)
143  float z = iter->first.get<0>() * cm;
144  float r = iter->first.get<1>() * cm;
145  float phi = iter->first.get<2>() * deg;
146  float Bz = iter->second.get<0>() * gauss;
147  float Br = iter->second.get<1>() * gauss;
148  float Bphi = iter->second.get<2>() * gauss;
149 
150  if (z > maxz_) maxz_ = z;
151  if (z < minz_) minz_ = z;
152 
153  // check for change in z value, when z changes we have a ton of updates to do
154  if (z != z_map_[iz])
155  {
156  ++iz;
157  ir = 0;
158  iphi = 0; // reset indices
159  }
160  else if (r != r_map_[ir])
161  { // check for change in r value
162  ++ir;
163  iphi = 0;
164  }
165  else if (phi != phi_map_[iphi])
166  { // change in phi value? (should be every time)
167  ++iphi;
168  }
169 
170  // shouldn't happen
171  if (iz > 0 && z < z_map_[iz - 1])
172  {
173  cout << "!!!!!!!!! Your map isn't ordered.... z: " << z << " zprev: " << z_map_[iz - 1] << endl;
174  }
175 
176  BFieldR_[iz][ir][iphi] = Br * magfield_rescale;
177  BFieldPHI_[iz][ir][iphi] = Bphi * magfield_rescale;
178  BFieldZ_[iz][ir][iphi] = Bz * magfield_rescale;
179 
180  // you can change this to check table values for correctness
181  // print_map prints the values in the root table, and the
182  // couts print the values entered into the vectors
183  if (fabs(z) < 10 && ir < 10 /*&& iphi==2*/ && verb_ > 3)
184  {
185  print_map(iter);
186 
187  cout << " B("
188  << r_map_[ir] << ", "
189  << phi_map_[iphi] << ", "
190  << z_map_[iz] << "): ("
191  << BFieldR_[iz][ir][iphi] << ", "
192  << BFieldPHI_[iz][ir][iphi] << ", "
193  << BFieldZ_[iz][ir][iphi] << ")" << endl;
194  }
195 
196  } // end loop over root field map file
197 
198  if (rootinput) rootinput->Close();
199 
200  cout << "\n ---> ... read file successfully "
201  << "\n ---> Z Boundaries ~ zlow, zhigh: "
202  << minz_ / cm << "," << maxz_ / cm << " cm " << endl;
203 
204  cout << "\n================= End Construct Mag Field ======================\n"
205  << endl;
206 }
207 
208 void PHField3DCylindrical::GetFieldValue(const double point[4], double *Bfield) const
209 {
210  if (verb_ > 2)
211  cout << "\nPHField3DCylindrical::GetFieldValue" << endl;
212  double x = point[0];
213  double y = point[1];
214  double z = point[2];
215  double r = sqrt(x * x + y * y);
216  double phi;
217  if (x == 0)
218  {
219  phi = atan2(y, 0.00000000001);
220  }
221  else
222  {
223  phi = atan2(y, x);
224  }
225  if (phi < 0) phi += 2 * M_PI; // normalize phi to be over the range [0,2*pi]
226 
227  // Check that the point is within the defined z region (check r in a second)
228  if ((z >= minz_) && (z <= maxz_))
229  {
230  double BFieldCyl[3];
231  double cylpoint[4] = {z, r, phi, 0};
232 
233  // take <z,r,phi> location and return a vector of <Bz, Br, Bphi>
234  GetFieldCyl(cylpoint, BFieldCyl);
235 
236  // X direction of B-field ( Bx = Br*cos(phi) - Bphi*sin(phi)
237  Bfield[0] = cos(phi) * BFieldCyl[1] - sin(phi) * BFieldCyl[2]; // unit vector transformations
238 
239  // Y direction of B-field ( By = Br*sin(phi) + Bphi*cos(phi)
240  Bfield[1] = sin(phi) * BFieldCyl[1] + cos(phi) * BFieldCyl[2];
241 
242  // Z direction of B-field
243  Bfield[2] = BFieldCyl[0];
244  }
245  else
246  { // x,y,z is outside of z range of the field map
247 
248  Bfield[0] = 0.0;
249  Bfield[1] = 0.0;
250  Bfield[2] = 0.0;
251  if (verb_ > 2)
252  cout << "!!!!!!!!!! Field point not in defined region (outside of z bounds)" << endl;
253  }
254 
255  if (verb_ > 2)
256  {
257  cout << "END PHField3DCylindrical::GetFieldValue\n"
258  << " ---> {Bx, By, Bz} : "
259  << "< " << Bfield[0] << ", " << Bfield[1] << ", " << Bfield[2] << " >" << endl;
260  }
261 
262  return;
263 }
264 
265 void PHField3DCylindrical::GetFieldCyl(const double CylPoint[4], double *BfieldCyl) const
266 {
267  float z = CylPoint[0];
268  float r = CylPoint[1];
269  float phi = CylPoint[2];
270 
271  BfieldCyl[0] = 0.0;
272  BfieldCyl[1] = 0.0;
273  BfieldCyl[2] = 0.0;
274 
275  if (verb_ > 2)
276  cout << "GetFieldCyl@ <z,r,phi>: {" << z << "," << r << "," << phi << "}" << endl;
277 
278  if (z < z_map_[0] || z > z_map_[z_map_.size() - 1])
279  {
280  if (verb_ > 2)
281  cout << "!!!! Point not in defined region (radius too large in specific z-plane)" << endl;
282  return;
283  }
284 
285  vector<float>::const_iterator ziter = upper_bound(z_map_.begin(), z_map_.end(), z);
286  unsigned int z_index0 = distance(z_map_.begin(), ziter) - 1;
287  unsigned int z_index1 = z_index0 + 1;
288 
289  vector<float>::const_iterator riter = upper_bound(r_map_.begin(), r_map_.end(), r);
290  unsigned int r_index0 = distance(r_map_.begin(), riter) - 1;
291  if (r_index0 >= r_map_.size())
292  {
293  if (verb_ > 2)
294  cout << "!!!! Point not in defined region (radius too large in specific z-plane)" << endl;
295  return;
296  }
297 
298  unsigned int r_index1 = r_index0 + 1;
299  if (r_index1 >= r_map_.size())
300  {
301  if (verb_ > 2)
302  cout << "!!!! Point not in defined region (radius too large in specific z-plane)" << endl;
303  return;
304  }
305 
306  vector<float>::const_iterator phiiter = upper_bound(phi_map_.begin(), phi_map_.end(), phi);
307  unsigned int phi_index0 = distance(phi_map_.begin(), phiiter) - 1;
308  unsigned int phi_index1 = phi_index0 + 1;
309  if (phi_index1 >= phi_map_.size())
310  phi_index1 = 0;
311 
312  double Br000 = BFieldR_[z_index0][r_index0][phi_index0];
313  double Br001 = BFieldR_[z_index0][r_index0][phi_index1];
314  double Br010 = BFieldR_[z_index0][r_index1][phi_index0];
315  double Br011 = BFieldR_[z_index0][r_index1][phi_index1];
316  double Br100 = BFieldR_[z_index1][r_index0][phi_index0];
317  double Br101 = BFieldR_[z_index1][r_index0][phi_index1];
318  double Br110 = BFieldR_[z_index1][r_index1][phi_index0];
319  double Br111 = BFieldR_[z_index1][r_index1][phi_index1];
320 
321  double Bphi000 = BFieldPHI_[z_index0][r_index0][phi_index0];
322  double Bphi001 = BFieldPHI_[z_index0][r_index0][phi_index1];
323  double Bphi010 = BFieldPHI_[z_index0][r_index1][phi_index0];
324  double Bphi011 = BFieldPHI_[z_index0][r_index1][phi_index1];
325  double Bphi100 = BFieldPHI_[z_index1][r_index0][phi_index0];
326  double Bphi101 = BFieldPHI_[z_index1][r_index0][phi_index1];
327  double Bphi110 = BFieldPHI_[z_index1][r_index1][phi_index0];
328  double Bphi111 = BFieldPHI_[z_index1][r_index1][phi_index1];
329 
330  double Bz000 = BFieldZ_[z_index0][r_index0][phi_index0];
331  double Bz001 = BFieldZ_[z_index0][r_index0][phi_index1];
332  double Bz100 = BFieldZ_[z_index1][r_index0][phi_index0];
333  double Bz101 = BFieldZ_[z_index1][r_index0][phi_index1];
334  double Bz010 = BFieldZ_[z_index0][r_index1][phi_index0];
335  double Bz110 = BFieldZ_[z_index1][r_index1][phi_index0];
336  double Bz011 = BFieldZ_[z_index0][r_index1][phi_index1];
337  double Bz111 = BFieldZ_[z_index1][r_index1][phi_index1];
338 
339  double zweight = z - z_map_[z_index0];
340  double zspacing = z_map_[z_index1] - z_map_[z_index0];
341  zweight /= zspacing;
342 
343  double rweight = r - r_map_[r_index0];
344  double rspacing = r_map_[r_index1] - r_map_[r_index0];
345  rweight /= rspacing;
346 
347  double phiweight = phi - phi_map_[phi_index0];
348  double phispacing = phi_map_[phi_index1] - phi_map_[phi_index0];
349  if (phi_index1 == 0)
350  phispacing += 2 * M_PI;
351  phiweight /= phispacing;
352 
353  // Z direction of B-field
354  BfieldCyl[0] =
355  (1 - zweight) * ((1 - rweight) * ((1 - phiweight) * Bz000 + phiweight * Bz001) +
356  rweight * ((1 - phiweight) * Bz010 + phiweight * Bz011)) +
357  zweight * ((1 - rweight) * ((1 - phiweight) * Bz100 + phiweight * Bz101) +
358  rweight * ((1 - phiweight) * Bz110 + phiweight * Bz111));
359 
360  // R direction of B-field
361  BfieldCyl[1] =
362  (1 - zweight) * ((1 - rweight) * ((1 - phiweight) * Br000 + phiweight * Br001) +
363  rweight * ((1 - phiweight) * Br010 + phiweight * Br011)) +
364  zweight * ((1 - rweight) * ((1 - phiweight) * Br100 + phiweight * Br101) +
365  rweight * ((1 - phiweight) * Br110 + phiweight * Br111));
366 
367  // PHI Direction of B-field
368  BfieldCyl[2] =
369  (1 - zweight) * ((1 - rweight) * ((1 - phiweight) * Bphi000 + phiweight * Bphi001) +
370  rweight * ((1 - phiweight) * Bphi010 + phiweight * Bphi011)) +
371  zweight * ((1 - rweight) * ((1 - phiweight) * Bphi100 + phiweight * Bphi101) +
372  rweight * ((1 - phiweight) * Bphi110 + phiweight * Bphi111));
373 
374  // cout << "wr: " << rweight << " wz: " << zweight << " wphi: " << phiweight << endl;
375  // cout << "Bz000: " << Bz000 << endl
376  // << "Bz001: " << Bz001 << endl
377  // << "Bz010: " << Bz010 << endl
378  // << "Bz011: " << Bz011 << endl
379  // << "Bz100: " << Bz100 << endl
380  // << "Bz101: " << Bz101 << endl
381  // << "Bz110: " << Bz110 << endl
382  // << "Bz111: " << Bz111 << endl
383  // << "Bz: " << BfieldCyl[0] << endl << endl;
384 
385  if (verb_ > 2)
386  {
387  cout << "End GFCyl Call: <bz,br,bphi> : {"
388  << BfieldCyl[0] / gauss << "," << BfieldCyl[1] / gauss << "," << BfieldCyl[2] / gauss << "}"
389  << endl;
390  }
391 
392  return;
393 }
394 
395 // a binary search algorithm that puts the location that "key" would be, into index...
396 // it returns true if key was found, and false if not.
397 bool PHField3DCylindrical::bin_search(const vector<float> &vec, unsigned start, unsigned end, const float &key, unsigned &index) const
398 {
399  // Termination condition: start index greater than end index
400  if (start > end)
401  {
402  index = start;
403  return false;
404  }
405 
406  // Find the middle element of the vector and use that for splitting
407  // the array into two pieces.
408  unsigned middle = start + ((end - start) / 2);
409 
410  if (vec[middle] == key)
411  {
412  index = middle;
413  return true;
414  }
415  else if (vec[middle] > key)
416  {
417  return bin_search(vec, start, middle - 1, key, index);
418  }
419  return bin_search(vec, middle + 1, end, key, index);
420 }
421 
422 // debug function to print key/value pairs in map
423 void PHField3DCylindrical::print_map(map<trio, trio>::iterator &it) const
424 {
425  cout << " Key: <"
426  << it->first.get<0>() * cm << ","
427  << it->first.get<1>() * cm << ","
428  << it->first.get<2>() * deg << ">"
429 
430  << " Value: <"
431  << it->second.get<0>() * gauss << ","
432  << it->second.get<1>() * gauss << ","
433  << it->second.get<2>() * gauss << ">\n";
434 }
std::vector< float > z_map_
std::vector< std::vector< std::vector< float > > > BFieldPHI_
std::vector< std::vector< std::vector< float > > > BFieldZ_
std::vector< float > r_map_
std::vector< std::vector< std::vector< float > > > BFieldR_
void GetFieldValue(const double Point[4], double *Bfield) const
void GetFieldCyl(const double CylPoint[4], double *Bfield) const
std::vector< float > phi_map_
PHField3DCylindrical(const std::string &filename, int verb=0, const float magfield_rescale=1.0)
transient DST object for field storage and access
Definition: PHField.h:14
unsigned verb_
Definition: PHField.h:34