REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
Public Member Functions | Private Member Functions | Private Attributes
TRestAxionMagneticField Class Reference

Detailed Description

A class to load magnetic field maps and evaluate the field on those maps including interpolation.

TRestAxionMagneticField is a class that allows to load externally defined magnetic fields, and create magnetic volume regions associated to those pre-generated definitions. Once the field maps have been loaded this class will be able to evaluate the field vector at any coordinate (x,y,z). If the coordinate (x,y,z) is ourside any defined region, the returned field will be (0,0,0).

TODO Description of magnetic field interpolation

RML definition

We can add any number of magnetic volumes inside the RML definition as shown in the following piece of code,

<addMagneticVolume file="magnetic.file" position="(30,0,0)cm" />
<addMagneticVolume file="magnetic.file" position="(-30,0,0)cm" />
A class to load magnetic field maps and evaluate the field on those maps including interpolation.

where we produce 2 magnetic regions, using the same magnetic map provided in file magnetic.file and shifted by x=-30cm and x=30cm. The following parameters might be defined at each addMagneticVolume entry.

All parameters are optional, and if not provided they will take their default values.

Note
If no magnetic field map is provided, i.e. we just want to define a constant magnetic field vector, we still need to provide a boundary box size and mesh size through the boundMax and meshSize parameters. We still have the possibility to define a constant magnetic field vector for that volume using the field parameter. In that case the method TRestMagneticField::IsFieldConstant will return true.

The magnetic field file format.

There are two possible formats to provide a magnetic field map.

A more detailed example

The following example shows different allowed volume definition entries.

<TRestAxionMagneticField name="bFieldBabyIAXO" title="First magnetic field definition"
verboseLevel="info" >
<!-- A volume from a text file centered at (0,0,0) -->
<addMagneticVolume fileName="Bykovskiy_201906.dat" position="(0,0,0)mm" />
<!-- A volume from binary file, include Bx=1T offset, and boundMax validation -->
<addMagneticVolume fileName="Bykovskiy_202004.bin" position="(800,800,800)mm"
field="(1,0,0)" boundMax="(350,350,4000)" />
<!-- A magnetic volume with a constant field definition (no field map) -->
<addMagneticVolume field="(0,0,3)T" position="(-800,-800,-800)mm"
boundMax="(100,200,300)" meshSize="(10,20,30)" />

Using this class

Once we have created an instance of this class we will be able to access the magnetic field directly, at any point by using trilinear interpolation.

For example, the following code will retrieve the field at a position found at a distance z=-2m at the magnet axis,

TRestAxionMagneticField *mag = new TRestAxionMagneticField("fields.rml", "babyIAXO" );
mag->GetMagneticField( 0, 0, -2000 );
TRestAxionMagneticField()
Default constructor.
TVector3 GetMagneticField(Double_t x, Double_t y, Double_t z)
It returns the magnetic field vector at x,y,z.

or in a similar way,

TRestAxionMagneticField *mag = new TRestAxionMagneticField("fields.rml", "babyIAXO" );
mag->GetMagneticField( TVector3(0, 0, -2000), true );

where the second argument, true will enable the warning message system.

The following code evaluates the transversal field component at 1m distance along the Z-axis, for a vector pointing in the direction (0,1,1).

TRestAxionMagneticField *mag = new TRestAxionMagneticField("fields.rml", "babyIAXO" );
mag->GetTransversalComponent( TVector3(0,0,1000), TVector3(0,1,1) );
Double_t GetTransversalComponent(TVector3 position, TVector3 direction)
It returns the intensity of the transversal magnetic field component for the defined propagation dire...

There are also geometric functions that allow to identify the boundaries of the magnetic volume. A test particle will penetrate in the bounding box and identify the moment where the field changes to a value different from (0,0,0) in order("fields.rml", "babyIAXO_2024"); to identify the entrance and exit point.

The following code will return two TVector3 with the magnetic volume entrance and exit coordinates, for a test particle being placed at x=10cm, y=10cm and z=-4m, pointing towards the magnetic field with a slight direction deviation from the magnet axis.

