16#ifndef ShipFieldMaker_H
17#define ShipFieldMaker_H
22#include "TVirtualMagField.h"
25#include "TG4VUserPostDetConstruction.h"
47 typedef std::map<TString, TVirtualMagField*>
SFMap;
66 fieldInfo(
const TString& volName,
const TString& fieldName, Double_t scale = 1.0) :
85 void defineUniform(
const TString& name,
const TVector3& BVector);
95 void defineConstant(
const TString& name,
const TVector2& xRange,
const TVector2& yRange,
96 const TVector2& zRange,
const TVector3& BVector);
109 void defineBell(
const TString& name, Double_t BPeak, Double_t zMiddle,
110 Int_t orient = 1, Double_t tubeR = 500.0,
111 Double_t xy = 0.0, Double_t z = 0.0, Double_t L = 0.0);
121 void defineFieldMap(
const TString& name,
const TString& mapFileName,
122 const TVector3& localCentre = TVector3(0.0, 0.0, 0.0),
123 const TVector3& localAngles = TVector3(0.0, 0.0, 0.0),
124 Bool_t useSymmetry = kFALSE);
134 const TVector3& translation,
135 const TVector3& eulerAngles = TVector3(0.0, 0.0, 0.0));
146 const TString& field2Name,
const TString& field3Name =
"",
147 const TString& field4Name =
"");
154 void defineComposite(
const TString& name, std::vector<TString> fieldNames);
163 void defineGlobalField(
const TString& field1Name,
const TString& field2Name =
"",
164 const TString& field3Name =
"",
const TString& field4Name =
"");
179 Double_t scale = 1.0);
188 Double_t scale = 1.0);
215 Bool_t
gotField(
const TString& name)
const;
222 TVirtualMagField*
getField(
const TString& name)
const;
230 void plotXYField(
const TVector3& xAxis,
const TVector3& yAxis,
const std::string& plotFile)
const;
238 void plotZXField(
const TVector3& zAxis,
const TVector3& xAxis,
const std::string& plotFile)
const;
246 void plotZYField(
const TVector3& zAxis,
const TVector3& yAxis,
const std::string& plotFile)
const;
255 void plotField(Int_t type,
const TVector3& xAxis,
const TVector3& yAxis,
256 const std::string& plotFile)
const;
258 void generateFieldMap(TString fileName,
const float step=2.5,
const float xRange=179,
const float yRange=317,
const float zRange=1515.5,
const float zShift=-4996);
351 void checkLocalFieldMap(TVirtualMagField*& localField,
const TString& volName, Double_t scale);
369 void findNode(TGeoVolume* aVolume,
const TString& volName);
Class that defines a magnetic field composed from many fields.
Creates various magnetic fields and assigns them to geometry regions.
void defineFieldMap(const TString &name, const TString &mapFileName, const TVector3 &localCentre=TVector3(0.0, 0.0, 0.0), const TVector3 &localAngles=TVector3(0.0, 0.0, 0.0), Bool_t useSymmetry=kFALSE)
std::vector< fieldInfo > localInfo_
Vector of fieldInfo for local fields.
Bool_t gotField(const TString &name) const
Check if we have a field stored with the given name name.
void findNode(TGeoVolume *aVolume, const TString &volName)
Update the current geometry node pointer that matches the required volume name.
void defineComposite(const TString &name, const TString &field1Name, const TString &field2Name, const TString &field3Name="", const TString &field4Name="")
Define a composite field from up to four fields.
ClassDef(ShipFieldMaker, 1)
Generate fieldMap csv file in the given region.
void plotField(Int_t type, const TVector3 &xAxis, const TVector3 &yAxis, const std::string &plotFile) const
Plot the magnetic field in "x-y" plane.
Bool_t gotNode_
Boolean to specify if we have found the volume node we need.
SFMap theFields_
The map storing all created magnetic fields.
TVirtualMagField * getField(const TString &name) const
Get the field stored with the given name name.
void plotZYField(const TVector3 &zAxis, const TVector3 &yAxis, const std::string &plotFile) const
Plot the magnetic field in z-y plane.
ShipCompField * globalField_
The global magnetic field.
std::vector< fieldInfo > regionInfo_
Vector of fieldInfo for regional fields.
std::vector< std::string > stringVect
Typedef of a vector of strings.
void getTransformation(const TString &volName, transformInfo &theInfo)
Get the transformation matrix for the volume position and orientation.
void readInputFile(const std::string &inputFile)
Read an input file to define fields and associated volumes.
TGeoNode * theNode_
The current volume node: used for finding volume transformations.
void plotXYField(const TVector3 &xAxis, const TVector3 &yAxis, const std::string &plotFile) const
Plot the magnetic field in x-y plane.
Double_t Tesla_
Double converting Tesla to kiloGauss (for VMC/FairRoot B field units)
TVirtualMagField * getVolumeField(const TString &volName) const
Get the magnetic field for the given volume.
ShipCompField * getGlobalField() const
Get the global magnetic field.
void defineGlobalField(const TString &field1Name, const TString &field2Name="", const TString &field3Name="", const TString &field4Name="")
Define a composite field from up to four fields.
void defineLocalField(const TString &volName, const TString &fieldName, Double_t scale=1.0)
Define a localised field and volume pairing.
std::map< TString, TVirtualMagField * > SFMap
Typedef for a TString-TVirtualMagField* map.
Bool_t verbose_
Verbose boolean.
void setAllRegionFields()
void defineUniform(const TString &name, const TVector3 &BVector)
Define a uniform field.
virtual void Construct()
Set-up all local and regional fields and assign them to volumes.
virtual ~ShipFieldMaker()
Destructor.
SFMap getAllFields() const
Get the map storing volume names and their associated local magnetic fields.
void defineBell(const TString &name, Double_t BPeak, Double_t zMiddle, Int_t orient=1, Double_t tubeR=500.0, Double_t xy=0.0, Double_t z=0.0, Double_t L=0.0)
Define a Bell field.
void defineFieldMapCopy(const TString &name, const TString &mapNameToCopy, const TVector3 &translation, const TVector3 &eulerAngles=TVector3(0.0, 0.0, 0.0))
Define a copy of a field map with a coordinate translation and optional rotation.
void defineConstant(const TString &name, const TVector2 &xRange, const TVector2 &yRange, const TVector2 &zRange, const TVector3 &BVector)
Define a constant field.
stringVect splitString(std::string &theString, std::string &splitter) const
Split a string.
void checkLocalFieldMap(TVirtualMagField *&localField, const TString &volName, Double_t scale)
Check if we have a local field map and store the volume global transformation.
void plotZXField(const TVector3 &zAxis, const TVector3 &xAxis, const std::string &plotFile) const
Plot the magnetic field in z-x plane.
void defineRegionField(const TString &volName, const TString &fieldName, Double_t scale=1.0)
Define a regional (local + global) field and volume pairing.
void generateFieldMap(TString fileName, const float step=2.5, const float xRange=179, const float yRange=317, const float zRange=1515.5, const float zShift=-4996)
Structure to hold volume name, field name and field scaling factor.
TString fieldName_
The name of the field.
Double_t scale_
The field scaling factor.
fieldInfo()
Default constructor.
fieldInfo(const TString &volName, const TString &fieldName, Double_t scale=1.0)
Constructor.
TString volName_
The name of the volume.