4.3.  Electromagnetic Field

4.3.1.  An Overview of Propagation in a Field

Geant4 is capable of describing and propagating in a variety of fields. Magnetic fields, electric fields, electromagnetic fields, and gravity fields, uniform or non-uniform, can specified for a Geant4 setup. The propagation of tracks inside them can be performed to a user-defined accuracy.

In order to propagate a track inside a field, the equation of motion of the particle in the field is integrated. In general, this is done using a Runge-Kutta method for the integration of ordinary differential equations. However, for specific cases where an analytical solution is known, it is possible to utilize this instead. Several Runge-Kutta methods are available, suitable for different conditions. In specific cases (such as a uniform field where the analytical solution is known) different solvers can also be used. In addition, when an approximate analytical solution is known, it is possible to utilize it in an iterative manner in order to converge to the solution to the precision required. This latter method is currently implemented and can be used particularly well for magnetic fields that are almost uniform.

Once a method is chosen that calculates the track's propagation in a specific field, the curved path is broken up into linear chord segments. These chord segments are determined so that they closely approximate the curved path. The chords are then used to interrogate the Navigator as to whether the track has crossed a volume boundary. Several parameters are available to adjust the accuracy of the integration and the subsequent interrogation of the model geometry.

How closely the set of chords approximates a curved trajectory is governed by a parameter called the miss distance (also called the chord distance ). This is an upper bound for the value of the sagitta - the distance between the 'real' curved trajectory and the approximate linear trajectory of the chord. By setting this parameter, the user can control the precision of the volume interrogation. Every attempt has been made to ensure that all volume interrogations will be made to an accuracy within this miss distance.

Miss Distance

Figure 4.6.  The curved trajectory will be approximated by chords, so that the maximum estimated distance between curve and chord is less than the the miss distance.


In addition to the miss distance there are two more parameters which the user can set in order to adjust the accuracy (and performance) of tracking in a field. In particular these parameters govern the accuracy of the intersection with a volume boundary and the accuracy of the integration of other steps. As such they play an important role for tracking.

The delta intersection parameter is the accuracy to which an intersection with a volume boundary is calculated. If a candidate boundary intersection is estimated to have a precision better than this, it is accepted. This parameter is especially important because it is used to limit a bias that our algorithm (for boundary crossing in a field) exhibits. This algorithm calculates the intersection with a volume boundary using a chord between two points on the curved particle trajectory. As such, the intersection point is always on the 'inside' of the curve. By setting a value for this parameter that is much smaller than some acceptable error, the user can limit the effect of this bias on, for example, the future estimation of the reconstructed particle momentum.

Figure: The distance between the calculated chord intersection point and a computed curve point.

Figure 4.7.  The distance between the calculated chord intersection point C and a computed curve point D is used to determine whether C is an accurate representation of the intersection of the curved path ADB with a volume boundary. Here CD is likely too large, and a new intersection on the chord AD will be calculated.


The delta one step parameter is the accuracy for the endpoint of 'ordinary' integration steps, those which do not intersect a volume boundary. This parameter is a limit on the estimated error of the endpoint of each physics step. It can be seen as akin to a statistical uncertainty and is not expected to contribute any systematic behavior to physical quantities. In contrast, the bias addressed by delta intersection is clearly correlated with potential systematic errors in the momentum of reconstructed tracks. Thus very strict limits on the intersection parameter should be used in tracking detectors or wherever the intersections are used to reconstruct a track's momentum.

Delta intersection and delta one step are parameters of the Field Manager; the user can set them according to the demands of his application. Because it is possible to use more than one field manager, different values can be set for different detector regions.

Note that reasonable values for the two parameters are strongly coupled: it does not make sense to request an accuracy of 1 nm for delta intersection and accept 100 &#956m for the delta one step error value. Nevertheless delta intersection is the more important of the two. It is recommended that these parameters should not differ significantly - certainly not by more than an order of magnitude.

4.3.2.  Practical Aspects

4.3.2.1.  Creating a Magnetic Field for a Detector

The simplest way to define a field for a detector involves the following steps:

  1. create a field:

        G4UniformMagField* magField
          = new G4UniformMagField(G4ThreeVector(0.,0.,fieldValue));
        

  2. set it as the default field:

        G4FieldManager* fieldMgr
          = G4TransportationManager::GetTransportationManager()
            ->GetFieldManager();
           fieldMgr->SetDetectorField(magField);
        

  3. create the objects which calculate the trajectory:

        fieldMgr->CreateChordFinder(magField);
        

To change the accuracy of volume intersection use the SetDeltaChord method:

    fieldMgr->GetChordFinder()->SetDeltaChord( G4double newValue);

Since 10.0 version, it is also possible to perform all three steps above at once using the G4GlobalMagFieldMessenger class:

    G4ThreeVector fieldValue = G4ThreeVector();
    fMagFieldMessenger = new G4GlobalMagFieldMessenger(fieldValue);
    fMagFieldMessenger->SetVerboseLevel(1);
    

The messenger creates the global uniform magnetic field, which is activated (set to the G4TransportationManager object) only when the fieldValue is non zero vector. The messenger class setter functions can be then used to change the field value (and activate or inactivate the field again) or the level of output messages. The messenger also takes care of deleting the field.

As its class name suggests, the messenger creates also UI commands which can be used to change the field value and the verbose level interactively or from a macro:

    /globalField/setValue vx vy vz unit
    /globalField/verbose level
    

4.3.2.2.  Creating a Field for a Part of the Volume Hierarchy

It is possible to create a field for a part of the detector. In particular it can describe the field (with pointer fEmField, for example) inside a logical volume and all its daughters. This can be done by simply creating a G4FieldManager and attaching it to a logical volume (with pointer, logicVolumeWithField, for example) or set of logical volumes.

  G4bool allLocal = true;       
  logicVolumeWithField->SetFieldManager(localFieldManager, allLocal);

Using the second parameter to SetFieldManager you choose whether daughter volumes of this logical volume will also be given this new field. If it has the value true, the field will be assigned also to its daughters, and all their sub-volumes. Else, if it is false, it will be copied only to those daughter volumes which do not have a field manager already.

4.3.2.3.  Creating an Electric or Electromagnetic Field

The design and implementation of the Field category allows and enables the use of an electric or combined electromagnetic field. These fields can also vary with time, as can magnetic fields.

Source listing Example 4.12 shows how to define a uniform electric field for the whole of a detector.

Example 4.12.  How to define a uniform electric field for the whole of a detector, extracted from example in examples/extended/field/field02 .

 // in the header file (or first)
  #include "G4EqMagElectricField.hh"
  #include "G4UniformElectricField.hh"

  ...
    G4ElectricField*        fEMfield;
    G4EqMagElectricField*   fEquation;
    G4MagIntegratorStepper* fStepper;
    G4FieldManager*         fFieldMgr;
    G4double                fMinStep ;
    G4ChordFinder*          fChordFinder ;

 // in the source file

  {
    fEMfield = new G4UniformElectricField(
                 G4ThreeVector(0.0,100000.0*kilovolt/cm,0.0));

    // Create an equation of motion for this field
    fEquation = new G4EqMagElectricField(fEMfield); 

    G4int nvar = 8;
    fStepper = new G4ClassicalRK4( fEquation, nvar );       

    // Get the global field manager 
    fFieldManager= G4TransportationManager::GetTransportationManager()->
         GetFieldManager();
    // Set this field to the global field manager 
    fFieldManager->SetDetectorField(fEMfield );

    fMinStep     = 0.010*mm ; // minimal step of 10 microns

    fIntgrDriver = new G4MagInt_Driver(fMinStep, 
                                     fStepper, 
                                     fStepper->GetNumberOfVariables() );

    fChordFinder = new G4ChordFinder(fIntgrDriver);
    fFieldManager->SetChordFinder( fChordFinder );

  }


An example with an electric field is examples/extended/field/field02, where the class F02ElectricFieldSetup demonstrates how to set these and other parameters, and how to choose different Integration Steppers. An example with a uniform gravity field (G4UniformGravityField) is examples/extended/field/field06.