TRestAxionMagneticField *mag = new TRestAxionMagneticField("fields.rml", "babyIAXO" );
mag->GetFieldBoundaries(TVector3(100,100, -4000) , TVector3(0.02, 0.03, 1))
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...

In the other hand, TRestAxionMagneticField::GetVolumeBoundaries will return the bounding box containing that magnetic field.

Remapping the field

It is possible to remap the magnetic field in order to increase the size of the grid cells in order to gain in performance. We may only increase the size of the cells by a multiple of the original mesh size.

Given the original mesh granularity by 10mmx10mmx50mm, we might increase the granularity of the first volume, index 0 to 50mmx50mmx200mm as follows:

TRestAxionMagneticField *mag = new TRestAxionMagneticField("fields.rml", "babyIAXO" );
mag->ReMap( 0, TVector3(50, 50, 200) );
void PrintMetadata()
Prints on screen the information about the metadata members of TRestAxionMagneticField.
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...

Visualizing the magnetic field

You may visualize the magnetic field profile along tracks towards a vanishing point by using the TRestAxionMagneticField::DrawTracks method, using:

TRestAxionMagneticField *field = new TRestAxionMagneticField( "fields.rml", "babyIAXO" );
field->DrawTracks( TVector3(0,0,8000), 100 );
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...

That will produce the following plot:

Tracks through the magnetic field volume and its corresponding T-field component

TODO Review and validate DrawHistogram drawing method and describe its use here.


RESTsoft - Software for Rare Event Searches with TPCs

History of developments:

2019-June: First concept and implementation of TRestAxionMagneticField class. Eve Pachoud

2020-April: Reviewing and validating TRestAxionMagneticField class. Javier Galan and Krešimir Jakovčić

Author
Eve Pachoud
Javier Galan javie.nosp@m.r.ga.nosp@m.lan@u.nosp@m.niza.nosp@m.r.es
Krešimir Jakovčić kjako.nosp@m.v@ir.nosp@m.b.hr

Definition at line 50 of file TRestAxionMagneticField.h.

#include <TRestAxionMagneticField.h>

Inheritance diagram for TRestAxionMagneticField:
TRestMetadata

Public Member Functions

Bool_t CheckOverlaps ()
 It will return true if the magnetic the regions overlap. More...
 
TCanvas * DrawComponents (Int_t volId=0)
 A method that creates a canvas where tracks traversing the magnetic volume are drawm together with their corresponding field intensity profile along the Z-axis. More...
 
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. More...
 
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 their corresponding field intensity profile along the Z-axis. More...
 
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 which may take values 0,1,2 for X,Y,Z components. The field profile will be constructed in the track defined between from and to positions given by argument. More...
 
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 that connects the 3-d coordinates from and to with respect to that line. More...
 
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 method TRestAxionMagneticField::GetVolumeBoudaries, but requiring that the in/out points are the first/last points where the transversal field intensity is not zero. More...
 
TVector3 GetMagneticField (Double_t x, Double_t y, Double_t z)
 It returns the magnetic field vector at x,y,z. More...
 
TVector3 GetMagneticField (TVector3 pos, Bool_t showWarning=true)
 It returns the magnetic field vector at TVector3(pos) using trilinear interpolation that is implemented following instructions given at https://en.wikipedia.org/wiki/Trilinear_interpolation. More...
 
unsigned int GetNumberOfVolumes ()
 The number of magnetic volumes loaded into the object. More...
 
Double_t GetPhotonAbsorptionLength (Double_t x, Double_t y, Double_t z, Double_t en)
 
Double_t GetPhotonAbsorptionLength (Int_t id, Double_t en)
 
Double_t GetPhotonAbsorptionLength (TVector3 pos, Double_t en)
 
Double_t GetPhotonMass (Double_t x, Double_t y, Double_t z, Double_t en)
 
Double_t GetPhotonMass (Int_t id, Double_t en)
 
Double_t GetPhotonMass (TVector3 pos, Double_t en)
 
TVector3 GetTrackDirection () const
 
Double_t GetTrackLength () const
 
