Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FieldParameters.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25
26/// \file G4FieldParameters.hh
27/// \brief Definition of the G4FieldParameters class
28///
29/// This code was initially developed in Geant4 VMC package
30/// (https://github.com/vmc-project)
31/// and adapted to Geant4.
32///
33/// \author I. Hrivnacova; IJCLab, Orsay
34
35#ifndef G4FIELDPARAMETERS_HH
36#define G4FIELDPARAMETERS_HH
37
38#include "G4MagneticField.hh"
39#include "globals.hh"
40
42
44
47
48/// The available fields in Geant4
50{
51 kMagnetic, ///< magnetic field
52 kElectroMagnetic, ///< electromagnetic field
53 kGravity ///< gravity field
54};
55
56/// The available equations of motion of a particle in a field
57/// in Geant4
59{
60 kEqMagnetic, ///< G4Mag_UsualEqRhs: the standard right-hand side for
61 ///< equation of motion.
62 kEqMagneticWithSpin,///< G4Mag_SpinEqRhs: the equation of motion for a particle
63 ///< with spin
64 ///< in a pure magnetic field
65 kEqElectroMagnetic, ///< G4EqMagElectricField: Equation of motion in a combined
66 ///< electric and magnetic field
67 kEqEMfieldWithSpin, ///< G4EqEMFieldWithSpin: Equation of motion for a
68 ///< particle with spin
69 ///< in a combined electric and magnetic field
70 kEqEMfieldWithEDM, ///< G4EqEMFieldWithEDM: Equation of motion in a combined
71 ///< electric and magnetic field, with spin tracking for
72 ///< both MDM and EDM terms
73 kUserEquation ///< User defined equation of motion
74};
75
76/// The available integrator of particle's equation of motion
77/// in Geant4
79{
80 // steppers with equation of motion of generic type (G4EquationOfMotion)
81 kCashKarpRKF45, ///< G4CashKarpRKF45
82 kClassicalRK4, ///< G4ClassicalRK4
83 kBogackiShampine23, ///< G4BogackiShampine23
84 kBogackiShampine45, ///< G4BogackiShampine45
85 kDormandPrince745, ///< G4DormandPrince745
86 kDormandPrinceRK56, ///< G4DormandPrinceRK56
87 kDormandPrinceRK78, ///< G4DormandPrinceRK78
88 kExplicitEuler, ///< G4ExplicitEuler
89 kImplicitEuler, ///< G4ImplicitEuler
90 kSimpleHeum, ///< G4SimpleHeum
91 kSimpleRunge, ///< G4SimpleRunge
92 kTsitourasRK45, ///< G4TsitourasRK45
93
94 // steppers with equation of motion of G4Mag_UsualEqRhs type
95 kConstRK4, ///< G4ConstRK4
96 kExactHelixStepper, ///< G4ExactHelixStepper
97 kHelixExplicitEuler, ///< G4HelixExplicitEuler
98 kHelixHeum, ///< G4HelixHeum
99 kHelixImplicitEuler, ///< G4HelixImplicitEuler
100 kHelixMixedStepper, ///< G4HelixMixedStepper
101 kHelixSimpleRunge, ///< G4HelixSimpleRunge
102 kNystromRK4, ///< G4NystromRK4
103 kRKG3Stepper, ///< G4RKG3_Stepper
104 kUserStepper, ///< User defined stepper
105
106 // FSAL steppers
107 kRK547FEq1, ///< G4RK547FEq1
108 kRK547FEq2, ///< G4RK547FEq2
109 kRK547FEq3 ///< G4RK547FEq3
110};
111
112/// \brief The magnetic field parameters
113///
114/// The class defines the type of equation of motion of a particle
115/// in a field and the integration method, as well as other accuracy
116/// parameters.
117///
118/// The default values correspond to the defaults set in Geant4
119/// (taken from Geant4 9.3 release.)
120/// As Geant4 classes to not provide access methods for these defaults,
121/// the defaults have to be checked with each new Geant4 release.
122///
123/// \author I. Hrivnacova; IJCLab, Orsay
124
126{
127 public:
128 /// Standard and default constructor
129 G4FieldParameters(const G4String& volumeName = "");
130 /// Destructor
132
133 // Methods
134 //
135
136 /// Return the field type as a string
137 static G4String FieldTypeName(G4FieldType field);
138 /// Return the equation type as a string
140 /// Return the stepper type as a string
141 static G4String StepperTypeName(G4StepperType stepper);
142 /// Return the field type for given field type name
143 static G4FieldType GetFieldType(const G4String& name);
144 /// Return the equation type for given equation type name
145 static G4EquationType GetEquationType(const G4String& name);
146 /// Return the stepper type for given stepper type name
147 static G4StepperType GetStepperType(const G4String& name);
148
149 /// Prints all customizable accuracy parameters
150 void PrintParameters() const;
151
152 // Set methods
153 //
154
155 /// Set type of field
156 void SetFieldType(G4FieldType field);
157 /// Set Type of equation of motion of a particle in a field
158 void SetEquationType(G4EquationType equation);
159 /// Type of integrator of particle's equation of motion
160 void SetStepperType(G4StepperType stepper);
161 /// Set user defined equation of motion
163 /// Set user defined integrator of particle's equation of motion
165
166 /// Set minimum step in G4ChordFinder
167 void SetMinimumStep(G4double value);
168 /// Set delta chord in G4ChordFinder
169 void SetDeltaChord(G4double value);
170 /// Set delta one step in global field manager
171 void SetDeltaOneStep(G4double value);
172 /// Set delta intersection in global field manager
173 void SetDeltaIntersection(G4double value);
174 /// Set minimum epsilon step in global field manager
176 /// Set maximum epsilon step in global field manager
178 /// Set the distance within which the field is considered constant
179 void SetConstDistance(G4double value);
180
181 // Get methods
182 //
183
184 // Get the name of associated volume, if local field
185 G4String GetVolumeName() const;
186
187 /// Get type of field
189 /// Get type of equation of motion of a particle in a field
191 /// Get rype of integrator of particle's equation of motion
193 /// Get user defined equation of motion
195 /// Get user defined integrator of particle's equation of motion
197
198 /// Get minimum step in G4ChordFinder
199 G4double GetMinimumStep() const;
200 /// Get delta chord in G4ChordFinder
201 G4double GetDeltaChord() const;
202 /// Get delta one step in global field manager
204 /// Get delta intersection in global field manager
206 /// Get minimum epsilon step in global field manager
208 /// Get maximum epsilon step in global field manager
210 /// Get the distance within which the field is considered constant
212
213 private:
214 /// Not implemented
215 G4FieldParameters(const G4FieldParameters& right) = delete;
216 /// Not implemented
217 G4FieldParameters& operator=(const G4FieldParameters& right) = delete;
218
219 // static data members
220 //
221 /// Default minimum step in G4ChordFinder
222 inline static const G4double fgkDefaultMinimumStep = 0.01 * CLHEP::mm;
223 /// Default delta chord in G4ChordFinder
224 inline static const G4double fgkDefaultDeltaChord = 0.25 * CLHEP::mm;
225 /// Default delta one step in global field manager
226 inline static const G4double fgkDefaultDeltaOneStep = 0.01 * CLHEP::mm;
227 /// Delta intersection in global field manager
228 inline static const G4double fgkDefaultDeltaIntersection = 0.001 * CLHEP::mm;
229 /// Default minimum epsilon step in global field manager
230 inline static const G4double fgkDefaultMinimumEpsilonStep = 5.0e-5;
231 /// Default maximum epsilon step in global field manager
232 inline static const G4double fgkDefaultMaximumEpsilonStep = 0.001;
233 /// Default constant distance
234 inline static const G4double fgkDefaultConstDistance = 0.;
235
236 // data members
237 //
238 /// Messenger for this class
239 G4FieldParametersMessenger* fMessenger = nullptr;
240
241 /// The name of associated volume, if local field
242 G4String fVolumeName;
243
244 /// Minimum step in G4ChordFinder
245 G4double fMinimumStep = fgkDefaultMinimumStep;
246 /// Delta chord in G4ChordFinder
247 G4double fDeltaChord = fgkDefaultDeltaChord;
248 /// Delta one step in global field manager
249 G4double fDeltaOneStep = fgkDefaultDeltaOneStep;
250 /// Delta intersection in global field manager
251 G4double fDeltaIntersection = fgkDefaultDeltaIntersection;
252 /// Minimum epsilon step in global field manager
253 G4double fMinimumEpsilonStep = fgkDefaultMinimumEpsilonStep;
254 /// Maximum epsilon step in global field manager
255 G4double fMaximumEpsilonStep = fgkDefaultMaximumEpsilonStep;
256
257 /// Type of field
258 G4FieldType fField = kMagnetic;
259
260 /// Type of equation of motion of a particle in a field
261 G4EquationType fEquation = kEqMagnetic;
262
263 /// Type of integrator of particle's equation of motion
265
266 /// User defined equation of motion
267 G4EquationOfMotion* fUserEquation = nullptr;
268
269 /// User defined integrator of particle's equation of motion
270 G4MagIntegratorStepper* fUserStepper = nullptr;
271
272 /// The distance within which the field is considered constant
273 G4double fConstDistance = fgkDefaultConstDistance;
274};
275
276// inline functions
277
278// Set type of field
280{
281 fField = field;
282}
283
284// Set the type of equation of motion of a particle in a field
286{
287 fEquation = equation;
288}
289
290// Set the type of integrator of particle's equation of motion
292{
293 fStepper = stepper;
294}
295
296// Set minimum step in G4ChordFinder
298{
299 fMinimumStep = value;
300}
301
302// Set delta chord in G4ChordFinder
304{
305 fDeltaChord = value;
306}
307
308// Set delta one step in global field manager
310{
311 fDeltaOneStep = value;
312}
313
314// Set delta intersection in global field manager
316{
317 fDeltaIntersection = value;
318}
319
320// Set minimum epsilon step in global field manager
322{
323 fMinimumEpsilonStep = value;
324}
325
326// Set maximum epsilon step in global field manager
328{
329 fMaximumEpsilonStep = value;
330}
331
332// Set the distance within which the field is considered constant
334{
335 fConstDistance = value;
336}
337
338// Return the name of associated volume, if local field
340{
341 return fVolumeName;
342}
343
344// Return the type of field
345inline G4FieldType G4FieldParameters::GetFieldType() const { return fField; }
346
347// Return the type of equation of motion of a particle in a field
349{
350 return fEquation;
351}
352
353// Return the type of integrator of particle's equation of motion
355{
356 return fStepper;
357}
358
359// Return the user defined equation of motion
361{
362 return fUserEquation;
363}
364
365// Return the user defined integrator of particle's equation of motion
367{
368 return fUserStepper;
369}
370
371// Return minimum step in G4ChordFinder
373{
374 return fMinimumStep;
375}
376
377// Return delta chord in G4ChordFinder
379{
380 return fDeltaChord;
381}
382
383// Return delta one step in global field manager
385{
386 return fDeltaOneStep;
387}
388
389// Return delta intersection in global field manager
391{
392 return fDeltaIntersection;
393}
394
395// Return minimum epsilon step in global field manager
397{
398 return fMinimumEpsilonStep;
399}
400
401// Return maximum epsilon step in global field manager
403{
404 return fMaximumEpsilonStep;
405}
406
407// Return the distance within which the field is considered constant
409{
410 return fConstDistance;
411}
412
413#endif // G4FIELDPARAMETERS_HH
G4EquationType
@ kEqMagneticWithSpin
@ kUserEquation
User defined equation of motion.
@ kEqEMfieldWithSpin
@ kEqMagnetic
@ kEqElectroMagnetic
@ kEqEMfieldWithEDM
G4FieldType
The available fields in Geant4.
@ kElectroMagnetic
electromagnetic field
@ kGravity
gravity field
@ kMagnetic
magnetic field
G4StepperType
@ kRKG3Stepper
G4RKG3_Stepper.
@ kRK547FEq2
G4RK547FEq2.
@ kHelixSimpleRunge
G4HelixSimpleRunge.
@ kNystromRK4
G4NystromRK4.
@ kDormandPrince745
G4DormandPrince745.
@ kCashKarpRKF45
G4CashKarpRKF45.
@ kDormandPrinceRK78
G4DormandPrinceRK78.
@ kSimpleRunge
G4SimpleRunge.
@ kHelixImplicitEuler
G4HelixImplicitEuler.
@ kConstRK4
G4ConstRK4.
@ kUserStepper
User defined stepper.
@ kSimpleHeum
G4SimpleHeum.
@ kHelixHeum
G4HelixHeum.
@ kHelixExplicitEuler
G4HelixExplicitEuler.
@ kDormandPrinceRK56
G4DormandPrinceRK56.
@ kTsitourasRK45
G4TsitourasRK45.
@ kImplicitEuler
G4ImplicitEuler.
@ kExactHelixStepper
G4ExactHelixStepper.
@ kHelixMixedStepper
G4HelixMixedStepper.
@ kBogackiShampine45
G4BogackiShampine45.
@ kExplicitEuler
G4ExplicitEuler.
@ kRK547FEq1
G4RK547FEq1.
@ kRK547FEq3
G4RK547FEq3.
@ kBogackiShampine23
G4BogackiShampine23.
@ kClassicalRK4
G4ClassicalRK4.
double G4double
Definition G4Types.hh:83
Messenger class that defines commands for TG4DetConstruction.
void SetEquationType(G4EquationType equation)
Set Type of equation of motion of a particle in a field.
G4double GetDeltaIntersection() const
Get delta intersection in global field manager.
G4String GetVolumeName() const
void SetFieldType(G4FieldType field)
Set type of field.
void SetUserStepper(G4MagIntegratorStepper *stepper)
Set user defined integrator of particle's equation of motion.
G4double GetMinimumEpsilonStep() const
Get minimum epsilon step in global field manager.
G4double GetDeltaChord() const
Get delta chord in G4ChordFinder.
G4MagIntegratorStepper * GetUserStepper() const
Get user defined integrator of particle's equation of motion.
void SetUserEquationOfMotion(G4EquationOfMotion *equation)
Set user defined equation of motion.
void SetConstDistance(G4double value)
Set the distance within which the field is considered constant.
G4double GetMinimumStep() const
Get minimum step in G4ChordFinder.
void SetMinimumStep(G4double value)
Set minimum step in G4ChordFinder.
void SetStepperType(G4StepperType stepper)
Type of integrator of particle's equation of motion.
void SetMinimumEpsilonStep(G4double value)
Set minimum epsilon step in global field manager.
G4double GetConstDistance() const
Get the distance within which the field is considered constant.
static G4String EquationTypeName(G4EquationType equation)
Return the equation type as a string.
void PrintParameters() const
Prints all customizable accuracy parameters.
G4EquationOfMotion * GetUserEquationOfMotion() const
Get user defined equation of motion.
G4FieldParameters(const G4String &volumeName="")
Standard and default constructor.
G4double GetDeltaOneStep() const
Get delta one step in global field manager.
void SetDeltaOneStep(G4double value)
Set delta one step in global field manager.
void SetDeltaIntersection(G4double value)
Set delta intersection in global field manager.
G4double GetMaximumEpsilonStep() const
Get maximum epsilon step in global field manager.
G4FieldType GetFieldType() const
Get type of field.
static G4String StepperTypeName(G4StepperType stepper)
Return the stepper type as a string.
void SetDeltaChord(G4double value)
Set delta chord in G4ChordFinder.
static G4String FieldTypeName(G4FieldType field)
Return the field type as a string.
G4StepperType GetStepperType() const
Get rype of integrator of particle's equation of motion.
~G4FieldParameters()
Destructor.
void SetMaximumEpsilonStep(G4double value)
Set maximum epsilon step in global field manager.
G4EquationType GetEquationType() const
Get type of equation of motion of a particle in a field.