Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FieldParameters Class Reference

The magnetic field parameters. More...

#include <G4FieldParameters.hh>

Public Member Functions

 G4FieldParameters (const G4String &volumeName="")
 Standard and default constructor.
 
 ~G4FieldParameters ()
 Destructor.
 
void PrintParameters () const
 Prints all customizable accuracy parameters.
 
void SetFieldType (G4FieldType field)
 Set type of field.
 
void SetEquationType (G4EquationType equation)
 Set Type of equation of motion of a particle in a field.
 
void SetStepperType (G4StepperType stepper)
 Type of integrator of particle's equation of motion.
 
void SetUserEquationOfMotion (G4EquationOfMotion *equation)
 Set user defined equation of motion.
 
void SetUserStepper (G4MagIntegratorStepper *stepper)
 Set user defined integrator of particle's equation of motion.
 
void SetMinimumStep (G4double value)
 Set minimum step in G4ChordFinder.
 
void SetDeltaChord (G4double value)
 Set delta chord in G4ChordFinder.
 
void SetDeltaOneStep (G4double value)
 Set delta one step in global field manager.
 
void SetDeltaIntersection (G4double value)
 Set delta intersection in global field manager.
 
void SetMinimumEpsilonStep (G4double value)
 Set minimum epsilon step in global field manager.
 
void SetMaximumEpsilonStep (G4double value)
 Set maximum epsilon step in global field manager.
 
void SetConstDistance (G4double value)
 Set the distance within which the field is considered constant.
 
G4String GetVolumeName () const
 
G4FieldType GetFieldType () const
 Get type of field.
 
G4EquationType GetEquationType () const
 Get type of equation of motion of a particle in a field.
 
G4StepperType GetStepperType () const
 Get rype of integrator of particle's equation of motion.
 
G4EquationOfMotionGetUserEquationOfMotion () const
 Get user defined equation of motion.
 
G4MagIntegratorStepperGetUserStepper () const
 Get user defined integrator of particle's equation of motion.
 
G4double GetMinimumStep () const
 Get minimum step in G4ChordFinder.
 
G4double GetDeltaChord () const
 Get delta chord in G4ChordFinder.
 
G4double GetDeltaOneStep () const
 Get delta one step in global field manager.
 
G4double GetDeltaIntersection () const
 Get delta intersection in global field manager.
 
G4double GetMinimumEpsilonStep () const
 Get minimum epsilon step in global field manager.
 
G4double GetMaximumEpsilonStep () const
 Get maximum epsilon step in global field manager.
 
G4double GetConstDistance () const
 Get the distance within which the field is considered constant.
 

Static Public Member Functions

static G4String FieldTypeName (G4FieldType field)
 Return the field type as a string.
 
static G4String EquationTypeName (G4EquationType equation)
 Return the equation type as a string.
 
static G4String StepperTypeName (G4StepperType stepper)
 Return the stepper type as a string.
 
static G4FieldType GetFieldType (const G4String &name)
 Return the field type for given field type name.
 
static G4EquationType GetEquationType (const G4String &name)
 Return the equation type for given equation type name.
 
static G4StepperType GetStepperType (const G4String &name)
 Return the stepper type for given stepper type name.
 

Detailed Description

The magnetic field parameters.

The class defines the type of equation of motion of a particle in a field and the integration method, as well as other accuracy parameters.

The default values correspond to the defaults set in Geant4 (taken from Geant4 9.3 release.) As Geant4 classes to not provide access methods for these defaults, the defaults have to be checked with each new Geant4 release.

Author
I. Hrivnacova; IJCLab, Orsay

Definition at line 125 of file G4FieldParameters.hh.

Constructor & Destructor Documentation

◆ G4FieldParameters()

G4FieldParameters::G4FieldParameters ( const G4String & volumeName = "")

Standard and default constructor.

Definition at line 230 of file G4FieldParameters.cc.

231 : fVolumeName(volumeName)
232{
233 // Standard constructor
234
235 fMessenger = new G4FieldParametersMessenger(this);
236}

◆ ~G4FieldParameters()

G4FieldParameters::~G4FieldParameters ( )

Destructor.

Definition at line 239 of file G4FieldParameters.cc.

240{
241 // Destructor
242
243 delete fMessenger;
244}

Member Function Documentation

◆ EquationTypeName()

G4String G4FieldParameters::EquationTypeName ( G4EquationType equation)
static

Return the equation type as a string.

Definition at line 67 of file G4FieldParameters.cc.

68{
69 // Return the equation type as a string
70
71 switch (equation) {
72 case kEqMagnetic:
73 return G4String("EqMagnetic");
75 return G4String("EqMagneticWithSpin");
77 return G4String("EqElectroMagnetic");
79 return G4String("EqEMfieldWithSpin");
81 return G4String("EqEMfieldWithEDM");
82 case kUserEquation:
83 return G4String("UserDefinedEq");
84 }
85
87 "G4FieldParameters::EquationTypeName:", "GeomFieldParameters0001",
88 FatalErrorInArgument, "Unknown equation value.");
89 return G4String();
90}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
@ kEqMagneticWithSpin
@ kUserEquation
User defined equation of motion.
@ kEqEMfieldWithSpin
@ kEqMagnetic
@ kEqElectroMagnetic
@ kEqEMfieldWithEDM

Referenced by G4FieldParametersMessenger::G4FieldParametersMessenger(), GetEquationType(), PrintParameters(), and G4FieldParametersMessenger::SetNewValue().

◆ FieldTypeName()

G4String G4FieldParameters::FieldTypeName ( G4FieldType field)
static

Return the field type as a string.

Definition at line 47 of file G4FieldParameters.cc.

48{
49 // Return the field type as a string
50
51 switch (field) {
52 case kMagnetic:
53 return G4String("Magnetic");
55 return G4String("ElectroMagnetic");
56 case kGravity:
57 return G4String("Gravity");
58 }
59
61 "G4FieldParameters::FieldTypeName:", "GeomFieldParameters0001",
62 FatalErrorInArgument, "Unknown field value.");
63 return G4String();
64}
@ kElectroMagnetic
electromagnetic field
@ kGravity
gravity field
@ kMagnetic
magnetic field

Referenced by G4FieldParametersMessenger::G4FieldParametersMessenger(), GetFieldType(), G4FieldSetup::PrintInfo(), PrintParameters(), and G4FieldParametersMessenger::SetNewValue().

◆ GetConstDistance()

G4double G4FieldParameters::GetConstDistance ( ) const
inline

Get the distance within which the field is considered constant.

Definition at line 408 of file G4FieldParameters.hh.

409{
410 return fConstDistance;
411}

◆ GetDeltaChord()

G4double G4FieldParameters::GetDeltaChord ( ) const
inline

Get delta chord in G4ChordFinder.

Definition at line 378 of file G4FieldParameters.hh.

379{
380 return fDeltaChord;
381}

◆ GetDeltaIntersection()

G4double G4FieldParameters::GetDeltaIntersection ( ) const
inline

Get delta intersection in global field manager.

Definition at line 390 of file G4FieldParameters.hh.

391{
392 return fDeltaIntersection;
393}

◆ GetDeltaOneStep()

G4double G4FieldParameters::GetDeltaOneStep ( ) const
inline

Get delta one step in global field manager.

Definition at line 384 of file G4FieldParameters.hh.

385{
386 return fDeltaOneStep;
387}

◆ GetEquationType() [1/2]

G4EquationType G4FieldParameters::GetEquationType ( ) const
inline

Get type of equation of motion of a particle in a field.

Definition at line 348 of file G4FieldParameters.hh.

349{
350 return fEquation;
351}

◆ GetEquationType() [2/2]

G4EquationType G4FieldParameters::GetEquationType ( const G4String & name)
static

Return the equation type for given equation type name.

Definition at line 172 of file G4FieldParameters.cc.