TVector3 GetTrackStart () const
 
Double_t GetTransversalComponent (TVector3 position, TVector3 direction)
 It returns the intensity of the transversal magnetic field component for the defined propagation direction and position given by argument. More...
 
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 given by argument. More...
 
Double_t GetTransversalComponentInParametricTrack (Double_t x)
 It will return the transversal magnetic field component evaluated at a parametric distance x (given by argument) for the track defined inside the class. The track will be defined by the data members fStartTrack and fEndTrack which should be initialized by external means by invoking the method SetTrack( position, direction );. More...
 
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 and to. More...
 
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. More...
 
TVector3 GetVolumeCenter (Int_t id)
 it returns the volume position (or center) for the given volume id. More...
 
Int_t GetVolumeIndex (TVector3 pos)
 It returns the corresponding volume index at the given position. If not found it will return -1. More...
 
TVector3 GetVolumePosition (Int_t id)
 it returns the volume position (or center) for the given volume id. More...
 
Bool_t IsFieldConstant (Int_t id)
 It returns true if no magnetic field map was loaded for that volume. More...
 
Bool_t IsInside (TVector3 pos)
 It returns true if the given position is found inside a magnetic volume. False otherwise. More...
 
void LoadMagneticVolumes ()
 It will load the magnetic field data from the data filenames specified at the RML definition. More...
 
void PrintMetadata ()
 Prints on screen the information about the metadata members of TRestAxionMagneticField. More...
 
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 provided by argument, and each dimension must be a factor of the present mesh size. More...
 
void SetTrack (const TVector3 &position, const TVector3 &direction)
 It initializes the field track boundaries data members of this class using a track position and direction so that these values can be used later on in parameterization. More...
 
 TRestAxionMagneticField ()
 Default constructor. More...
 
 TRestAxionMagneticField (const char *cfgFileName, std::string name="")
 Constructor loading data from a config file. More...
 
 ~TRestAxionMagneticField ()
 Default destructor. More...
 

Private Member Functions

Bool_t FieldLoaded ()
 This private method returns true if the magnetic field volumes loaded are the same as the volumes defined. More...
 
MagneticFieldVolume * GetMagneticVolume (Int_t id)
 It returns a pointer to the corresponding magnetic volume id. More...
 
TVector3 GetMagneticVolumeNode (size_t id, TVector3 pos)
 It returns the corresponding mesh node in the magnetic volume. More...
 
void InitFromConfigFile ()
 Initialization of TRestAxionMagnetic field members through a RML file. More...
 
void Initialize ()
 Initialization of TRestAxionMagneticField members. More...
 
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 using its corresponding mesh. More...
 

Private Attributes

std::vector< TVector3 > fBoundMax
 A vector to store the maximum bounding box values. More...
 
TCanvas * fCanvas
 A canvas to insert the histogram drawing. More...
 
std::vector< TVector3 > fConstantField
 A constant field component that will be added to the field map. More...
 
std::vector< std::string > fFileNames
 The name of the filenames containing the field data. More...
 
TH2D * fHisto
 A helper histogram to plot the field. More...
 
std::vector< MagneticFieldVolume > fMagneticFieldVolumes
 A magnetic field volume structure to store field data and mesh. More...
 
std::vector< TVector3 > fMeshSize
 
std::vector< TString > fMeshType
 The type of the mesh used (default is cylindrical) More...
 
std::vector< TVector3 > fPositions
 The absolute position of each of the magnetic volumes defined in this class. More...
 
TVector3 fReMap = TVector3(0, 0, 0)
 A vector that defines the new mesh cell volume. It will re-scale the original fMeshSize. More...
 
TVector3 fTrackDirection = TVector3(0, 0, 0)
 The track direction used to parameterize the field along a track. More...
 
Double_t fTrackLength = 0
 The total length of the track which defines the limit for field parameterization. More...
 
TVector3 fTrackStart = TVector3(0, 0, 0)
 The start track position used to parameterize the field along a track. More...
 

Additional Inherited Members

Constructor & Destructor Documentation