The user can also create their own type of field, inheriting from G4VField, and an associated Equation of Motion class (inheriting from G4EqRhs) to simulate other types of fields.

4.3.2.4.  Choosing a Stepper

Runge-Kutta integration is used to compute the motion of a charged track in a general field. There are many general steppers from which to choose, of low and high order, and specialized steppers for pure magnetic fields. By default, Geant4 uses the classical fourth-order Runge-Kutta stepper, which is general purpose and robust. If the field is known to have specific properties, lower or higher order steppers can be used to obtain the same quality results using fewer computing cycles.

In particular, if the field is calculated from a field map, a lower order stepper is recommended. The less smooth the field is, the lower the order of the stepper that should be used. The choice of lower order steppers includes the third order stepper G4SimpleHeum, the second order G4ImplicitEuler and G4SimpleRunge, and the first order G4ExplicitEuler. A first order stepper would be useful only for very rough fields. For somewhat smooth fields (intermediate), the choice between second and third order steppers should be made by trial and error. Trying a few different types of steppers for a particular field or application is suggested if maximum performance is a goal.

The choice of stepper depends on the type of field: magnetic or general. A general field can be an electric or electromagnetic field, it can be a magnetic field or a user-defined field (which requires a user-defined equation of motion.) For a general field several steppers are available as alternatives to the default (G4ClassicalRK4):

   G4int nvar = 8;  // To integrate time & energy 
                    //    in addition to position, momentum
   G4EqMagElectricField* fEquation= new G4EqMagElectricField(fEMfield); 

   fStepper = new G4SimpleHeum( fEquation, nvar );         
            //  3rd  order, a good alternative to ClassicalRK
   fStepper = new G4SimpleRunge( fEquation, nvar ); 
            //  2nd  order, for less smooth fields
   fStepper = new G4CashKarpRKF45( fEquation );     
            // 4/5th order for very smooth fields 

Specialized steppers for pure magnetic fields are also available. They take into account the fact that a local trajectory in a slowly varying field will not vary significantly from a helix. Combining this in with a variation the Runge-Kutta method can provide higher accuracy at lower computational cost when large steps are possible.

   G4Mag_UsualEqRhs* 
      fEquation = new G4Mag_UsualEqRhs(fMagneticField); 
   fStepper = new G4HelixImplicitEuler( fEquation ); 
      // Note that for magnetic field that do not vary with time,
      //  the default number of variables suffices.
 
      // or ..
     fStepper = new G4HelixExplicitEuler( fEquation ); 
     fStepper = new G4HelixSimpleRunge( fEquation );   

A new stepper for propagation in magnetic field is available in release 9.3. Choosing the G4NystromRK4 stepper provides accuracy near that of G4ClassicalRK4 (4th order) with a significantly reduced cost in field evaluation. Using a novel analytical expression for estimating the error of a proposed step and the Nystrom reuse of the mid-point field value, it requires only 2 additional field evaluations per attempted step, in place of 10 field evaluations of ClassicalRK4 (which uses the general midpoint method for estimating the step error.)

   G4Mag_UsualEqRhs* 
      pMagFldEquation = new G4Mag_UsualEqRhs(fMagneticField); 
   fStepper = new G4NystromRK4( pMagFldEquation ); 

It is proposed as an alternative stepper in the case of a pure magnetic field. It is not applicable for the simulation of electric or full electromagnetic or other types of field. For a pure magnetic field, results should be fully compatible with the results of ClassicalRK4 in nearly all cases. ( The only potential exceptions are large steps for tracks with small momenta - which cannot be integrated well by any RK method except the Helical extended methods.)

You can choose an alternative stepper either when the field manager is constructed or later. At the construction of the ChordFinder it is an optional argument:

 G4ChordFinder( G4MagneticField* itsMagField,
                G4double         stepMinimum = 1.0e-2 * mm,
                G4MagIntegratorStepper* pItsStepper = 0 );

