Garfield++ 5.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
wire2d.cc
Go to the documentation of this file.
1/**
2 * wire2d.cc
3 *
4 * Demonstrates importing of a 2D Elmer finite element
5 * field map.
6 *
7*/
8#include <iostream>
9#include <cmath>
10
11#include <TCanvas.h>
12#include <TApplication.h>
13#include <TFile.h>
14
17#include "Garfield/Sensor.hh"
18#include "Garfield/ViewField.hh"
21#include "Garfield/Random.hh"
23
24using namespace Garfield;
25
26int main(int argc, char* argv[]) {
27
28 TApplication app("app", &argc, argv);
29
30 // Set relevant parameters.
31 // Wire radius
32 const double rwire = 1.0;
33 // X-width of drift simulation will cover between +/- axis_x
34 const double axis_x = 5;
35 // Y-width of drift simulation will cover between +/- axis_y
36 const double axis_y = 5;
37 const double axis_z = 5;
38
39 // Define the medium (Ar/CO2 70:30).
40 MediumMagboltz gas("ar", 70., "co2", 30.);
41 // Set the temperature (K).
42 gas.SetTemperature(293.15);
43 // Set the pressure (Torr).
44 gas.SetPressure(740.);
45
46 // Import an Elmer-created field map.
48 "wire2d/mesh.header", "wire2d/mesh.elements", "wire2d/mesh.nodes",
49 "wire2d/dielectrics.dat", "wire2d/wire2d.result", "cm");
50 elm.SetGas(&gas);
51 elm.SetRangeZ(-5.,5.);
52
53 // Set up a sensor object.
54 Sensor sensor;
55 sensor.AddComponent(&elm);
56 sensor.SetArea(-axis_x, -axis_y, -axis_z, axis_x, axis_y, axis_z);
57
58 // Create an avalanche object
59 AvalancheMicroscopic aval(&sensor);
60
61 // Set up the object for drift line visualization.
62 ViewDrift viewDrift;
63 viewDrift.SetArea(-axis_x, -axis_y, -axis_z, axis_x, axis_y, axis_z);
64 aval.EnablePlotting(&viewDrift, 100);
65
66 // Set the electron start parameters.
67 const double zi = 1.0;
68 double ri = rwire + 2.0;
69 double thetai = RndmUniform() * TwoPi;
70 double xi = ri * cos(thetai);
71 double yi = ri * sin(thetai);
72 // Calculate the avalanche.
73 std::cout << "Avalanche of a single electron starting from ("
74 << xi << ", " << yi << ", " << zi << ")..." << std::endl;
75 aval.AvalancheElectron(xi, yi, zi, 0., 0., 0., 0., 0.);
76
77 std::cout << "... avalanche complete with "
78 << aval.GetNumberOfElectronEndpoints() << " electron tracks.\n";
79
80 // Plot the geometry, field and drift lines.
81 TCanvas* cGeom = new TCanvas("geom", "Geometry/Avalanche/Fields");
82 cGeom->SetLeftMargin(0.14);
83 const bool plotContours = false;
84 if (plotContours) {
85 ViewField* vf = new ViewField();
86 vf->SetSensor(&sensor);
87 vf->SetCanvas(cGeom);
88 vf->SetArea(-axis_x, -axis_y, axis_x, axis_y);
89 vf->SetNumberOfContours(40);
90 vf->SetNumberOfSamples2d(30, 30);
91 vf->SetPlane(0, 0, 1, 0, 0, 0);
92 vf->PlotContour("v");
93 }
94
95 // Set up the object for FE mesh visualization.
96 ViewFEMesh vFE;
97 vFE.SetArea(-axis_x, -axis_z, -axis_y, axis_x, axis_z, axis_y);
98 vFE.SetCanvas(cGeom);
99 vFE.SetComponent(&elm);
100 vFE.SetPlane(0, 0, 1, 0, 0, 0);
101 vFE.SetFillMesh(true);
102 vFE.SetColor(1, kGray);
103 if (!plotContours) {
104 vFE.EnableAxes();
105 vFE.SetXaxisTitle("x (cm)");
106 vFE.SetYaxisTitle("y (cm)");
107 }
108 vFE.SetViewDrift(&viewDrift);
109 vFE.Plot();
110 app.Run();
111 return 0;
112}
Calculate electron drift lines and avalanches using microscopic tracking.
bool AvalancheElectron(const double x, const double y, const double z, const double t, const double e, const double dx=0., const double dy=0., const double dz=0.)
Calculate an avalanche initiated by a given electron.
void EnablePlotting(ViewDrift *view, const size_t nColl=100)
Switch on drift line plotting.
Component for importing two-dimensional field maps computed by Elmer.
void SetRangeZ(const double zmin, const double zmax)
void SetTemperature(const double t)
Set the temperature [K].
Definition Medium.cc:72
void SetPressure(const double p)
Definition Medium.cc:82
void AddComponent(Component *comp)
Add a component.
Definition Sensor.cc:355
bool SetArea(const bool verbose=false)
Set the user area to the default.
Definition Sensor.cc:187
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
virtual void SetPlane(const double fx, const double fy, const double fz, const double x0, const double y0, const double z0)
Definition ViewBase.cc:142
Visualize drift lines and tracks.
Definition ViewDrift.hh:19
Draw the mesh of a field-map component.
Definition ViewFEMesh.hh:22
bool Plot(const bool twod=true)
Plot method to be called by user.
Definition ViewFEMesh.cc:60
void SetYaxisTitle(const std::string &ytitle)
Definition ViewFEMesh.hh:42
void SetViewDrift(ViewDrift *vd)
Set the optional associated ViewDrift.
Definition ViewFEMesh.hh:64
void SetXaxisTitle(const std::string &xtitle)
Definition ViewFEMesh.hh:41
void SetPlane(const double fx, const double fy, const double fz, const double x0, const double y0, const double z0) override
void SetComponent(ComponentFieldMap *cmp)
Set the component from which to retrieve the mesh and field.
Definition ViewFEMesh.cc:49
void SetColor(int matID, int colorID)
Definition ViewFEMesh.hh:58
void SetFillMesh(const bool f)
Element fill switch; 2D only, set false for wireframe mesh.
Definition ViewFEMesh.hh:50
Visualize the potential or electric field of a component or sensor.
Definition ViewField.hh:15
void SetNumberOfSamples2d(const unsigned int nx, const unsigned int ny)
Set the number of points used for drawing 2D functions.
Definition ViewField.cc:123
void PlotContour(const std::string &option="v")
Definition ViewField.cc:129
void SetNumberOfContours(const unsigned int n)
Set the number of contour levels.
Definition ViewField.cc:115
void SetSensor(Sensor *s)
Set the sensor for which to plot the field.
Definition ViewField.cc:73
double RndmUniform()
Draw a random number uniformly distributed in the range [0, 1).
Definition Random.hh:14
int main(int argc, char *argv[])
Definition wire2d.cc:26