Geant4 VMC Version 6.6
Loading...
Searching...
No Matches
TG4FieldParameters.cxx
Go to the documentation of this file.
1//------------------------------------------------
2// The Geant4 Virtual Monte Carlo package
3// Copyright (C) 2007 - 2014 Ivana Hrivnacova
4// All rights reserved.
5//
6// For the licensing terms see geant4_vmc/LICENSE.
7// Contact: root-vmc@cern.ch
8//-------------------------------------------------
9
14
15#include "TG4FieldParameters.h"
17#include "TG4G3Units.h"
18#include "TG4Globals.h"
19
20#include <TVirtualMC.h>
21#include <TVirtualMCApplication.h>
22
23// Moved after Root includes to avoid shadowed variables
24// generated from short units names
25#include <G4SystemOfUnits.hh>
26
27const G4double TG4FieldParameters::fgkDefaultStepMinimum = 0.01 * mm;
28const G4double TG4FieldParameters::fgkDefaultDeltaChord = 0.25 * mm;
29const G4double TG4FieldParameters::fgkDefaultDeltaOneStep = 0.01 * mm;
30const G4double TG4FieldParameters::fgkDefaultDeltaIntersection = 0.001 * mm;
34
35//
36// static methods
37//
38
39//_____________________________________________________________________________
41{
43
44 switch (field) {
45 case kMagnetic:
46 return G4String("Magnetic");
48 return G4String("ElectroMagnetic");
49 case kGravity:
50 return G4String("Gravity");
51 }
52
54 "TG4FieldParameters", "FieldTypeName:", "Unknown field value.");
55 return G4String();
56}
57
58//_____________________________________________________________________________
60{
62
63 switch (equation) {
64 case kMagUsualEqRhs:
65 return G4String("MagUsualEqRhs");
66 case kMagSpinEqRhs:
67 return G4String("MagSpinEqRhs");
68 case kEqMagElectric:
69 return G4String("EqMagElectric");
71 return G4String("EqEMFieldWithSpin");
73 return G4String("EqEMFieldWithEDM");
74 case kUserEquation:
75 return G4String("UserDefinedEq");
76 }
77
79 "TG4FieldParameters", "EquationTypeName:", "Unknown equation value.");
80 return G4String();
81}
82
83//_____________________________________________________________________________
85{
87
88 switch (stepper) {
90 return G4String("BogackiShampine23");
92 return G4String("BogackiShampine45");
93 case kCashKarpRKF45:
94 return G4String("CashKarpRKF45");
95 case kClassicalRK4:
96 return G4String("ClassicalRK4");
98 return G4String("DormandPrince745");
100 return G4String("DormandPrinceRK56");
102 return G4String("DormandPrinceRK78");
103 case kExplicitEuler:
104 return G4String("ExplicitEuler");
105 case kImplicitEuler:
106 return G4String("ImplicitEuler");
107 case kSimpleHeum:
108 return G4String("SimpleHeum");
109 case kSimpleRunge:
110 return G4String("SimpleRunge");
111 case kConstRK4:
112 return G4String("ConstRK4");
114 return G4String("ExactHelixStepper");
116 return G4String("HelixExplicitEuler");
117 case kHelixHeum:
118 return G4String("HelixHeum");
120 return G4String("HelixImplicitEuler");
122 return G4String("HelixMixedStepper");
124 return G4String("HelixSimpleRunge");
125 case kNystromRK4:
126 return G4String("NystromRK4");
127 case kRKG3Stepper:
128 return G4String("RKG3_Stepper");
129 case kTsitourasRK45:
130 return G4String("TsitourasRK45");
131 case kUserStepper:
132 return G4String("UserDefinedStepper");
133 case kRK547FEq1:
134 return G4String("RK547FEq1");
135 case kRK547FEq2:
136 return G4String("RK547FEq2");
137 case kRK547FEq3:
138 return G4String("RK547FEq3");
139 }
140
142 "TG4FieldParameters", "StepperTypeName:", "Unknown stepper value.");
143 return G4String();
144}
145
146//_____________________________________________________________________________
148{
150
151 if (name == FieldTypeName(kMagnetic)) return kMagnetic;
153 if (name == FieldTypeName(kGravity)) return kGravity;
154
156 "TG4FieldParameters", "GetFieldType:", "Unknown field name.");
157 return kMagnetic;
158}
159
160//_____________________________________________________________________________
162{
164
166 if (name == EquationTypeName(kMagSpinEqRhs)) return kMagSpinEqRhs;
170 if (name == EquationTypeName(kUserEquation)) return kUserEquation;
171
173 "TG4FieldParameters", "GetEquationType:", "Unknown equation name.");
174 return kMagUsualEqRhs;
175}
176
177//_____________________________________________________________________________
179{
183 if (name == StepperTypeName(kCashKarpRKF45)) return kCashKarpRKF45;
184 if (name == StepperTypeName(kClassicalRK4)) return kClassicalRK4;
188 if (name == StepperTypeName(kExplicitEuler)) return kExplicitEuler;
189 if (name == StepperTypeName(kImplicitEuler)) return kImplicitEuler;
190 if (name == StepperTypeName(kSimpleHeum)) return kSimpleHeum;
191 if (name == StepperTypeName(kSimpleRunge)) return kSimpleRunge;
192 if (name == StepperTypeName(kConstRK4)) return kConstRK4;
195 if (name == StepperTypeName(kHelixHeum)) return kHelixHeum;
199 if (name == StepperTypeName(kNystromRK4)) return kNystromRK4;
200 if (name == StepperTypeName(kRKG3Stepper)) return kRKG3Stepper;
201 if (name == StepperTypeName(kRK547FEq1)) return kRK547FEq1;
202 if (name == StepperTypeName(kRK547FEq2)) return kRK547FEq2;
203 if (name == StepperTypeName(kRK547FEq3)) return kRK547FEq3;
204 if (name == StepperTypeName(kTsitourasRK45)) return kTsitourasRK45;
205 if (name == StepperTypeName(kUserStepper)) return kUserStepper;
206
208 "TG4FieldParameters", "GetStepperType:", "Unknown stepper name.");
209 return kClassicalRK4;
210}
211
212//
213// ctors, dtor
214//
215
216//_____________________________________________________________________________
217TG4FieldParameters::TG4FieldParameters(const G4String& volumeName)
218 : fMessenger(0),
219 fVolumeName(volumeName),
220 fStepMinimum(fgkDefaultStepMinimum),
221 fDeltaChord(fgkDefaultDeltaChord),
222 fDeltaOneStep(fgkDefaultDeltaOneStep),
223 fDeltaIntersection(fgkDefaultDeltaIntersection),
224 fMinimumEpsilonStep(fgkDefaultMinimumEpsilonStep),
225 fMaximumEpsilonStep(fgkDefaultMaximumEpsilonStep),
226 fField(kMagnetic),
227 fEquation(kMagUsualEqRhs),
228 fStepper(kDormandPrince745),
229 fUserEquation(0),
230 fUserStepper(0),
231 fConstDistance(0),
232 fIsMonopole(false)
233{
235
237}
238
239//_____________________________________________________________________________
246
247//
248// public methods
249//
250
251//_____________________________________________________________________________
253{
255
256 G4cout << "Magnetic field parameters: " << G4endl;
257 if (fVolumeName.size()) {
258 G4cout << " volume name = " << fVolumeName << G4endl;
259 }
260 G4cout << " field type = " << FieldTypeName(fField) << G4endl
261 << " equation type = " << EquationTypeName(fEquation) << G4endl
262 << " stepper type = " << StepperTypeName(fStepper) << G4endl
263 << " minStep = " << fStepMinimum << " mm" << G4endl
264 << " constDistance = " << fConstDistance << " mm" << G4endl
265 << " isMonopole = " << std::boolalpha << fIsMonopole << G4endl
266 << " deltaChord = " << fDeltaChord << " mm" << G4endl
267 << " deltaOneStep = " << fDeltaOneStep << " mm" << G4endl
268 << " deltaIntersection = " << fDeltaIntersection << " mm" << G4endl
269 << " epsMin = " << fMinimumEpsilonStep << G4endl
270 << " epsMax= " << fMaximumEpsilonStep << G4endl;
271}
272
273//_____________________________________________________________________________
281
282//_____________________________________________________________________________
283void TG4FieldParameters::SetUserStepper(G4MagIntegratorStepper* stepper)
284{
286
287 fUserStepper = stepper;
289}
Definition of the TG4FieldParametersMessenger class.
Definition of the TG4FieldParameters class.
@ 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
Definition of the TG4G3Units class.
Definition of the TG4Globals class and basic container types.
Messenger class that defines commands for TG4DetConstruction.
G4double fStepMinimum
Minimum step in G4ChordFinder.
G4String fVolumeName
The name of associated volume, if local field.
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)
G4EquationOfMotion * fUserEquation
User defined equation of motion.
G4double fDeltaChord
Delta chord in G4ChordFinder.
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.
G4double fDeltaOneStep
Delta one step in global field manager.
TG4FieldParameters(const G4String &volumeName="")
static const G4double fgkDefaultStepMinimum
Default minimum step in G4ChordFinder.
StepperType fStepper
Type of integrator of particle's equation of motion.
static G4String StepperTypeName(StepperType stepper)
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.
static G4String FieldTypeName(FieldType field)
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.
static void Exception(const TString &className, const TString &methodName, const TString &text)