File indexing completed on 2025-08-06 08:17:58
0001 #include "PHField3DCylindrical.h"
0002
0003 #include <TDirectory.h> // for TDirectory, gDirectory
0004 #include <TFile.h>
0005 #include <TNtuple.h>
0006
0007 #include <Geant4/G4SystemOfUnits.hh>
0008
0009 #include <algorithm>
0010 #include <cassert>
0011 #include <cmath>
0012 #include <cstdlib>
0013 #include <iostream>
0014 #include <iterator>
0015 #include <set>
0016 #include <utility>
0017
0018 PHField3DCylindrical::PHField3DCylindrical(const std::string &filename, const int verb, const float magfield_rescale)
0019 : PHField(verb)
0020 {
0021 std::cout << "\n================ Begin Construct Mag Field =====================" << std::endl;
0022 std::cout << "\n-----------------------------------------------------------"
0023 << "\n Magnetic field Module - Verbosity:" << Verbosity()
0024 << "\n-----------------------------------------------------------";
0025
0026
0027 TFile *rootinput = TFile::Open(filename.c_str());
0028 if (!rootinput)
0029 {
0030 std::cout << "\n could not open " << filename << " exiting now" << std::endl;
0031 exit(1);
0032 }
0033 std::cout << "\n ---> "
0034 "Reading the field grid from "
0035 << filename << " ... " << std::endl;
0036 rootinput->cd();
0037
0038
0039 TNtuple *field_map = (TNtuple *) gDirectory->Get("map");
0040 Float_t ROOT_Z, ROOT_R, ROOT_PHI;
0041 Float_t ROOT_BZ, ROOT_BR, ROOT_BPHI;
0042 field_map->SetBranchAddress("z", &ROOT_Z);
0043 field_map->SetBranchAddress("r", &ROOT_R);
0044 field_map->SetBranchAddress("phi", &ROOT_PHI);
0045 field_map->SetBranchAddress("bz", &ROOT_BZ);
0046 field_map->SetBranchAddress("br", &ROOT_BR);
0047 field_map->SetBranchAddress("bphi", &ROOT_BPHI);
0048
0049
0050 int nz, nr, nphi;
0051 nz = field_map->GetEntries("z>-1e6");
0052 nr = field_map->GetEntries("r>-1e6");
0053 nphi = field_map->GetEntries("phi>-1e6");
0054 static const int NENTRIES = field_map->GetEntries();
0055
0056
0057 std::cout << " ---> The field grid contained " << NENTRIES << " entries" << std::endl;
0058 if (Verbosity() > 0)
0059 {
0060 std::cout << "\n NENTRIES should be the same as the following values:"
0061 << "\n [ Number of values r,phi,z: "
0062 << nr << " " << nphi << " " << nz << " ]! " << std::endl;
0063 }
0064
0065 if (nz != nr || nz != nphi || nr != nphi)
0066 {
0067 std::cout << "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!"
0068 << "\n The file you entered is not a \"table\" of values"
0069 << "\n Something very likely went oh so wrong"
0070 << "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!" << std::endl;
0071 }
0072
0073
0074 std::set<float> z_set, r_set, phi_set;
0075
0076
0077
0078
0079
0080
0081 if (Verbosity() > 0)
0082 {
0083 std::cout << " --> Sorting Entries..." << std::endl;
0084 }
0085 std::map<trio, trio> sorted_map;
0086 for (int i = 0; i < field_map->GetEntries(); i++)
0087 {
0088 field_map->GetEntry(i);
0089 trio coord_key(ROOT_Z, ROOT_R, ROOT_PHI);
0090 trio field_val(ROOT_BZ, ROOT_BR, ROOT_BPHI);
0091 sorted_map[coord_key] = field_val;
0092
0093 z_set.insert(ROOT_Z * cm);
0094 r_set.insert(ROOT_R * cm);
0095 phi_set.insert(ROOT_PHI * deg);
0096 }
0097
0098
0099 if (Verbosity() > 4)
0100 {
0101 std::map<trio, trio>::iterator it = sorted_map.begin();
0102 print_map(it);
0103 float last_z = std::get<0>(it->first);
0104 for (it = sorted_map.begin(); it != sorted_map.end(); ++it)
0105 {
0106 if (std::get<0>(it->first) != last_z)
0107 {
0108 last_z = std::get<0>(it->first);
0109 print_map(it);
0110 }
0111 }
0112 }
0113
0114 if (Verbosity() > 0)
0115 {
0116 std::cout << " --> Putting entries into containers... " << std::endl;
0117 }
0118
0119
0120 minz_ = *(z_set.begin());
0121 std::set<float>::iterator ziter = z_set.end();
0122 --ziter;
0123 maxz_ = *ziter;
0124
0125
0126 nz = z_set.size();
0127 nr = r_set.size();
0128 nphi = phi_set.size();
0129 r_map_.resize(nr, 0);
0130 z_map_.resize(nz, 0);
0131 phi_map_.resize(nphi, 0);
0132
0133 std::copy(z_set.begin(), z_set.end(), z_map_.begin());
0134 std::copy(phi_set.begin(), phi_set.end(), phi_map_.begin());
0135 std::copy(r_set.begin(), r_set.end(), r_map_.begin());
0136
0137
0138 BFieldR_.resize(nz, std::vector<std::vector<float> >(nr, std::vector<float>(nphi, 0)));
0139 BFieldPHI_.resize(nz, std::vector<std::vector<float> >(nr, std::vector<float>(nphi, 0)));
0140 BFieldZ_.resize(nz, std::vector<std::vector<float> >(nr, std::vector<float>(nphi, 0)));
0141
0142
0143 unsigned int ir = 0, iphi = 0, iz = 0;
0144 std::map<trio, trio>::iterator iter = sorted_map.begin();
0145 for (; iter != sorted_map.end(); ++iter)
0146 {
0147
0148 float z = std::get<0>(iter->first) * cm;
0149 float r = std::get<1>(iter->first) * cm;
0150 float phi = std::get<2>(iter->first) * deg;
0151 float Bz = std::get<0>(iter->second) * gauss;
0152 float Br = std::get<1>(iter->second) * gauss;
0153 float Bphi = std::get<2>(iter->second) * gauss;
0154
0155 if (z > maxz_)
0156 {
0157 maxz_ = z;
0158 }
0159 if (z < minz_)
0160 {
0161 minz_ = z;
0162 }
0163
0164
0165 if (z != z_map_[iz])
0166 {
0167 ++iz;
0168 ir = 0;
0169 iphi = 0;
0170 }
0171 else if (r != r_map_[ir])
0172 {
0173 ++ir;
0174 iphi = 0;
0175 }
0176 else if (phi != phi_map_[iphi])
0177 {
0178 ++iphi;
0179 }
0180
0181
0182 if (iz > 0 && z < z_map_[iz - 1])
0183 {
0184 std::cout << "!!!!!!!!! Your map isn't ordered.... z: " << z << " zprev: " << z_map_[iz - 1] << std::endl;
0185 }
0186
0187 BFieldR_[iz][ir][iphi] = Br * magfield_rescale;
0188 BFieldPHI_[iz][ir][iphi] = Bphi * magfield_rescale;
0189 BFieldZ_[iz][ir][iphi] = Bz * magfield_rescale;
0190
0191
0192
0193
0194 if (std::fabs(z) < 10 && ir < 10 && Verbosity() > 3)
0195 {
0196 print_map(iter);
0197
0198 std::cout << " B("
0199 << r_map_[ir] << ", "
0200 << phi_map_[iphi] << ", "
0201 << z_map_[iz] << "): ("
0202 << BFieldR_[iz][ir][iphi] << ", "
0203 << BFieldPHI_[iz][ir][iphi] << ", "
0204 << BFieldZ_[iz][ir][iphi] << ")" << std::endl;
0205 }
0206
0207 }
0208
0209 rootinput->Close();
0210
0211 std::cout << "\n ---> ... read file successfully "
0212 << "\n ---> Z Boundaries ~ zlow, zhigh: "
0213 << minz_ / cm << "," << maxz_ / cm << " cm " << std::endl;
0214
0215 std::cout << "\n================= End Construct Mag Field ======================\n"
0216 << std::endl;
0217 }
0218
0219 void PHField3DCylindrical::GetFieldValue(const double point[4], double *Bfield) const
0220 {
0221 if (Verbosity() > 2)
0222 {
0223 std::cout << "\nPHField3DCylindrical::GetFieldValue" << std::endl;
0224 }
0225 double x = point[0];
0226 double y = point[1];
0227 double z = point[2];
0228 double r = sqrt(x * x + y * y);
0229 double phi;
0230 if (x == 0)
0231 {
0232 phi = atan2(y, 0.00000000001);
0233 }
0234 else
0235 {
0236 phi = atan2(y, x);
0237 }
0238 if (phi < 0)
0239 {
0240 phi += 2 * M_PI;
0241 }
0242
0243
0244 if ((z >= minz_) && (z <= maxz_))
0245 {
0246 double BFieldCyl[3];
0247 double cylpoint[4] = {z, r, phi, 0};
0248
0249
0250 GetFieldCyl(cylpoint, BFieldCyl);
0251
0252
0253 Bfield[0] = cos(phi) * BFieldCyl[1] - sin(phi) * BFieldCyl[2];
0254
0255
0256 Bfield[1] = sin(phi) * BFieldCyl[1] + cos(phi) * BFieldCyl[2];
0257
0258
0259 Bfield[2] = BFieldCyl[0];
0260 }
0261 else
0262 {
0263
0264 Bfield[0] = 0.0;
0265 Bfield[1] = 0.0;
0266 Bfield[2] = 0.0;
0267 if (Verbosity() > 2)
0268 {
0269 std::cout << "!!!!!!!!!! Field point not in defined region (outside of z bounds)" << std::endl;
0270 }
0271 }
0272
0273 if (Verbosity() > 2)
0274 {
0275 std::cout << "END PHField3DCylindrical::GetFieldValue\n"
0276 << " ---> {Bx, By, Bz} : "
0277 << "< " << Bfield[0] << ", " << Bfield[1] << ", " << Bfield[2] << " >" << std::endl;
0278 }
0279
0280 return;
0281 }
0282
0283 void PHField3DCylindrical::GetFieldCyl(const double CylPoint[4], double *BfieldCyl) const
0284 {
0285 float z = CylPoint[0];
0286 float r = CylPoint[1];
0287 float phi = CylPoint[2];
0288
0289 BfieldCyl[0] = 0.0;
0290 BfieldCyl[1] = 0.0;
0291 BfieldCyl[2] = 0.0;
0292
0293 if (Verbosity() > 2)
0294 {
0295 std::cout << "GetFieldCyl@ <z,r,phi>: {" << z << "," << r << "," << phi << "}" << std::endl;
0296 }
0297
0298 if (z <= z_map_[0] || z >= z_map_[z_map_.size() - 1])
0299 {
0300 if (Verbosity() > 2)
0301 {
0302 std::cout << "!!!! Point not in defined region (|z| too large)" << std::endl;
0303 }
0304 return;
0305 }
0306 if (r < r_map_[0])
0307 {
0308 r = r_map_[0];
0309 if (Verbosity() > 2)
0310 {
0311 std::cout << "!!!! Point not in defined region (radius too small in specific z-plane). Use min radius" << std::endl;
0312 }
0313
0314 }
0315 if (r > r_map_[r_map_.size() - 1])
0316 {
0317 if (Verbosity() > 2)
0318 {
0319 std::cout << "!!!! Point not in defined region (radius too large in specific z-plane)" << std::endl;
0320 }
0321 return;
0322 }
0323
0324 std::vector<float>::const_iterator ziter = upper_bound(z_map_.begin(), z_map_.end(), z);
0325 int z_index0 = distance(z_map_.begin(), ziter) - 1;
0326 int z_index1 = z_index0 + 1;
0327
0328 assert(z_index0 >= 0);
0329 assert(z_index1 >= 0);
0330 assert(z_index0 < (int) z_map_.size());
0331 assert(z_index1 < (int) z_map_.size());
0332
0333 std::vector<float>::const_iterator riter = upper_bound(r_map_.begin(), r_map_.end(), r);
0334 int r_index0 = distance(r_map_.begin(), riter) - 1;
0335 if (r_index0 >= (int) r_map_.size())
0336 {
0337 if (Verbosity() > 2)
0338 {
0339 std::cout << "!!!! Point not in defined region (radius too large in specific z-plane)" << std::endl;
0340 }
0341 return;
0342 }
0343
0344 int r_index1 = r_index0 + 1;
0345 if (r_index1 >= (int) r_map_.size())
0346 {
0347 if (Verbosity() > 2)
0348 {
0349 std::cout << "!!!! Point not in defined region (radius too large in specific z-plane)" << std::endl;
0350 }
0351 return;
0352 }
0353
0354 assert(r_index0 >= 0);
0355 assert(r_index1 >= 0);
0356
0357 std::vector<float>::const_iterator phiiter = upper_bound(phi_map_.begin(), phi_map_.end(), phi);
0358 int phi_index0 = distance(phi_map_.begin(), phiiter) - 1;
0359 int phi_index1 = phi_index0 + 1;
0360 if (phi_index1 >= (int) phi_map_.size())
0361 {
0362 phi_index1 = 0;
0363 }
0364
0365 assert(phi_index0 >= 0);
0366 assert(phi_index0 < (int) phi_map_.size());
0367 assert(phi_index1 >= 0);
0368
0369 double Br000 = BFieldR_[z_index0][r_index0][phi_index0];
0370 double Br001 = BFieldR_[z_index0][r_index0][phi_index1];
0371 double Br010 = BFieldR_[z_index0][r_index1][phi_index0];
0372 double Br011 = BFieldR_[z_index0][r_index1][phi_index1];
0373 double Br100 = BFieldR_[z_index1][r_index0][phi_index0];
0374 double Br101 = BFieldR_[z_index1][r_index0][phi_index1];
0375 double Br110 = BFieldR_[z_index1][r_index1][phi_index0];
0376 double Br111 = BFieldR_[z_index1][r_index1][phi_index1];
0377
0378 double Bphi000 = BFieldPHI_[z_index0][r_index0][phi_index0];
0379 double Bphi001 = BFieldPHI_[z_index0][r_index0][phi_index1];
0380 double Bphi010 = BFieldPHI_[z_index0][r_index1][phi_index0];
0381 double Bphi011 = BFieldPHI_[z_index0][r_index1][phi_index1];
0382 double Bphi100 = BFieldPHI_[z_index1][r_index0][phi_index0];
0383 double Bphi101 = BFieldPHI_[z_index1][r_index0][phi_index1];
0384 double Bphi110 = BFieldPHI_[z_index1][r_index1][phi_index0];
0385 double Bphi111 = BFieldPHI_[z_index1][r_index1][phi_index1];
0386
0387 double Bz000 = BFieldZ_[z_index0][r_index0][phi_index0];
0388 double Bz001 = BFieldZ_[z_index0][r_index0][phi_index1];
0389 double Bz100 = BFieldZ_[z_index1][r_index0][phi_index0];
0390 double Bz101 = BFieldZ_[z_index1][r_index0][phi_index1];
0391 double Bz010 = BFieldZ_[z_index0][r_index1][phi_index0];
0392 double Bz110 = BFieldZ_[z_index1][r_index1][phi_index0];
0393 double Bz011 = BFieldZ_[z_index0][r_index1][phi_index1];
0394 double Bz111 = BFieldZ_[z_index1][r_index1][phi_index1];
0395
0396 double zweight = z - z_map_[z_index0];
0397 double zspacing = z_map_[z_index1] - z_map_[z_index0];
0398 zweight /= zspacing;
0399
0400 double rweight = r - r_map_[r_index0];
0401 double rspacing = r_map_[r_index1] - r_map_[r_index0];
0402 rweight /= rspacing;
0403
0404 double phiweight = phi - phi_map_[phi_index0];
0405 double phispacing = phi_map_[phi_index1] - phi_map_[phi_index0];
0406 if (phi_index1 == 0)
0407 {
0408 phispacing += 2 * M_PI;
0409 }
0410 phiweight /= phispacing;
0411
0412
0413 BfieldCyl[0] =
0414 (1 - zweight) * ((1 - rweight) * ((1 - phiweight) * Bz000 + phiweight * Bz001) +
0415 rweight * ((1 - phiweight) * Bz010 + phiweight * Bz011)) +
0416 zweight * ((1 - rweight) * ((1 - phiweight) * Bz100 + phiweight * Bz101) +
0417 rweight * ((1 - phiweight) * Bz110 + phiweight * Bz111));
0418
0419
0420 BfieldCyl[1] =
0421 (1 - zweight) * ((1 - rweight) * ((1 - phiweight) * Br000 + phiweight * Br001) +
0422 rweight * ((1 - phiweight) * Br010 + phiweight * Br011)) +
0423 zweight * ((1 - rweight) * ((1 - phiweight) * Br100 + phiweight * Br101) +
0424 rweight * ((1 - phiweight) * Br110 + phiweight * Br111));
0425
0426
0427 BfieldCyl[2] =
0428 (1 - zweight) * ((1 - rweight) * ((1 - phiweight) * Bphi000 + phiweight * Bphi001) +
0429 rweight * ((1 - phiweight) * Bphi010 + phiweight * Bphi011)) +
0430 zweight * ((1 - rweight) * ((1 - phiweight) * Bphi100 + phiweight * Bphi101) +
0431 rweight * ((1 - phiweight) * Bphi110 + phiweight * Bphi111));
0432
0433
0434
0435
0436
0437
0438
0439
0440
0441
0442
0443
0444 if (Verbosity() > 2)
0445 {
0446 std::cout << "End GFCyl Call: <bz,br,bphi> : {"
0447 << BfieldCyl[0] / gauss << "," << BfieldCyl[1] / gauss << "," << BfieldCyl[2] / gauss << "}"
0448 << std::endl;
0449 }
0450
0451 return;
0452 }
0453
0454
0455
0456 bool PHField3DCylindrical::bin_search(const std::vector<float> &vec, unsigned start, unsigned end, const float &key, unsigned &index) const
0457 {
0458
0459 if (start > end)
0460 {
0461 index = start;
0462 return false;
0463 }
0464
0465
0466
0467 unsigned middle = start + ((end - start) / 2);
0468
0469 if (vec[middle] == key)
0470 {
0471 index = middle;
0472 return true;
0473 }
0474 else if (vec[middle] > key)
0475 {
0476 return bin_search(vec, start, middle - 1, key, index);
0477 }
0478 return bin_search(vec, middle + 1, end, key, index);
0479 }
0480
0481
0482 void PHField3DCylindrical::print_map(std::map<trio, trio>::iterator &it) const
0483 {
0484 std::cout << " Key: <"
0485 << std::get<0>(it->first) * cm << ","
0486 << std::get<1>(it->first) * cm << ","
0487 << std::get<2>(it->first) * deg << ">"
0488
0489 << " Value: <"
0490 << std::get<0>(it->second) * gauss << ","
0491 << std::get<1>(it->second) * gauss << ","
0492 << std::get<2>(it->second) * gauss << ">" << std::endl;
0493 }