QwGeant4
QweakSimMagneticField Class Reference

#include <QweakSimMagneticField.hh>

Inherits G4MagneticField.

+ Collaboration diagram for QweakSimMagneticField:

Public Member Functions

 QweakSimMagneticField ()
 Default constructor. More...
 
virtual ~QweakSimMagneticField ()
 Virtual destructor. More...
 
void GetFieldValue (const G4double point[4], G4double *field) const
 Get the field value. More...
 
void PrintFieldValue (const G4double point[4]) const
 Print the field value. More...
 
void SetScaleFactor (const double scalefactor)
 Set the scale factor. More...
 
double GetScaleFactor () const
 Get the scale factor. More...
 
void SetReferenceCurrent (const double referencecurrent)
 Set the reference current. More...
 
double GetReferenceCurrent () const
 Get the reference current. More...
 
void SetActualCurrent (const double actualcurrent)
 Set the actual current. More...
 
double GetActualCurrent () const
 Get the actual current. More...
 
void SetRotation (const double rotation)
 Set the field rotation around z. More...
 
double GetRotation () const
 Get the field rotation around z. More...
 
void SetTranslation (const double translation)
 Set the field translation along z. More...
 
double GetTranslation () const
 Get the field translation along z. More...
 
void SetTranslationVector (const G4ThreeVector translationvector)
 Set the field translation vector. More...
 
G4ThreeVector GetTranslationVector () const
 Get the field translation vector. More...
 
void SetMinimum (const double value, const int i)
 Set the minimum. More...
 
void SetMaximum (const double value, const int i)
 Set the maximum. More...
 
void SetStep (const double value, const int i)
 Set the step. More...
 
void SetWrap (const int value, const int i)
 Set the wrap. More...
 
void ReadFieldMap (const G4String &filename)
 Read the field map. More...
 
void ReadFieldMapText (const G4String &filename)
 
void ReadFieldMapBinary (const G4String &filename)
 
void TestFieldMap (const G4double point[4], const G4double exact[3]) const
 Test the field map. More...
 

Private Types

typedef float value_t
 Field map storage data type. More...
 

Private Attributes

QweakSimMagneticFieldMessengerfMagneticFieldMessenger
 Messenger. More...
 
QweakSimFieldMap< value_t,
value_n > * 
fField
 
double fRotation
 Field rotation and translation with respect to the z axis. More...
 
double fRotationCos
 
double fRotationSin
 
double fTranslation
 
G4ThreeVector fTranslationVector
 
double fScaleFactor
 Field map reference and actual currents. More...
 
double fActualCurrent
 
double fReferenceCurrent
 
std::vector< double > fMin
 Field map grid definitions. More...
 
std::vector< double > fMax
 
std::vector< double > fStep
 
std::vector< size_t > fWrap
 

Static Private Attributes

static const unsigned int value_n = 3
 Number of field component in the map (x,y,z,r,phi) More...
 

Detailed Description

Definition at line 24 of file QweakSimMagneticField.hh.

Member Typedef Documentation

typedef float QweakSimMagneticField::value_t
private

Field map storage data type.

Definition at line 32 of file QweakSimMagneticField.hh.

Constructor & Destructor Documentation

QweakSimMagneticField::QweakSimMagneticField ( )

Default constructor.

Definition at line 48 of file QweakSimMagneticField.cc.

References fActualCurrent, fField, fMagneticFieldMessenger, fMax, fMin, fReferenceCurrent, fScaleFactor, fStep, fWrap, SetRotation(), SetTranslation(), and SetTranslationVector().

