Calculation of induced currents
Weighting potentials/fields
The steps for preparing the weighting potential and/or weighting field depend on the type of Component
.
Analytic fields
When setting up the cell, each element (wire, plane, tube) is given a label, e.g. in case of a plane:
ComponentAnalyticField cmp;
cmp.AddPlaneY(0., 0., "readout");
If this label is a non-empty string, then the class will prepare the weighting potential and weighting field for this element. If elements share the same label, they are treated as one electrode, i.e. the calculated signal will be the sum of the signals induced on the individual elements.
ComponentAnalyticField
can calculate weighting fields and potentials for planes, wires, strips, and pixels.
Finite-element field maps
- Compute a weighting field map using the same mesh as your main field map. See the Ansys123 tutorial for a recipe.
-
After initialising your field map class (
ComponentAnsys123
),ComponentAnsys123 cmp; cmp.Initialise("EFIELD.lis", "NLIST.lis", "PRNSOL.lis", "MPLIST.lis");
load the weighting field map using
cmp.SetWeightingField("WSOL.lis", "readout");
where "WSOL.lis" is the name of the Ansys file (you can name it as you like, of course) that contains the nodal solutions for the weighting field configuration, and "readout" is a label you give to the electrode.
The same recipe applies to Elmer, Comsol, and CST field maps.
TCAD field maps
ComponentTcad2d
and ComponentTcad3d
require two field maps (.dat
files) to define the weighting potential and field. The first one corresponds to the "drift" field. The second one corresponds to the solution for a configuration where the potential on the electrode to be read out is increased by a small voltage ΔV. The two solutions have to be calculated using the same mesh.
After importing the mesh and drift electric field,
ComponentTcad2d cmp;
cmp.Initialise("sensor_des.grd", "sensor_des.dat");
set up the weighting field using
cmp.SetWeightingField("sensor_des.dat", "sensor1_des.dat", 1., "readout");
The third parameter in SetWeightingField
corresponds to the extra voltage ΔV (1 V in this example) applied to the electrode to be read out in the second field map. The last argument is the identifier to be given to the electrode.
neBEM
When defining the detector geometry, we need to assign a "label" to the solids for which we would like to calculate the weighting field. As an example, let us consider a simple parallel-plate geometry.
SolidBox box1(0, 0, -0.5, 5, 5, 0);
SolidBox box2(0, 0, 0.5, 5, 5, 0);
box1.SetBoundaryPotential(0.);
box2.SetBoundaryPotential(1000.);
box1.SetLabel("readout");
neBEM will then compute the weighting potential and field corresponding to a given label by setting potential of all elements associated to solids with this label to 1 V (and grounding all other conducting elements).
Sensor
- Tell the
Sensor
that you want to calculate the signal on the electrode labelled "readout" using the weighting field with the same name.sensor.AddElectrode(&cmp, "readout");
-
Setup the binning for the signal calculation, e.g.
sensor.SetTimeWindow(0., 1., 100);
The first parameter is the lower time limit, the second one is the bin width (in ns), the third one is the number of bins. The above setting would thus request the signal to be recorded from 0 to 100 ns in steps of 1 ns. You'll probably have to adapt it.
After that the Sensor
accumulates the signals of all drift lines which are calculated. The induced current can be calculated either using the weighting field or the weighting potential, with the latter being the default. To use the weighting field instead, add a line
aval.UseWeightingPotential(false);
To tell the sensor that the signal from a particular particle track, photon conversion etc. is finished and that the drift lines calculated from now on belong to a new track call
sensor.ClearSignal();
In this case the previous signal is deleted. Alternatively, you can call
sensor.NewSignal();
To retrieve the signal, you can call
sensor.GetSignal("readout", i);
i
is the index of the time bin. The value is in fC / ns.
To make a plot, you can use the class ViewSignal
:
ViewSignal signalView;
signalView.SetSensor(&sensor);
signalView.PlotSignal("readout");