Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
example.cpp File Reference
#include <iostream>
#include <cmath>
#include <TApplication.h>
#include <TCanvas.h>
#include "Garfield/ViewField.hh"
#include "Garfield/ViewDrift.hh"
#include "Garfield/ViewSignal.hh"
#include "Garfield/MediumMagboltz.hh"
#include "Garfield/ComponentAnalyticField.hh"
#include "Garfield/Sensor.hh"
#include "Garfield/TrackHeed.hh"
#include "Garfield/AvalancheMicroscopic.hh"
#include "Garfield/Plotting.hh"

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 19 of file example.cpp.

19 {
20
21 TApplication app("app", &argc, argv);
23
24 // Set the gas mixture.
26 gas.SetComposition("ar", 90., "c4h10", 10.);
27
28 // Build the Micromegas geometry.
29 // We use separate components for drift and amplification regions.
32 cmpD.SetMedium(&gas);
33 cmpA.SetMedium(&gas);
34 // Set the magnetic field [Tesla].
35 constexpr double bfield = 2.;
36 cmpD.SetMagneticField(0, 0, bfield);
37 cmpA.SetMagneticField(0, 0, bfield);
38
39 // 100 um amplification gap.
40 constexpr double yMesh = 0.01;
41 constexpr double vMesh = -500.;
42 cmpA.AddPlaneY(yMesh, vMesh, "mesh");
43 cmpA.AddPlaneY(0., 0., "bottom");
44
45 // 3 mm drift gap.
46 constexpr double yDrift = yMesh + 0.3;
47 constexpr double vDrift = -3000.;
48 cmpD.AddPlaneY(yDrift, vDrift, "top");
49 cmpD.AddPlaneY(yMesh, vMesh, "mesh");
50
51 // We place three strips along z direction on the anode plane
52 // in order to be able to resolve the avalanche signal.
53 // Strip pitch [cm].
54 constexpr double pitch = 0.07;
55 // Distance between two strips [cm].
56 constexpr double interstrip = 0.01;
57 // Half-width of a strip.
58 constexpr double hw = 0.5 * (pitch - interstrip);
59 const double xStrip1 = hw;
60 cmpA.AddStripOnPlaneY('z', 0., xStrip1 - hw, xStrip1 + hw, "strip1");
61 const double xStrip2 = xStrip1 + pitch;
62 cmpA.AddStripOnPlaneY('z', 0., xStrip2 - hw, xStrip2 + hw, "strip2");
63 const double xStrip3 = xStrip2 + pitch;
64 cmpA.AddStripOnPlaneY('z', 0., xStrip3 - hw, xStrip3 + hw, "strip3");
65
66 // We want to calculate the signals induced on the strips.
67 cmpA.AddReadout("strip1");
68 cmpA.AddReadout("strip2");
69 cmpA.AddReadout("strip3");
70
71 // Assemble a sensor.
72 Sensor sensor;
73 sensor.AddComponent(&cmpD);
74 sensor.AddComponent(&cmpA);
75 sensor.SetArea(0., 0., -1., 5 * pitch, yDrift, 1.);
76
77 // Request signal calculation for the strip electrodes.
78 sensor.AddElectrode(&cmpA, "strip1");
79 sensor.AddElectrode(&cmpA, "strip2");
80 sensor.AddElectrode(&cmpA, "strip3");
81
82 // Set the time window [ns] for the signal calculation.
83 const double tMin = 0.;
84 const double tMax = 100.;
85 const double tStep = 0.05;
86 const int nTimeBins = int((tMax - tMin) / tStep);
87 sensor.SetTimeWindow(0., tStep, nTimeBins);
88
89 // We use microscopic tracking for simulating the electron avalanche.
91 aval.SetSensor(&sensor);
92 // Switch on signal calculation.
95
96 // Simulate an ionizing particle (negative pion) using Heed.
97 TrackHeed track;
98 track.SetParticle("pi-");
99 constexpr double momentum = 300.e6; // [eV / c]
100 track.SetMomentum(momentum);
101 track.SetSensor(&sensor);
102 track.EnableMagneticField();
103 // track.EnableElectricField();
104
105 // Construct a viewer to visualise the drift lines.
106 ViewDrift driftView;
107 track.EnablePlotting(&driftView);
108 aval.EnablePlotting(&driftView);
109 aval.EnableExcitationMarkers(false);
110 aval.EnableIonisationMarkers(false);
111
112 // Set the starting point of the incoming ionising track.
113 double xt = 0.05; // [cm]
114 double yt = yDrift;
115 double zt = 0.0;
116
117 // Now simulate a track.
118 track.NewTrack(xt, yt, zt, 0, 0, -1, 0);
119 // Loop over the clusters.
120 double xc, yc, zc, tc, ec, extra;
121 int nc;
122 while (track.GetCluster(xc, yc, zc, tc, nc, ec, extra)) {
123 for (int j = 0; j < nc; ++j) {
124 double xe, ye, ze, te, ee, dxe, dye, dze;
125 track.GetElectron(j, xe, ye, ze, te, ee, dxe, dye, dze);
126 // Simulate the drift/avalanche of this electron.
127 aval.AvalancheElectron(xe, ye, ze, te, 0.1, dxe, dye, dze);
128 }
129 }
130
131 // Calculate the accumulated charge on each strip.
132 sensor.IntegrateSignals();
133 const double q1 = sensor.GetSignal("strip1", nTimeBins - 1);
134 const double q2 = sensor.GetSignal("strip2", nTimeBins - 1);
135 const double q3 = sensor.GetSignal("strip3", nTimeBins - 1);
136
137 // Now calculate the reconstructed pion position in the Micromegas
138 // using the weighted average of the strip signals.
139 const double sum = q1 + q2 + q3;
140 double mean = (q1 * xStrip1 + q2 * xStrip2 + q3 * xStrip3) / sum;
141 const double res = fabs(xt - mean) * 10000.0;
142 std::cout << "---------------------------\n";
143 std::cout << "XMean: " << mean << " cm, XTrue: " << xt << " cm\n";
144 std::cout << "XTrue - XMean: " << res << " um\n";
145 std::cout << "---------------------------\n";
146
147 // Plotting
148 // This canvas will be used to display the drift lines and the field
149 gStyle->SetPadRightMargin(0.15);
150 TCanvas* c1 = new TCanvas("c1", "", 1400, 600);
151 c1->Divide(2, 1);
152
153 // Plot the drift lines.
154 driftView.SetArea(0.0, 0.0, 0.2, yDrift);
155 driftView.SetClusterMarkerSize(0.1);
156 driftView.SetCollisionMarkerSize(0.5);
157 driftView.SetCanvas((TPad*)c1->cd(1));
158 driftView.Plot(true);
159 // Plot the potential.
160 ViewField fieldView;
161 fieldView.SetSensor(&sensor);
162 fieldView.SetCanvas((TPad*)c1->cd(2));
163 fieldView.Plot("v", "cont1z");
164
165 c1->SaveAs("plot.pdf");
166 app.Run(true);
167 return 0;
168
169}
Calculate electron drift lines and avalanches using microscopic tracking.
void EnablePlotting(ViewDrift *view)
Switch on drift line plotting.
void EnableIonisationMarkers(const bool on=true)
Draw a marker at every ionising collision or not.
void SetSensor(Sensor *sensor)
Set the sensor.
void EnableExcitationMarkers(const bool on=true)
Draw a marker at every excitation or not.
void EnableSignalCalculation(const bool on=true)
Switch on calculation of induced currents (default: off).
void EnableMagneticField(const bool on=true)
Enable magnetic field in stepping algorithm (default: off).
bool AvalancheElectron(const double x0, const double y0, const double z0, const double t0, const double e0, const double dx0=0., const double dy0=0., const double dz0=0.)
Calculate an avalanche initiated by a given electron.
void SetMedium(Medium *medium)
Set the medium inside the cell.
void AddPlaneY(const double y, const double voltage, const std::string &label)
Add a plane at constant y.
void AddStripOnPlaneY(const char direction, const double y, const double smin, const double smax, const std::string &label, const double gap=-1.)
Add a strip in the x or z direction on an existing plane at constant y.
void AddReadout(const std::string &label)
Setup the weighting field for a given group of wires or planes.
void SetMagneticField(const double bx, const double by, const double bz)
Set a constant magnetic field.
Definition: Component.cc:94
bool SetComposition(const std::string &gas1, const double f1=1., const std::string &gas2="", const double f2=0., const std::string &gas3="", const double f3=0., const std::string &gas4="", const double f4=0., const std::string &gas5="", const double f5=0., const std::string &gas6="", const double f6=0.)
Set the gas mixture.
Definition: MediumGas.cc:135
void SetDefaultStyle()
Apply the default Garfield ROOT style.
void SetTimeWindow(const double tstart, const double tstep, const unsigned int nsteps)
Definition: Sensor.cc:820
void AddElectrode(Component *comp, const std::string &label)
Add an electrode.
Definition: Sensor.cc:349
void AddComponent(Component *comp)
Add a component.
Definition: Sensor.cc:316
bool SetArea()
Set the user area to the default.
Definition: Sensor.cc:183
double GetSignal(const std::string &label, const unsigned int bin)
Retrieve the total signal for a given electrode and time bin.
Definition: Sensor.cc:912
bool IntegrateSignals()
Replace the signals on all electrodes curve by their integrals.
Definition: Sensor.cc:1234
Generate tracks using Heed++.
Definition: TrackHeed.hh:35
void EnableMagneticField()
Take the magnetic field into account in the stepping algorithm.
Definition: TrackHeed.cc:678
bool GetCluster(double &xcls, double &ycls, double &zcls, double &tcls, int &n, double &e, double &extra) override
Definition: TrackHeed.cc:222
bool NewTrack(const double x0, const double y0, const double z0, const double t0, const double dx0, const double dy0, const double dz0) override
Definition: TrackHeed.cc:59
bool GetElectron(const unsigned int i, double &x, double &y, double &z, double &t, double &e, double &dx, double &dy, double &dz)
Definition: TrackHeed.cc:354
void EnablePlotting(ViewDrift *viewer)
Switch on plotting.
Definition: Track.cc:175
void SetSensor(Sensor *s)
Set the sensor through which to transport the particle.
Definition: Track.cc:166
void SetMomentum(const double p)
Set the particle momentum.
Definition: Track.cc:140
virtual void SetParticle(const std::string &part)
Definition: Track.cc:15
void SetCanvas(TPad *pad)
Set the canvas to be painted on.
Definition: ViewBase.hh:28
void SetArea(const double xmin, const double ymin, const double xmax, const double ymax)
Definition: ViewBase.cc:109
Visualize drift lines and tracks.
Definition: ViewDrift.hh:18
void Plot(const bool twod=false, const bool axis=true)
Draw the drift lines.
Definition: ViewDrift.cc:185
void SetClusterMarkerSize(const double size)
Set the size of the cluster markers (see TAttMarker).
Definition: ViewDrift.cc:38
void SetCollisionMarkerSize(const double size)
Set the size of the collision markers (see TAttMarker).
Definition: ViewDrift.cc:46
Visualize the potential or electric field of a component or sensor.
Definition: ViewField.hh:15
void Plot(const std::string &option="v", const std::string &drawopt="arr")
Definition: ViewField.cc:126
void SetSensor(Sensor *s)
Set the sensor for which to plot the field.
Definition: ViewField.cc:72
PlottingEngine plottingEngine
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615