◆ TRestAxionMagneticField() [1/2]

TRestAxionMagneticField::TRestAxionMagneticField ( )

Default constructor.

Definition at line 248 of file TRestAxionMagneticField.cxx.

◆ TRestAxionMagneticField() [2/2]

TRestAxionMagneticField::TRestAxionMagneticField ( const char *  cfgFileName,
std::string  name = "" 
)

Constructor loading data from a config file.

If no configuration path is defined using TRestMetadata::SetConfigFilePath the path to the config file must be specified using full path, absolute or relative.

The default behaviour is that the config file must be specified with full path, absolute or relative.

Parameters
cfgFileNameA const char* giving the path to an RML file.
nameThe name of the specific metadata. It will be used to find the corresponding TRestAxionMagneticField section inside the RML.

Definition at line 264 of file TRestAxionMagneticField.cxx.

◆ ~TRestAxionMagneticField()

TRestAxionMagneticField::~TRestAxionMagneticField ( )

Default destructor.

Definition at line 278 of file TRestAxionMagneticField.cxx.

Member Function Documentation

◆ CheckOverlaps()

Bool_t TRestAxionMagneticField::CheckOverlaps ( )

It will return true if the magnetic the regions overlap.

Definition at line 1443 of file TRestAxionMagneticField.cxx.

◆ DrawComponents()

TCanvas * TRestAxionMagneticField::DrawComponents ( Int_t  volId = 0)

A method that creates a canvas where tracks traversing the magnetic volume are drawm together with their corresponding field intensity profile along the Z-axis.

Definition at line 591 of file TRestAxionMagneticField.cxx.

◆ DrawHistogram()

TCanvas * TRestAxionMagneticField::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.

This method is used to create various 2D visualisations of the Bx, By and Bz magnetic field components for a given magnetic field region.

The following input parameters can be specified :

  • projection : it specifies the plane in the magnetic field region from which the plotted values of the magnetic field component were taken. The allowed values are: "XY", "XZ" and "YZ" for X-Y, X-Z and Y-Z plane, respectively.
  • Bcomp : it specifies which component of the magnetic field is plotted. The allowed values are: "X", "Y" and "Z" for Bx, By and Bz components, respectively.
  • volIndex : it specifies the index of the magnetic field region/volume for which the plot is shown.
  • step : it specifies the step/bin size for the 2D histogram that is shown. If this parameter is not specified, or if it has value < = 0, the values are taken from the corresponding mesh size values stored in data member fMeshSize.
  • style : it specifies the plotting style for the histogram. It can correspond to any draw option for 2D histograms in ROOT. The list of options can be found at: https://root.cern.ch/root/htmldoc/guides/users-guide/Histograms.html#drawing-histograms The default value is "COLZ0" which draws a box for each cell in the histogram with a color scale varying with values. The other useful option is "SURF3Z" which draws a surface plot with a coloured contour view on the top.
  • depth : it specifies the position of the plane in the magnetic field region from which the plotted values of the magnetic field component were taken. If this parameter is not specified, the histogram is plotted for the plane that goes through the middle of the region.

The following example will plot the values of the Bx component of the magnetic field taken in the X-Y plane positioned at z=5000mm in the magnetic field region with index 0. The step/bin size is 50mm, and the plotting style is "SURF3Z".

field->DrawHistogram("XY","X",0,50.0,"SURF3Z",5000.0)

where field is a pointer to the TRestAxionMagneticField object that describes the magnetic field.

Definition at line 337 of file TRestAxionMagneticField.cxx.

◆ DrawTracks()

TCanvas * TRestAxionMagneticField::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 their corresponding field intensity profile along the Z-axis.

Definition at line 650 of file TRestAxionMagneticField.cxx.

◆ FieldLoaded()

Bool_t TRestAxionMagneticField::FieldLoaded ( )
inlineprivate

This private method returns true if the magnetic field volumes loaded are the same as the volumes defined.

Definition at line 102 of file TRestAxionMagneticField.h.

◆ GetComponentAlongPath()

