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

G4FieldParameters 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. More...

#include <G4FieldParameters.hh>

Public Member Functions

 G4FieldParameters (const G4String &volumeName="")
 ~G4FieldParameters ()
 G4FieldParameters (const G4FieldParameters &right)=delete
G4FieldParametersoperator= (const G4FieldParameters &right)=delete
void PrintParameters () const
void SetFieldType (G4FieldType field)
void SetEquationType (G4EquationType equation)
void SetStepperType (G4StepperType stepper)
void SetUserEquationOfMotion (G4EquationOfMotion *equation)
void SetUserStepper (G4MagIntegratorStepper *stepper)
void SetMinimumStep (G4double value)
void SetDeltaChord (G4double value)
void SetDeltaOneStep (G4double value)
void SetDeltaIntersection (G4double value)
void SetMinimumEpsilonStep (G4double value)
void SetMaximumEpsilonStep (G4double value)
void SetConstDistance (G4double value)
const G4StringGetVolumeName () const
const G4FieldTypeGetFieldType () const
const G4EquationTypeGetEquationType () const
const G4StepperTypeGetStepperType () const
G4EquationOfMotionGetUserEquationOfMotion () const
G4MagIntegratorStepperGetUserStepper () const
G4double GetMinimumStep () const
G4double GetDeltaChord () const
G4double GetDeltaOneStep () const
G4double GetDeltaIntersection () const
G4double GetMinimumEpsilonStep () const
G4double GetMaximumEpsilonStep () const
G4double GetConstDistance () const

Static Public Member Functions

static G4String FieldTypeName (G4FieldType field)
static G4String EquationTypeName (G4EquationType equation)
static G4String StepperTypeName (G4StepperType stepper)
static G4FieldType GetFieldType (const G4String &name)
static G4EquationType GetEquationType (const G4String &name)
static G4StepperType GetStepperType (const G4String &name)

Detailed Description

G4FieldParameters 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.

Definition at line 167 of file G4FieldParameters.hh.

Constructor & Destructor Documentation

◆ G4FieldParameters() [1/2]

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

Constructor for G4FieldParameters.

Parameters
[in]volumeNameThe volume name where field is applied.

Definition at line 255 of file G4FieldParameters.cc.

256 : fVolumeName(volumeName)
257{
258 // Standard constructor
259
260 fMessenger = new G4FieldParametersMessenger(this);
261}

Referenced by G4FieldParameters(), and operator=().

◆ ~G4FieldParameters()

G4FieldParameters::~G4FieldParameters ( )

Destructor.

Definition at line 264 of file G4FieldParameters.cc.

265{
266 // Destructor
267
268 delete fMessenger;
269}

◆ G4FieldParameters() [2/2]

G4FieldParameters::G4FieldParameters ( const G4FieldParameters & right)
delete

Copy constructor and assignment operator not allowed.

Member Function Documentation

◆ EquationTypeName()

G4String G4FieldParameters::EquationTypeName ( G4EquationType equation)
static

Returns the equation type as a string.

Definition at line 65 of file G4FieldParameters.cc.

66{
67 // Return the equation type as a string
68
69 switch (equation)
70 {
71 case kEqMagnetic:
72 return {"EqMagnetic"};
74 return {"EqMagneticWithSpin"};
76 return {"EqElectroMagnetic"};
78 return {"EqEMfieldWithSpin"};
80 return {"EqEMfieldWithEDM"};
81 case kEqGravity:
82 return {"EqGravity"};
83 case kEqMonopole:
84 return {"EqMonopole"};
85 case kEqReplate:
86 return {"EqReplate"};
87 case kUserEquation:
88 return {"UserDefinedEq"};
89 }
90
92 "G4FieldParameters::EquationTypeName:", "GeomFieldParameters0001",
93 FatalErrorInArgument, "Unknown equation value.");
94 return {};
95}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
@ kEqMagneticWithSpin
@ kUserEquation
User defined equation of motion.
@ kEqEMfieldWithSpin
@ kEqMagnetic
@ kEqGravity
@ kEqReplate
@ kEqMonopole
@ kEqElectroMagnetic
@ kEqEMfieldWithEDM

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

◆ FieldTypeName()

G4String G4FieldParameters::FieldTypeName ( G4FieldType field)
static

Returns the field type as a string.

Definition at line 42 of file G4FieldParameters.cc.

43{
44 // Return the field type as a string
45
46 switch (field)
47 {
48 case kMagnetic:
49 return {"Magnetic"};
51 return {"ElectroMagnetic"};
52 case kGravity:
53 return {"Gravity"};
54 case kUserFieldType:
55 return {"UserDefinedField"};
56 }
57
59 "G4FieldParameters::FieldTypeName:", "GeomFieldParameters0001",
60 FatalErrorInArgument, "Unknown field value.");
61 return {};
62}
@ kElectroMagnetic
electromagnetic field
@ kUserFieldType
User defined field type.
@ kGravity
gravity field
@ kMagnetic
magnetic field

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

