File indexing completed on 2025-08-02 08:28:36
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035 #ifndef G4FIELDPARAMETERS_HH
0036 #define G4FIELDPARAMETERS_HH
0037
0038 #include "G4MagneticField.hh"
0039 #include "globals.hh"
0040
0041 #include <CLHEP/Units/SystemOfUnits.h>
0042
0043 class G4FieldParametersMessenger;
0044
0045 class G4EquationOfMotion;
0046 class G4MagIntegratorStepper;
0047
0048
0049 enum G4FieldType
0050 {
0051 kMagnetic,
0052 kElectroMagnetic,
0053 kGravity
0054 };
0055
0056
0057
0058 enum G4EquationType
0059 {
0060 kEqMagnetic,
0061
0062 kEqMagneticWithSpin,
0063
0064
0065 kEqElectroMagnetic,
0066
0067 kEqEMfieldWithSpin,
0068
0069
0070 kEqEMfieldWithEDM,
0071
0072
0073 kUserEquation
0074 };
0075
0076
0077
0078 enum G4StepperType
0079 {
0080
0081 kCashKarpRKF45,
0082 kClassicalRK4,
0083 kBogackiShampine23,
0084 kBogackiShampine45,
0085 kDormandPrince745,
0086 kDormandPrinceRK56,
0087 kDormandPrinceRK78,
0088 kExplicitEuler,
0089 kImplicitEuler,
0090 kSimpleHeum,
0091 kSimpleRunge,
0092 kTsitourasRK45,
0093
0094
0095 kConstRK4,
0096 kExactHelixStepper,
0097 kHelixExplicitEuler,
0098 kHelixHeum,
0099 kHelixImplicitEuler,
0100 kHelixMixedStepper,
0101 kHelixSimpleRunge,
0102 kNystromRK4,
0103 kRKG3Stepper,
0104 kUserStepper,
0105
0106
0107 kRK547FEq1,
0108 kRK547FEq2,
0109 kRK547FEq3
0110 };
0111
0112
0113
0114
0115
0116
0117
0118
0119
0120
0121
0122
0123
0124
0125 class G4FieldParameters
0126 {
0127 public:
0128
0129 G4FieldParameters(const G4String& volumeName = "");
0130
0131 ~G4FieldParameters();
0132
0133
0134
0135
0136
0137 static G4String FieldTypeName(G4FieldType field);
0138
0139 static G4String EquationTypeName(G4EquationType equation);
0140
0141 static G4String StepperTypeName(G4StepperType stepper);
0142
0143 static G4FieldType GetFieldType(const G4String& name);
0144
0145 static G4EquationType GetEquationType(const G4String& name);
0146
0147 static G4StepperType GetStepperType(const G4String& name);
0148
0149
0150 void PrintParameters() const;
0151
0152
0153
0154
0155
0156 void SetFieldType(G4FieldType field);
0157
0158 void SetEquationType(G4EquationType equation);
0159
0160 void SetStepperType(G4StepperType stepper);
0161
0162 void SetUserEquationOfMotion(G4EquationOfMotion* equation);
0163
0164 void SetUserStepper(G4MagIntegratorStepper* stepper);
0165
0166
0167 void SetMinimumStep(G4double value);
0168
0169 void SetDeltaChord(G4double value);
0170
0171 void SetDeltaOneStep(G4double value);
0172
0173 void SetDeltaIntersection(G4double value);
0174
0175 void SetMinimumEpsilonStep(G4double value);
0176
0177 void SetMaximumEpsilonStep(G4double value);
0178
0179 void SetConstDistance(G4double value);
0180
0181
0182
0183
0184
0185 G4String GetVolumeName() const;
0186
0187
0188 G4FieldType GetFieldType() const;
0189
0190 G4EquationType GetEquationType() const;
0191
0192 G4StepperType GetStepperType() const;
0193
0194 G4EquationOfMotion* GetUserEquationOfMotion() const;
0195
0196 G4MagIntegratorStepper* GetUserStepper() const;
0197
0198
0199 G4double GetMinimumStep() const;
0200
0201 G4double GetDeltaChord() const;
0202
0203 G4double GetDeltaOneStep() const;
0204
0205 G4double GetDeltaIntersection() const;
0206
0207 G4double GetMinimumEpsilonStep() const;
0208
0209 G4double GetMaximumEpsilonStep() const;
0210
0211 G4double GetConstDistance() const;
0212
0213 private:
0214
0215 G4FieldParameters(const G4FieldParameters& right) = delete;
0216
0217 G4FieldParameters& operator=(const G4FieldParameters& right) = delete;
0218
0219
0220
0221
0222 inline static const G4double fgkDefaultMinimumStep = 0.01 * CLHEP::mm;
0223
0224 inline static const G4double fgkDefaultDeltaChord = 0.25 * CLHEP::mm;
0225
0226 inline static const G4double fgkDefaultDeltaOneStep = 0.01 * CLHEP::mm;
0227
0228 inline static const G4double fgkDefaultDeltaIntersection = 0.001 * CLHEP::mm;
0229
0230 inline static const G4double fgkDefaultMinimumEpsilonStep = 5.0e-5;
0231
0232 inline static const G4double fgkDefaultMaximumEpsilonStep = 0.001;
0233
0234 inline static const G4double fgkDefaultConstDistance = 0.;
0235
0236
0237
0238
0239 G4FieldParametersMessenger* fMessenger = nullptr;
0240
0241
0242 G4String fVolumeName;
0243
0244
0245 G4double fMinimumStep = fgkDefaultMinimumStep;
0246
0247 G4double fDeltaChord = fgkDefaultDeltaChord;
0248
0249 G4double fDeltaOneStep = fgkDefaultDeltaOneStep;
0250
0251 G4double fDeltaIntersection = fgkDefaultDeltaIntersection;
0252
0253 G4double fMinimumEpsilonStep = fgkDefaultMinimumEpsilonStep;
0254
0255 G4double fMaximumEpsilonStep = fgkDefaultMaximumEpsilonStep;
0256
0257
0258 G4FieldType fField = kMagnetic;
0259
0260
0261 G4EquationType fEquation = kEqMagnetic;
0262
0263
0264 G4StepperType fStepper = kDormandPrince745;
0265
0266
0267 G4EquationOfMotion* fUserEquation = nullptr;
0268
0269
0270 G4MagIntegratorStepper* fUserStepper = nullptr;
0271
0272
0273 G4double fConstDistance = fgkDefaultConstDistance;
0274 };
0275
0276
0277
0278
0279 inline void G4FieldParameters::SetFieldType(G4FieldType field)
0280 {
0281 fField = field;
0282 }
0283
0284
0285 inline void G4FieldParameters::SetEquationType(G4EquationType equation)
0286 {
0287 fEquation = equation;
0288 }
0289
0290
0291 inline void G4FieldParameters::SetStepperType(G4StepperType stepper)
0292 {
0293 fStepper = stepper;
0294 }
0295
0296
0297 inline void G4FieldParameters::SetMinimumStep(G4double value)
0298 {
0299 fMinimumStep = value;
0300 }
0301
0302
0303 inline void G4FieldParameters::SetDeltaChord(G4double value)
0304 {
0305 fDeltaChord = value;
0306 }
0307
0308
0309 inline void G4FieldParameters::SetDeltaOneStep(G4double value)
0310 {
0311 fDeltaOneStep = value;
0312 }
0313
0314
0315 inline void G4FieldParameters::SetDeltaIntersection(G4double value)
0316 {
0317 fDeltaIntersection = value;
0318 }
0319
0320
0321 inline void G4FieldParameters::SetMinimumEpsilonStep(G4double value)
0322 {
0323 fMinimumEpsilonStep = value;
0324 }
0325
0326
0327 inline void G4FieldParameters::SetMaximumEpsilonStep(G4double value)
0328 {
0329 fMaximumEpsilonStep = value;
0330 }
0331
0332
0333 inline void G4FieldParameters::SetConstDistance(G4double value)
0334 {
0335 fConstDistance = value;
0336 }
0337
0338
0339 inline G4String G4FieldParameters::GetVolumeName() const
0340 {
0341 return fVolumeName;
0342 }
0343
0344
0345 inline G4FieldType G4FieldParameters::GetFieldType() const { return fField; }
0346
0347
0348 inline G4EquationType G4FieldParameters::GetEquationType() const
0349 {
0350 return fEquation;
0351 }
0352
0353
0354 inline G4StepperType G4FieldParameters::GetStepperType() const
0355 {
0356 return fStepper;
0357 }
0358
0359
0360 inline G4EquationOfMotion* G4FieldParameters::GetUserEquationOfMotion() const
0361 {
0362 return fUserEquation;
0363 }
0364
0365
0366 inline G4MagIntegratorStepper* G4FieldParameters::GetUserStepper() const
0367 {
0368 return fUserStepper;
0369 }
0370
0371
0372 inline G4double G4FieldParameters::GetMinimumStep() const
0373 {
0374 return fMinimumStep;
0375 }
0376
0377
0378 inline G4double G4FieldParameters::GetDeltaChord() const
0379 {
0380 return fDeltaChord;
0381 }
0382
0383
0384 inline G4double G4FieldParameters::GetDeltaOneStep() const
0385 {
0386 return fDeltaOneStep;
0387 }
0388
0389
0390 inline G4double G4FieldParameters::GetDeltaIntersection() const
0391 {
0392 return fDeltaIntersection;
0393 }
0394
0395
0396 inline G4double G4FieldParameters::GetMinimumEpsilonStep() const
0397 {
0398 return fMinimumEpsilonStep;
0399 }
0400
0401
0402 inline G4double G4FieldParameters::GetMaximumEpsilonStep() const
0403 {
0404 return fMaximumEpsilonStep;
0405 }
0406
0407
0408 inline G4double G4FieldParameters::GetConstDistance() const
0409 {
0410 return fConstDistance;
0411 }
0412
0413 #endif