std::vector< Double_t > TRestAxionMagneticField::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 which may take values 0,1,2 for X,Y,Z components. The field profile will be constructed in the track defined between from and to positions given by argument.

The differential element dl is by default 1mm, but it can be modified through the third argument of this function.

The maximum number of divisions (unlimited by default) of the output vector can be fixed by the forth argument. In that case, the differential element dl length might be increased to fullfil such condition.

Definition at line 1288 of file TRestAxionMagneticField.cxx.

◆ GetFieldAverageTransverseVector()

TVector3 TRestAxionMagneticField::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 that connects the 3-d coordinates from and to with respect to that line.

The differential element dl defines the step, and it is by default 10mm, but it can be modified through the third argument of this function.

The maximum number of divisions (unlimited by default) can be fixed by the forth argument. In that case, the differential element dl length might be increased to fullfil such condition.

Definition at line 1385 of file TRestAxionMagneticField.cxx.

◆ GetFieldBoundaries()

std::vector< TVector3 > TRestAxionMagneticField::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 method TRestAxionMagneticField::GetVolumeBoudaries, but requiring that the in/out points are the first/last points where the transversal field intensity is not zero.

If no precision is given, the mesh size of the corresponding volume will be used as reference. The precision will be meshSize/2.

If no intersection is found the returned std::vector will be empty.

Definition at line 1497 of file TRestAxionMagneticField.cxx.

◆ GetMagneticField() [1/2]

TVector3 TRestAxionMagneticField::GetMagneticField ( Double_t  x,
Double_t  y,
Double_t  z 
)

It returns the magnetic field vector at x,y,z.

Definition at line 1006 of file TRestAxionMagneticField.cxx.

◆ GetMagneticField() [2/2]

TVector3 TRestAxionMagneticField::GetMagneticField ( TVector3  pos,
Bool_t  showWarning = true 
)

It returns the magnetic field vector at TVector3(pos) using trilinear interpolation that is implemented following instructions given at https://en.wikipedia.org/wiki/Trilinear_interpolation.

The warning in case the evaluated point is found outside any volume might be disabled using the showWarning argument.

Definition at line 1018 of file TRestAxionMagneticField.cxx.

◆ GetMagneticVolume()

MagneticFieldVolume * TRestAxionMagneticField::GetMagneticVolume ( Int_t  id)
inlineprivate

It returns a pointer to the corresponding magnetic volume id.

Definition at line 105 of file TRestAxionMagneticField.h.

◆ GetMagneticVolumeNode()

TVector3 TRestAxionMagneticField::GetMagneticVolumeNode ( size_t  id,
TVector3  pos 
)
private

It returns the corresponding mesh node in the magnetic volume.

The corresponging node to x,y,z is the bottom, down, left node in the cell volume defined by 8-nodes.

This method will be made private, no reason to use it outside this class.

Definition at line 1433 of file TRestAxionMagneticField.cxx.

◆ GetNumberOfVolumes()

unsigned int TRestAxionMagneticField::GetNumberOfVolumes ( )
inline

The number of magnetic volumes loaded into the object.

Definition at line 133 of file TRestAxionMagneticField.h.

◆ GetTrackDirection()

TVector3 TRestAxionMagneticField::GetTrackDirection ( ) const
inline

Definition at line 124 of file TRestAxionMagneticField.h.

◆ GetTrackLength()

Double_t TRestAxionMagneticField::GetTrackLength ( ) const
inline

Definition at line 122 of file TRestAxionMagneticField.h.

◆ GetTrackStart()

TVector3 TRestAxionMagneticField::GetTrackStart ( ) const
inline

Definition at line 123 of file TRestAxionMagneticField.h.

◆ GetTransversalComponent()

Double_t TRestAxionMagneticField::GetTransversalComponent ( TVector3  position,
TVector3  direction 
)

It returns the intensity of the transversal magnetic field component for the defined propagation direction and position given by argument.

Definition at line 1236 of file TRestAxionMagneticField.cxx.

◆ GetTransversalComponentAlongPath()

std::vector< Double_t > TRestAxionMagneticField::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 given by argument.

