3 #include <TDirectory.h>
7 #include <Geant4/G4SystemOfUnits.hh>
21 std::cout <<
"\n================ Begin Construct Mag Field =====================" << std::endl;
22 std::cout <<
"\n-----------------------------------------------------------"
23 <<
"\n Magnetic field Module - Verbosity:" <<
Verbosity()
24 <<
"\n-----------------------------------------------------------";
27 TFile *rootinput = TFile::Open(filename.c_str());
30 std::cout <<
"\n could not open " << filename <<
" exiting now" << std::endl;
33 std::cout <<
"\n ---> "
34 "Reading the field grid from "
35 << filename <<
" ... " << std::endl;
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);
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();
57 std::cout <<
" ---> The field grid contained " << NENTRIES <<
" entries" << std::endl;
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;
65 if (nz != nr || nz != nphi || nr != nphi)
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;
74 std::set<float> z_set, r_set, phi_set;
83 std::cout <<
" --> Sorting Entries..." << std::endl;
85 std::map<trio, trio> sorted_map;
86 for (
int i = 0;
i < field_map->GetEntries();
i++)
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;
93 z_set.insert(ROOT_Z *
cm);
94 r_set.insert(ROOT_R * cm);
95 phi_set.insert(ROOT_PHI *
deg);
101 std::map<trio, trio>::iterator
it = sorted_map.begin();
103 float last_z = std::get<0>(it->first);
104 for (it = sorted_map.begin(); it != sorted_map.end(); ++
it)
106 if (std::get<0>(it->first) != last_z)
108 last_z = std::get<0>(it->first);
116 std::cout <<
" --> Putting entries into containers... " << std::endl;
120 minz_ = *(z_set.begin());
121 std::set<float>::iterator ziter = z_set.end();
128 nphi = phi_set.size();
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());
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)));
143 unsigned int ir = 0, iphi = 0, iz = 0;
144 std::map<trio, trio>::iterator iter = sorted_map.begin();
145 for (; iter != sorted_map.end(); ++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;
182 if (iz > 0 && z <
z_map_[iz - 1])
184 std::cout <<
"!!!!!!!!! Your map isn't ordered.... z: " << z <<
" zprev: " <<
z_map_[iz - 1] << std::endl;
194 if (std::fabs(z) < 10 && ir < 10 &&
Verbosity() > 3)
204 <<
BFieldZ_[iz][ir][iphi] <<
")" << std::endl;
211 std::cout <<
"\n ---> ... read file successfully "
212 <<
"\n ---> Z Boundaries ~ zlow, zhigh: "
215 std::cout <<
"\n================= End Construct Mag Field ======================\n"
223 std::cout <<
"\nPHField3DCylindrical::GetFieldValue" << std::endl;
228 double r = sqrt(x * x + y * y);
232 phi = atan2(y, 0.00000000001);
247 double cylpoint[4] = {
z,
r,
phi, 0};
253 Bfield[0] = cos(phi) * BFieldCyl[1] - sin(phi) * BFieldCyl[2];
256 Bfield[1] = sin(phi) * BFieldCyl[1] + cos(phi) * BFieldCyl[2];
259 Bfield[2] = BFieldCyl[0];
269 std::cout <<
"!!!!!!!!!! Field point not in defined region (outside of z bounds)" << std::endl;
275 std::cout <<
"END PHField3DCylindrical::GetFieldValue\n"
276 <<
" ---> {Bx, By, Bz} : "
277 <<
"< " << Bfield[0] <<
", " << Bfield[1] <<
", " << Bfield[2] <<
" >" << std::endl;
285 float z = CylPoint[0];
286 float r = CylPoint[1];
287 float phi = CylPoint[2];
295 std::cout <<
"GetFieldCyl@ <z,r,phi>: {" << z <<
"," << r <<
"," << phi <<
"}" << std::endl;
302 std::cout <<
"!!!! Point not in defined region (|z| too large)" << std::endl;
311 std::cout <<
"!!!! Point not in defined region (radius too small in specific z-plane). Use min radius" << std::endl;
319 std::cout <<
"!!!! Point not in defined region (radius too large in specific z-plane)" << std::endl;
324 std::vector<float>::const_iterator ziter = upper_bound(
z_map_.begin(),
z_map_.end(),
z);
326 int z_index1 = z_index0 + 1;
333 std::vector<float>::const_iterator riter = upper_bound(
r_map_.begin(),
r_map_.end(),
r);
335 if (r_index0 >= (
int)
r_map_.size())
339 std::cout <<
"!!!! Point not in defined region (radius too large in specific z-plane)" << std::endl;
344 int r_index1 = r_index0 + 1;
345 if (r_index1 >= (
int)
r_map_.size())
349 std::cout <<
"!!!! Point not in defined region (radius too large in specific z-plane)" << std::endl;
357 std::vector<float>::const_iterator phiiter = upper_bound(
phi_map_.begin(),
phi_map_.end(),
phi);
359 int phi_index1 = phi_index0 + 1;
360 if (phi_index1 >= (
int)
phi_map_.size())
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];
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];
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];
396 double zweight = z -
z_map_[z_index0];
397 double zspacing = z_map_[z_index1] - z_map_[z_index0];
400 double rweight = r -
r_map_[r_index0];
401 double rspacing = r_map_[r_index1] - r_map_[r_index0];
404 double phiweight = phi -
phi_map_[phi_index0];
405 double phispacing = phi_map_[phi_index1] - phi_map_[phi_index0];
408 phispacing += 2 * M_PI;
410 phiweight /= phispacing;
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));
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));
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));
446 std::cout <<
"End GFCyl Call: <bz,br,bphi> : {"
447 << BfieldCyl[0] /
gauss <<
"," << BfieldCyl[1] /
gauss <<
"," << BfieldCyl[2] /
gauss <<
"}"
467 unsigned middle = start + ((end -
start) / 2);
469 if (vec[middle] == key)
474 else if (vec[middle] > key)
476 return bin_search(vec, start, middle - 1, key, index);
478 return bin_search(vec, middle + 1, end, key, index);
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 <<
">"
490 << std::get<0>(it->second) *
gauss <<
","
491 << std::get<1>(it->second) *
gauss <<
","
492 << std::get<2>(it->second) *
gauss <<
">" << std::endl;