11 using namespace CLHEP;
16 cout <<
"\n================ Begin Construct Mag Field =====================" << endl;
17 cout <<
"\n-----------------------------------------------------------"
18 <<
"\n Magnetic field Module - Verbosity:" <<
verb_
19 <<
"\n-----------------------------------------------------------";
22 TFile *rootinput = TFile::Open(filename.c_str());
25 cout <<
"\n could not open " << filename <<
" exiting now" << endl;
29 "Reading the field grid from "
30 << filename <<
" ... " << endl;
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);
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();
52 cout <<
" ---> The field grid contained " << NENTRIES <<
" entries" << endl;
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;
60 if (nz != nr || nz != nphi || nr != nphi)
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;
69 std::set<float> z_set, r_set, phi_set;
78 cout <<
" --> Sorting Entries..." << endl;
80 std::map<trio, trio> sorted_map;
81 for (
int i = 0; i < field_map->GetEntries(); i++)
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;
88 z_set.insert(ROOT_Z * cm);
89 r_set.insert(ROOT_R * cm);
90 phi_set.insert(ROOT_PHI * deg);
96 map<trio, trio>::iterator it = sorted_map.begin();
98 float last_z = it->first.get<0>();
99 for (it = sorted_map.begin(); it != sorted_map.end(); ++it)
101 if (it->first.get<0>() != last_z)
103 last_z = it->first.get<0>();
111 cout <<
" --> Putting entries into containers... " << endl;
115 minz_ = *(z_set.begin());
116 std::set<float>::iterator ziter = z_set.end();
123 nphi = phi_set.size();
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());
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)));
138 unsigned int ir = 0, iphi = 0, iz = 0;
139 map<trio, trio>::iterator iter = sorted_map.begin();
140 for (; iter != sorted_map.end(); ++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;
171 if (iz > 0 && z <
z_map_[iz - 1])
173 cout <<
"!!!!!!!!! Your map isn't ordered.... z: " << z <<
" zprev: " <<
z_map_[iz - 1] << endl;
176 BFieldR_[iz][ir][iphi] = Br * magfield_rescale;
177 BFieldPHI_[iz][ir][iphi] = Bphi * magfield_rescale;
178 BFieldZ_[iz][ir][iphi] = Bz * magfield_rescale;
183 if (fabs(z) < 10 && ir < 10 &&
verb_ > 3)
193 <<
BFieldZ_[iz][ir][iphi] <<
")" << endl;
198 if (rootinput) rootinput->Close();
200 cout <<
"\n ---> ... read file successfully "
201 <<
"\n ---> Z Boundaries ~ zlow, zhigh: "
202 <<
minz_ / cm <<
"," <<
maxz_ / cm <<
" cm " << endl;
204 cout <<
"\n================= End Construct Mag Field ======================\n"
211 cout <<
"\nPHField3DCylindrical::GetFieldValue" << endl;
215 double r = sqrt(x * x + y * y);
219 phi = atan2(y, 0.00000000001);
225 if (phi < 0) phi += 2 * M_PI;
231 double cylpoint[4] = {z, r, phi, 0};
237 Bfield[0] = cos(phi) * BFieldCyl[1] - sin(phi) * BFieldCyl[2];
240 Bfield[1] = sin(phi) * BFieldCyl[1] + cos(phi) * BFieldCyl[2];
243 Bfield[2] = BFieldCyl[0];
252 cout <<
"!!!!!!!!!! Field point not in defined region (outside of z bounds)" << endl;
257 cout <<
"END PHField3DCylindrical::GetFieldValue\n"
258 <<
" ---> {Bx, By, Bz} : "
259 <<
"< " << Bfield[0] <<
", " << Bfield[1] <<
", " << Bfield[2] <<
" >" << endl;
267 float z = CylPoint[0];
268 float r = CylPoint[1];
269 float phi = CylPoint[2];
276 cout <<
"GetFieldCyl@ <z,r,phi>: {" << z <<
"," << r <<
"," << phi <<
"}" << endl;
281 cout <<
"!!!! Point not in defined region (radius too large in specific z-plane)" << endl;
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;
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())
294 cout <<
"!!!! Point not in defined region (radius too large in specific z-plane)" << endl;
298 unsigned int r_index1 = r_index0 + 1;
299 if (r_index1 >=
r_map_.size())
302 cout <<
"!!!! Point not in defined region (radius too large in specific z-plane)" << endl;
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;
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];
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];
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];
339 double zweight = z -
z_map_[z_index0];
343 double rweight = r -
r_map_[r_index0];
347 double phiweight = phi -
phi_map_[phi_index0];
350 phispacing += 2 * M_PI;
351 phiweight /= phispacing;
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));
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));
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));
387 cout <<
"End GFCyl Call: <bz,br,bphi> : {"
388 << BfieldCyl[0] / gauss <<
"," << BfieldCyl[1] / gauss <<
"," << BfieldCyl[2] / gauss <<
"}"
397 bool PHField3DCylindrical::bin_search(
const vector<float> &vec,
unsigned start,
unsigned end,
const float &key,
unsigned &index)
const
408 unsigned middle = start + ((end - start) / 2);
410 if (vec[middle] == key)
415 else if (vec[middle] > key)
417 return bin_search(vec, start, middle - 1, key, index);
419 return bin_search(vec, middle + 1, end, key, index);
423 void PHField3DCylindrical::print_map(map<trio, trio>::iterator &it)
const
426 << it->first.get<0>() * cm <<
","
427 << it->first.get<1>() * cm <<
","
428 << it->first.get<2>() * deg <<
">"
431 << it->second.get<0>() * gauss <<
","
432 << it->second.get<1>() * gauss <<
","
433 << it->second.get<2>() * gauss <<
">\n";
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