The differential element dl is by default 1mm, but it can be modified through the third argument of this function.

The maximum number of divisions (unlimited by default) of the output vector can be fixed by the forth argument. In that case, the differential element dl length might be increased to fullfil such condition.

Definition at line 1251 of file TRestAxionMagneticField.cxx.

◆ GetTransversalComponentInParametricTrack()

Double_t TRestAxionMagneticField::GetTransversalComponentInParametricTrack ( Double_t  x)

It will return the transversal magnetic field component evaluated at a parametric distance x (given by argument) for the track defined inside the class. The track will be defined by the data members fStartTrack and fEndTrack which should be initialized by external means by invoking the method SetTrack( position, direction );.

This method will be used by the integration method

Definition at line 1346 of file TRestAxionMagneticField.cxx.

◆ GetTransversalFieldAverage()

Double_t TRestAxionMagneticField::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 and to.

The differential element dl defines the integration step, and it is by default 1mm, but it can be modified through the third argument of this function.

The maximum number of divisions of the output vector can be fixed by the forth argument. In that case, the differential element dl length might be increased to fullfil such condition.

Definition at line 1362 of file TRestAxionMagneticField.cxx.

◆ GetVolumeBoundaries()

std::vector< TVector3 > TRestAxionMagneticField::GetVolumeBoundaries ( TVector3  pos,
TVector3  dir,
Int_t  id = 0 
)

Finds the in/out particle trajectory boundaries for a particular magnetic region bounding box.

This method checks if the trajectory defined by the position pos and direction dir passes through the magnetic field region/volume id given. If two such points (entry point and exit point) are found, their coordinates are returned. In the example shown in Fig. 1 from TRestAxionFieldPropagationProcess these points are: IN 1 and OUT 1 for the region #1 and IN2 and OUT 2 for the region #2.

If no intersection is found, or the particle is not moving towards the volume, the returned std::vector will be empty.

Definition at line 1477 of file TRestAxionMagneticField.cxx.

◆ GetVolumeCenter()

TVector3 TRestAxionMagneticField::GetVolumeCenter ( Int_t  id)

it returns the volume position (or center) for the given volume id.

Definition at line 1216 of file TRestAxionMagneticField.cxx.

◆ GetVolumeIndex()

Int_t TRestAxionMagneticField::GetVolumeIndex ( TVector3  pos)

It returns the corresponding volume index at the given position. If not found it will return -1.

Definition at line 1196 of file TRestAxionMagneticField.cxx.

◆ GetVolumePosition()

TVector3 TRestAxionMagneticField::GetVolumePosition ( Int_t  id)

it returns the volume position (or center) for the given volume id.

Definition at line 1221 of file TRestAxionMagneticField.cxx.

◆ InitFromConfigFile()

void TRestAxionMagneticField::InitFromConfigFile ( )
privatevirtual

Initialization of TRestAxionMagnetic field members through a RML file.

Reimplemented from TRestMetadata.

Definition at line 1535 of file TRestAxionMagneticField.cxx.

◆ Initialize()

void TRestAxionMagneticField::Initialize ( )
privatevirtual

Initialization of TRestAxionMagneticField members.

Reimplemented from TRestMetadata.

Definition at line 285 of file TRestAxionMagneticField.cxx.

◆ IsFieldConstant()

Bool_t TRestAxionMagneticField::IsFieldConstant ( Int_t  id)
inline

It returns true if no magnetic field map was loaded for that volume.

Definition at line 127 of file TRestAxionMagneticField.h.

◆ IsInside()

Bool_t TRestAxionMagneticField::IsInside ( TVector3  pos)

It returns true if the given position is found inside a magnetic volume. False otherwise.

Definition at line 1208 of file TRestAxionMagneticField.cxx.

◆ LoadMagneticFieldData()

void TRestAxionMagneticField::LoadMagneticFieldData ( MagneticFieldVolume &  mVol,
std::vector< std::vector< Float_t > >  data 
)
private

A method to help loading magnetic field data, as x,y,z,Bx,By,Bz into a magnetic volume definition using its corresponding mesh.

