REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestAxionMagneticField.h
1/*************************************************************************
2 * This file is part of the REST software framework. *
3 * *
4 * Copyright (C) 2016 GIFNA/TREX (University of Zaragoza) *
5 * For more information see http://gifna.unizar.es/trex *
6 * *
7 * REST is free software: you can redistribute it and/or modify *
8 * it under the terms of the GNU General Public License as published by *
9 * the Free Software Foundation, either version 3 of the License, or *
10 * (at your option) any later version. *
11 * *
12 * REST is distributed in the hope that it will be useful, *
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15 * GNU General Public License for more details. *
16 * *
17 * You should have a copy of the GNU General Public License along with *
18 * REST in $REST_PATH/LICENSE. *
19 * If not, see http://www.gnu.org/licenses/. *
20 * For the list of contributors see $REST_PATH/CREDITS. *
21 *************************************************************************/
22
23#ifndef _TRestAxionMagneticField
24#define _TRestAxionMagneticField
25
26#include <TRestMetadata.h>
27
28#include <iostream>
29
30#include "TCanvas.h"
31#include "TH2D.h"
32#include "TRestMesh.h"
33#include "TVector3.h"
34#include "TVectorD.h"
35
38struct MagneticFieldVolume {
40 TVector3 position;
41
43 TRestMesh mesh;
44
46 std::vector<std::vector<std::vector<TVector3>>> field;
47};
48
51 private:
53 std::vector<std::string> fFileNames; //<
54
56 std::vector<TVector3> fPositions; //<
57
59 std::vector<TVector3> fConstantField; //<
60
63 std::vector<TVector3> fMeshSize; //<
64
66 std::vector<TString> fMeshType; //<
67
69 std::vector<TVector3> fBoundMax; //<
70
72 TVector3 fReMap = TVector3(0, 0, 0); //<
73
75 std::vector<MagneticFieldVolume> fMagneticFieldVolumes;
76
78 TVector3 fTrackStart = TVector3(0, 0, 0);
79
81 TVector3 fTrackDirection = TVector3(0, 0, 0);
82
84 Double_t fTrackLength = 0;
85
87 TH2D* fHisto;
88
90 TCanvas* fCanvas;
91
92 void Initialize();
93
94 void InitFromConfigFile();
95
96 void LoadMagneticFieldData(MagneticFieldVolume& mVol, std::vector<std::vector<Float_t>> data);
97
98 TVector3 GetMagneticVolumeNode(size_t id, TVector3 pos);
99
102 Bool_t FieldLoaded() { return GetNumberOfVolumes() == fMagneticFieldVolumes.size(); }
103
105 MagneticFieldVolume* GetMagneticVolume(Int_t id) {
107 if (id >= 0 && fMagneticFieldVolumes.size() > (unsigned int)id)
108 return &fMagneticFieldVolumes[id];
109 else {
110 RESTError << "TRestAxionMagneticField::GetMagneticVolume. Id outside limits!" << RESTendl;
111 return NULL;
112 }
113 }
114
115 public:
116 void LoadMagneticVolumes();
117
118 void ReMap(const size_t& n, const TVector3& newMapSize);
119
120 void SetTrack(const TVector3& position, const TVector3& direction);
121
122 Double_t GetTrackLength() const { return fTrackLength; }
123 TVector3 GetTrackStart() const { return fTrackStart; }
124 TVector3 GetTrackDirection() const { return fTrackDirection; }
125
127 Bool_t IsFieldConstant(Int_t id) {
128 if (GetMagneticVolume(id)) return GetMagneticVolume(id)->field.size() == 0;
129 return true;
130 }
131
133 unsigned int GetNumberOfVolumes() { return fPositions.size(); }
134
135 Bool_t CheckOverlaps();
136
137 std::vector<TVector3> GetVolumeBoundaries(TVector3 pos, TVector3 dir, Int_t id = 0);
138 std::vector<TVector3> GetFieldBoundaries(TVector3 pos, TVector3 dir, Double_t precision = 0,
139 Int_t id = 0);
140
141 TVector3 GetMagneticField(Double_t x, Double_t y, Double_t z);
142 TVector3 GetMagneticField(TVector3 pos, Bool_t showWarning = true);
143
144 Double_t GetPhotonMass(Double_t x, Double_t y, Double_t z, Double_t en);
145 Double_t GetPhotonMass(TVector3 pos, Double_t en);
146 Double_t GetPhotonMass(Int_t id, Double_t en);
147
148 Double_t GetPhotonAbsorptionLength(Double_t x, Double_t y, Double_t z, Double_t en);
149 Double_t GetPhotonAbsorptionLength(TVector3 pos, Double_t en);
150 Double_t GetPhotonAbsorptionLength(Int_t id, Double_t en);
151
152 Int_t GetVolumeIndex(TVector3 pos);
153 Bool_t IsInside(TVector3 pos);
154
155 TVector3 GetVolumePosition(Int_t id);
156 TVector3 GetVolumeCenter(Int_t id);
157
158 Double_t GetTransversalComponent(TVector3 position, TVector3 direction);
159 Double_t GetTransversalComponentInParametricTrack(Double_t x);
160
161 std::vector<Double_t> GetTransversalComponentAlongPath(TVector3 from, TVector3 to, Double_t dl = 1.,
162 Int_t Nmax = 0);
163
164 std::vector<Double_t> GetComponentAlongPath(Int_t axis, TVector3 from, TVector3 to, Double_t dl = 1.,
165 Int_t Nmax = 0);
166
167 Double_t GetTransversalFieldAverage(TVector3 from, TVector3 to, Double_t dl = 1., Int_t Nmax = 0);
168
169 TVector3 GetFieldAverageTransverseVector(TVector3 from, TVector3 to, Double_t dl = 10., Int_t Nmax = 0);
170
171 TCanvas* DrawHistogram(TString projection, TString Bcomp, Int_t volIndex = -1, Double_t step = -1,
172 TString style = "COLZ0", Double_t depth = -100010.0);
173
174 TCanvas* DrawTracks(TVector3 vanishingPoint, Int_t divisions, Int_t volId = 0);
175
176 TCanvas* DrawComponents(Int_t volId = 0);
177
178 void PrintMetadata();
179
181 TRestAxionMagneticField(const char* cfgFileName, std::string name = "");
183
184 ClassDef(TRestAxionMagneticField, 3);
185};
186#endif
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.
Definition: TRestMesh.h:38
A base class for any REST metadata class.
Definition: TRestMetadata.h:74
endl_t RESTendl
Termination flag object for TRestStringOutput.