Geant4 VMC Version 6.6
Loading...
Searching...
No Matches
TG4FieldParameters.h
Go to the documentation of this file.
1#ifndef TG4_FIELD_PARAMETERS_H
2#define TG4_FIELD_PARAMETERS_H
3
4//-------------------------------------------------
5// The Geant4 Virtual Monte Carlo package
6// Copyright (C) 2007 - 2014 Ivana Hrivnacova
7// All rights reserved.
8//
9// For the licensing terms see geant4_vmc/LICENSE.
10// Contact: root-vmc@cern.ch
11//-------------------------------------------------
12
17
18#include <G4MagneticField.hh>
19#include <globals.hh>
20
22
24class G4MagIntegratorStepper;
25
33
54
90
105
107{
108 public:
109 TG4FieldParameters(const G4String& volumeName = "");
111
112 // methods
113 void PrintParameters() const;
114
115 static G4String FieldTypeName(FieldType field);
116 static G4String EquationTypeName(EquationType equation);
117 static G4String StepperTypeName(StepperType stepper);
118 static FieldType GetFieldType(const G4String& name);
119 static EquationType GetEquationType(const G4String& name);
120 static StepperType GetStepperType(const G4String& name);
121
122 // set methods
123 void SetFieldType(FieldType field);
124 void SetEquationType(EquationType equation);
125 void SetStepperType(StepperType stepper);
127 void SetUserStepper(G4MagIntegratorStepper* stepper);
128
129 void SetStepMinimum(G4double value);
130 void SetDeltaChord(G4double value);
131 void SetDeltaOneStep(G4double value);
132 void SetDeltaIntersection(G4double value);
133 void SetMinimumEpsilonStep(G4double value);
134 void SetMaximumEpsilonStep(G4double value);
135 void SetConstDistance(G4double value);
136 void SetIsMonopole(G4bool isMonopole);
137
138 // get methods
139 G4String GetVolumeName() const;
140
141 FieldType GetFieldType() const;
145 G4MagIntegratorStepper* GetUserStepper() const;
146
147 G4double GetStepMinimum() const;
148 G4double GetDeltaChord() const;
149 G4double GetDeltaOneStep() const;
150 G4double GetDeltaIntersection() const;
151 G4double GetMinimumEpsilonStep() const;
152 G4double GetMaximumEpsilonStep() const;
153 G4double GetConstDistance() const;
154 G4bool GetIsMonopole() const;
155
156 private:
157 // static data members
158 //
160 static const G4double fgkDefaultStepMinimum;
162 static const G4double fgkDefaultDeltaChord;
164 static const G4double fgkDefaultDeltaOneStep;
166 static const G4double fgkDefaultDeltaIntersection;
168 static const G4double fgkDefaultMinimumEpsilonStep;
170 static const G4double fgkDefaultMaximumEpsilonStep;
172 static const G4double fgkDefaultConstDistance;
173
174 // data members
175 //
178
180 G4String fVolumeName;
181
183 G4double fStepMinimum;
185 G4double fDeltaChord;
194
197
200
203
206
208 G4MagIntegratorStepper* fUserStepper;
209
212
216};
217
218// inline functions
219
222{
223 fField = field;
224}
225
228{
229 fEquation = equation;
230}
231
234{
235 fStepper = stepper;
236}
237
239inline void TG4FieldParameters::SetStepMinimum(G4double value)
240{
241 fStepMinimum = value;
242}
243
245inline void TG4FieldParameters::SetDeltaChord(G4double value)
246{
247 fDeltaChord = value;
248}
249
251inline void TG4FieldParameters::SetDeltaOneStep(G4double value)
252{
253 fDeltaOneStep = value;
254}
255
258{
259 fDeltaIntersection = value;
260}
261
264{
265 fMinimumEpsilonStep = value;
266}
267
270{
271 fMaximumEpsilonStep = value;
272}
273
275inline void TG4FieldParameters::SetConstDistance(G4double value)
276{
277 fConstDistance = value;
278}
279
282inline void TG4FieldParameters::SetIsMonopole(G4bool isMonopole)
283{
284 fIsMonopole = isMonopole;
285}
286
289{
290 return fVolumeName;
291}
292
295
298{
299 return fEquation;
300}
301
304{
305 return fStepper;
306}
307
313
315inline G4MagIntegratorStepper* TG4FieldParameters::GetUserStepper() const
316{
317 return fUserStepper;
318}
319
322{
323 return fStepMinimum;
324}
325
328{
329 return fDeltaChord;
330}
331
334{
335 return fDeltaOneStep;
336}
337
340{
341 return fDeltaIntersection;
342}
343
346{
347 return fMinimumEpsilonStep;
348}
349
352{
353 return fMaximumEpsilonStep;
354}
355
358{
359 return fConstDistance;
360}
361
364inline G4bool TG4FieldParameters::GetIsMonopole() const { return fIsMonopole; }
365
366#endif // TG4_FIELD_PARAMETERS_H
@ 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.
@ kMagUsualEqRhs
@ kUserEquation
User defined equation of motion.
@ kEqEMFieldWithSpin
@ kMagSpinEqRhs
@ kEqMagElectric
@ kEqEMFieldWithEDM
FieldType
The available fields in Geant4.
@ kElectroMagnetic
electromagnetic field
@ kGravity
gravity field
@ kMagnetic
magnetic field
Messenger class that defines commands for TG4DetConstruction.
The magnetic field parameters.
G4double GetMaximumEpsilonStep() const
Return maximum epsilon step in global field manager.
G4EquationOfMotion * GetUserEquationOfMotion() const
Return the user defined equation of motion.
G4double fStepMinimum
Minimum step in G4ChordFinder.
G4double GetMinimumEpsilonStep() const
Return minimum epsilon step in global field manager.
void SetEquationType(EquationType equation)
Set the type of equation of motion of a particle in a field.
G4double GetDeltaOneStep() const
Return delta one step in global field manager.
G4String fVolumeName
The name of associated volume, if local field.
void SetConstDistance(G4double value)
Set the distance within which the field is considered constant.
G4double GetDeltaIntersection() const
Return delta intersection in global field manager.
G4MagIntegratorStepper * fUserStepper
User defined integrator of particle's equation of motion.
static const G4double fgkDefaultDeltaChord
Default delta chord in G4ChordFinder.
static const G4double fgkDefaultMaximumEpsilonStep
Default maximum epsilon step in global field manager.
EquationType fEquation
Type of equation of motion of a particle in a field.
G4double fMaximumEpsilonStep
Maximum epsilon step in global field manager.
void SetUserEquationOfMotion(G4EquationOfMotion *equation)
void SetStepperType(StepperType stepper)
Set the type of integrator of particle's equation of motion.
G4EquationOfMotion * fUserEquation
User defined equation of motion.
G4double fDeltaChord
Delta chord in G4ChordFinder.
void SetIsMonopole(G4bool isMonopole)
G4double GetStepMinimum() const
Return minimum step in G4ChordFinder.
G4String GetVolumeName() const
Return the name of associated volume, if local field.
void SetMaximumEpsilonStep(G4double value)
Set maximum epsilon step in global field manager.
FieldType fField
Type of field.
static const G4double fgkDefaultDeltaIntersection
Delta intersection in global field manager.
static G4String EquationTypeName(EquationType equation)
G4double fMinimumEpsilonStep
Minimum epsilon step in global field manager.
void SetUserStepper(G4MagIntegratorStepper *stepper)
EquationType GetEquationType() const
Return the type of equation of motion of a particle in a field.
FieldType GetFieldType() const
Return the type of field.
void SetStepMinimum(G4double value)
Set minimum step in G4ChordFinder.
G4double fDeltaOneStep
Delta one step in global field manager.
TG4FieldParameters(const G4String &volumeName="")
static const G4double fgkDefaultStepMinimum
Default minimum step in G4ChordFinder.
void SetDeltaChord(G4double value)
Set delta chord in G4ChordFinder.
G4double GetDeltaChord() const
Return delta chord in G4ChordFinder.
G4double GetConstDistance() const
Return the distance within which the field is considered constant.
StepperType fStepper
Type of integrator of particle's equation of motion.
G4bool GetIsMonopole() const
static G4String StepperTypeName(StepperType stepper)
void SetDeltaIntersection(G4double value)
Set delta intersection in global field manager.
void SetFieldType(FieldType field)
Set type of field.
TG4FieldParametersMessenger * fMessenger
Messenger for this class.
StepperType GetStepperType() const
Return the type of integrator of particle's equation of motion.
static const G4double fgkDefaultConstDistance
Default constant distance.
G4double fConstDistance
The distance within which the field is considered constant.
void SetDeltaOneStep(G4double value)
Set delta one step in global field manager.
void SetMinimumEpsilonStep(G4double value)
Set minimum epsilon step in global field manager.
static G4String FieldTypeName(FieldType field)
G4MagIntegratorStepper * GetUserStepper() const
Return the user defined integrator of particle's equation of motion.
static const G4double fgkDefaultDeltaOneStep
Default delta one step in global field manager.
static const G4double fgkDefaultMinimumEpsilonStep
Default minimum epsilon step in global field manager.
G4double fDeltaIntersection
Delta intersection in global field manager.