235#include "TRestAxionMagneticField.h"
239#include "TRestPhysics.h"
266 RESTDebug <<
"Entering TRestAxionMagneticField constructor( cfgFileName, name )" <<
RESTendl;
279 RESTDebug <<
"Entering ... TRestAxionMagneticField() destructor." <<
RESTendl;
338 Double_t step, TString style, Double_t depth) {
339 Double_t step_x, step_y, step_z;
352 if (volIndex < 0) volIndex = 0;
355 RESTError << volIndex <<
" corresponds to none volume index " <<
RESTendl;
357 RESTError <<
"Setting volIndex to the first volume" <<
RESTendl;
366 step_x = step_y = step_z = step;
371 if (!(projection ==
"XY" || projection ==
"XZ" || projection ==
"YZ")) {
372 RESTError <<
"You entered : " << projection <<
" as a projection but you have to choose XY, XZ or YZ"
381 Double_t halfSizeX = vol->mesh.GetNetSizeX() / 2.;
382 Double_t halfSizeY = vol->mesh.GetNetSizeY() / 2.;
383 Double_t halfSizeZ = vol->mesh.GetNetSizeZ() / 2.;
385 Double_t xMin = centerX - halfSizeX;
386 Double_t yMin = centerY - halfSizeY;
387 Double_t zMin = centerZ - halfSizeZ;
389 Double_t xMax = centerX + halfSizeX;
390 Double_t yMax = centerY + halfSizeY;
391 Double_t zMax = centerZ + halfSizeZ;
393 Int_t nBinsX = (xMax - xMin) / step_x;
394 Int_t nBinsY = (yMax - yMin) / step_y;
395 Int_t nBinsZ = (zMax - zMin) / step_z;
397 Double_t x = -1, y = -1, z = -1;
401 if (projection ==
"XY") {
402 fCanvas =
new TCanvas(
"fCanvas",
"");
403 fHisto =
new TH2D(
"",
"", nBinsX, xMin, xMax, nBinsY, yMin, yMax);
405 if (depth < -100000.0)
406 z = (zMin + zMax) / 2.0;
407 else if ((depth >= zMin) && (depth <= zMax))
410 RESTError <<
"You entered depth = " << depth <<
", but you have to choose depth between " << zMin
414 for (Int_t i = 0; i < nBinsX; i++) {
416 for (Int_t j = 0; j < nBinsY; j++) {
427 RESTError <<
"You entered : " << Bcomp
428 <<
" as a B component but you have to choose X, Y or Z" <<
RESTendl;
438 fHisto->SetBit(TH1::kNoStats);
439 fHisto->GetXaxis()->SetTitle(
"x (mm)");
440 fHisto->GetYaxis()->SetTitle(
"y (mm)");
443 TString title = Form(
"B_{x} against x and y @ z = %.1f", z);
448 TString title = Form(
"B_{y} against x and y @ z = %.1f", z);
453 TString title = Form(
"B_{z} against x and y @ z = %.1f", z);
462 if (projection ==
"XZ") {
463 TCanvas*
fCanvas =
new TCanvas(
"fCanvas",
"");
464 fHisto =
new TH2D(
"",
"", nBinsX, xMin, xMax, nBinsZ, zMin, zMax);
466 if (depth < -100000.0)
467 y = (yMin + yMax) / 2.0;
468 else if ((depth >= yMin) && (depth <= yMax))
471 RESTError <<
"You entered depth = " << depth <<
", but you have to choose depth between " << yMin
475 for (Int_t i = 0; i < nBinsX; i++) {
477 for (Int_t j = 0; j < nBinsZ; j++) {
488 RESTError <<
"You entered : " << Bcomp
489 <<
" as a B component but you have to choose X, Y or Z" <<
RESTendl;
499 fHisto->SetBit(TH1::kNoStats);
500 fHisto->GetXaxis()->SetTitle(
"x (mm)");
501 fHisto->GetYaxis()->SetTitle(
"z (mm)");
504 TString title = Form(
"B_{x} against x and z @ y = %.1f", y);
509 TString title = Form(
"B_{y} against x and z @ y = %.1f", y);
514 TString title = Form(
"B_{z} against x and z @ y = %.1f", y);
523 if (projection ==
"YZ") {
524 TCanvas*
fCanvas =
new TCanvas(
"fCanvas",
"");
525 fHisto =
new TH2D(
"",
"", nBinsY, yMin, yMax, nBinsZ, zMin, zMax);
527 if (depth < -100000.0)
528 x = (xMin + xMax) / 2.0;
529 else if ((depth >= xMin) && (depth <= xMax))
532 RESTError <<
"You entered depth = " << depth <<
", but you have to choose depth between " << xMin
536 for (Int_t i = 0; i < nBinsY; i++) {
538 for (Int_t j = 0; j < nBinsZ; j++) {
549 RESTError <<
"You entered : " << Bcomp
550 <<
" as a B component but you have to choose X, Y or Z" <<
RESTendl;
560 fHisto->SetBit(TH1::kNoStats);
561 fHisto->GetXaxis()->SetTitle(
"y (mm)");
562 fHisto->GetYaxis()->SetTitle(
"z (mm)");
565 TString title = Form(
"B_{x} against y and z @ x = %.1f", x);
570 TString title = Form(
"B_{y} against y and z @ x = %.1f", x);
575 TString title = Form(
"B_{z} against y and z @ x = %.1f", x);
592 Int_t divisions = 100;
597 fCanvas =
new TCanvas(
"fCanvas",
"", 1600, 600);
599 TPad* pad1 =
new TPad(
"pad1",
"This is pad1", 0.01, 0.02, 0.99, 0.97);
607 TVector3 position(0, 0, genPositionZ);
608 TVector3 direction(0, 0, 1);
611 RESTDebug <<
"Start moving along" <<
RESTendl;
612 RESTDebug <<
"++++++++++++++++++" <<
RESTendl;
614 TGraph* fieldGr =
new TGraph();
616 Double_t delta =
fBoundMax[volId][2] * 2. / divisions;
621 position =
MoveToPlane(position, direction, TVector3(0, 0, 1), posAlongAxis);
624 fieldGr->SetPoint(fieldGr->GetN(), posZ, fieldY);
629 fieldGr->SetLineWidth(3);
630 fieldGr->SetLineColor(38 + n);
631 fieldGr->GetXaxis()->SetLimits(genPositionZ - 500, finalPositionZ + 500);
632 fieldGr->GetHistogram()->SetMaximum(2.5);
633 fieldGr->GetHistogram()->SetMinimum(0);
634 fieldGr->GetXaxis()->SetTitle(
"Z [mm]");
635 fieldGr->GetXaxis()->SetTitleSize(0.05);
636 fieldGr->GetXaxis()->SetLabelSize(0.05);
637 fieldGr->GetYaxis()->SetTitle(
"B [T]");
638 fieldGr->GetYaxis()->SetTitleOffset(1.3);
639 fieldGr->GetYaxis()->SetTitleSize(0.05);
640 fieldGr->GetYaxis()->SetLabelSize(0.05);
655 fCanvas =
new TCanvas(
"fCanvas",
"", 1600, 600);
657 TPad* pad1 =
new TPad(
"pad1",
"This is pad1", 0.01, 0.02, 0.99, 0.97);
663 Double_t genSizeY =
fBoundMax[volId].Y() * 3.;
668 TGraph* bBox =
new TGraph();
680 RESTDebug <<
"Gen position : " << genPositionZ <<
RESTendl;
682 bBox->GetXaxis()->SetLimits(genPositionZ - 500, finalPositionZ + 500);
683 bBox->GetHistogram()->SetMaximum(genSizeY + 100);
684 bBox->GetHistogram()->SetMinimum(-genSizeY - 100);
686 bBox->GetXaxis()->SetTitle(
"Z [mm]");
687 bBox->GetXaxis()->SetTitleSize(0.05);
688 bBox->GetXaxis()->SetLabelSize(0.05);
689 bBox->GetYaxis()->SetTitle(
"Y [mm]");
690 bBox->GetYaxis()->SetTitleOffset(1.3);
691 bBox->GetYaxis()->SetTitleSize(0.05);
692 bBox->GetYaxis()->SetLabelSize(0.05);
693 bBox->SetLineWidth(2);
697 for (Double_t y = genSizeY; y >= -genSizeY; y -=
fBoundMax[volId].Y()) {
698 TVector3 position(0, y, genPositionZ);
699 TVector3 direction = (vanishingPoint - position).Unit();
702 TGraph* grB =
new TGraph();
703 grB->SetPoint(0, trackBounds[0].Z(), trackBounds[0].Y());
704 grB->SetPoint(1, trackBounds[1].Z(), trackBounds[1].Y());
711 RESTDebug <<
"Start moving along" <<
RESTendl;
712 RESTDebug <<
"++++++++++++++++++" <<
RESTendl;
714 TGraph* fieldGr =
new TGraph();
716 Double_t delta =
fBoundMax[volId][2] * 2. / divisions;
719 fieldGr->SetPoint(fieldGr->GetN(), posZ, field);
724 position =
MoveToPlane(position, direction, TVector3(0, 0, 1), posAlongAxis);
727 fieldGr->SetPoint(fieldGr->GetN(), posZ, field);
733 position =
MoveToPlane(position, direction, TVector3(0, 0, 1), posAlongAxis);
736 fieldGr->SetPoint(fieldGr->GetN(), posZ, field2);
739 fieldGr->SetLineWidth(3);
740 fieldGr->SetLineColor(38 + n);
741 fieldGr->GetXaxis()->SetLimits(genPositionZ - 500, finalPositionZ + 500);
742 fieldGr->GetHistogram()->SetMaximum(2.5);
743 fieldGr->GetHistogram()->SetMinimum(0);
744 fieldGr->GetXaxis()->SetTitle(
"Z [mm]");
745 fieldGr->GetXaxis()->SetTitleSize(0.05);
746 fieldGr->GetXaxis()->SetLabelSize(0.05);
747 fieldGr->GetYaxis()->SetTitle(
"B [T]");
748 fieldGr->GetYaxis()->SetTitleOffset(1.3);
749 fieldGr->GetYaxis()->SetTitleSize(0.05);
750 fieldGr->GetYaxis()->SetLabelSize(0.05);
757 TGraph* gr =
new TGraph();
758 gr->SetPoint(0, genPositionZ, y);
759 position =
MoveToPlane(position, direction, TVector3(0, 0, 1), TVector3(0, 0, finalPositionZ));
760 gr->SetPoint(1, finalPositionZ, position.Y());
764 grB->SetLineColor(38 + n);
766 grB->SetLineWidth(3);
781 std::vector<std::vector<Float_t>> data) {
782 mVol.field.resize(mVol.mesh.GetNodesX());
783 for (
unsigned int n = 0; n < mVol.field.size(); n++) {
784 mVol.field[n].resize(mVol.mesh.GetNodesY());
785 for (
unsigned int m = 0; m < mVol.field[n].size(); m++)
786 mVol.field[n][m].resize(mVol.mesh.GetNodesZ());
789 RESTDebug <<
"TRestAxionMagneticField::LoadMagneticFieldData. Printing first 5 data rows" <<
RESTendl;
790 for (
unsigned int n = 0; n < data.size(); n++) {
794 Int_t nX = mVol.mesh.GetNodeX((Int_t)(data[n][0] + mVol.mesh.GetNetSizeX() / 2.),
true);
795 Int_t nY = mVol.mesh.GetNodeY((Int_t)(data[n][1] + mVol.mesh.GetNetSizeY() / 2.),
true);
796 Int_t nZ = mVol.mesh.GetNodeZ((Int_t)(data[n][2] + mVol.mesh.GetNetSizeZ() / 2.),
true);
799 RESTDebug <<
"X: " << data[n][0] <<
" Y: " << data[n][1] <<
" Z: " << data[n][2] <<
RESTendl;
800 RESTDebug <<
"absX: " << data[n][0] + mVol.position.X()
801 <<
" absY: " << data[n][1] + mVol.position.Y()
802 <<
" absZ: " << data[n][2] + mVol.position.Z() <<
RESTendl;
803 RESTDebug <<
"nX: " << nX <<
" nY: " << nY <<
" nZ: " << nZ <<
RESTendl;
804 RESTDebug <<
"Bx: " << data[n][3] <<
" By: " << data[n][4] <<
" Bz: " << data[n][5] <<
RESTendl;
808 if (mVol.field[nX][nY][nZ] != TVector3(0.0, 0.0, 0.0)) {
809 RESTWarning <<
"X: " << data[n][0] <<
" Y: " << data[n][1] <<
" Z: " << data[n][2] <<
RESTendl;
810 RESTWarning <<
"nX: " << nX <<
" nY: " << nY <<
" nZ: " << nZ <<
RESTendl;
811 RESTWarning <<
"WARNING: field[nX][nY][nZ] element not equal to initial value (0, 0, 0) !!"
813 RESTWarning <<
"It has value: "
814 <<
"mVol.field[" << nX <<
"][" << nY <<
"][" << nZ <<
"] = ("
815 << mVol.field[nX][nY][nZ].X() <<
" , " << mVol.field[nX][nY][nZ].Y() <<
" , "
816 << mVol.field[nX][nY][nZ].Z() <<
")" <<
RESTendl;
817 RESTWarning <<
"Values to write: "
818 <<
"Bx: " << data[n][3] <<
" By: " << data[n][4] <<
" Bz: " << data[n][5] <<
RESTendl
821 this->
SetError(
"There was a problem assigning the field matrix!");
825 mVol.field[nX][nY][nZ] = TVector3(data[n][3], data[n][4], data[n][5]);
835 for (
unsigned int n = 0; n <
fPositions.size(); n++) {
838 RESTDebug <<
"Full path : " << fullPathName <<
RESTendl;
840 if (
fFileNames[n] !=
"none" && fullPathName ==
"") {
841 RESTError <<
"TRestAxionMagneticField::LoadMagneticVolumes. File " <<
fFileNames[n]
844 <<
"REST will look for this file at any path given by <searchPath at globals definitions"
849 std::vector<std::vector<Float_t>> fieldData;
851 if (fullPathName.find(
".dat") != string::npos) {
852 RESTDebug <<
"Reading ASCII format" <<
RESTendl;
854 RESTError <<
"Problem reading file : " << fullPathName <<
RESTendl;
858 if (fullPathName.find(
".bin") != string::npos) {
859 RESTDebug <<
"Reading binary format" <<
RESTendl;
861 RESTError <<
"Problem reading file : " << fullPathName <<
RESTendl;
867 RESTError <<
"Filename : " << fullPathName <<
RESTendl;
868 RESTError <<
"File format not recognized!" <<
RESTendl;
872 if (
fFileNames[n] !=
"none" && fieldData.size() < 2) {
873 RESTError <<
"Field data size is no more than 2 grid points!" <<
RESTendl;
874 RESTError <<
"Filename : " << fullPathName <<
RESTendl;
875 RESTError <<
"Probably something went wrong loading the file" <<
RESTendl;
884 if (fieldData.size() > 0) {
885 RESTDebug <<
"Reading max boundary values" <<
RESTendl;
891 if (
fBoundMax[n] != TVector3(xMax, yMax, zMax)) {
892 RESTWarning <<
"Volume : " << n <<
RESTendl;
893 RESTWarning <<
"boundMax was defined in RML but does not match the field map boundaries!"
895 RESTWarning <<
"Max. Field map boundaries : (" << xMax <<
", " << yMax <<
", " << zMax
900 RESTDebug <<
"Reading min boundary values" <<
RESTendl;
906 if (-
fBoundMax[n] != TVector3(xMin, yMin, zMin)) {
907 RESTWarning <<
"Volume : " << n <<
RESTendl;
908 RESTWarning <<
"boundMax was defined in RML but does not match the field map boundaries"
910 RESTWarning <<
"Min. Field map boundaries : (" << xMin <<
", " << yMin <<
", " << zMin
914 fBoundMax[n] = TVector3(xMax, yMax, zMax);
916 RESTDebug <<
"Reading mesh size" <<
RESTendl;
922 if (
fMeshSize[n] != TVector3(meshSizeX, meshSizeY, meshSizeZ)) {
923 RESTWarning <<
"Volume : " << n <<
RESTendl;
925 <<
"MeshSize was defined in RML but does not match the mesh size deduced from "
928 RESTWarning <<
"Mesh size : (" << meshSizeX <<
", " << meshSizeY <<
", " << meshSizeZ
932 fMeshSize[n] = TVector3(meshSizeX, meshSizeY, meshSizeZ);
936 RESTDebug <<
"Reading magnetic field map" <<
RESTendl;
937 RESTDebug <<
"--------------------------" <<
RESTendl;
939 RESTDebug <<
"Full path : " << fullPathName <<
RESTendl;
941 RESTDebug <<
"Boundaries" <<
RESTendl;
942 RESTDebug <<
"xMin: " << xMin <<
" yMin: " << yMin <<
" zMin: " << zMin <<
RESTendl;
943 RESTDebug <<
"xMax: " << xMax <<
" yMax: " << yMax <<
" zMax: " << zMax <<
RESTendl;
944 RESTDebug <<
"Mesh size" <<
RESTendl;
946 RESTDebug <<
"sX: " << meshSizeX <<
" sY: " << meshSizeY <<
" sZ: " << meshSizeZ <<
RESTendl;
948 if (fieldData.size() > 4) {
949 RESTDebug <<
"Printing beginning of magnetic file table : " << fieldData.size() <<
RESTendl;
952 RESTDebug <<
"The data file contains no field map" <<
RESTendl;
958 Int_t nx = (Int_t)(2 * xMax / meshSizeX) + 1;
959 Int_t ny = (Int_t)(2 * yMax / meshSizeY) + 1;
960 Int_t nz = (Int_t)(2 * zMax / meshSizeZ) + 1;
966 restMesh.
SetSize(2 * xMax, 2 * yMax, 2 * zMax);
974 MagneticFieldVolume mVolume;
976 mVolume.mesh = restMesh;
981 RESTError <<
"The bounding box was not defined for volume " << n <<
"!" <<
RESTendl;
982 RESTError <<
"Please review RML configuration for TRestAxionMagneticField" <<
RESTendl;
984 }
else if (
fMeshSize[n] == TVector3(0, 0, 0)) {
985 RESTError <<
"The mesh grid size was not defined for volume " << n <<
"!" <<
RESTendl;
986 RESTError <<
"Please review RML configuration for TRestAxionMagneticField" <<
RESTendl;
993 RESTError <<
"TRestAxionMagneticField::LoadMagneticVolumes. Volumes overlap!" <<
RESTendl;
1000 RESTDebug <<
"Finished loading magnetic volumes" <<
RESTendl;
1023 RESTWarning <<
"TRestAxionMagneticField::GetMagneticField position (" << pos.X() <<
", "
1024 << pos.Y() <<
", " << pos.Z() <<
") is outside any volume" <<
RESTendl;
1025 return TVector3(0, 0, 0);
1029 Int_t nX = node.X();
1030 Int_t nY = node.Y();
1031 Int_t nZ = node.Z();
1033 Int_t nX_1, nY_1, nZ_1;
1063 xd = (pos.X() - x0) / (x1 - x0);
1064 if ((xd < -0.00001) || (xd > 1.00001))
1065 RESTWarning <<
"TRestAxionMagneticField::GetMagneticField Error: xd NOT between 0 and 1"
1074 yd = (pos.Y() - y0) / (y1 - y0);
1075 if ((yd < -0.00001) || (yd > 1.00001))
1076 RESTWarning <<
"TRestAxionMagneticField::GetMagneticField Error: yd NOT between 0 and 1"
1085 zd = (pos.Z() - z0) / (z1 - z0);
1086 if ((zd < -0.00001) || (zd > 1.00001))
1087 RESTWarning <<
"TRestAxionMagneticField::GetMagneticField Error: zd NOT between 0 and 1"
1091 TVector3 C00 = C000 * (1.0 - xd) + C100 * xd;
1092 TVector3 C01 = C001 * (1.0 - xd) + C101 * xd;
1093 TVector3 C10 = C010 * (1.0 - xd) + C110 * xd;
1094 TVector3 C11 = C011 * (1.0 - xd) + C111 * xd;
1097 TVector3 C0 = C00 * (1.0 - yd) + C10 * yd;
1098 TVector3 C1 = C01 * (1.0 - yd) + C11 * yd;
1101 TVector3 C = C0 * (1.0 - zd) + C1 * zd;
1103 RESTDebug <<
"position = (" << pos.X() <<
", " << pos.Y() <<
", " << pos.Z() <<
") ";
1104 RESTDebug <<
"nX = " << nX <<
" nY = " << nY <<
" nZ = " << nZ <<
" nX_1 = " << nX_1
1106 RESTDebug <<
"C000 = (" << C000.X() <<
", " << C000.Y() <<
", " << C000.Z() <<
")" <<
RESTendl
1108 RESTDebug <<
"C100 = (" << C100.X() <<
", " << C100.Y() <<
", " << C100.Z() <<
")" <<
RESTendl
1110 RESTDebug <<
"C010 = (" << C010.X() <<
", " << C010.Y() <<
", " << C010.Z() <<
")" <<
RESTendl
1112 RESTDebug <<
"C110 = (" << C110.X() <<
", " << C110.Y() <<
", " << C110.Z() <<
")" <<
RESTendl
1114 RESTDebug <<
"C001 = (" << C001.X() <<
", " << C001.Y() <<
", " << C001.Z() <<
")" <<
RESTendl
1116 RESTDebug <<
"C101 = (" << C101.X() <<
", " << C101.Y() <<
", " << C101.Z() <<
")" <<
RESTendl
1118 RESTDebug <<
"C011 = (" << C011.X() <<
", " << C011.Y() <<
", " << C011.Z() <<
")" <<
RESTendl
1120 RESTDebug <<
"C111 = (" << C111.X() <<
", " << C111.Y() <<
", " << C111.Z() <<
")" <<
RESTendl
1122 RESTDebug <<
" -------------------------------------------------------" <<
RESTendl;
1123 RESTDebug <<
"C00 = (" << C00.X() <<
", " << C00.Y() <<
", " << C00.Z() <<
")" <<
RESTendl
1125 RESTDebug <<
"C01 = (" << C01.X() <<
", " << C01.Y() <<
", " << C01.Z() <<
")" <<
RESTendl
1127 RESTDebug <<
"C10 = (" << C10.X() <<
", " << C10.Y() <<
", " << C10.Z() <<
")" <<
RESTendl
1129 RESTDebug <<
"C11 = (" << C11.X() <<
", " << C11.Y() <<
", " << C11.Z() <<
")" <<
RESTendl
1131 RESTDebug <<
" -------------------------------------------------------" <<
RESTendl;
1132 RESTDebug <<
"C0 = (" << C0.X() <<
", " << C0.Y() <<
", " << C0.Z() <<
")" <<
RESTendl <<
RESTendl;
1133 RESTDebug <<
"C1 = (" << C1.X() <<
", " << C1.Y() <<
", " << C1.Z() <<
")" <<
RESTendl <<
RESTendl;
1134 RESTDebug <<
" -------------------------------------------------------" <<
RESTendl;
1135 RESTDebug <<
"C = (" << C.X() <<
", " << C.Y() <<
", " << C.Z() <<
")" <<
RESTendl <<
RESTendl;
1147 if (newMapSize.X() == 0 || newMapSize.Y() == 0 || newMapSize.Z() == 0) {
1148 RESTError <<
"TRestAxionMagneticField::ReMap. The mesh granularity cannot be 0" <<
RESTendl;
1149 RESTError <<
"Remapping will not take effect" <<
RESTendl;
1153 Double_t remainder = std::fmod(newMapSize.X(),
fMeshSize[n].X()) +
1154 std::fmod(newMapSize.Y(),
fMeshSize[n].Y()) +
1155 std::fmod(newMapSize.Z(),
fMeshSize[n].Z());
1156 if (remainder != 0) {
1157 RESTError <<
"TRestAxionMagneticField::ReMap. The field cannot be remapped." <<
RESTendl;
1158 RESTError <<
"The new mesh granularity must be a multiple of the existing granularity." <<
RESTendl;
1159 RESTError <<
"Present mesh size : (" <<
fMeshSize[n].X() <<
", " <<
fMeshSize[n].Y() <<
", "
1161 RESTError <<
"Requested mesh size : (" << newMapSize.X() <<
", " << newMapSize.Y() <<
", "
1162 << newMapSize.Z() <<
")" <<
RESTendl;
1163 RESTError <<
"Remapping will not take effect" <<
RESTendl;
1167 Int_t scaleX = (Int_t)(newMapSize.X() /
fMeshSize[n].X());
1168 Int_t scaleY = (Int_t)(newMapSize.Y() /
fMeshSize[n].Y());
1169 Int_t scaleZ = (Int_t)(newMapSize.Z() /
fMeshSize[n].Z());
1175 for (Int_t nx = 0; nx < newNodesX; nx++)
1176 for (Int_t ny = 0; ny < newNodesY; ny++)
1177 for (Int_t nz = 0; nz < newNodesZ; nz++)
1225 RESTWarning <<
"TRestAxionMagneticField::GetVolumePosition. Id : " <<
id <<
" out of range!"
1228 return TVector3(0, 0, 0);
1252 Double_t dl, Int_t Nmax) {
1253 Double_t length = (to - from).Mag();
1257 if (length / dl > Nmax) {
1258 diff = length / Nmax;
1259 RESTWarning <<
"TRestAxionMagneticField::GetTransversalComponentAlongPath. Nmax reached!"
1261 RESTWarning <<
"Nmax = " << Nmax <<
RESTendl;
1262 RESTWarning <<
"Adjusting differential step to : " << diff <<
" mm" <<
RESTendl;
1266 TVector3 direction = (to - from).Unit();
1268 std::vector<Double_t> Bt;
1269 for (Double_t d = 0; d < length; d += diff) {
1289 Double_t dl, Int_t Nmax) {
1290 std::vector<Double_t> Bt;
1291 if (axis != 0 && axis != 1 && axis != 2) {
1292 RESTError <<
"TRestAxionMagneticField::GetComponentAlongPath. Axis must take values 0,1 or 2"
1296 Double_t length = (to - from).Mag();
1300 if (length / dl > Nmax) {
1301 diff = length / Nmax;
1302 RESTWarning <<
"TRestAxionMagneticField::GetComponentAlongPath. Nmax reached!" <<
RESTendl;
1303 RESTWarning <<
"Nmax = " << Nmax <<
RESTendl;
1304 RESTWarning <<
"Adjusting differential step to : " << diff <<
" mm" <<
RESTendl;
1308 TVector3 direction = (to - from).Unit();
1310 for (Double_t d = 0; d < length; d += diff) {
1326 if (trackBounds.size() != 2) {
1334 fTrackLength = (trackBounds[1] - trackBounds[0]).Mag() - 1;
1366 for (
auto& b : Bt) Bavg += b;
1368 if (Bt.size() > 0)
return Bavg / Bt.size();
1370 RESTError <<
"TRestAxionMagneticField::GetTransversalFieldAverage. Lenght is zero!" <<
RESTendl;
1387 Double_t length = (to - from).Mag();
1391 if (length / dl > Nmax) {
1392 diff = length / Nmax;
1393 RESTWarning <<
"TRestAxionMagneticField::GetFieldAverageTransverseVector Nmax reached!"
1395 RESTWarning <<
"Nmax = " << Nmax <<
RESTendl;
1396 RESTWarning <<
"Adjusting differential step to : " << diff <<
" mm" <<
RESTendl;
1400 TVector3 direction = (to - from).Unit();
1401 TVector3 Bavg = TVector3(0.0, 0.0, 0.0);
1402 TVector3 BTavg = TVector3(0.0, 0.0, 0.0);
1403 Int_t numberofpoints = 0;
1405 for (Double_t d = 0; d <= length; d += diff) {
1407 numberofpoints = numberofpoints + 1;
1410 if ((length > 0) && (numberofpoints > 0)) {
1411 Bavg = Bavg * (1.0 / numberofpoints);
1413 Bavg - (Bavg * direction) *
1415 RESTDebug <<
"B average vector = (" << Bavg.x() <<
", " << Bavg.y() <<
", " << Bavg.z() <<
")"
1417 RESTDebug <<
"Transverse B average vector = (" << BTavg.x() <<
", " << BTavg.y() <<
", " << BTavg.z()
1421 RESTError <<
"TRestAxionMagneticField::GetTransversalFieldAverage. Lenght is zero!" <<
RESTendl;
1422 return TVector3(0.0, 0.0, 0.0);
1437 return TVector3(nx, ny, nz);
1444 RESTDebug <<
"Checking overlaps" <<
RESTendl;
1447 if (m == n)
continue;
1448 RESTDebug <<
"Volume : " << m <<
RESTendl;
1451 RESTDebug <<
"Relative bottom vertex : (" << b.X() <<
", " << b.Y() <<
", " << b.Z() <<
")"
1456 RESTDebug <<
"Relative top vertex : (" << t.X() <<
", " << t.Y() <<
", " << t.Z() <<
")"
1480 std::vector<TVector3> boundaries;
1482 if (vol) boundaries = vol->mesh.GetTrackBoundaries(pos, dir);
1498 Double_t precision, Int_t
id) {
1500 if (volumeBoundaries.size() != 2)
return volumeBoundaries;
1505 if (!vol)
return volumeBoundaries;
1509 TVector3 unit = dir.Unit();
1510 std::vector<TVector3> fieldBoundaries;
1512 TVector3 in = volumeBoundaries[0];
1515 if (((volumeBoundaries[1] - in) * dir) > 0)
1516 fieldBoundaries.push_back(in);
1518 return fieldBoundaries;
1520 TVector3 out = volumeBoundaries[1];
1522 (((out - in) * dir) > 0))
1524 if ((((volumeBoundaries[0] - out) * dir) < 0) && (((out - in) * dir) > 0))
1525 fieldBoundaries.push_back(out);
1527 return fieldBoundaries;
1529 return fieldBoundaries;
1538 auto magVolumeDef =
GetElement(
"addMagneticVolume");
1539 while (magVolumeDef) {
1541 if (filename ==
"Not defined")
1546 TVector3 position = Get3DVectorParameterWithUnits(
"position", magVolumeDef);
1547 if (position == TVector3(-1, -1, -1))
1552 TVector3 field = Get3DVectorParameterWithUnits(
"field", magVolumeDef);
1553 if (field == TVector3(-1, -1, -1))
1558 TVector3 boundMax = Get3DVectorParameterWithUnits(
"boundMax", magVolumeDef);
1559 if (boundMax == TVector3(-1, -1, -1))
1564 TVector3 meshSize = Get3DVectorParameterWithUnits(
"meshSize", magVolumeDef);
1565 if (meshSize == TVector3(-1, -1, -1))
1571 if (type ==
"NO_SUCH_PARA")
1578 RESTWarning <<
"Mesh type is cylinder. But X and Y inside boundMax are not the same!" <<
RESTendl;
1579 RESTWarning <<
"Making second bound component Y equal to the X bound component!" <<
RESTendl;
1583 RESTDebug <<
"Reading new magnetic volume" <<
RESTendl;
1585 RESTDebug <<
"Filename : " << filename <<
RESTendl;
1586 RESTDebug <<
"Position: ( " << position.X() <<
", " << position.Y() <<
", " << position.Z() <<
") mm"
1588 RESTDebug <<
"Field: ( " << field.X() <<
", " << field.Y() <<
", " << field.Z() <<
") T" <<
RESTendl;
1589 RESTDebug <<
"Max bounding box ( " << boundMax.X() <<
", " << boundMax.Y() <<
", " << boundMax.Z()
1591 RESTDebug <<
"Mesh size ( " << meshSize.X() <<
", " << meshSize.Y() <<
", " << meshSize.Z() <<
")"
1609 if (
fReMap != TVector3(0, 0, 0)) {
1610 RESTMetadata <<
"Fields original mesh size has been remapped" <<
RESTendl;
1615 RESTMetadata <<
" ------------------------------------------------ " <<
RESTendl;
1617 if (p > 0) RESTMetadata <<
" ------------------------------------------------ " <<
RESTendl;
1627 Double_t xMin = centerX - halfSizeX;
1628 Double_t yMin = centerY - halfSizeY;
1629 Double_t zMin = centerZ - halfSizeZ;
1631 Double_t xMax = centerX + halfSizeX;
1632 Double_t yMax = centerY + halfSizeY;
1633 Double_t zMax = centerZ + halfSizeZ;
1635 RESTMetadata <<
"* Volume " << p <<
" centered at (" << centerX <<
"," << centerY <<
"," << centerZ
1637 RESTMetadata <<
" - Grid mesh element size. X: " <<
fMeshSize[p].X() <<
"mm "
1644 RESTMetadata <<
" - Bounds : " <<
RESTendl;
1645 RESTMetadata <<
" xmin : " << xMin <<
" mm , xmax : " << xMax <<
" mm" <<
RESTendl;
1646 RESTMetadata <<
" ymin : " << yMin <<
" mm, ymax : " << yMax <<
" mm" <<
RESTendl;
1647 RESTMetadata <<
" zmin : " << zMin <<
" mm, zmax : " << zMax <<
" mm" <<
RESTendl;
1650 RESTMetadata <<
"+++++++++++++++++++++++++++++++++++++++++++++++++" <<
RESTendl;
A class to load magnetic field maps and evaluate the field on those maps including interpolation.
TVector3 GetVolumeCenter(Int_t id)
it returns the volume position (or center) for the given volume id.
void PrintMetadata()
Prints on screen the information about the metadata members of TRestAxionMagneticField.
TVector3 GetVolumePosition(Int_t id)
it returns the volume position (or center) for the given volume id.
TVector3 fTrackDirection
The track direction used to parameterize the field along a track.
void InitFromConfigFile()
Initialization of TRestAxionMagnetic field members through a RML file.
Double_t fTrackLength
The total length of the track which defines the limit for field parameterization.
void LoadMagneticVolumes()
It will load the magnetic field data from the data filenames specified at the RML definition.
std::vector< TVector3 > fConstantField
A constant field component that will be added to the field map.
std::vector< TVector3 > fPositions
The absolute position of each of the magnetic volumes defined in this class.
void SetTrack(const TVector3 &position, const TVector3 &direction)
It initializes the field track boundaries data members of this class using a track position and direc...
std::vector< Double_t > GetComponentAlongPath(Int_t axis, TVector3 from, TVector3 to, Double_t dl=1., Int_t Nmax=0)
It returns a vector describing the magnetic field profile component requested using the axis argument...
TVector3 fReMap
A vector that defines the new mesh cell volume. It will re-scale the original fMeshSize.
std::vector< TVector3 > GetFieldBoundaries(TVector3 pos, TVector3 dir, Double_t precision=0, Int_t id=0)
Finds the in/out particle trajectory boundaries for a particular magnetic region, similar to the meth...
unsigned int GetNumberOfVolumes()
The number of magnetic volumes loaded into the object.
TCanvas * DrawComponents(Int_t volId=0)
A method that creates a canvas where tracks traversing the magnetic volume are drawm together with th...
Bool_t CheckOverlaps()
It will return true if the magnetic the regions overlap.
TVector3 GetFieldAverageTransverseVector(TVector3 from, TVector3 to, Double_t dl=10., Int_t Nmax=0)
It returns the transverse component of the average magnetic field vector calculated along the line th...
TCanvas * DrawTracks(TVector3 vanishingPoint, Int_t divisions, Int_t volId=0)
A method that creates a canvas where tracks traversing the magnetic volume are drawm together with th...
TRestAxionMagneticField()
Default constructor.
std::vector< TVector3 > GetVolumeBoundaries(TVector3 pos, TVector3 dir, Int_t id=0)
Finds the in/out particle trajectory boundaries for a particular magnetic region bounding box.
void ReMap(const size_t &n, const TVector3 &newMapSize)
It allows to remap the magnetic field to a larger mesh size. The new mesh size granularity must be pr...
Int_t GetVolumeIndex(TVector3 pos)
It returns the corresponding volume index at the given position. If not found it will return -1.
std::vector< TVector3 > fBoundMax
A vector to store the maximum bounding box values.
TCanvas * DrawHistogram(TString projection, TString Bcomp, Int_t volIndex=-1, Double_t step=-1, TString style="COLZ0", Double_t depth=-100010.0)
A method that creates a canvas where magnetic field map is drawn.
Bool_t IsFieldConstant(Int_t id)
It returns true if no magnetic field map was loaded for that volume.
MagneticFieldVolume * GetMagneticVolume(Int_t id)
It returns a pointer to the corresponding magnetic volume id.
TVector3 fTrackStart
The start track position used to parameterize the field along a track.
TVector3 GetMagneticVolumeNode(size_t id, TVector3 pos)
It returns the corresponding mesh node in the magnetic volume.
TVector3 GetMagneticField(Double_t x, Double_t y, Double_t z)
It returns the magnetic field vector at x,y,z.
std::vector< TVector3 > fMeshSize
std::vector< MagneticFieldVolume > fMagneticFieldVolumes
A magnetic field volume structure to store field data and mesh.
std::vector< Double_t > GetTransversalComponentAlongPath(TVector3 from, TVector3 to, Double_t dl=1., Int_t Nmax=0)
It returns a vector describing the transversal magnetic field component between from and to positions...
std::vector< TString > fMeshType
The type of the mesh used (default is cylindrical)
void LoadMagneticFieldData(MagneticFieldVolume &mVol, std::vector< std::vector< Float_t > > data)
A method to help loading magnetic field data, as x,y,z,Bx,By,Bz into a magnetic volume definition usi...
Double_t GetTransversalComponent(TVector3 position, TVector3 direction)
It returns the intensity of the transversal magnetic field component for the defined propagation dire...
std::vector< std::string > fFileNames
The name of the filenames containing the field data.
TH2D * fHisto
A helper histogram to plot the field.
TCanvas * fCanvas
A canvas to insert the histogram drawing.
~TRestAxionMagneticField()
Default destructor.
Bool_t IsInside(TVector3 pos)
It returns true if the given position is found inside a magnetic volume. False otherwise.
void Initialize()
Initialization of TRestAxionMagneticField members.
Double_t GetTransversalFieldAverage(TVector3 from, TVector3 to, Double_t dl=1., Int_t Nmax=0)
It returns the average of the transversal magnetic field intensity between the 3-d coordinates from a...
Double_t GetTransversalComponentInParametricTrack(Double_t x)
It will return the transversal magnetic field component evaluated at a parametric distance x (given b...
Bool_t FieldLoaded()
This private method returns true if the magnetic field volumes loaded are the same as the volumes def...
A basic class inheriting from TObject to help creating a node grid definition.
void SetSize(Double_t sX, Double_t sY, Double_t sZ)
Sets the origin of the bounding-box and initializes the nodes vector to zero.
void SetOrigin(Double_t oX, Double_t oY, Double_t oZ)
Sets the origin of the bounding-box and initializes the nodes vector to zero.
void SetCylindrical(Bool_t v)
Sets the coordinate system to cylindrical.
void SetNodes(Int_t nX, Int_t nY, Int_t nZ)
Sets the number of nodes and initializes the nodes vector to zero.
@ REST_Extreme
show everything
@ REST_Info
+show most of the information for each steps
@ REST_Debug
+show the defined debug messages
This namespace serves to define physics constants and other basic physical operations.
TVector3 MoveToPlane(const TVector3 &pos, const TVector3 &dir, const TVector3 &n, const TVector3 &a)
This method will translate the vector with direction dir starting at position pos to the plane define...
TVector3 MoveByDistanceFast(const TVector3 &pos, const TVector3 &dir, Double_t d)
This method transports a position pos by a distance d in the direction defined by dir....
Int_t GetChar(std::string hint="Press a KEY to continue ...")
Helps to pause the program, printing a message before pausing.