Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
diode.cpp File Reference
#include <iostream>
#include <fstream>
#include <cmath>
#include <TCanvas.h>
#include <TROOT.h>
#include <TApplication.h>
#include <TSystem.h>
#include "Garfield/MediumSilicon.hh"
#include "Garfield/ComponentTcad3d.hh"
#include "Garfield/Sensor.hh"
#include "Garfield/TrackHeed.hh"
#include "Garfield/AvalancheMC.hh"
#include "Garfield/ViewDrift.hh"
#include "Garfield/ViewField.hh"
#include "Garfield/FundamentalConstants.hh"
#include "Garfield/Random.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 22 of file diode.cpp.

22 {
23
24 TApplication app("app", &argc, argv);
25
26 // Create a drift medium.
28 // Set the temperature [K].
29 si.SetTemperature(266.15);
30
31 // Make a component with three-dimensional TCAD field map.
33 // Load the mesh (.grd file) and electric field profile (.dat).
34 cmp.Initialise("diode_msh.grd", "diode_des.dat");
35
36 // Associate the TCAD regions with the silicon medium.
37 auto nRegions = cmp.GetNumberOfRegions();
38 for (size_t i = 0; i < nRegions; ++i) {
39 cmp.SetMedium(i, &si);
40 }
41
42 // Make a sensor.
43 Sensor sensor;
44 sensor.AddComponent(&cmp);
45 sensor.SetArea();
46
47 // Plot the electrostatic potential.
48 ViewField fieldView;
49 fieldView.SetSensor(&sensor);
50 fieldView.SetPlaneXZ();
51 fieldView.PlotContour();
52
53 // Use MC transport for drifting the electrons and holes.
54 AvalancheMC drift;
55 drift.SetDistanceSteps(1.e-4);
56 drift.SetSensor(&sensor);
57
58 // Visualize the drift lines.
59 ViewDrift driftView;
60 drift.EnablePlotting(&driftView);
61
62 TrackHeed track;
63 track.SetSensor(&sensor);
64 track.SetParticle("pi");
65 track.SetMomentum(3.e9);
66
67 // Sample the inclination and starting point of the track.
68 const double theta = (RndmUniform() - 0.5) * 20. * DegreeToRad;
69 const double x0 = (RndmUniform() - 0.5) * 40.e-4;
70 const double y0 = (RndmUniform() - 0.5) * 40.e-4;
71 track.NewTrack(x0, y0, 0., 0., sin(theta), 0., cos(theta));
72 // Retrieve the clusters along the track.
73 double xc, yc, zc, tc, ec, dummy;
74 int ne;
75 while (track.GetCluster(xc, yc, zc, tc, ne, ec, dummy)) {
76 // Retrieve the electrons in the cluster.
77 for (int i = 0; i < ne; ++i) {
78 double xe, ye, ze, te, ee, dxe, dye, dze;
79 track.GetElectron(i, xe, ye, ze, te, ee, dxe, dye, dze);
80 // Simulate and plot only a small fraction of the drift lines.
81 constexpr double fPlot = 0.01;
82 if (RndmUniform() > fPlot) continue;
83 drift.DriftElectron(xe, ye, ze, te);
84 drift.DriftHole(xe, ye, ze, te);
85 }
86 }
87 driftView.SetPlaneXZ();
88 driftView.Plot2d(true);
89 app.Run(true);
90}
void EnablePlotting(ViewDrift *view)
Switch on drift line plotting.
Definition: AvalancheMC.cc:38
bool DriftElectron(const double x0, const double y0, const double z0, const double t0)
Simulate the drift line of an electron from a given starting point.
Definition: AvalancheMC.cc:186
void SetSensor(Sensor *s)
Set the sensor.
Definition: AvalancheMC.cc:29
void SetDistanceSteps(const double d=0.001)
Use fixed distance steps (default 10 um).
Definition: AvalancheMC.cc:62
bool DriftHole(const double x0, const double y0, const double z0, const double t0)
Simulate the drift line of a hole from a given starting point.
Definition: AvalancheMC.cc:205
Interpolation in a three-dimensional field map created by Sentaurus Device.
size_t GetNumberOfRegions() const
Get the number of regions in the device.
bool Initialise(const std::string &gridfilename, const std::string &datafilename)
void SetMedium(const size_t ireg, Medium *m)
Set the medium to be associated to a given region.
Solid crystalline silicon
void SetTemperature(const double t)
Set the temperature [K].
Definition: Medium.cc:71
void AddComponent(Component *comp)
Add a component.
Definition: Sensor.cc:316
bool SetArea()
Set the user area to the default.
Definition: Sensor.cc:183
Generate tracks using Heed++.
Definition: TrackHeed.hh:35
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 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 SetPlaneXZ()
Set the viewing plane to x-z.
Definition: ViewBase.cc:278
Visualize drift lines and tracks.
Definition: ViewDrift.hh:18
void Plot2d(const bool axis)
Make a 2D plot of the drift lines in the current viewing plane.
Definition: ViewDrift.cc:193
Visualize the potential or electric field of a component or sensor.
Definition: ViewField.hh:15
void PlotContour(const std::string &option="v")
Definition: ViewField.cc:122
void SetSensor(Sensor *s)
Set the sensor for which to plot the field.
Definition: ViewField.cc:72
double RndmUniform()
Draw a random number uniformly distributed in the range [0, 1).
Definition: Random.hh:14
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:432
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:384