◆ GetConstDistance()

G4double G4FieldParameters::GetConstDistance ( ) const

Gets the distance within which the field is considered constant.

◆ GetDeltaChord()

G4double G4FieldParameters::GetDeltaChord ( ) const

Gets the delta chord in G4ChordFinder.

◆ GetDeltaIntersection()

G4double G4FieldParameters::GetDeltaIntersection ( ) const

Gets the delta intersection in global field manager.

◆ GetDeltaOneStep()

G4double G4FieldParameters::GetDeltaOneStep ( ) const

Gets the delta one step in global field manager.

◆ GetEquationType() [1/2]

const G4EquationType & G4FieldParameters::GetEquationType ( ) const

Gets the type of equation of motion of a particle in a field.

◆ GetEquationType() [2/2]

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

Returns the equation type for given equation type name.

Definition at line 189 of file G4FieldParameters.cc.

190{
191 // Return the equation type for given equation type name
192
193 if (name == EquationTypeName(kEqMagnetic)) { return kEqMagnetic; }
198 if (name == EquationTypeName(kEqGravity)) { return kEqGravity; }
199 if (name == EquationTypeName(kEqMonopole)) { return kEqMonopole; }
200 if (name == EquationTypeName(kEqReplate)) { return kEqReplate; }
201 if (name == EquationTypeName(kUserEquation)) { return kUserEquation; }
202
204 "G4FieldParameters::GetEquationType:", "GeomFieldParameters0001",
205 FatalErrorInArgument, "Unknown equation name.");
206 return kEqMagnetic;
207}
static G4String EquationTypeName(G4EquationType equation)

◆ GetFieldType() [1/2]

const G4FieldType & G4FieldParameters::GetFieldType ( ) const

Gets the type of field.

◆ GetFieldType() [2/2]

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

Returns the field type for given field type name.

Definition at line 173 of file G4FieldParameters.cc.

174{
175 // Return the field type for given field type name
176
177 if (name == FieldTypeName(kMagnetic)) { return kMagnetic; }
178 if (name == FieldTypeName(kElectroMagnetic)) { return kElectroMagnetic; }
179 if (name == FieldTypeName(kGravity)) { return kGravity; }
180 if (name == FieldTypeName(kUserFieldType)) { return kUserFieldType; }
181
183 "G4FieldParameters::GetFieldType:", "GeomFieldParameters0001",
184 FatalErrorInArgument, "Unknown field name.");
185 return kMagnetic;
186}
static G4String FieldTypeName(G4FieldType field)

◆ GetMaximumEpsilonStep()

G4double G4FieldParameters::GetMaximumEpsilonStep ( ) const

Gets the maximum epsilon step in global field manager.

◆ GetMinimumEpsilonStep()

G4double G4FieldParameters::GetMinimumEpsilonStep ( ) const

Gets the minimum epsilon step in global field manager.

◆ GetMinimumStep()

G4double G4FieldParameters::GetMinimumStep ( ) const

Gets the minimum step in G4ChordFinder.

◆ GetStepperType() [1/2]

const G4StepperType & G4FieldParameters::GetStepperType ( ) const

Gets the type of integrator of particle's equation of motion.

◆ GetStepperType() [2/2]

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

Returns the stepper type for given stepper type name.

Definition at line 210 of file G4FieldParameters.cc.

