Class Reference for E1039 Core & Analysis Software
PHField2D.cc
Go to the documentation of this file.
1 
2 #include "PHField2D.h"
3 
4 //root framework
5 #include <TFile.h>
6 #include <TNtuple.h>
7 
8 #include <set>
9 #include <iostream>
10 #include <algorithm>
11 
12 using namespace std;
13 using namespace CLHEP;// units
14 
15 PHField2D::PHField2D( const string &filename, const int verb, const float magfield_rescale) :
16  PHField(verb)
17 {
18  r_index0_cache = 0;
19  r_index1_cache = 0;
20  z_index0_cache = 0;
21  z_index1_cache = 0;
22 
23  if (verb_ > 0)
24  cout << " ------------- PHField2D::PHField2D() ------------------" << endl;
25 
26  // open file
27  TFile *rootinput = TFile::Open(filename.c_str());
28  if (!rootinput)
29  {
30  cout << " could not open " << filename << " exiting now" << endl;
31  exit(1);
32  }
33  if (verb_ > 0) cout << " Field grid file: " << filename << endl;
34  rootinput->cd();
35 
36  Float_t ROOT_Z, ROOT_R;
37  Float_t ROOT_BZ, ROOT_BR;
38  // get root NTuple objects
39  TNtuple *field_map = (TNtuple*)gDirectory->Get("fieldmap");
40  if (field_map)
41  {
42  magfield_unit = tesla;
43  }
44  else
45  {
46  field_map = (TNtuple*)gDirectory->Get("map");
47  if (! field_map)
48  {
49  cout << "PHField2D: could not locate ntuple of name map or fieldmap, exiting now" << endl;
50  exit(1);
51  }
52  magfield_unit = gauss;
53  }
54  field_map->SetBranchAddress("z", &ROOT_Z);
55  field_map->SetBranchAddress("r", &ROOT_R);
56  field_map->SetBranchAddress("bz", &ROOT_BZ);
57  field_map->SetBranchAddress("br", &ROOT_BR);
58 
59  // get the number of entries in the tree
60  int nz, nr;
61  nz = field_map->GetEntries("z>-1e6");
62  nr = field_map->GetEntries("r>-1e6");
63  static const int NENTRIES = field_map->GetEntries();
64 
65  // run checks on entries
66  if (verb_ > 0) cout << " The field grid contained " << NENTRIES << " entries" << endl;
67  if ( verb_ > 1 )
68  {
69  cout << "\n NENTRIES should be the same as the following values:"
70  << "\n [ Number of values r,z: "
71  << nr << " " << nz << " ]! " << endl;
72  }
73 
74  if (nz != nr)
75  {
76  cout << "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!"
77  << "\n The file you entered is not a \"table\" of values"
78  << "\n Something very likely went oh so wrong"
79  << "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!" << endl;
80  }
81 
82  // Keep track of the unique z, r, phi values in the grid using sets
83  std::set<float> z_set, r_set;
84 
85  // Sort the entries to get rid of any stupid ordering problems...
86 
87  // We copy the TNtuple into a std::map (which is always sorted)
88  // using a 3-tuple of (z, r, phi) so it is sorted in z, then r, then
89  // phi.
90  if ( verb_ > 1 )
91  {
92  cout << " --> Sorting Entries..." << endl;
93  }
94  std::map<trio, trio> sorted_map;
95  for ( int i = 0; i < field_map->GetEntries(); i++)
96  {
97  field_map->GetEntry(i);
98  trio coord_key( ROOT_Z, ROOT_R );
99  trio field_val( ROOT_BZ, ROOT_BR );
100  sorted_map[coord_key] = field_val;
101 
102  z_set.insert(ROOT_Z*cm);
103  r_set.insert(ROOT_R*cm);
104  }
105 
106  // couts for assurance
107  if (verb_ > 4)
108  {
109  map<trio, trio>::iterator it = sorted_map.begin();
110  print_map(it);
111  float last_z = it->first.get<0>();
112  for ( it = sorted_map.begin(); it != sorted_map.end(); ++it )
113  {
114  if ( it->first.get<0>() != last_z )
115  {
116  last_z = it->first.get<0>();
117  print_map( it );
118  }
119  }
120  }
121 
122  if (verb_ > 1)
123  {
124  cout << " --> Putting entries into containers... " << endl;
125  }
126 
127  // grab the minimum and maximum z values
128  minz_ = *(z_set.begin());
129  std::set<float>::iterator ziter = z_set.end();
130  --ziter;
131  maxz_ = *ziter;
132 
133  // initialize maps
134  nz = z_set.size();
135  nr = r_set.size();
136  r_map_.resize( nr, 0 );
137  z_map_.resize( nz, 0 );
138 
139  std::copy(z_set.begin(), z_set.end(), z_map_.begin());
140  std::copy(r_set.begin(), r_set.end(), r_map_.begin());
141 
142  // initialize the field map vectors to the correct sizes
143  BFieldR_.resize( nz, vector<float>(nr, 0) );
144  BFieldZ_.resize( nz, vector<float>(nr, 0) );
145 
146  // all of this assumes that z_prev < z , i.e. the table is ordered (as of right now)
147  unsigned int ir = 0, iz = 0; // useful indexes to keep track of
148  map<trio, trio>::iterator iter = sorted_map.begin();
149  for ( ; iter != sorted_map.end(); ++iter )
150  {
151 
152  // equivalent to ->GetEntry(iter)
153  float z = iter-> first.get<0>() * cm;
154  float r = iter-> first.get<1>() * cm;
155  float Bz = iter->second.get<0>() * magfield_unit;
156  float Br = iter->second.get<1>() * magfield_unit;
157 
158  if ( z > maxz_ ) maxz_ = z;
159  if ( z < minz_ ) minz_ = z;
160 
161  // check for change in z value, when z changes we have a ton of updates to do
162  if ( z != z_map_[iz] )
163  {
164  ++iz;
165  ir = 0;
166 
167  }
168  else if ( r != r_map_[ir] ) // check for change in r value
169  {
170  ++ir;
171  }
172 
173  // shouldn't happen
174  if ( iz > 0 && z < z_map_[iz-1] )
175  {
176  cout << "!!!!!!!!! Your map isn't ordered.... z: " << z << " zprev: " << z_map_[iz-1] << endl;
177  }
178 
179  BFieldR_[iz][ir] = Br * magfield_rescale;
180  BFieldZ_[iz][ir] = Bz * magfield_rescale;
181 
182  // you can change this to check table values for correctness
183  // print_map prints the values in the root table, and the
184  // couts print the values entered into the vectors
185  if ( fabs(z) < 10 && ir < 10 /*&& iphi==2*/ && verb_ > 3)
186  {
187  print_map(iter);
188 
189  cout << " B("
190  << r_map_[ir] << ", "
191  << z_map_[iz] << "): ("
192  << BFieldR_[iz][ir] << ", "
193  << BFieldZ_[iz][ir] << ")" << endl;
194  }
195 
196  } // end loop over root field map file
197 
198  if (rootinput)
199  rootinput->Close();
200 
201  if (verb_ > 0) cout << " Mag field z boundaries (min,max): (" << minz_ / cm << ", " << maxz_ / cm << ") cm" << endl;
202  if (verb_ > 0) cout << " Mag field r max boundary: " << r_map_.back()/ cm << " cm" << endl;
203 
204  if (verb_ > 0)
205  cout << " -----------------------------------------------------------" << endl;
206 }
207 
208 void PHField2D::GetFieldValue(const double point[4], double *Bfield ) const
209 {
210  if (verb_ > 2)
211  cout << "\nPHField2D::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  phi = atan2(y, x);
218  if ( phi < 0 ) phi += 2 * M_PI; // normalize phi to be over the range [0,2*pi]
219 
220  // Check that the point is within the defined z region (check r in a second)
221  if ( (z >= minz_) && (z <= maxz_) )
222  {
223 
224  double BFieldCyl[3];
225  double cylpoint[4] = { z, r, phi, 0 };
226 
227  // take <z,r,phi> location and return a vector of <Bz, Br, Bphi>
228  GetFieldCyl( cylpoint, BFieldCyl );
229 
230  // X direction of B-field ( Bx = Br*cos(phi) - Bphi*sin(phi)
231  Bfield[0] = cos(phi) * BFieldCyl[1] - sin(phi) * BFieldCyl[2]; // unit vector transformations
232 
233  // Y direction of B-field ( By = Br*sin(phi) + Bphi*cos(phi)
234  Bfield[1] = sin(phi) * BFieldCyl[1] + cos(phi) * BFieldCyl[2];
235 
236  // Z direction of B-field
237  Bfield[2] = BFieldCyl[0];
238 
239  }
240  else // x,y,z is outside of z range of the field map
241  {
242 
243  Bfield[0] = 0.0;
244  Bfield[1] = 0.0;
245  Bfield[2] = 0.0;
246  if ( verb_ > 2 )
247  cout << "!!!!!!!!!! Field point not in defined region (outside of z bounds)" << endl;
248  }
249 
250  if (verb_ > 2)
251  {
252  cout << "END PHField2D::GetFieldValue\n"
253  << " ---> {Bx, By, Bz} : "
254  << "< " << Bfield[0] << ", " << Bfield[1] << ", " << Bfield[2] << " >" << endl;
255  }
256 
257  return;
258 }
259 
260 void PHField2D::GetFieldCyl( const double CylPoint[4], double *BfieldCyl ) const
261 {
262  float z = CylPoint[0];
263  float r = CylPoint[1];
264 
265  BfieldCyl[0] = 0.0;
266  BfieldCyl[1] = 0.0;
267  BfieldCyl[2] = 0.0;
268 
269  if ( verb_ > 2 )
270  cout << "GetFieldCyl@ <z,r>: {" << z << "," << r << "}" << endl;
271 
272  if (z < z_map_[0] || z > z_map_[z_map_.size()-1])
273  {
274  if ( verb_ > 2 )
275  cout << "!!!! Point not in defined region (radius too large in specific z-plane)" << endl;
276  return;
277  }
278 
279  // since GEANT4 looks up the field ~95% of the time in the same voxel
280  // between subsequent calls, we can save on the expense of the upper_bound
281  // lookup (~10-15% of central event run time) with some caching between calls
282 
283  unsigned int r_index0 = r_index0_cache;
284  unsigned int r_index1 = r_index1_cache;
285 
286  if (!((r > r_map_[r_index0])&&(r < r_map_[r_index1]))) {
287 
288  // if miss cached r values, search through the lookup table
289  vector<float>::const_iterator riter = upper_bound(r_map_.begin(), r_map_.end(), r);
290  r_index0 = distance(r_map_.begin(), riter) - 1;
291  if (r_index0 >= r_map_.size()) {
292  if ( verb_ > 2 )
293  cout << "!!!! Point not in defined region (radius too large in specific z-plane)" << endl;
294  return;
295  }
296 
297  r_index1 = r_index0 + 1;
298  if (r_index1 >= r_map_.size()) {
299  if ( verb_ > 2 )
300  cout << "!!!! Point not in defined region (radius too large in specific z-plane)" << endl;
301  return;
302  }
303 
304  // update cache
305  r_index0_cache = r_index0;
306  r_index1_cache = r_index1;
307  }
308 
309  unsigned int z_index0 = z_index0_cache;
310  unsigned int z_index1 = z_index1_cache;
311 
312  if (!((z > z_map_[z_index0])&&(z < z_map_[z_index1]))) {
313 
314  // if miss cached z values, search through the lookup table
315  vector<float>::const_iterator ziter = upper_bound(z_map_.begin(), z_map_.end(), z);
316  z_index0 = distance(z_map_.begin(), ziter) - 1;
317  z_index1 = z_index0 + 1;
318  if (z_index1 >= z_map_.size()) {
319  if ( verb_ > 2 )
320  cout << "!!!! Point not in defined region (z too large in specific r-plane)" << endl;
321  return;
322  }
323 
324  // update cache
325  z_index0_cache = z_index0;
326  z_index1_cache = z_index1;
327  }
328 
329  double Br000 = BFieldR_[z_index0][r_index0];
330  double Br010 = BFieldR_[z_index0][r_index1];
331  double Br100 = BFieldR_[z_index1][r_index0];
332  double Br110 = BFieldR_[z_index1][r_index1];
333 
334  double Bz000 = BFieldZ_[z_index0][r_index0];
335  double Bz100 = BFieldZ_[z_index1][r_index0];
336  double Bz010 = BFieldZ_[z_index0][r_index1];
337  double Bz110 = BFieldZ_[z_index1][r_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  // Z direction of B-field
348  BfieldCyl[0] =
349  (1 - zweight) * ((1 - rweight) * Bz000 +
350  rweight * Bz010) +
351  zweight * ((1 - rweight) * Bz100 +
352  rweight * Bz110);
353 
354  // R direction of B-field
355  BfieldCyl[1] =
356  (1 - zweight) * ((1 - rweight) * Br000 +
357  rweight * Br010) +
358  zweight * ((1 - rweight) * Br100 +
359  rweight * Br110);
360 
361  // PHI Direction of B-field
362  BfieldCyl[2] = 0;
363 
364  if ( verb_ > 2 )
365  {
366  cout << "End GFCyl Call: <bz,br,bphi> : {"
367  << BfieldCyl[0] / gauss << "," << BfieldCyl[1] / gauss << "," << BfieldCyl[2] / gauss << "}"
368  << endl;
369  }
370 
371  return;
372 }
373 
374 // debug function to print key/value pairs in map
375 void PHField2D::print_map( map<trio, trio>::iterator& it ) const
376 {
377 
378  cout << " Key: <"
379  << it->first.get<0>() / cm << ","
380  << it->first.get<1>() / cm << ">"
381 
382  << " Value: <"
383  << it->second.get<0>() / magfield_unit << ","
384  << it->second.get<1>() / magfield_unit << ">\n";
385 
386 }
387 
std::vector< std::vector< float > > BFieldR_
Definition: PHField2D.h:32
void GetFieldCyl(const double CylPoint[4], double *Bfield) const
Definition: PHField2D.cc:260
PHField2D(const std::string &filename, const int verb=0, const float magfield_rescale=1.0)
Definition: PHField2D.cc:15
std::vector< float > r_map_
Definition: PHField2D.h:37
std::vector< float > z_map_
Definition: PHField2D.h:36
double magfield_unit
Definition: PHField2D.h:41
float maxz_
Definition: PHField2D.h:40
std::vector< std::vector< float > > BFieldZ_
Definition: PHField2D.h:31
float minz_
Definition: PHField2D.h:40
void GetFieldValue(const double Point[4], double *Bfield) const
Definition: PHField2D.cc:208
transient DST object for field storage and access
Definition: PHField.h:14
unsigned verb_
Definition: PHField.h:34