173{
174 // Return the equation type for given equation type name
175
176 if (name == EquationTypeName(kEqMagnetic)) return kEqMagnetic;
181 if (name == EquationTypeName(kUserEquation)) return kUserEquation;
182
184 "G4FieldParameters::GetEquationType:", "GeomFieldParameters0001",
185 FatalErrorInArgument, "Unknown equation name.");
186 return kEqMagnetic;
187}
static G4String EquationTypeName(G4EquationType equation)
Return the equation type as a string.

◆ GetFieldType() [1/2]

G4FieldType G4FieldParameters::GetFieldType ( ) const
inline

Get type of field.

Definition at line 345 of file G4FieldParameters.hh.

345{ return fField; }

◆ GetFieldType() [2/2]

G4FieldType G4FieldParameters::GetFieldType ( const G4String & name)
static

Return the field type for given field type name.

Definition at line 157 of file G4FieldParameters.cc.

158{
159 // Return the field type for given field type name
160
161 if (name == FieldTypeName(kMagnetic)) return kMagnetic;
163 if (name == FieldTypeName(kGravity)) return kGravity;
164
166 "G4FieldParameters::GetFieldType:", "GeomFieldParameters0001",
167 FatalErrorInArgument, "Unknown field name.");
168 return kMagnetic;
169}
static G4String FieldTypeName(G4FieldType field)
Return the field type as a string.

◆ GetMaximumEpsilonStep()

G4double G4FieldParameters::GetMaximumEpsilonStep ( ) const
inline

Get maximum epsilon step in global field manager.

Definition at line 402 of file G4FieldParameters.hh.

403{
404 return fMaximumEpsilonStep;
405}

◆ GetMinimumEpsilonStep()

G4double G4FieldParameters::GetMinimumEpsilonStep ( ) const
inline

Get minimum epsilon step in global field manager.

Definition at line 396 of file G4FieldParameters.hh.

397{
398 return fMinimumEpsilonStep;
399}

◆ GetMinimumStep()

G4double G4FieldParameters::GetMinimumStep ( ) const
inline

Get minimum step in G4ChordFinder.

Definition at line 372 of file G4FieldParameters.hh.

373{
374 return fMinimumStep;
375}

◆ GetStepperType() [1/2]

G4StepperType G4FieldParameters::GetStepperType ( ) const
inline

Get rype of integrator of particle's equation of motion.

Definition at line 354 of file G4FieldParameters.hh.

355{
356 return fStepper;
357}

◆ GetStepperType() [2/2]

G4StepperType G4FieldParameters::GetStepperType ( const G4String & name)
static

Return the stepper type for given stepper type name.

Definition at line 190 of file G4FieldParameters.cc.

191{
192 // Return the stepper type for given stepper type name
195 if (name == StepperTypeName(kCashKarpRKF45)) return kCashKarpRKF45;
196 if (name == StepperTypeName(kClassicalRK4)) return kClassicalRK4;
200 if (name == StepperTypeName(kExplicitEuler)) return kExplicitEuler;
201 if (name == StepperTypeName(kImplicitEuler)) return kImplicitEuler;
202 if (name == StepperTypeName(kSimpleHeum)) return kSimpleHeum;
203 if (name == StepperTypeName(kSimpleRunge)) return kSimpleRunge;
204 if (name == StepperTypeName(kConstRK4)) return kConstRK4;
207 if (name == StepperTypeName(kHelixHeum)) return kHelixHeum;
211 if (name == StepperTypeName(kNystromRK4)) return kNystromRK4;
212 if (name == StepperTypeName(kRKG3Stepper)) return kRKG3Stepper;
213 if (name == StepperTypeName(kRK547FEq1)) return kRK547FEq1;
214 if (name == StepperTypeName(kRK547FEq2)) return kRK547FEq2;
215 if (name == StepperTypeName(kRK547FEq3)) return kRK547FEq3;
216 if (name == StepperTypeName(kTsitourasRK45)) return kTsitourasRK45;
217 if (name == StepperTypeName(kUserStepper)) return kUserStepper;
218
220 "G4FieldParameters::GetStepperType:", "GeomFieldParameters0001",
221 FatalErrorInArgument, "Unknown stepper name.");
222 return kClassicalRK4;
223}
@ 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.
static G4String StepperTypeName(G4StepperType stepper)
Return the stepper type as a string.

◆ GetUserEquationOfMotion()

G4EquationOfMotion * G4FieldParameters::GetUserEquationOfMotion ( ) const
inline

Get user defined equation of motion.

Definition at line 360 of file G4FieldParameters.hh.

361{
362 return fUserEquation;
363}

◆ GetUserStepper()

G4MagIntegratorStepper * G4FieldParameters::GetUserStepper ( ) const
inline

Get user defined integrator of particle's equation of motion.

Definition at line 366 of file G4FieldParameters.hh.

367{
368 return fUserStepper;
369}

◆ GetVolumeName()

G4String G4FieldParameters::GetVolumeName ( ) const
inline

Definition at line 339 of file G4FieldParameters.hh.

340{
341 return fVolumeName;
342}

Referenced by G4FieldParametersMessenger::G4FieldParametersMessenger(), and G4FieldBuilder::GetFieldParameters().

◆ PrintParameters()

void G4FieldParameters::PrintParameters ( ) const

Prints all customizable accuracy parameters.

Definition at line 251 of file G4FieldParameters.cc.

252{
253 // Prints all customizable accuracy parameters
254
255 G4cout << "Magnetic field parameters: " << G4endl;
256 if (fVolumeName.size()) {
257 G4cout << " volume name = " << fVolumeName << G4endl;
258 }
259 G4cout << " field type = " << FieldTypeName(fField) << G4endl
260 << " equation type = " << EquationTypeName(fEquation) << G4endl
261 << " stepper type = " << StepperTypeName(fStepper) << G4endl
262 << " minStep = " << fMinimumStep << " mm" << G4endl
263 << " constDistance = " << fConstDistance << " mm" << G4endl
264 << " deltaChord = " << fDeltaChord << " mm" << G4endl
265 << " deltaOneStep = " << fDeltaOneStep << " mm" << G4endl
266 << " deltaIntersection = " << fDeltaIntersection << " mm" << G4endl
267 << " epsMin = " << fMinimumEpsilonStep << G4endl
268 << " epsMax= " << fMaximumEpsilonStep << G4endl;
269}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout

◆ SetConstDistance()

void G4FieldParameters::SetConstDistance ( G4double value)
inline

Set the distance within which the field is considered constant.

Definition at line 333 of file G4FieldParameters.hh.

334{
335 fConstDistance = value;
336}

◆ SetDeltaChord()

void G4FieldParameters::SetDeltaChord ( G4double value)
inline

Set delta chord in G4ChordFinder.

Definition at line 303 of file G4FieldParameters.hh.

304{
305 fDeltaChord = value;
306}

◆ SetDeltaIntersection()

void G4FieldParameters::SetDeltaIntersection ( G4double value)
inline

Set delta intersection in global field manager.

Definition at line 315 of file G4FieldParameters.hh.

316{
317 fDeltaIntersection = value;
318}

◆ SetDeltaOneStep()

void G4FieldParameters::SetDeltaOneStep ( G4double value)
inline

Set delta one step in global field manager.

Definition at line 309 of file G4FieldParameters.hh.

310{
311 fDeltaOneStep = value;
312}

◆ SetEquationType()

void G4FieldParameters::SetEquationType ( G4EquationType equation)
inline

Set Type of equation of motion of a particle in a field.

Definition at line 285 of file G4FieldParameters.hh.

286{
287 fEquation = equation;
288}

◆ SetFieldType()

void G4FieldParameters::SetFieldType ( G4FieldType field)
inline

Set type of field.

Definition at line 279 of file G4FieldParameters.hh.

280{
281 fField = field;
282}

◆ SetMaximumEpsilonStep()

void G4FieldParameters::SetMaximumEpsilonStep ( G4double value)
inline

Set maximum epsilon step in global field manager.

Definition at line 327 of file G4FieldParameters.hh.