This method will be made private since it will only be used internally.

Definition at line 780 of file TRestAxionMagneticField.cxx.

◆ LoadMagneticVolumes()

void TRestAxionMagneticField::LoadMagneticVolumes ( )

It will load the magnetic field data from the data filenames specified at the RML definition.

This method will be made private since it will only be used internally.

Definition at line 834 of file TRestAxionMagneticField.cxx.

◆ PrintMetadata()

void TRestAxionMagneticField::PrintMetadata ( )
virtual

Prints on screen the information about the metadata members of TRestAxionMagneticField.

Reimplemented from TRestMetadata.

Definition at line 1606 of file TRestAxionMagneticField.cxx.

◆ ReMap()

void TRestAxionMagneticField::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 provided by argument, and each dimension must be a factor of the present mesh size.

Definition at line 1146 of file TRestAxionMagneticField.cxx.

◆ SetTrack()

void TRestAxionMagneticField::SetTrack ( const TVector3 &  position,
const TVector3 &  direction 
)

It initializes the field track boundaries data members of this class using a track position and direction so that these values can be used later on in parameterization.

Definition at line 1323 of file TRestAxionMagneticField.cxx.

Field Documentation

◆ fBoundMax

std::vector<TVector3> TRestAxionMagneticField::fBoundMax
private

A vector to store the maximum bounding box values.

Definition at line 69 of file TRestAxionMagneticField.h.

◆ fCanvas

TCanvas* TRestAxionMagneticField::fCanvas
private

A canvas to insert the histogram drawing.

Definition at line 90 of file TRestAxionMagneticField.h.

◆ fConstantField

std::vector<TVector3> TRestAxionMagneticField::fConstantField
private

A constant field component that will be added to the field map.

Definition at line 59 of file TRestAxionMagneticField.h.

◆ fFileNames

std::vector<std::string> TRestAxionMagneticField::fFileNames
private

The name of the filenames containing the field data.

Definition at line 53 of file TRestAxionMagneticField.h.

◆ fHisto

TH2D* TRestAxionMagneticField::fHisto
private

A helper histogram to plot the field.

Definition at line 87 of file TRestAxionMagneticField.h.

◆ fMagneticFieldVolumes

std::vector<MagneticFieldVolume> TRestAxionMagneticField::fMagneticFieldVolumes
private

A magnetic field volume structure to store field data and mesh.

Definition at line 75 of file TRestAxionMagneticField.h.

◆ fMeshSize

std::vector<TVector3> TRestAxionMagneticField::fMeshSize
private

The size of a grid element from the mesh in mm. Initially, it must be the same as the binary input data

Definition at line 63 of file TRestAxionMagneticField.h.

◆ fMeshType

std::vector<TString> TRestAxionMagneticField::fMeshType
private

The type of the mesh used (default is cylindrical)

Definition at line 66 of file TRestAxionMagneticField.h.

◆ fPositions

std::vector<TVector3> TRestAxionMagneticField::fPositions
private

The absolute position of each of the magnetic volumes defined in this class.

Definition at line 56 of file TRestAxionMagneticField.h.

◆ fReMap

TVector3 TRestAxionMagneticField::fReMap = TVector3(0, 0, 0)
private

A vector that defines the new mesh cell volume. It will re-scale the original fMeshSize.

Definition at line 72 of file TRestAxionMagneticField.h.

◆ fTrackDirection

TVector3 TRestAxionMagneticField::fTrackDirection = TVector3(0, 0, 0)
private

The track direction used to parameterize the field along a track.

Definition at line 81 of file TRestAxionMagneticField.h.

◆ fTrackLength

Double_t TRestAxionMagneticField::fTrackLength = 0
private

The total length of the track which defines the limit for field parameterization.

Definition at line 84 of file TRestAxionMagneticField.h.

◆ fTrackStart

TVector3 TRestAxionMagneticField::fTrackStart = TVector3(0, 0, 0)
private

The start track position used to parameterize the field along a track.

Definition at line 78 of file TRestAxionMagneticField.h.


The documentation for this class was generated from the following files: