3.2.  Electromagnetic Fields

3.2.1.  Creating a New Type of Field

Geant4 currently handles magnetic and electric fields and, in future releases, will handle combined electromagnetic fields. Fields due to other forces, not yet included in Geant4, can be provided by describing the new field and the force it exerts on a particle passing through it. For the time being, all fields must be time-independent. This restriction may be lifted in the future.

In order to accommodate a new type of field, two classes must be created: a field type and a class that determines the force. The Geant4 system must then be informed of the new field.

A new Field class

A new type of Field class may be created by inheriting from G4Field

      class NewField : public G4Field
      {
         public: 
            void  GetFieldValue( const  double Point[3],
                                        double *pField )=0;
      }

and deciding how many components your field will have, and what each component represents. For example, three components are required to describe a vector field while only one component is required to describe a scalar field.

If you want your field to be a combination of different fields, you must choose your convention for which field goes first, which second etc. For example, to define an electromagnetic field we follow the convention that components 0,1 and 2 refer to the magnetic field and components 3, 4 and 5 refer to the electric field.

By leaving the GetFieldValue method pure virtual, you force those users who want to describe their field to create a class that implements it for their detector's instance of this field. So documenting what each component means is required, to give them the necessary information.

For example someone can describe DetectorAbc's field by creating a class DetectorAbcField, that derives from your NewField

class DetectorAbcField : public NewField 
{
  public:
     void  MyFieldGradient::GetFieldValue( const double Point[3], 
                                           double *pField );
}

They then implement the function GetFieldValue

      void  MyFieldGradient::GetFieldValue( const  double Point[3],
                                            double *pField )
      {
         // We expect pField to point to pField[9]; 
         // This & the order of the components of pField is your own 
         // convention
       
         // We calculate the value of pField at Point ...
      }

A new Equation of Motion for the new Field

Once you have created a new type of field, you must create an Equation of Motion for this Field. This is required in order to obtain the force that a particle feels.

To do this you must inherit from G4Mag_EqRhs and create your own equation of motion that understands your field. In it you must implement the virtual function EvaluateRhsGivenB. Given the value of the field, this function calculates the value of the generalised force. This is the only function that a subclass must define.

     virtual void EvaluateRhsGivenB( const  G4double y[],
                              const  G4double B[3],
                                     G4double dydx[] ) const = 0;

In particular, the derivative vector dydx is a vector with six components. The first three are the derivative of the position with respect to the curve length. Thus they should set equal to the normalised velocity, which is components 3, 4 and 5 of y.

  (dydx[0], dydx[1], dydx[2]) = (y[3], y[4], y[5])

The next three components are the derivatives of the velocity vector with respect to the path length. So you should write the "force" components for

      dydx[3], dydx[4] and dydx[5]

for your field.

Get a G4FieldManager to use your field

In order to inform the Geant4 system that you want it to use your field as the global field, you must do the following steps:

  1. Create a Stepper of your choice:

           
          yourStepper = new G4ClassicalRK( yourEquationOfMotion );
                    // or if your field is not smooth eg
                    //     new G4ImplicitEuler( yourEquationOfMotion );
        

  2. Create a chord finder that uses your Field and Stepper. You must also give it a minimum step size, below which it does not make sense to attempt to integrate:

          yourChordFinder= new G4ChordFinder( yourField,
                                     yourMininumStep, // say 0.01*mm
                                     yourStepper );
        

  3. Next create a G4FieldManager and give it that chord finder,

          yourFieldManager= new G4FieldManager();
          yourFieldManager.SetChordFinder(yourChordFinder);
        

  4. Finally we tell the Geometry that this FieldManager is responsible for creating a field for the detector.

          G4TransportationManager::GetTransportationManager()
                                -> SetFieldManager( yourFieldManager );
        

Changes for non-electromagnetic fields

If the field you are interested in simulating is not electromagnetic, another minor modification may be required. The transportation currently chooses whether to propagate a particle in a field or rectilinearly based on whether the particle is charged or not. If your field affects non-charged particles, you must inherit from the G4Transportation and re-implement the part of GetAlongStepPhysicalInteractionLength that decides whether the particles is affected by your force.

In particular the relevant section of code does the following:

  //  Does the particle have an (EM) field force exerting upon it?
  //
  if( (particleCharge!=0.0) ){
     
     fieldExertsForce= this->DoesGlobalFieldExist();
     // Future: will/can also check whether current volume's field is Zero or
     //  set by the user (in the logical volume) to be zero.
  }

and you want it to ask whether it feels your force. If, for the sake of an example, you wanted to see the effects of gravity on a heavy hypothetical particle, you could say

  
  //  Does the particle have my field's force exerted on it?
  //
  if (particle->GetName() == "VeryHeavyWIMP") {
     fieldExertsForce= this->DoesGlobalFieldExist();  // For gravity
  }

After doing all these steps, you will be able to see the effects of your force on a particle's motion.