328{
329 fMaximumEpsilonStep = value;
330}

◆ SetMinimumEpsilonStep()

void G4FieldParameters::SetMinimumEpsilonStep ( G4double value)
inline

Set minimum epsilon step in global field manager.

Definition at line 321 of file G4FieldParameters.hh.

322{
323 fMinimumEpsilonStep = value;
324}

◆ SetMinimumStep()

void G4FieldParameters::SetMinimumStep ( G4double value)
inline

Set minimum step in G4ChordFinder.

Definition at line 297 of file G4FieldParameters.hh.

298{
299 fMinimumStep = value;
300}

◆ SetStepperType()

void G4FieldParameters::SetStepperType ( G4StepperType stepper)
inline

Type of integrator of particle's equation of motion.

Definition at line 291 of file G4FieldParameters.hh.

292{
293 fStepper = stepper;
294}

◆ SetUserEquationOfMotion()

void G4FieldParameters::SetUserEquationOfMotion ( G4EquationOfMotion * equation)

Set user defined equation of motion.

Definition at line 272 of file G4FieldParameters.cc.

273{
274 // Set user defined equation of motion
275
276 fUserEquation = equation;
277 fEquation = kUserEquation;
278}

Referenced by G4FieldBuilder::SetUserEquationOfMotion().

◆ SetUserStepper()

void G4FieldParameters::SetUserStepper ( G4MagIntegratorStepper * stepper)

Set user defined integrator of particle's equation of motion.

Definition at line 281 of file G4FieldParameters.cc.

282{
283 // Set user defined integrator of particle's equation of motion
284
285 fUserStepper = stepper;
286 fStepper = kUserStepper;
287}

Referenced by G4FieldBuilder::SetUserStepper().

◆ StepperTypeName()

G4String G4FieldParameters::StepperTypeName ( G4StepperType stepper)
static

Return the stepper type as a string.

Definition at line 93 of file G4FieldParameters.cc.

94{
95 // Return the stepper type as a string
96
97 switch (stepper) {
99 return G4String("BogackiShampine23");
101 return G4String("BogackiShampine45");
102 case kCashKarpRKF45:
103 return G4String("CashKarpRKF45");
104 case kClassicalRK4:
105 return G4String("ClassicalRK4");
107 return G4String("DormandPrince745");
109 return G4String("DormandPrinceRK56");
111 return G4String("DormandPrinceRK78");
112 case kExplicitEuler:
113 return G4String("ExplicitEuler");
114 case kImplicitEuler:
115 return G4String("ImplicitEuler");
116 case kSimpleHeum:
117 return G4String("SimpleHeum");
118 case kSimpleRunge:
119 return G4String("SimpleRunge");
120 case kConstRK4:
121 return G4String("ConstRK4");
123 return G4String("ExactHelixStepper");
125 return G4String("HelixExplicitEuler");
126 case kHelixHeum:
127 return G4String("HelixHeum");
129 return G4String("HelixImplicitEuler");
131 return G4String("HelixMixedStepper");
133 return G4String("HelixSimpleRunge");
134 case kNystromRK4:
135 return G4String("NystromRK4");
136 case kRKG3Stepper:
137 return G4String("RKG3_Stepper");
138 case kTsitourasRK45:
139 return G4String("TsitourasRK45");
140 case kUserStepper:
141 return G4String("UserDefinedStepper");
142 case kRK547FEq1:
143 return G4String("RK547FEq1");
144 case kRK547FEq2:
145 return G4String("RK547FEq2");
146 case kRK547FEq3:
147 return G4String("RK547FEq3");
148 }
149
151 "G4FieldParameters::StepperTypeName:", "GeomFieldParameters0001",
152 FatalErrorInArgument, "Unknown stepper value.");
153 return G4String();
154}

Referenced by G4FieldParametersMessenger::G4FieldParametersMessenger(), GetStepperType(), G4FieldSetup::PrintInfo(), PrintParameters(), and G4FieldParametersMessenger::SetNewValue().


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