SND@LHC Software
Loading...
Searching...
No Matches
ShipFieldMaker.h
Go to the documentation of this file.
1
16#ifndef ShipFieldMaker_H
17#define ShipFieldMaker_H
18
19#include "ShipCompField.h"
20
21#include "TString.h"
22#include "TVirtualMagField.h"
23#include "TVector2.h"
24#include "TVector3.h"
25#include "TG4VUserPostDetConstruction.h"
26
27#include <map>
28#include <string>
29#include <vector>
30
31class TGeoMatrix;
32class TGeoNode;
33class TGeoVolume;
34
35class ShipFieldMaker : public TG4VUserPostDetConstruction
36{
37
38 public:
39
41 ShipFieldMaker(Bool_t verbose = kFALSE);
42
44 virtual ~ShipFieldMaker();
45
47 typedef std::map<TString, TVirtualMagField*> SFMap;
48
50 typedef std::vector<std::string> stringVect;
51
53 struct fieldInfo {
54
56 TString volName_;
58 TString fieldName_;
60 Double_t scale_;
61
63 fieldInfo() : volName_(""), fieldName_(""), scale_(1.0) {};
64
66 fieldInfo(const TString& volName, const TString& fieldName, Double_t scale = 1.0) :
67 volName_(volName), fieldName_(fieldName), scale_(scale) {};
68
69 };
70
72 virtual void Construct();
73
75
78 void readInputFile(const std::string& inputFile);
79
81
85 void defineUniform(const TString& name, const TVector3& BVector);
86
88
95 void defineConstant(const TString& name, const TVector2& xRange, const TVector2& yRange,
96 const TVector2& zRange, const TVector3& BVector);
97
99
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);
112
113 // ! Define a field map
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);
125
127
133 void defineFieldMapCopy(const TString& name, const TString& mapNameToCopy,
134 const TVector3& translation,
135 const TVector3& eulerAngles = TVector3(0.0, 0.0, 0.0));
136
138
145 void defineComposite(const TString& name, const TString& field1Name,
146 const TString& field2Name, const TString& field3Name = "",
147 const TString& field4Name = "");
148
150
154 void defineComposite(const TString& name, std::vector<TString> fieldNames);
155
157
163 void defineGlobalField(const TString& field1Name, const TString& field2Name = "",
164 const TString& field3Name = "", const TString& field4Name = "");
165
167
170 void defineGlobalField(std::vector<TString> fieldNames);
171
173
178 void defineRegionField(const TString& volName, const TString& fieldName,
179 Double_t scale = 1.0);
180
182
187 void defineLocalField(const TString& volName, const TString& fieldName,
188 Double_t scale = 1.0);
189
190
192
196
198
201 SFMap getAllFields() const {return theFields_;}
202
204
208 TVirtualMagField* getVolumeField(const TString& volName) const;
209
211
215 Bool_t gotField(const TString& name) const;
216
218
222 TVirtualMagField* getField(const TString& name) const;
223
225
230 void plotXYField(const TVector3& xAxis, const TVector3& yAxis, const std::string& plotFile) const;
231
233
238 void plotZXField(const TVector3& zAxis, const TVector3& xAxis, const std::string& plotFile) const;
239
241
246 void plotZYField(const TVector3& zAxis, const TVector3& yAxis, const std::string& plotFile) const;
247
249
255 void plotField(Int_t type, const TVector3& xAxis, const TVector3& yAxis,
256 const std::string& plotFile) const;
257
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);
260
263
264
265 protected:
266
269
271 Double_t x0_;
273 Double_t y0_;
275 Double_t z0_;
276
278 Double_t phi_;
280 Double_t theta_;
282 Double_t psi_;
283
284 };
285
287
290 void defineUniform(const stringVect& inputLine);
291
293
296 void defineConstant(const stringVect& inputLine);
297
299
302 void defineBell(const stringVect& inputLine);
303
305
309 void defineFieldMap(const stringVect& inputLine, Bool_t useSymmetry = kFALSE);
310
312
315 void defineFieldMapCopy(const stringVect& inputLine);
316
318
321 void defineComposite(const stringVect& inputLine);
322
324
327 void defineGlobalField(const stringVect& inputLine);
328
330
333 void defineRegionField(const stringVect& inputLine);
334
335 // ! Setup all of the regional fields. Called by Construct()
336 void setAllRegionFields();
337
338
340
343 void defineLocalField(const stringVect& inputLine);
344
346
351 void checkLocalFieldMap(TVirtualMagField*& localField, const TString& volName, Double_t scale);
352
353 // ! Setup all of the local fields. Called by Construct()
354 void setAllLocalFields();
355
356
358
362 void getTransformation(const TString& volName, transformInfo& theInfo);
363
365
369 void findNode(TGeoVolume* aVolume, const TString& volName);
370
371 private:
372
375
378
380 std::vector<fieldInfo> regionInfo_;
381
383 std::vector<fieldInfo> localInfo_;
384
386 Bool_t verbose_;
387
389 Double_t Tesla_;
390
392 TGeoNode* theNode_;
393
395 Bool_t gotNode_;
396
398
403 stringVect splitString(std::string& theString, std::string& splitter) const;
404
405};
406
407#endif
408
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 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.
Structure to hold transformation information.
Double_t y0_
The y translation displacement.
Double_t theta_
The second Euler rotation angle (about new X' axis)
Double_t z0_
The z translation displacement.
Double_t x0_
The x translation displacement.
Double_t phi_
The first Euler rotation angle (about Z axis)
Double_t psi_
The third Euler rotation angle (about new Z' axis)