211{
212 // Return the stepper type for given stepper type name
215 if (name == StepperTypeName(kCashKarpRKF45)) { return kCashKarpRKF45; }
216 if (name == StepperTypeName(kClassicalRK4)) { return kClassicalRK4; }
217 if (name == StepperTypeName(kDoLoMcPriRK34)) { return kDoLoMcPriRK34; }
218 if (name == StepperTypeName(kDormandPrince745)) { return kDormandPrince745; }
221 if (name == StepperTypeName(kExplicitEuler)) { return kExplicitEuler; }
222 if (name == StepperTypeName(kImplicitEuler)) { return kImplicitEuler; }
223 if (name == StepperTypeName(kSimpleHeum)) { return kSimpleHeum; }
224 if (name == StepperTypeName(kSimpleRunge)) { return kSimpleRunge; }
225 if (name == StepperTypeName(kConstRK4)) { return kConstRK4; }
228 if (name == StepperTypeName(kHelixHeum)) { return kHelixHeum; }
231 if (name == StepperTypeName(kHelixSimpleRunge)) { return kHelixSimpleRunge; }
232 if (name == StepperTypeName(kNystromRK4)) { return kNystromRK4; }
233 if (name == StepperTypeName(kQSStepper)) { return kQSStepper; }
234 if (name == StepperTypeName(kRKG3Stepper)) { return kRKG3Stepper; }
235 if (name == StepperTypeName(kRK547FEq1)) { return kRK547FEq1; }
236 if (name == StepperTypeName(kRK547FEq2)) { return kRK547FEq2; }
237 if (name == StepperTypeName(kRK547FEq3)) { return kRK547FEq3; }
238 if (name == StepperTypeName(kTsitourasRK45)) { return kTsitourasRK45; }
239 if (name == StepperTypeName(kTCashKarpRKF45)) { return kTCashKarpRKF45; }
240 if (name == StepperTypeName(kTDormandPrince45)) { return kTDormandPrince45; }
241 if (name == StepperTypeName(kTMagErrorStepper)) { return kTMagErrorStepper; }
242 if (name == StepperTypeName(kUserStepper)) { return kUserStepper; }
243
245 "G4FieldParameters::GetStepperType:", "GeomFieldParameters0001",
246 FatalErrorInArgument, "Unknown stepper name.");
247 return kClassicalRK4;
248}
@ kRKG3Stepper
G4RKG3_Stepper.
@ kRK547FEq2
G4RK547FEq2.
@ kHelixSimpleRunge
G4HelixSimpleRunge.
@ kNystromRK4
G4NystromRK4.
@ kDormandPrince745
G4DormandPrince745.
@ kCashKarpRKF45
G4CashKarpRKF45.
@ kDormandPrinceRK78
G4DormandPrinceRK78.
@ kQSStepper
G4QSStepper.
@ kSimpleRunge
G4SimpleRunge.
@ kHelixImplicitEuler
G4HelixImplicitEuler.
@ kConstRK4
G4ConstRK4.
@ kUserStepper
User defined stepper.
@ kDoLoMcPriRK34
G4DoLoMcPriRK34.
@ kSimpleHeum
G4SimpleHeum.
@ kHelixHeum
G4HelixHeum.
@ kHelixExplicitEuler
G4HelixExplicitEuler.
@ kDormandPrinceRK56
G4DormandPrinceRK56.
@ kTDormandPrince45
G4TDormandPrince45.
@ kTsitourasRK45
G4TsitourasRK45.
@ kImplicitEuler
G4ImplicitEuler.
@ kExactHelixStepper
G4ExactHelixStepper.
@ kHelixMixedStepper
G4HelixMixedStepper.
@ kBogackiShampine45
G4BogackiShampine45.
@ kExplicitEuler
G4ExplicitEuler.
@ kTCashKarpRKF45
G4TCashKarpRKF45.
@ kRK547FEq1
G4RK547FEq1.
@ kRK547FEq3
G4RK547FEq3.
@ kBogackiShampine23
G4BogackiShampine23.
@ kClassicalRK4
G4ClassicalRK4.
@ kTMagErrorStepper
G4TMagErrorStepper.
static G4String StepperTypeName(G4StepperType stepper)

◆ GetUserEquationOfMotion()

G4EquationOfMotion * G4FieldParameters::GetUserEquationOfMotion ( ) const

Gets the user defined equation of motion.

◆ GetUserStepper()

G4MagIntegratorStepper * G4FieldParameters::GetUserStepper ( ) const

Gets the user defined integrator of particle's equation of motion.

◆ GetVolumeName()

const G4String & G4FieldParameters::GetVolumeName ( ) const

Gets the name of associated volume, if local field.

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

◆ operator=()

G4FieldParameters & G4FieldParameters::operator= ( const G4FieldParameters & right)
delete

◆ PrintParameters()

void G4FieldParameters::PrintParameters ( ) const

Prints all customisable accuracy parameters.

Definition at line 276 of file G4FieldParameters.cc.

277{
278 // Prints all customizable accuracy parameters
279
280 G4cout << "Magnetic field parameters: " << G4endl;
281 if (!fVolumeName.empty())
282 {
283 G4cout << " volume name = " << fVolumeName << G4endl;
284 }
285 G4cout << " field type = " << FieldTypeName(fField) << G4endl
286 << " equation type = " << EquationTypeName(fEquation) << G4endl
287 << " stepper type = " << StepperTypeName(fStepper) << G4endl
288 << " minStep = " << fMinimumStep << " mm" << G4endl
289 << " constDistance = " << fConstDistance << " mm" << G4endl
290 << " deltaChord = " << fDeltaChord << " mm" << G4endl
291 << " deltaOneStep = " << fDeltaOneStep << " mm" << G4endl
292 << " deltaIntersection = " << fDeltaIntersection << " mm" << G4endl
293 << " epsMin = " << fMinimumEpsilonStep << G4endl
294 << " epsMax= " << fMaximumEpsilonStep << G4endl;
295}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout

◆ SetConstDistance()

void G4FieldParameters::SetConstDistance ( G4double value)

Sets the distance within which the field is considered constant.

◆ SetDeltaChord()

void G4FieldParameters::SetDeltaChord ( G4double value)

Sets the delta chord in G4ChordFinder.

◆ SetDeltaIntersection()

void G4FieldParameters::SetDeltaIntersection ( G4double value)

Sets the delta intersection in global field manager.

◆ SetDeltaOneStep()

void G4FieldParameters::SetDeltaOneStep ( G4double value)

Sets the delta one step in global field manager.

◆ SetEquationType()

void G4FieldParameters::SetEquationType ( G4EquationType equation)

Sets the type of equation of motion of a particle in a field.

◆ SetFieldType()

void G4FieldParameters::SetFieldType ( G4FieldType field)

Sets the type of field.

◆ SetMaximumEpsilonStep()

void G4FieldParameters::SetMaximumEpsilonStep ( G4double value)

Sets the maximum epsilon step in global field manager.

◆ SetMinimumEpsilonStep()

void G4FieldParameters::SetMinimumEpsilonStep ( G4double value)

Sets the minimum epsilon step in global field manager.

◆ SetMinimumStep()

void G4FieldParameters::SetMinimumStep ( G4double value)

Sets the minimum step in G4ChordFinder.

◆ SetStepperType()

void G4FieldParameters::SetStepperType ( G4StepperType stepper)

Sets the type of integrator of particle's equation of motion.

◆ SetUserEquationOfMotion()

void G4FieldParameters::SetUserEquationOfMotion ( G4EquationOfMotion * equation)

Sets the user defined equation of motion.

Definition at line 298 of file G4FieldParameters.cc.

299{
300 // Set user defined equation of motion
301
302 fUserEquation = equation;
303 fEquation = kUserEquation;
304}

Referenced by G4FieldBuilder::SetUserEquationOfMotion().

◆ SetUserStepper()

void G4FieldParameters::SetUserStepper ( G4MagIntegratorStepper * stepper)

Sets the user defined integrator of particle's equation of motion.

Definition at line 307 of file G4FieldParameters.cc.

308{
309 // Set user defined integrator of particle's equation of motion
310
311 fUserStepper = stepper;
312 fStepper = kUserStepper;
313}

Referenced by G4FieldBuilder::SetUserStepper().

◆ StepperTypeName()

G4String G4FieldParameters::StepperTypeName ( G4StepperType stepper)
static

Returns the stepper type as a string.

Definition at line 98 of file G4FieldParameters.cc.

99{
100 // Return the stepper type as a string
101
102 switch (stepper)
103 {
105 return {"BogackiShampine23"};
107 return {"BogackiShampine45"};
108 case kCashKarpRKF45:
109 return {"CashKarpRKF45"};
110 case kClassicalRK4:
111 return {"ClassicalRK4"};
112 case kDoLoMcPriRK34:
113 return {"DoLoMcPriRK34"};
115 return {"DormandPrince745"};
117 return {"DormandPrinceRK56"};
119 return {"DormandPrinceRK78"};
120 case kExplicitEuler:
121 return {"ExplicitEuler"};
122 case kImplicitEuler:
123 return {"ImplicitEuler"};
124 case kSimpleHeum:
125 return {"SimpleHeum"};
126 case kSimpleRunge:
127 return {"SimpleRunge"};
128 case kConstRK4:
129 return {"ConstRK4"};
131 return {"ExactHelixStepper"};
133 return {"HelixExplicitEuler"};
134 case kHelixHeum:
135 return {"HelixHeum"};
137 return {"HelixImplicitEuler"};
139 return {"HelixMixedStepper"};
141 return {"HelixSimpleRunge"};
142 case kNystromRK4:
143 return {"NystromRK4"};
144 case kRKG3Stepper:
145 return {"RKG3_Stepper"};
146 case kTsitourasRK45:
147 return {"TsitourasRK45"};
148 case kUserStepper:
149 return {"UserDefinedStepper"};
150 case kRK547FEq1:
151 return {"RK547FEq1"};
152 case kRK547FEq2:
153 return {"RK547FEq2"};
154 case kRK547FEq3:
155 return {"RK547FEq3"};
156 case kTCashKarpRKF45:
157 return {"TCashKarpRKF45"};
159 return {"TDormandPrince45"};
161 return {"TMagErrorStepper"};
162 case kQSStepper:
163 return {"QSStepper"};
164 }
165
167 "G4FieldParameters::StepperTypeName:", "GeomFieldParameters0001",
168 FatalErrorInArgument, "Unknown stepper value.");
169 return {};
170}

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


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