Skip to content

Gas files

As a prerequisite to simulating electron drift lines in gases using AvalancheMC or DriftLineRKF, one needs to know the relevant macroscopic transport parameters (drift velocity, diffusion coefficients etc.) as functions of the electric and magnetic field. The Magboltz program computes these properties using a Monte Carlo simulation based on the microscopic electron-atom/molecule scattering cross-sections. As the calculation can be quite time consuming, it is common practice to run Magboltz for a set of electric and magnetic fields, and save a table of the resulting transport properties to a file.

Generating a gas table using Magboltz

There are two options for specifying the set of electric fields, magnetic fields, and angles between E and B to be included in the table.

One option is to set the lower and upper limits of the field and angular ranges and the number of points.

MediumMagboltz gas;
const size_t nE = 20;
const double emin =    100.;
const double emax = 100000.;
// Flag to request logarithmic spacing.
constexpr bool useLog = true;
const size_t nB = 3;
// Range of magnetic fields [Tesla]
const double bmin = 0.;
const double bmax = 2.;
const size_t nA = 4;
// Range of angles [rad]
const double amin = 0.;
const double amax = HalfPi;
gas.SetFieldGrid(emin, emax, nE, useLog, bmin, bmax, nB, amin, amax, nA);

For the electric field, one has the choice between linear and logarithmic spacing of the grid points. Magnetic fields and angles will be evenly spaced. The parameters for magnetic field and angle are optional. The default is one magnetic field (B = 0).

Alternatively, one can set the list of electric fields, magnetic fields, and angles explicitly.

MediumMagboltz gas;
std::vector<double> efields = {100., 200., 500., 1000., 2000., 5000., 10000.};
std::vector<double> bfields = {0., 1., 2.};
std::vector<double> angles = {0., 30., 60., 90.};
gas.SetFieldGrid(efields, bfields, angles);

The field values have to be strictly positive and in increasing order. The angles are in rad and have to be in the range [0, π/2]. Make sure that the all fields present in the detector are covered by the table. Try to be reasonable when selecting the number of electric fields at which you run Magboltz. On the one hand the number of points must be large enough so as to represent any fine structure that may be present in the transport properties. On the other hand, an excessively large number of points may lead to oscillations in the interpolation.

Temperature and pressure should also be specified. The default values are 293.15 K and 760 Torr.