49 {
50  G4cout << G4endl << "###### Calling QweakSimMagneticField::QweakSimMagneticField() " << G4endl << G4endl;
51 
52  // No field map
53  fField = 0;
54 
55  // Messenger
57 
58  // Default settings
59  SetRotation(0.0);
60  SetTranslation(0.0);
61  SetTranslationVector(G4ThreeVector(0.0,0.0,0.0));
62 
63  // Default field map grid (peiqing_2007.dat)
64  fMin.push_back(-250.0); fMax.push_back(250.0); fStep.push_back(2.0); fWrap.push_back(0);
65  fMin.push_back(2.0); fMax.push_back(300.0); fStep.push_back(2.0); fWrap.push_back(0);
66  fMin.push_back(-0.5 * degree/radian); fMax.push_back(360.5 * degree/radian); fStep.push_back(1.0 * degree/radian); fWrap.push_back(2);
67 
68  // Set current and scale factors
69  G4double current = 8921.0;
70  fActualCurrent = current;
71  fReferenceCurrent = current;
72  fScaleFactor = 1.0;
73 
74  G4cout << G4endl << "###### Leaving QweakSimMagneticField::QweakSimMagneticField " << G4endl << G4endl;
75 }
QweakSimFieldMap< value_t, value_n > * fField
void SetTranslationVector(const G4ThreeVector translationvector)
Set the field translation vector.
void SetTranslation(const double translation)
Set the field translation along z.
QweakSimMagneticFieldMessenger * fMagneticFieldMessenger
Messenger.
std::vector< size_t > fWrap
std::vector< double > fMax
std::vector< double > fMin
Field map grid definitions.
void SetRotation(const double rotation)
Set the field rotation around z.
double fScaleFactor
Field map reference and actual currents.
std::vector< double > fStep
Scans input file for /MagnetField/xyz commands.

+ Here is the call graph for this function:

QweakSimMagneticField::~QweakSimMagneticField ( )
virtual

Virtual destructor.

Definition at line 78 of file QweakSimMagneticField.cc.

References fField, and fMagneticFieldMessenger.

79 {
80  G4cout << G4endl << "###### Calling QweakSimMagneticField::~QweakSimMagneticField() " << G4endl << G4endl;
81 
82  // Delete existing magnetic field
83  if (fField) delete fField;
84 
85  // Delete messenger
87 
88  G4cout << G4endl << "###### Leaving QweakSimMagneticField::~QweakSimMagneticField() " << G4endl << G4endl;
89 }
QweakSimFieldMap< value_t, value_n > * fField
QweakSimMagneticFieldMessenger * fMagneticFieldMessenger
Messenger.

Member Function Documentation

double QweakSimMagneticField::GetActualCurrent ( ) const
inline

Get the actual current.

Definition at line 78 of file QweakSimMagneticField.hh.

References fActualCurrent.

78 { return fActualCurrent; };
void QweakSimMagneticField::GetFieldValue ( const G4double  point[4],
G4double *  field 
) const

Get the field value.

Definition at line 262 of file QweakSimMagneticField.cc.

References fField, fScaleFactor, and fTranslationVector.

Referenced by PrintFieldValue(), and TestFieldMap().

263 {
264  // Shift by translation vector
265  G4double x = point[0] - fTranslationVector[0];
266  G4double y = point[1] - fTranslationVector[1];
267  G4double z = point[2] - fTranslationVector[2];
268 
269  // Convert from cartesian to cylindrical coordinates
270  G4double r = sqrt(x*x + y*y);
271  G4double phi = atan2(y, x);
272  if (phi < 0.0) phi += 360.0 * degree;
273 
274  // Retrieve field value
275  field[0] = 0.0;
276  field[1] = 0.0;
277  field[2] = 0.0;
278  if (fField) {
279  // Get field value
280  G4double coord[3] = {z/cm, r/cm, phi/radian};
281  fField->GetValue(coord,field);
282 
283  // Units and scale factor
284  field[0] *= fScaleFactor * tesla;
285  field[1] *= fScaleFactor * tesla;
286  field[2] *= fScaleFactor * tesla;
287  }
288 }
QweakSimFieldMap< value_t, value_n > * fField
double fScaleFactor
Field map reference and actual currents.

+ Here is the caller graph for this function:

double QweakSimMagneticField::GetReferenceCurrent ( ) const
inline

Get the reference current.

Definition at line 68 of file QweakSimMagneticField.hh.

References fReferenceCurrent.

double QweakSimMagneticField::GetRotation ( ) const
inline

Get the field rotation around z.

Definition at line 89 of file QweakSimMagneticField.hh.

References fRotation.

89 { return fRotation; };
double fRotation
Field rotation and translation with respect to the z axis.
double QweakSimMagneticField::GetScaleFactor ( ) const
inline

Get the scale factor.

Definition at line 58 of file QweakSimMagneticField.hh.

References fScaleFactor.

58 { return fScaleFactor; };
double fScaleFactor
Field map reference and actual currents.
double QweakSimMagneticField::GetTranslation ( ) const
inline