To change the stepper at a later time use

 pChordFinder->GetIntegrationDriver()
             ->RenewStepperAndAdjust( newStepper );

4.3.2.5.  How to Adjust the Accuracy of Propagation

In order to obtain a particular accuracy in tracking particles through an electromagnetic field, it is necessary to adjust the parameters of the field propagation module. In the following section, some of these additional parameters are discussed.

When integration is used to calculate the trajectory, it is necessary to determine an acceptable level of numerical imprecision in order to get performant simulation with acceptable errors. The parameters in Geant4 tell the field module what level of integration inaccuracy is acceptable.

In all quantities which are integrated (position, momentum, energy) there will be errors. Here, however, we focus on the error in two key quantities: the position and the momentum. (The error in the energy will come from the momentum integration).

Three parameters exist which are relevant to the integration accuracy. DeltaOneStep is a distance and is roughly the position error which is acceptable in an integration step. Since many integration steps may be required for a single physics step, DeltaOneStep should be a fraction of the average physics step size. The next two parameters impose a further limit on the relative error of the position/momentum inaccuracy. EpsilonMin and EpsilonMax impose a minimum and maximum on this relative error - and take precedence over DeltaOneStep. (Note: if you set EpsilonMin=EpsilonMax=your-value, then all steps will be made to this relative precision.

Example 4.13.  How to set accuracy parameters for the 'global' field of the setup.

  G4FieldManager *globalFieldManager;

  G4TransportationManager *transportMgr= 
       G4TransportationManager::GetTransportationManager();

  globalFieldManager = transportMgr->GetFieldManager();
                            // Relative accuracy values:
  G4double minEps= 1.0e-5;  //   Minimum & value for smallest steps
  G4double maxEps= 1.0e-4;  //   Maximum & value for largest steps

  globalFieldManager->SetMinimumEpsilonStep( minEps );
  globalFieldManager->SetMaximumEpsilonStep( maxEps );
  globalFieldManager->SetDeltaOneStep( 0.5e-3 * mm );  // 0.5 micrometer

  G4cout << "EpsilonStep: set min= " << minEps << " max= " << maxEps << G4endl;


We note that the relevant parameters above limit the inaccuracy in each step. The final inaccuracy due to the full trajectory will accumulate!

The exact point at which a track crosses a boundary is also calculated with finite accuracy. To limit this inaccuracy, a parameter called DeltaIntersection is used. This is a maximum for the inaccuracy of a single boundary crossing. Thus the accuracy of the position of the track after a number of boundary crossings is directly proportional to the number of boundaries.

4.3.2.6.  Reducing the number of field calls to speed-up simulation

An additional method to reduce the number of field evaluations is to utilise the new class G4CachedMagneticField class. It is applicable only for pure magnetic fields which do not vary with time.

  G4MagneticField * pMagField;  // Your field - Defined elsewhere

  G4double                 distanceConst = 2.5 * millimeter; 
  G4MagneticField * pCachedMagField= new G4CachedMagneticField(  pMagField,  distanceConst); 

4.3.2.7.  Choosing different accuracies for the same volume

It is possible to create a FieldManager which has different properties for particles of different momenta (or depending on other parameters of a track). This is useful, for example, in obtaining high accuracy for 'important' tracks (e.g. muons) and accept less accuracy in tracking others (e.g. electrons). To use this, you must create your own field manager which uses the method

  void ConfigureForTrack( const G4Track * );

to configure itself using the parameters of the current track. An example of this will be available in examples/extended/field05.

4.3.2.8.  Parameters that must scale with problem size

The default settings of this module are for problems with the physical size of a typical high energy physics setup, that is, distances smaller than about one kilometer. A few parameters are necessary to carry this information to the magnetic field module, and must typically be rescaled for problems of vastly different sizes in order to get reasonable performance and robustness. Two of these parameters are the maximum acceptable step and the minimum step size.

The maximum acceptable step should be set to a distance larger than the biggest reasonable step. If the apparatus in a setup has a diameter of two meters, a likely maximum acceptable steplength would be 10 meters. A particle could then take large spiral steps, but would not attempt to take, for example, a 1000-meter-long step in the case of a very low-density material. Similarly, for problems of a planetary scale, such as the earth with its radius of roughly 6400 km, a maximum acceptabe steplength of a few times this value would be reasonable.

An upper limit for the size of a step is a parameter of G4PropagatorInField, and can be set by calling its SetLargestAcceptableStep method.

The minimum step size is used during integration to limit the amount of work in difficult cases. It is possible that strong fields or integration problems can force the integrator to try very small steps; this parameter stops them from becoming unnecessarily small.

Trial steps smaller than this parameter will be treated with less accuracy, and may even be ignored, depending on the situation.

The minimum step size is a parameter of the MagInt_Driver, but can be set in the contstructor of G4ChordFinder, as in the source listing above.

4.3.2.9.  Known Issues

Currently it is computationally expensive to change the miss distance to very small values, as it causes tracks to be limited to curved sections whose 'bend' is smaller than this value. (The bend is the distance of the mid-point from the chord between endpoints.) For tracks with small curvature (typically low momentum particles in strong fields) this can cause a large number of steps

  • even in areas where there are no volumes to intersect (something that is expected to be addressed in future development, in which the safety will be utilized to partially alleviate this limitation)

  • especially in a region near a volume boundary (in which case it is necessary in order to discover whether a track might intersect a volume for only a short distance.)

Requiring such precision at the intersection is clearly expensive, and new development would be necessary to minimize the expense.

By contrast, changing the intersection parameter is less computationally expensive. It causes further calculation for only a fraction of the steps, in particular those that intersect a volume boundary.

4.3.3.  Spin Tracking

The effects of a particle's motion on the precession of its spin angular momentum in slowly varying external fields are simulated. The relativistic equation of motion for spin is known as the BMT equation. The equation demonstrates a remarkable property; in a purely magnetic field, in vacuum, and neglecting small anomalous magnetic moments, the particle's spin precesses in such a manner that the longitudinal polarization remains a constant, whatever the motion of the particle. But when the particle interacts with electric fields of the medium and multiple scatters, the spin, which is related to the particle's magnetic moment, does not participate, and the need thus arises to propagate it independent of the momentum vector. In the case of a polarized muon beam, for example, it is important to predict the muon's spin direction at decay-time in order to simulate the decay electron (Michel) distribution correctly.

In order to track the spin of a particle in a magnetic field, you need to code the following:

  1. in your DetectorConstruction

    #include "G4Mag_SpinEqRhs.hh"
    
      G4Mag_EqRhs* fEquation = new G4Mag_SpinEqRhs(magField);
    
      G4MagIntegratorStepper* pStepper = new G4ClassicalRK4(fEquation,12);
                                                           notice the 12 
        

  2. in your PrimaryGenerator

      particleGun->SetParticlePolarization(G4ThreeVector p)
        

    for example:

      particleGun->
      SetParticlePolarization(-(particleGun->GetParticleMomentumDirection()));
    
      // or
      particleGun->
      SetParticlePolarization(particleGun->GetParticleMomentumDirection()
                                             .cross(G4ThreeVector(0.,1.,0.)));
        

    where you set the initial spin direction.

While the G4Mag_SpinEqRhs class constructor

  G4Mag_SpinEqRhs::G4Mag_SpinEqRhs( G4MagneticField* MagField )
     : G4Mag_EqRhs( MagField ) 
  {
      anomaly = 1.165923e-3;
  }

sets the muon anomaly by default, the class also comes with the public method:

  inline void SetAnomaly(G4double a) { anomaly = a; }

with which you can set the magnetic anomaly to any value you require.

The code has been rewritten (in Release 9.5) such that field tracking of the spin can now be done for charged and neutral particles with a magnetic moment, for example spin tracking of ultra cold neutrons. This requires the user to set EnableUseMagneticMoment, a method of the G4Transportation process. The force resulting from the term, μ⋅∇B, is not yet implemented in Geant4 (for example, simulated trajectory of a neutral hydrogen atom trapped by its magnetic moment in a gradient B-field.)