Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-08-02 08:28:36

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 
0026 /// \file G4FieldParameters.hh
0027 /// \brief Definition of the G4FieldParameters class
0028 ///
0029 /// This code was initially developed in Geant4 VMC package
0030 /// (https://github.com/vmc-project)
0031 /// and adapted to Geant4.
0032 ///
0033 /// \author I. Hrivnacova; IJCLab, Orsay
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 /// The available fields in Geant4
0049 enum G4FieldType
0050 {
0051   kMagnetic,        ///< magnetic field
0052   kElectroMagnetic, ///< electromagnetic field
0053   kGravity          ///< gravity field
0054 };
0055 
0056 /// The available equations of motion of a particle in a field
0057 /// in Geant4
0058 enum G4EquationType
0059 {
0060   kEqMagnetic,        ///< G4Mag_UsualEqRhs: the standard right-hand side for
0061                       ///< equation of motion.
0062   kEqMagneticWithSpin,///< G4Mag_SpinEqRhs: the equation of motion for a particle
0063                       ///< with spin
0064                       ///< in a pure magnetic field
0065   kEqElectroMagnetic, ///< G4EqMagElectricField: Equation of motion in a combined
0066                       ///< electric and magnetic field
0067   kEqEMfieldWithSpin, ///< G4EqEMFieldWithSpin: Equation of motion for a
0068                       ///< particle with spin
0069                       ///< in a combined electric and magnetic field
0070   kEqEMfieldWithEDM,  ///< G4EqEMFieldWithEDM: Equation of motion in a combined
0071                       ///< electric and magnetic field, with spin tracking for
0072                       ///< both MDM and EDM terms
0073   kUserEquation       ///< User defined equation of motion
0074 };
0075 
0076 /// The available integrator of particle's equation of motion
0077 /// in Geant4
0078 enum G4StepperType
0079 {
0080   // steppers with equation of motion of generic type (G4EquationOfMotion)
0081   kCashKarpRKF45,     ///< G4CashKarpRKF45
0082   kClassicalRK4,      ///< G4ClassicalRK4
0083   kBogackiShampine23, ///< G4BogackiShampine23
0084   kBogackiShampine45, ///< G4BogackiShampine45
0085   kDormandPrince745,  ///< G4DormandPrince745
0086   kDormandPrinceRK56, ///< G4DormandPrinceRK56
0087   kDormandPrinceRK78, ///< G4DormandPrinceRK78
0088   kExplicitEuler,     ///< G4ExplicitEuler
0089   kImplicitEuler,     ///< G4ImplicitEuler
0090   kSimpleHeum,        ///< G4SimpleHeum
0091   kSimpleRunge,       ///< G4SimpleRunge
0092   kTsitourasRK45,     ///< G4TsitourasRK45
0093 
0094   // steppers with equation of motion of G4Mag_UsualEqRhs type
0095   kConstRK4,           ///< G4ConstRK4
0096   kExactHelixStepper,  ///< G4ExactHelixStepper
0097   kHelixExplicitEuler, ///< G4HelixExplicitEuler
0098   kHelixHeum,          ///< G4HelixHeum
0099   kHelixImplicitEuler, ///< G4HelixImplicitEuler
0100   kHelixMixedStepper,  ///< G4HelixMixedStepper
0101   kHelixSimpleRunge,   ///< G4HelixSimpleRunge
0102   kNystromRK4,         ///< G4NystromRK4
0103   kRKG3Stepper,        ///< G4RKG3_Stepper
0104   kUserStepper,        ///< User defined stepper
0105 
0106   // FSAL steppers
0107   kRK547FEq1, ///< G4RK547FEq1
0108   kRK547FEq2, ///< G4RK547FEq2
0109   kRK547FEq3  ///< G4RK547FEq3
0110 };
0111 
0112 /// \brief The magnetic field parameters
0113 ///
0114 /// The class defines the type of equation of motion of a particle
0115 /// in a field and the integration method, as well as other accuracy
0116 /// parameters.
0117 ///
0118 /// The default values correspond to the defaults set in Geant4
0119 /// (taken from Geant4 9.3 release.)
0120 /// As Geant4 classes to not provide access methods for these defaults,
0121 /// the defaults have to be checked with each new Geant4 release.
0122 ///
0123 /// \author I. Hrivnacova; IJCLab, Orsay
0124 
0125 class G4FieldParameters
0126 {
0127  public:
0128   /// Standard and default constructor
0129   G4FieldParameters(const G4String& volumeName = "");
0130   /// Destructor
0131   ~G4FieldParameters();
0132 
0133   // Methods
0134   //
0135 
0136   /// Return the field type as a string
0137   static G4String FieldTypeName(G4FieldType field);
0138   /// Return the equation type as a string
0139   static G4String EquationTypeName(G4EquationType equation);
0140   /// Return the stepper type as a string
0141   static G4String StepperTypeName(G4StepperType stepper);
0142   /// Return the field type for given field type name
0143   static G4FieldType GetFieldType(const G4String& name);
0144   /// Return the equation type for given equation type name
0145   static G4EquationType GetEquationType(const G4String& name);
0146   /// Return the stepper type for given stepper type name
0147   static G4StepperType GetStepperType(const G4String& name);
0148 
0149   /// Prints all customizable accuracy parameters
0150   void PrintParameters() const;
0151 
0152   // Set methods
0153   //
0154 
0155   /// Set type of field
0156   void SetFieldType(G4FieldType field);
0157   /// Set Type of equation of motion of a particle in a field
0158   void SetEquationType(G4EquationType equation);
0159   /// Type of integrator of particle's equation of motion
0160   void SetStepperType(G4StepperType stepper);
0161   /// Set user defined equation of motion
0162   void SetUserEquationOfMotion(G4EquationOfMotion* equation);
0163   /// Set user defined integrator of particle's equation of motion
0164   void SetUserStepper(G4MagIntegratorStepper* stepper);
0165 
0166   /// Set minimum step in G4ChordFinder
0167   void SetMinimumStep(G4double value);
0168   /// Set delta chord in G4ChordFinder
0169   void SetDeltaChord(G4double value);
0170   /// Set delta one step in global field manager
0171   void SetDeltaOneStep(G4double value);
0172   /// Set delta intersection in global field manager
0173   void SetDeltaIntersection(G4double value);
0174   /// Set minimum epsilon step in global field manager
0175   void SetMinimumEpsilonStep(G4double value);
0176   /// Set maximum epsilon step in global field manager
0177   void SetMaximumEpsilonStep(G4double value);
0178   /// Set the distance within which the field is considered constant
0179   void SetConstDistance(G4double value);
0180 
0181   // Get methods
0182   //
0183 
0184   // Get the name of associated volume, if local field
0185   G4String GetVolumeName() const;
0186 
0187   /// Get type of field
0188   G4FieldType GetFieldType() const;
0189   /// Get type of equation of motion of a particle in a field
0190   G4EquationType GetEquationType() const;
0191   /// Get rype of integrator of particle's equation of motion
0192   G4StepperType GetStepperType() const;
0193   /// Get user defined equation of motion
0194   G4EquationOfMotion* GetUserEquationOfMotion() const;
0195   /// Get user defined integrator of particle's equation of motion
0196   G4MagIntegratorStepper* GetUserStepper() const;
0197 
0198   /// Get minimum step in G4ChordFinder
0199   G4double GetMinimumStep() const;
0200   /// Get delta chord in G4ChordFinder
0201   G4double GetDeltaChord() const;
0202   /// Get delta one step in global field manager
0203   G4double GetDeltaOneStep() const;
0204   /// Get delta intersection in global field manager
0205   G4double GetDeltaIntersection() const;
0206   /// Get minimum epsilon step in global field manager
0207   G4double GetMinimumEpsilonStep() const;
0208   /// Get maximum epsilon step in global field manager
0209   G4double GetMaximumEpsilonStep() const;
0210   /// Get the distance within which the field is considered constant
0211   G4double GetConstDistance() const;
0212 
0213  private:
0214   /// Not implemented
0215   G4FieldParameters(const G4FieldParameters& right) = delete;
0216   /// Not implemented
0217   G4FieldParameters& operator=(const G4FieldParameters& right) = delete;
0218 
0219   // static data members
0220   //
0221   /// Default minimum step in G4ChordFinder
0222   inline static const G4double fgkDefaultMinimumStep  = 0.01 * CLHEP::mm;
0223   /// Default delta chord in G4ChordFinder
0224   inline static const G4double fgkDefaultDeltaChord = 0.25 * CLHEP::mm;
0225   /// Default delta one step in global field manager
0226   inline static const G4double fgkDefaultDeltaOneStep = 0.01 * CLHEP::mm;
0227   /// Delta intersection in global field manager
0228   inline static const G4double fgkDefaultDeltaIntersection = 0.001 * CLHEP::mm;
0229   /// Default minimum epsilon step in global field manager
0230   inline static const G4double fgkDefaultMinimumEpsilonStep = 5.0e-5;
0231   /// Default maximum epsilon step in global field manager
0232   inline static const G4double fgkDefaultMaximumEpsilonStep = 0.001;
0233   /// Default constant distance
0234   inline static const G4double fgkDefaultConstDistance = 0.;
0235 
0236   // data members
0237   //
0238   /// Messenger for this class
0239   G4FieldParametersMessenger* fMessenger = nullptr;
0240 
0241   /// The name of associated volume, if local field
0242   G4String fVolumeName;
0243 
0244   /// Minimum step in G4ChordFinder
0245   G4double fMinimumStep = fgkDefaultMinimumStep;
0246   /// Delta chord in G4ChordFinder
0247   G4double fDeltaChord = fgkDefaultDeltaChord;
0248   /// Delta one step in global field manager
0249   G4double fDeltaOneStep = fgkDefaultDeltaOneStep;
0250   /// Delta intersection in global field manager
0251   G4double fDeltaIntersection = fgkDefaultDeltaIntersection;
0252   /// Minimum epsilon step in global field manager
0253   G4double fMinimumEpsilonStep = fgkDefaultMinimumEpsilonStep;
0254   /// Maximum epsilon step in global field manager
0255   G4double fMaximumEpsilonStep = fgkDefaultMaximumEpsilonStep;
0256 
0257   /// Type of field
0258   G4FieldType fField = kMagnetic;
0259 
0260   /// Type of equation of motion of a particle in a field
0261   G4EquationType fEquation = kEqMagnetic;
0262 
0263   /// Type of integrator of particle's equation of motion
0264   G4StepperType fStepper = kDormandPrince745;
0265 
0266   /// User defined equation of motion
0267   G4EquationOfMotion* fUserEquation = nullptr;
0268 
0269   /// User defined integrator of particle's equation of motion
0270   G4MagIntegratorStepper* fUserStepper = nullptr;
0271 
0272   /// The distance within which the field is considered constant
0273   G4double fConstDistance = fgkDefaultConstDistance;
0274 };
0275 
0276 // inline functions
0277 
0278 // Set type of field
0279 inline void G4FieldParameters::SetFieldType(G4FieldType field)
0280 {
0281   fField = field;
0282 }
0283 
0284 // Set the type of equation of motion of a particle in a field
0285 inline void G4FieldParameters::SetEquationType(G4EquationType equation)
0286 {
0287   fEquation = equation;
0288 }
0289 
0290 // Set the type of integrator of particle's equation of motion
0291 inline void G4FieldParameters::SetStepperType(G4StepperType stepper)
0292 {
0293   fStepper = stepper;
0294 }
0295 
0296 // Set minimum step in G4ChordFinder
0297 inline void G4FieldParameters::SetMinimumStep(G4double value)
0298 {
0299   fMinimumStep = value;
0300 }
0301 
0302 // Set delta chord in G4ChordFinder
0303 inline void G4FieldParameters::SetDeltaChord(G4double value)
0304 {
0305   fDeltaChord = value;
0306 }
0307 
0308 // Set delta one step in global field manager
0309 inline void G4FieldParameters::SetDeltaOneStep(G4double value)
0310 {
0311   fDeltaOneStep = value;
0312 }
0313 
0314 // Set delta intersection in global field manager
0315 inline void G4FieldParameters::SetDeltaIntersection(G4double value)
0316 {
0317   fDeltaIntersection = value;
0318 }
0319 
0320 // Set minimum epsilon step in global field manager
0321 inline void G4FieldParameters::SetMinimumEpsilonStep(G4double value)
0322 {
0323   fMinimumEpsilonStep = value;
0324 }
0325 
0326 // Set maximum epsilon step in global field manager
0327 inline void G4FieldParameters::SetMaximumEpsilonStep(G4double value)
0328 {
0329   fMaximumEpsilonStep = value;
0330 }
0331 
0332 // Set the distance within which the field is considered constant
0333 inline void G4FieldParameters::SetConstDistance(G4double value)
0334 {
0335   fConstDistance = value;
0336 }
0337 
0338 // Return the name of associated volume, if local field
0339 inline G4String G4FieldParameters::GetVolumeName() const
0340 {
0341   return fVolumeName;
0342 }
0343 
0344 // Return the type of field
0345 inline G4FieldType G4FieldParameters::GetFieldType() const { return fField; }
0346 
0347 // Return the type of equation of motion of a particle in a field
0348 inline G4EquationType G4FieldParameters::GetEquationType() const
0349 {
0350   return fEquation;
0351 }
0352 
0353 // Return the type of integrator of particle's equation of motion
0354 inline G4StepperType G4FieldParameters::GetStepperType() const
0355 {
0356   return fStepper;
0357 }
0358 
0359 // Return the user defined equation of motion
0360 inline G4EquationOfMotion* G4FieldParameters::GetUserEquationOfMotion() const
0361 {
0362   return fUserEquation;
0363 }
0364 
0365 // Return the user defined integrator of particle's equation of motion
0366 inline G4MagIntegratorStepper* G4FieldParameters::GetUserStepper() const
0367 {
0368   return fUserStepper;
0369 }
0370 
0371 // Return minimum step in G4ChordFinder
0372 inline G4double G4FieldParameters::GetMinimumStep() const
0373 {
0374   return fMinimumStep;
0375 }
0376 
0377 // Return delta chord in G4ChordFinder
0378 inline G4double G4FieldParameters::GetDeltaChord() const
0379 {
0380   return fDeltaChord;
0381 }
0382 
0383 // Return delta one step in global field manager
0384 inline G4double G4FieldParameters::GetDeltaOneStep() const
0385 {
0386   return fDeltaOneStep;
0387 }
0388 
0389 // Return delta intersection in global field manager
0390 inline G4double G4FieldParameters::GetDeltaIntersection() const
0391 {
0392   return fDeltaIntersection;
0393 }
0394 
0395 // Return minimum epsilon step in global field manager
0396 inline G4double G4FieldParameters::GetMinimumEpsilonStep() const
0397 {
0398   return fMinimumEpsilonStep;
0399 }
0400 
0401 // Return maximum epsilon step in global field manager
0402 inline G4double G4FieldParameters::GetMaximumEpsilonStep() const
0403 {
0404   return fMaximumEpsilonStep;
0405 }
0406 
0407 // Return the distance within which the field is considered constant
0408 inline G4double G4FieldParameters::GetConstDistance() const
0409 {
0410   return fConstDistance;
0411 }
0412 
0413 #endif // G4FIELDPARAMETERS_HH