Get the field translation along z.

Definition at line 97 of file QweakSimMagneticField.hh.

References fTranslation.

97 { return fTranslation; };
G4ThreeVector QweakSimMagneticField::GetTranslationVector ( ) const
inline

Get the field translation vector.

Definition at line 105 of file QweakSimMagneticField.hh.

References fTranslationVector.

105 { return fTranslationVector; };
void QweakSimMagneticField::PrintFieldValue ( const G4double  point[4]) const

Print the field value.

Definition at line 243 of file QweakSimMagneticField.cc.

References GetFieldValue().

Referenced by QweakSimMagneticFieldMessenger::SetNewValue().

244 {
245  G4double field[3];
246  GetFieldValue(point,field);
247 
248  // Calculate field value
249  G4cout << "Calculating test field value at cartesian position "
250  << "(" << point[0]/cm << "," << point[1]/cm << "," << point[2]/cm << ") cm "
251  << G4endl;
252  G4double value[3] = {0.0, 0.0, 0.0};
253  GetFieldValue(point,value);
254 
255  // Output
256  G4cout << " Value = "
257  << "(" << value[0]/kilogauss << ", " << value[1]/kilogauss << ", " << value[2]/kilogauss << ") kG."
258  << G4endl;
259 }
void GetFieldValue(const G4double point[4], G4double *field) const
Get the field value.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QweakSimMagneticField::ReadFieldMap ( const G4String &  filename)

Read the field map.

Definition at line 92 of file QweakSimMagneticField.cc.

References ReadFieldMapBinary(), ReadFieldMapText(), and TestFieldMap().

Referenced by QweakSimMagneticFieldMessenger::SetNewValue().

93 {
94  G4cout << G4endl << "###### Calling QweakSimMagneticField::ReadFieldMap() " << G4endl << G4endl;
95 
96  // Depending on form of filename, read regular/binary field map
97  if (filename.find(".dat") != std::string::npos) {
98  ReadFieldMapText(filename);
99  }
100  if (filename.find(".bin") != std::string::npos) {
101  ReadFieldMapBinary(filename);
102  }
103 
104  // Test new field map
105  G4double point1[4] = {100.0 * cm, 0.0 * cm, 100.0 * cm, 0.0};
106  G4double exact1[3] = {-0.0499845 * kilogauss, 3.28516 * kilogauss, -0.0112704 * kilogauss};
107  TestFieldMap(point1,exact1);
108  G4double point2[4] = {0.0 * cm, 100.0 * cm, 100.0 * cm, 0.0};
109  G4double exact2[3] = {-3.28505 * kilogauss, -0.0030823 * kilogauss, 0.0014373 * kilogauss};
110  TestFieldMap(point2,exact2);
111  G4double point3[4] = {100.0 * cm, 100.0 * cm, 100.0 * cm, 0.0};
112  G4double exact3[3] = {-0.838069 * kilogauss, 0.836675 * kilogauss, 0.00004997 * kilogauss};
113  TestFieldMap(point3,exact3);
114 
115  G4cout << G4endl << "###### Leaving QweakSimMagneticField::ReadFieldMap() " << G4endl << G4endl;
116 }
void ReadFieldMapText(const G4String &filename)
void ReadFieldMapBinary(const G4String &filename)
void TestFieldMap(const G4double point[4], const G4double exact[3]) const
Test the field map.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QweakSimMagneticField::ReadFieldMapBinary ( const G4String &  filename)

Definition at line 156 of file QweakSimMagneticField.cc.

References fField.

Referenced by ReadFieldMap().

157 {
158  // Delete old field map
159  if (fField) delete fField;
160 
161  // Create new field map
164 
165  // Set debug level to something reasonable
167 }
QweakSimFieldMap< value_t, value_n > * fField
A multi-dimensional grid of values with interpolation methods.

+ Here is the caller graph for this function:

void QweakSimMagneticField::ReadFieldMapText ( const G4String &  filename)

Definition at line 170 of file QweakSimMagneticField.cc.

References fField, fMax, fMin, fRotation, fRotationCos, fRotationSin, fStep, fTranslation, and fWrap.

Referenced by ReadFieldMap().

