Garfield++ 4.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"
19#include "Garfield/Plotting.hh"
22#include "Garfield/Random.hh"
24
25using namespace Garfield;
26
27int main(int argc, char* argv[]) {
28
29 TApplication app("app", &argc, argv);
30
31 // Set relevant parameters.
32 // Wire radius
33 const double rwire = 1.0;
34 // X-width of drift simulation will cover between +/- axis_x
35 const double axis_x = 5;
36 // Y-width of drift simulation will cover between +/- axis_y
37 const double axis_y = 5;
38 const double axis_z = 5;
39
40
41 // Define the medium.
42 MediumMagboltz* gas = new MediumMagboltz();
43 // Set the temperature (K)
44 gas->SetTemperature(293.15);
45 // Set the pressure (Torr)
46 gas->SetPressure(740.);
47 // Allow for drifting in this medium
48 gas->EnableDrift();
49 // Specify the gas mixture (Ar/CO2 70:30)
50 gas->SetComposition("ar", 70., "co2", 30.);
51
52 // Import an Elmer-created field map.
54 "wire2d/mesh.header", "wire2d/mesh.elements", "wire2d/mesh.nodes",
55 "wire2d/dielectrics.dat", "wire2d/wire2d.result", "cm");
56 elm->SetMedium(0, gas);
57 elm->SetRangeZ(-5.,5.);
58
59 // Set up a sensor object.
60 Sensor* sensor = new Sensor();
61 sensor->AddComponent(elm);
62 sensor->SetArea(-axis_x, -axis_y, -axis_z, axis_x, axis_y, axis_z);
63
64 // Create an avalanche object
66 aval->SetSensor(sensor);
67 aval->SetCollisionSteps(100);
68
69 // Set up the object for drift line visualization.
70 ViewDrift* viewDrift = new ViewDrift();
71 viewDrift->SetArea(-axis_x, -axis_y, -axis_z, axis_x, axis_y, axis_z);
72 aval->EnablePlotting(viewDrift);
73
74 // Set the electron start parameters.
75 const double zi = 1.0;
76 double ri = rwire + 2.0;
77 double thetai = RndmUniform() * TwoPi;
78 double xi = ri * cos(thetai);
79 double yi = ri * sin(thetai);
80 // Calculate the avalanche.
81 std::cout << "Avalanche of a single electron starting from (" << xi << ", " << yi << ", " << zi << ")..." << std::endl;
82 aval->AvalancheElectron(xi, yi, zi, 0., 0., 0., 0., 0.);
83
84 std::cout << "... avalanche complete with "
85 << aval->GetNumberOfElectronEndpoints() << " electron tracks.\n";
86
87 // Plot the geometry, field and drift lines.
88 TCanvas* cGeom = new TCanvas("geom", "Geometry/Avalanche/Fields");
89 cGeom->SetLeftMargin(0.14);
90 const bool plotContours = false;
91 if (plotContours) {
92 ViewField* vf = new ViewField();
93 vf->SetSensor(sensor);
94 vf->SetCanvas(cGeom);
95 vf->SetArea(-axis_x, -axis_y, axis_x, axis_y);
96 vf->SetNumberOfContours(40);
97 vf->SetNumberOfSamples2d(30, 30);
98 vf->SetPlane(0, 0, 1, 0, 0, 0);
99 vf->PlotContour("v");
100 }
101
102 // Set up the object for FE mesh visualization.
103 ViewFEMesh* vFE = new ViewFEMesh();
104 vFE->SetArea(-axis_x, -axis_z, -axis_y, axis_x, axis_z, axis_y);
105 vFE->SetCanvas(cGeom);
106 vFE->SetComponent(elm);
107 vFE->SetPlane(0, 0, 1, 0, 0, 0);
108 vFE->SetFillMesh(true);
109 vFE->SetColor(1, kGray);
110 if (!plotContours) {
111 vFE->EnableAxes();
112 vFE->SetXaxisTitle("x (cm)");
113 vFE->SetYaxisTitle("y (cm)");
114 vFE->SetViewDrift(viewDrift);
115 vFE->Plot();
116 }
117
118 app.Run(kTRUE);
119
120 return 0;
121}
Calculate electron drift lines and avalanches using microscopic tracking.
void EnablePlotting(ViewDrift *view)
Switch on drift line plotting.
void SetCollisionSteps(const unsigned int n)
Set number of collisions to be skipped for plotting.
void SetSensor(Sensor *sensor)
Set the sensor.
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.
unsigned int GetNumberOfElectronEndpoints() const
Component for importing field maps computed by Elmer.
void SetRangeZ(const double zmin, const double zmax)
void SetMedium(const unsigned int imat, Medium *medium)
Associate a field map material with a Medium class.
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 SetTemperature(const double t)
Set the temperature [K].
Definition: Medium.cc:71
virtual void EnableDrift(const bool on=true)
Switch electron/ion/hole on/off.
Definition: Medium.hh:67
void SetPressure(const double p)
Definition: Medium.cc:81
void AddComponent(Component *comp)
Add a component.
Definition: Sensor.cc:316
bool SetArea()
Set the user area to the default.
Definition: Sensor.cc:183
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:18
Draw the mesh of a field-map component.
Definition: ViewFEMesh.hh:21
void SetYaxisTitle(const std::string &ytitle)
Definition: ViewFEMesh.hh:41
void SetViewDrift(ViewDrift *vd)
Set the optional associated ViewDrift.
Definition: ViewFEMesh.hh:63
void SetXaxisTitle(const std::string &xtitle)
Definition: ViewFEMesh.hh:40
void SetPlane(const double fx, const double fy, const double fz, const double x0, const double y0, const double z0) override
Definition: ViewFEMesh.cc:202
bool Plot()
Plot method to be called by user.
Definition: ViewFEMesh.cc:32
void SetComponent(ComponentFieldMap *cmp)
Set the component from which to retrieve the mesh and field.
Definition: ViewFEMesh.cc:21
void SetColor(int matID, int colorID)
Definition: ViewFEMesh.hh:57
void SetFillMesh(const bool f)
Element fill switch; 2D only, set false for wireframe mesh.
Definition: ViewFEMesh.hh:49
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:116
void PlotContour(const std::string &option="v")
Definition: ViewField.cc:122
void SetNumberOfContours(const unsigned int n)
Set the number of contour levels.
Definition: ViewField.cc:108
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
int main(int argc, char *argv[])
Definition: wire2d.cc:27