FairShip
Loading...
Searching...
No Matches
ShipFieldMaker.h
Go to the documentation of this file.
1// SPDX-License-Identifier: LGPL-3.0-or-later
2// SPDX-FileCopyrightText: Copyright CERN for the benefit of the SHiP
3// Collaboration
4
20#ifndef FIELD_SHIPFIELDMAKER_H_
21#define FIELD_SHIPFIELDMAKER_H_
22
23#include <map>
24#include <string>
25#include <vector>
26
27#include "ShipCompField.h"
28#include "TG4VUserPostDetConstruction.h"
29#include "TString.h"
30#include "TVector2.h"
31#include "TVector3.h"
32#include "TVirtualMagField.h"
33
34class TGeoMatrix;
35class TGeoNode;
36class TGeoVolume;
37
38class ShipFieldMaker : public TG4VUserPostDetConstruction {
39 public:
41 explicit 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 {
55 TString volName_;
57 TString fieldName_;
59 Double_t scale_;
60
62 fieldInfo() : volName_(""), fieldName_(""), scale_(1.0) {};
63
65 fieldInfo(const TString& volName, const TString& fieldName,
66 Double_t scale = 1.0)
67 : volName_(volName), fieldName_(fieldName), scale_(scale) {};
68 };
69
71 virtual void Construct();
72
74
77 void readInputFile(const std::string& inputFile);
78
80
84 void defineUniform(const TString& name, const TVector3& BVector);
85
87
94 void defineConstant(const TString& name, const TVector2& xRange,
95 const TVector2& yRange, const TVector2& zRange,
96 const TVector3& BVector);
97
99
110 void defineBell(const TString& name, Double_t BPeak, Double_t zMiddle,
111 Int_t orient = 1, Double_t tubeR = 500.0, Double_t xy = 0.0,
112 Double_t z = 0.0, Double_t L = 0.0);
113
114 // ! Define a field map
126 void defineFieldMap(const TString& name, const TString& mapFileName,
127 const TVector3& localCentre = TVector3(0.0, 0.0, 0.0),
128 const TVector3& localAngles = TVector3(0.0, 0.0, 0.0),
129 Bool_t useSymmetry = kFALSE);
130
133
139 void defineFieldMapCopy(const TString& name, const TString& mapNameToCopy,
140 const TVector3& translation,
141 const TVector3& eulerAngles = TVector3(0.0, 0.0,
142 0.0));
143
145
152 void defineComposite(const TString& name, const TString& field1Name,
153 const TString& field2Name,
154 const TString& field3Name = "",
155 const TString& field4Name = "");
156
158
162 void defineComposite(const TString& name, std::vector<TString> fieldNames);
163
165
171 void defineGlobalField(const TString& field1Name,
172 const TString& field2Name = "",
173 const TString& field3Name = "",
174 const TString& field4Name = "");
175
177
180 void defineGlobalField(std::vector<TString> fieldNames);
181
183
188 void defineRegionField(const TString& volName, const TString& fieldName,
189 Double_t scale = 1.0);
190
192
197 void defineLocalField(const TString& volName, const TString& fieldName,
198 Double_t scale = 1.0);
199
201
205
208
212 SFMap getAllFields() const { return theFields_; }
213
215
219 TVirtualMagField* getVolumeField(const TString& volName) const;
220
222
227 Bool_t gotField(const TString& name) const;
228
230
234 TVirtualMagField* getField(const TString& name) const;
235
237
245 void plotXYField(const TVector3& xAxis, const TVector3& yAxis,
246 const std::string& plotFile) const;
247
249
257 void plotZXField(const TVector3& zAxis, const TVector3& xAxis,
258 const std::string& plotFile) const;
259
261
269 void plotZYField(const TVector3& zAxis, const TVector3& yAxis,
270 const std::string& plotFile) const;
271
273
282 void plotField(Int_t type, const TVector3& xAxis, const TVector3& yAxis,
283 const std::string& plotFile) const;
284
285 void generateFieldMap(TString fileName, const float step = 2.5,
286 const float xRange = 179, const float yRange = 317,
287 const float zRange = 1515.5,
288 const float zShift = -4996);
290
293
294 protected:
298 Double_t x0_;
300 Double_t y0_;
302 Double_t z0_;
303
305 Double_t phi_;
307 Double_t theta_;
309 Double_t psi_;
310 };
311
313
316 void defineUniform(const stringVect& inputLine);
317
319
322 void defineConstant(const stringVect& inputLine);
323
325
328 void defineBell(const stringVect& inputLine);
329
331
335 void defineFieldMap(const stringVect& inputLine, Bool_t useSymmetry = kFALSE);
336
339
342 void defineFieldMapCopy(const stringVect& inputLine);
343
345
348 void defineComposite(const stringVect& inputLine);
349
351
354 void defineGlobalField(const stringVect& inputLine);
355
358
361 void defineRegionField(const stringVect& inputLine);
362
363 // ! Setup all of the regional fields. Called by Construct()
364 void setAllRegionFields();
365
367
370 void defineLocalField(const stringVect& inputLine);
371
374
381 void checkLocalFieldMap(TVirtualMagField*& localField, const TString& volName,
382 Double_t scale);
383
384 // ! Setup all of the local fields. Called by Construct()
385 void setAllLocalFields();
386
388
392 void getTransformation(const TString& volName, transformInfo& theInfo);
393
396
400 void findNode(TGeoVolume* aVolume, const TString& volName);
401
402 private:
405
408
410 std::vector<fieldInfo> regionInfo_;
411
413 std::vector<fieldInfo> localInfo_;
414
416 Bool_t verbose_;
417
419 Double_t Tesla_;
420
422 TGeoNode* theNode_;
423
425 Bool_t gotNode_;
426
428
433 stringVect splitString(std::string& theString, std::string& splitter) const;
434};
435
436#endif // FIELD_SHIPFIELDMAKER_H_
Class that defines a magnetic field composed from many fields.
Definition: ShipCompField.h:19
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)
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
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))
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)
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)