171 {
172  // Delete old field map
173  if (fField) delete fField;
174 
175  // Create new field map
178 
179  // Fill vector of grid min, max, and step size
180  fField->SetDimensions(fMin.size());
181  fField->SetMinimumMaximumStep(fMin,fMax,fStep);
182  fField->SetWrapCoordinate(fWrap);
183 
184  // Debug output
185  G4cout << "Reading field grid: ";
186 
187  // Open the field map file
188  std::ifstream inputfile;
189  inputfile.open(filename);
190  while (inputfile.good()) {
191  // Progress bar
192  if (fField->GetCurrentEntries() % (fField->GetMaximumEntries() / 10) == 0) {
193  int pct = fField->GetCurrentEntries() / (fField->GetMaximumEntries() / 100);
194  G4cout << pct << "%" << std::flush;
195  }
196  if (fField->GetCurrentEntries() % (fField->GetMaximumEntries() / 10) != 0
197  && fField->GetCurrentEntries() % (fField->GetMaximumEntries() / 40) == 0) {
198  G4cout << "." << std::flush;
199  }
200 
201  // Read variables
202  G4double r = 0.0, z = 0.0, phi = 0.0;
203  value_t b[3] = {0.0, 0.0, 0.0};
204  inputfile >> r >> z >> phi >> b[0] >> b[1] >> b[2];
205 
206  // Transform units into Geant4 system of units
207  r *= cm;
208  z *= cm;
209  phi *= degree;
210  b[0] *= kilogauss;
211  b[1] *= kilogauss;
212  b[2] *= kilogauss;
213 
214  // Translation
215  z -= fTranslation;
216 
217  // Rotation
218  phi -= fRotation;
219  if (phi < 0.0 * degree) phi += 360.0 * degree;
220  if (phi > 360.0 * degree) phi -= 360.0 * degree;
221  value_t b_rot[2];
222  b_rot[0] = b[0] * fRotationCos + b[1] * fRotationSin;
223  b_rot[1] = -b[0] * fRotationSin + b[1] * fRotationCos;
224  b[0] = b_rot[0]; b[1] = b_rot[1];
225 
226  // Transform units into fieldmap system of units
227  double coord[3] = { z/cm, r/cm, phi/radian };
228  value_t field[3] = { value_t(b[0]/tesla), value_t(b[1]/tesla), value_t(b[2]/tesla) };
229 
230  // Store into fieldmap
231  fField->Set(coord,field);
232  }
233  G4cout << G4endl;
234 
235  // Close file
236  inputfile.close();
237 
238  // Set debug level to something reasonable
240 }
QweakSimFieldMap< value_t, value_n > * fField
std::vector< size_t > fWrap
A multi-dimensional grid of values with interpolation methods.
std::vector< double > fMax
double fRotation
Field rotation and translation with respect to the z axis.
std::vector< double > fMin
Field map grid definitions.
std::vector< double > fStep
float value_t
Field map storage data type.

+ Here is the caller graph for this function:

void QweakSimMagneticField::SetActualCurrent ( const double  actualcurrent)
inline

Set the actual current.

Definition at line 71 of file QweakSimMagneticField.hh.

References fActualCurrent, fReferenceCurrent, and SetScaleFactor().

Referenced by QweakSimMagneticFieldMessenger::SetNewValue().

71  {
72  fActualCurrent = actualcurrent;
73  G4cout << "Magnetic field actual current is now: " << fActualCurrent << G4endl;
74  if (fReferenceCurrent != 0)
76  };
void SetScaleFactor(const double scalefactor)
Set the scale factor.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QweakSimMagneticField::SetMaximum ( const double  value,
const int  i 
)
inline

Set the maximum.

Definition at line 114 of file QweakSimMagneticField.hh.

References fMax.

Referenced by QweakSimMagneticFieldMessenger::SetNewValue().

114  {
115  fMax.at(i) = value;
116  G4cout << "Magnetic field grid max " << i << " is " << fMax.at(i) << " cm/radian" << G4endl;
117  };
std::vector< double > fMax

+ Here is the caller graph for this function:

void QweakSimMagneticField::SetMinimum ( const double  value,
const int  i 
)
inline

Set the minimum.

Definition at line 109 of file QweakSimMagneticField.hh.

References fMin.

Referenced by QweakSimMagneticFieldMessenger::SetNewValue().

