Skip to content

Modelling a drift tube

In this example, we consider a drift tube with an outer diameter of 15 mm and a wire diameter of 50 μm, similar to the ATLAS small-diameter muon drift tubes (sMDTs). We will calculate the signal induced on the wire by the passage of a charged particle.

Gas table

As a first step, we have to prepare a table of transport parameters (drift velocity, diffusion coefficients, Townsend and attachment coefficients) as a function of the electric field (and, in general, also the magnetic field as well as the angle between the two fields). In this example we will use a gas mixture of 93% argon and 7% CO2 at a pressure of 3 atm and room temperature.

MediumMagboltz gas("ar", 93., "co2", 7.);
// Set temperature [K] and pressure [Torr].
gas.SetPressure(3 * 760.);
gas.SetTemperature(293.15);

We also have to specify the number of electric fields to be included in the table and the electric field range to be covered. Here we use 20 field points between 100 V/cm and 100 kV/cm with logarithmic spacing.

gas.SetFieldGrid(100., 100000., 20, true);

Now we run Magboltz to generate the gas table for this grid. As an input parameter we have to specify the number of collisions (in multiples of 107) over which the electron is traced by Magboltz.

gas.GenerateGasTable(10);

The calculation will take a while, don't panic. After the calculation is done, we save the gas table to a file for later use.

gas.WriteGasFile("ar_93_co2_7.gas");

Once we have saved the transport parameters to file we can skip the steps above, and simply import the table in our program using

gas.LoadGasFile("ar_93_co2_7.gas");

In order to make sure the calculation of the gas table was successful, it is a good idea to plot, for instance, the drift velocity as a function of the electric field.

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

Drift velocity

Elecric field

For calculating the electric field inside the tube, we use the class ComponentAnalyticField which can handle (two-dimensional) arrangements of wires, planes and tubes. First we set the medium associated with the active region of the tube.

ComponentAnalyticField cmp;
cmp.SetMedium(&gas);

Then we add the elements defining the electric field: a wire (which we label "s") and a tube.

// Wire radius [cm]
const double rWire = 25.e-4;
// Outer radius of the tube [cm]
const double rTube = 0.71;
// Voltages
const double vWire = 2730.;
const double vTube =    0.;
// Add the wire in the centre.
cmp.AddWire(0., 0., 2 * rWire, vWire, "s");
// Add the tube.
cmp.AddTube(rTube, vTube, 0);

Finally we assemble a Sensor object which acts as an interface to the transport classes discussed below.

Sensor sensor(&cmp);

Drift lines from a charged-particle track

We use Heed (see this tutorial for a more detailed description) to simulate the ionization produced by a charged particle crossing the tube (a 170 GeV muon in our example).

TrackHeed track(&sensor);
track.SetParticle("muon");
track.SetEnergy(170.e9);

The drift lines of the electrons created along the track are calculated using Runge-Kutta-Fehlberg (RKF) integration, implemented in the class DriftLineRKF. This method uses the previously computed tables of transport parameters to calculate drift lines and multiplication.

DriftLineRKF drift(&sensor);

If we want to simulate the ion tail (i. e. the drift lines of the ions that are produced in the avalanche), we need to import a table of ion mobilities

gas.LoadIonMobility("IonMobility_Ar+_Ar.txt");

and call

drift.EnableIonTail();

before simulating the drift lines.

Let us consider in this example a track that passes at a distance of 3 mm from the wire centre. After simulating the passage of the charged particle, we loop over the "clusters" (i. e. the ionizing collisions of the primary particle) along the track and calculate a drift line for each electron produced in the cluster.

const double rTrack = 0.3;
const double x0 = rTrack;
const double y0 = -sqrt(rTube * rTube - rTrack * rTrack);
track.NewTrack(x0, y0, 0, 0, 0, 1, 0);
// Loop over the clusters along the track.
for (const auto& cluster : track.GetClusters()) {
  // Loop over the electrons in the cluster.
  for (const auto& electron : cluster.electrons) {
    drift.DriftElectron(electron.x, electron.y, electron.z, electron.t);
  }
}

As a check whether the simulation is doing something sensible, it can be useful to visualize the drift lines. Before simulating the charged particle track and the electron drift lines, we have to instruct TrackHeed and DriftLineRKF to pass the coordinates of the clusters and the points along the drift line to a ViewDrift object which then takes care of plotting them.

// Create a canvas.
TCanvas* cD = new TCanvas("cD", "", 600, 600);
ViewDrift driftView;
driftView.SetCanvas(cD);
drift.EnablePlotting(&driftView);
track.EnablePlotting(&driftView);

We use the class ViewCell to draw the outline of the tube and the position of the wire on the same plot as the drift lines.

ViewCell cellView;
cellView.SetCanvas(cD);
cellView.SetComponent(&cmp);

After we've simulated all drift lines from a charged particle track, we create a plot using

cellView.Plot2d();
constexpr bool twod = true;
constexpr bool drawaxis = false;
driftView.Plot(twod, drawaxis);

Drift lines

Signals

As discussed in more detail in this tutorial, we need to instruct the Sensor object to use ComponentAnalyticField for calculating the weighting field of the wire (which we had labelled "s").

// Request signal calculation for the electrode named "s",
// using the weighting field provided by the Component object cmp.
sensor.AddElectrode(&cmp, "s");

In addition, we have to set the time interval within which the signal is recorded and the granularity (bin width). In this example, we use 1000 bins with a width of 0.5 ns.

// Set the signal time window.
const double tstep = 0.5;
const double tmin = -0.5 * tstep;
const unsigned int nbins = 1000;
sensor.SetTimeWindow(tmin, tstep, nbins);

Using the class ViewSignal, we plot the current induced on the wire by the drifting electrons.

ViewSignal signalView(&sensor);
signalView.PlotSignal("s");

To calculate the signal at the output of the front-end electronics, we convolute the induced current pulse with the so-called transfer function (or delta-response function) of the front-end circuit. In this example, we use a table read in from a text file. We store the values of the transfer function and the corresponding times (in ns) in separate vectors, and pass them to our Sensor object.

std::vector times;
std::vector values;
// Read data from file.
// ...
sensor.SetTransferFunction(times, values);

To perform the convolution after we've calculated all drift lines from a track, we add a line

sensor.ConvoluteSignals();

to our program.

Source files

Example applications for generating and importing gas files can be found in the folder Examples/GasFile of the source tree. An example for simulating electron drift lines from a charged particle track and the corresponding induced signal is available in the folder Examples/DriftTube.