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