109  {
110  fMin.at(i) = value;
111  G4cout << "Magnetic field grid min " << i << " is " << fMin.at(i) << " cm/radian" << G4endl;
112  };
std::vector< double > fMin
Field map grid definitions.

+ Here is the caller graph for this function:

void QweakSimMagneticField::SetReferenceCurrent ( const double  referencecurrent)
inline

Set the reference current.

Definition at line 61 of file QweakSimMagneticField.hh.

References fActualCurrent, fReferenceCurrent, and SetScaleFactor().

Referenced by QweakSimMagneticFieldMessenger::SetNewValue().

61  {
62  fReferenceCurrent = referencecurrent;
63  G4cout << "Magnetic field reference current is now: " << fReferenceCurrent << G4endl;
64  if (fReferenceCurrent != 0)
66  };
void SetScaleFactor(const double scalefactor)
Set the scale factor.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QweakSimMagneticField::SetRotation ( const double  rotation)
inline

Set the field rotation around z.

Definition at line 82 of file QweakSimMagneticField.hh.

References fRotation, fRotationCos, and fRotationSin.

Referenced by QweakSimMagneticField(), and QweakSimMagneticFieldMessenger::SetNewValue().

82  {
83  fRotation = rotation;
84  fRotationCos = cos(rotation);
85  fRotationSin = sin(rotation);
86  G4cout << "Magnetic field rotation is now: " << fRotation/degree << " degree" << G4endl;
87  };
double fRotation
Field rotation and translation with respect to the z axis.

+ Here is the caller graph for this function:

void QweakSimMagneticField::SetScaleFactor ( const double  scalefactor)
inline

Set the scale factor.

Definition at line 53 of file QweakSimMagneticField.hh.

References fScaleFactor.

Referenced by SetActualCurrent(), QweakSimMagneticFieldMessenger::SetNewValue(), and SetReferenceCurrent().

53  {
54  fScaleFactor = scalefactor;
55  G4cout << "Magnetic field scale factor is now: " << fScaleFactor << G4endl;
56  };
double fScaleFactor
Field map reference and actual currents.

+ Here is the caller graph for this function:

void QweakSimMagneticField::SetStep ( const double  value,
const int  i 
)
inline

Set the step.

Definition at line 119 of file QweakSimMagneticField.hh.

References fStep.

Referenced by QweakSimMagneticFieldMessenger::SetNewValue().

119  {
120  fStep.at(i) = value;
121  G4cout << "Magnetic field grid step " << i << " is " << fStep.at(i) << " cm/radian" << G4endl;
122  };
std::vector< double > fStep

+ Here is the caller graph for this function:

void QweakSimMagneticField::SetTranslation ( const double  translation)
inline

Set the field translation along z.

Definition at line 92 of file QweakSimMagneticField.hh.

References fTranslation.

Referenced by QweakSimMagneticField(), and QweakSimMagneticFieldMessenger::SetNewValue().

92  {
93  fTranslation = translation;
94  G4cout << "Magnetic field translation is now: " << fTranslation/cm << " cm" << G4endl;
95  };

+ Here is the caller graph for this function:

void QweakSimMagneticField::SetTranslationVector ( const G4ThreeVector  translationvector)
inline

Set the field translation vector.

Definition at line 100 of file QweakSimMagneticField.hh.

References fTranslationVector.

Referenced by QweakSimMagneticField(), and QweakSimMagneticFieldMessenger::SetNewValue().

100  {
101  fTranslationVector = translationvector;
102  G4cout << "Magnetic field translation is now: " << fTranslationVector/cm << " cm" << G4endl;
103  };

+ Here is the caller graph for this function:

void QweakSimMagneticField::SetWrap ( const int  value,
const int  i 
)
inline

Set the wrap.

Definition at line 124 of file QweakSimMagneticField.hh.

References fWrap.

Referenced by QweakSimMagneticFieldMessenger::SetNewValue().

124  {
125  fWrap.at(i) = value;
126  G4cout << "Magnetic field grid wrap " << i << " is " << fWrap.at(i) << G4endl;
127  };
std::vector< size_t > fWrap

+ Here is the caller graph for this function:

void QweakSimMagneticField::TestFieldMap ( const G4double  point[4],
const G4double  exact[3] 
) const

Test the field map.

Definition at line 119 of file QweakSimMagneticField.cc.