MediumMagboltz gas;
gas.SetComposition("Ar", 70., "CO2", 30.);
// Set the temperature [K].
gas.SetTemperature(293.15);
// Set the pressure [Torr].
gas.SetPressure(3 * 760.)`

We are now ready to run Magboltz.

const int ncoll = 10;
// Run Magboltz to generate the gas table.
gas.GenerateGasTable(ncoll);

The function GenerateGasTable loops over all combinations of E, B and angle, invokes Magboltz and retrieves the drift velocity vector, the longitudinal and transverse diffusion coefficients, the Townsend and attachment coefficients as well as the excitation and ionisation rates for the various gas molecules.

The argument of the function GenerateGasTable specifies the number of collisions, in multiples of 107, to be used to compute the transport properties for each (E, B, angle) combination. Low values (e.g. 2 - 5) are suitable for mixtures than contain substantial amounts of quenchers. Pure noble gases call for larger values (at least 10). The statistical accuracy of the drift velocity calculation improves with the square root of this parameter, while the CPU time consumption increases linearly.

Exporting and importing a gas table

After having generated a gas table using Magboltz, we can write it to file.

gas.WriteGasFile("ar_70_co2_30.gas");

In addition to the electron transport property tables, the gas file also contains the pressure and the temperature. For compatibility with classic (Fortran) Garfield, it also includes the ion mobility and fields related to primary ionisation (cluster size and density, Heed initialisation information etc.). These are not used by Garfield++ however. The gas files use a compact format and are not intended to serve as easily legible tables, nor should they be edited ''by hand''.

The dataset stored in a gas file can be imported using

MediumMagboltz gas;
gas.LoadGasFile("ar_70_co2_30.gas");

This will set the gas composition, pressure and temperature to the ones given in the gas file. It is technically possible to change the pressure after loading a gas file. The transport properties will then be scaled according to simple scaling laws. It is not recommended, however, to rely on these scaling laws since these are very approximate.

The method PrintGas prints out all available information on the gas.

gas.PrintGas();

Merging gas files

Suppose we have two gas files for the same gas mixture, one for zero magnetic field and one for B = 2 T. Using the method MergeGasFile we can merge these two datasets.

// Load the first gas table.
gas.LoadGasFile("ar_70_co2_30_0T.gas");
// In case of overlaps, use the values from the file being added.
constexpr bool replaceOld = true;
// Add the second gas table.
gas.MergeGasFile("ar_70_co2_30_2T.gas", replaceOld);
// Save the merged table.
gas.WriteGasFile("ar_70_co2_30_merged.gas");

This only works if the data in the two files have at least two of the following three axes in common: the set of electric field values, the set of angles between E and B and the set of magnetic field values. This condition is trivially met if they both concern a zero magnetic field situation. If precisely two of the three axes coincide, then the merged data will contain only the transport properties which are found in both datasets. If all three axes coincide, then the merged data will contain all the transport properties that are found in either the current data or the file to be merged in.

If a transport parameter for a given (E, B, angle) combination is found in both datasets, either the value from the existing ("old") dataset will be kept (if the flag replaceOld is set to false) or the value from the file to be merged in will be used (replaceOld set to true).

Plotting transport properties

The transport properties stored in the gas table can be plotted using the class ViewMedium.

MediumMagboltz gas;
gas.LoadGasFile("ar_70_co2_30.gas");

ViewMedium view(&gas);
view.PlotElectronVelocity('e');

Similar functions exist for the other transport parameters. The first argument indicates the variable to be used for the x-axis of the plot. Valid values are 'e' (electric field), 'b' (magnetic field), and 'a' (angle). By default, the plot is made for B = 0. To plot e.g. the drift velocity at B = 2 T, one would need to replace the above lines of code by

ViewMedium view(&gas);
view.SetMagneticField(2.);
view.PlotElectronVelocity('e');

If the magnetic field is non-zero, the drift velocity plot will include three curves: the component parallel to the electric field, the component parallel to E×B, and the component parallel to the orthogonal part of B.

It is possible to superimpose curves from different datasets on the same plot. An example is given in the snippet below.

MediumMagboltz gas1;
gas1.LoadGasFile("ar_80_co2_20.gas");
MediumMagboltz gas2;
gas2.LoadGasFile("ar_70_co2_30.gas");
ViewMedium view;
view.SetMedium(&gas1);
view.PlotElectronVelocity('e');
std::vector<std::string> labels = {"Ar/CO_{2} 80/20", "Ar/CO_{2} 70/30"};
view.SetLabels(labels);
view.SetMedium(&gas2);
view.PlotElectronVelocity('e', true);

One can of course also retrieve the transport properties directly from MediumMagboltz (without using ViewMedium) and plot them using another tool. The function ElectronVelocity returns the drift velocity at a given electric and magnetic field.

// Get the drift velocity at E = 10 kV/cm, B = 1 T,  
// with E and B at an angle of 30 degree.
double emag = 10.e3;
double bmag = 1.;
double theta = 30. * DegreeToRad;
double vx, vy, vz;
gas.ElectronVelocity(emag, 0, 0, bmag * cos(theta), bmag * sin(theta), 0,
                     vx, vy, vz);
// Convert from cm/ns to cm/us.
vx *= 1.e3;
vy *= 1.e3;
vz *= 1.e3;
std::printf("Velocity along E:   %10.3f cm/us\n", vx);
std::printf("Velocity along Bt:  %10.3f cm/us\n", vy);
std::printf("Velocity along ExB: %10.3f cm/us\n", vz); 

This function will interpolate between the drift velocity values at the (E, B, angle) points in the table. Similar functions are available for the other transport parameters.

One can also retrieve directly the values at the (E, B, angle) points in the table, as illustrated in the following snippet of code.

std::vector<double> efields;
std::vector<double> bfields;
std::vector<double> angles;
gas.GetFieldGrid(efields, bfields, angles);
const auto nE = efields.size();
const auto nB = bfields.size();
const auto nA = angles.size();
for (size_t j = 0; j < nB; ++j) {
  for (size_t k = 0; k < nA; ++k) {
    std::cout << "B = " << bfields[j] << " T, theta = "
              << angles[k] * RadToDegree << " degree\n";
    std::cout << "   E [V/cm]     vE [cm/us]    alpha [1/cm]\n";
    for (size_t i = 0; i < nE; ++i) {
       double ve = 0.;
       gas.GetElectronVelocityE(i, j, k, ve);
       // Convert from cm/ns to cm/us.
       ve *= 1.e3;
       double alpha = 0.;
       gas.GetElectronTownsend(i, j, k, alpha);
       alpha = exp(alpha);
       std::printf("%10.3f    %10.3f    %10.3f\n", efields[i], ve, alpha);
    }
  }
}

For the Townsend coefficient α and the attachment coefficient η, the logarithms lnα, lnη and not the actual values are stored in the tables, hence the line alpha = exp(alpha); in the above example.

Examples

Some simple example applications are available in the directory Examples/GasFile of the Garfield++ source tree.