References GetFieldValue().

Referenced by ReadFieldMap().

120 {
121  // Calculate field value
122  G4cout << "Calculating test field value at cartesian position "
123  << "(" << point[0]/cm << "," << point[1]/cm << "," << point[2]/cm << ") cm "
124  << G4endl;
125  G4double value[3] = {0.0, 0.0, 0.0};
126  GetFieldValue(point,value);
127 
128  // Calculate difference
129  G4double diff[3] = {0.0, 0.0, 0.0};
130  G4double norm = 0.0;
131  for (size_t i = 0; i < 3; i++) {
132  diff[i] = value[i] - exact[i];
133  norm += diff[i] * diff[i];
134  }
135  norm = sqrt(norm);
136 
137  // Output
138  G4cout << " Value = "
139  << "(" << value[0]/kilogauss << ", " << value[1]/kilogauss << ", " << value[2]/kilogauss << ") kG."
140  << G4endl;
141  G4cout << " Exact = "
142  << "(" << exact[0]/kilogauss << ", " << exact[1]/kilogauss << ", " << exact[2]/kilogauss << ") kG."
143  << G4endl;
144  G4cout << " Diff = "
145  << "(" << diff[0]/kilogauss << ", " << diff[1]/kilogauss << ", " << diff[2]/kilogauss << ") kG."
146  << G4endl;
147 
148  // Test difference (0.1 kilogauss is of the order of 2%)
149  if (norm > 0.1 * kilogauss) {
150  G4cerr << "Magnetic field is different from expected value by " << norm/kilogauss << " kG."
151  << G4endl;
152  }
153 }
void GetFieldValue(const G4double point[4], G4double *field) const
Get the field value.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

Field Documentation

double QweakSimMagneticField::fActualCurrent
private
QweakSimMagneticFieldMessenger* QweakSimMagneticField::fMagneticFieldMessenger
private

Messenger.

Definition at line 35 of file QweakSimMagneticField.hh.

Referenced by QweakSimMagneticField(), and ~QweakSimMagneticField().

std::vector<double> QweakSimMagneticField::fMax
private

Definition at line 157 of file QweakSimMagneticField.hh.

Referenced by QweakSimMagneticField(), ReadFieldMapText(), and SetMaximum().

std::vector<double> QweakSimMagneticField::fMin
private

Field map grid definitions.

Definition at line 157 of file QweakSimMagneticField.hh.

Referenced by QweakSimMagneticField(), ReadFieldMapText(), and SetMinimum().

double QweakSimMagneticField::fReferenceCurrent
private
double QweakSimMagneticField::fRotation
private

Field rotation and translation with respect to the z axis.

Definition at line 147 of file QweakSimMagneticField.hh.

Referenced by GetRotation(), ReadFieldMapText(), and SetRotation().

double QweakSimMagneticField::fRotationCos
private

Definition at line 147 of file QweakSimMagneticField.hh.

Referenced by ReadFieldMapText(), and SetRotation().

double QweakSimMagneticField::fRotationSin
private

Definition at line 147 of file QweakSimMagneticField.hh.

Referenced by ReadFieldMapText(), and SetRotation().

double QweakSimMagneticField::fScaleFactor
private

Field map reference and actual currents.

Definition at line 152 of file QweakSimMagneticField.hh.

Referenced by GetFieldValue(), GetScaleFactor(), QweakSimMagneticField(), and SetScaleFactor().

std::vector<double> QweakSimMagneticField::fStep
private

Definition at line 157 of file QweakSimMagneticField.hh.

Referenced by QweakSimMagneticField(), ReadFieldMapText(), and SetStep().

double QweakSimMagneticField::fTranslation
private

Definition at line 148 of file QweakSimMagneticField.hh.

Referenced by GetTranslation(), ReadFieldMapText(), and SetTranslation().

G4ThreeVector QweakSimMagneticField::fTranslationVector
private
std::vector<size_t> QweakSimMagneticField::fWrap
private

Definition at line 158 of file QweakSimMagneticField.hh.

Referenced by QweakSimMagneticField(), ReadFieldMapText(), and SetWrap().

const unsigned int QweakSimMagneticField::value_n = 3
staticprivate

Number of field component in the map (x,y,z,r,phi)

Definition at line 29 of file QweakSimMagneticField.hh.


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