Garfield++ 3.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
ComponentNeBem2d.hh
Go to the documentation of this file.
1#ifndef G_COMPONENT_NEBEM_2D_H
2#define G_COMPONENT_NEBEM_2D_H
3
4#include "ComponentBase.hh"
5
6namespace Garfield {
7
8/// Two-dimensional implementation of the nearly exact Boundary %Element Method.
9
11 public:
12 /// Constructor
14 /// Destructor
16
17 Medium* GetMedium(const double x, const double y, const double z) override;
18
19 void ElectricField(const double x, const double y, const double z, double& ex,
20 double& ey, double& ez, Medium*& m, int& status) override;
21 void ElectricField(const double x, const double y, const double z, double& ex,
22 double& ey, double& ez, double& v, Medium*& m,
23 int& status) override;
24 bool GetVoltageRange(double& vmin, double& vmax) override;
25
26 bool GetBoundingBox(double& xmin, double& ymin, double& zmin,
27 double& xmax, double& ymax, double& zmax) override;
28
29 bool IsWireCrossed(const double x0, const double y0, const double z0,
30 const double x1, const double y1, const double z1,
31 double& xc, double& yc, double& zc, const bool centre,
32 double& rc) override;
33 bool IsInTrapRadius(const double q0, const double x0, const double y0,
34 const double z0, double& xw, double& yx,
35 double& rw) override;
36
37 /// Set the "background" medium.
38 void SetMedium(Medium* medium) { m_medium = medium; }
39
40 /** Add a conducting straight-line segment.
41 * \param x0,y0,x1,y1 coordinates of start and end point.
42 * \param v applied potential.
43 * \param ndiv number of elements in which to split the segment.
44 */
45 bool AddSegment(const double x0, const double y0, const double x1,
46 const double y1, const double v, const int ndiv = -1);
47 /** Add a wire.
48 * \param x,y centre of the wire.
49 * \param d wire diameter.
50 * \param v applied potential.
51 * \param ntrap multiple of the wire radius within which a particle is
52 considered to be trapped by the wire.
53 */
54 bool AddWire(const double x, const double y, const double d, const double v,
55 const int ntrap = 5);
56 /** Add a region bounded by a polygon.
57 * \param xp,yp x/y-coordinates of the vertices of the polygon.
58 * \param medium pointer to the medium associated to the region.
59 * \param bctype 1: fixed voltage, 4: dielectric-dielectric interface.
60 * \param v applied potential.
61 * \param ndiv number of elements on each edge segment.
62 */
63 bool AddRegion(const std::vector<double>& xp,
64 const std::vector<double>& yp, Medium* medium,
65 const unsigned int bctype = 4, const double v = 0.,
66 const int ndiv = -1);
67 /// Set the extent of the drift region along z.
68 void SetRangeZ(const double zmin, const double zmax);
69
70 /// Discretise the geometry and compute the solution.
71 bool Initialise();
72
73 /// Set the default number of elements per segment.
74 void SetNumberOfDivisions(const unsigned int ndiv);
75
76 void SetNumberOfCollocationPoints(const unsigned int ncoll);
77 void EnableAutoResizing(const bool on = true) { m_autoSize = on; }
78 void EnableRandomCollocation(const bool on = true) {
79 m_randomCollocation = on;
80 }
81 void SetMaxNumberOfIterations(const unsigned int niter);
82
83 /// Return the number of regions.
84 unsigned int GetNumberOfRegions() const { return m_regions.size(); }
85 /// Return the properties of a given region.
86 bool GetRegion(const unsigned int i,
87 std::vector<double>& xv, std::vector<double>& yv,
88 Medium*& medium, unsigned int& bctype, double& v);
89 /// Return the number of conducting straight-line segments.
90 unsigned int GetNumberOfSegments() const { return m_segments.size(); }
91 /// Return the coordinates and voltage of a given straight-line segment.
92 bool GetSegment(const unsigned int i, double& x0, double& y0,
93 double& x1, double& x2, double& v) const;
94 /// Return the number of wires.
95 unsigned int GetNumberOfWires() const { return m_wires.size(); }
96 /// Return the coordinates, diameter, potential and charge of a given wire.
97 bool GetWire(const unsigned int i, double& x, double& y, double& d,
98 double& v, double& q) const;
99 /// Return the number of boundary elements.
100 unsigned int GetNumberOfElements() const { return m_elements.size(); }
101 /// Return the coordinates and charge of a given boundary element.
102 bool GetElement(const unsigned int i, double& x0, double& y0,
103 double& x1, double& y1, double& q) const;
104 private:
105 static const double InvEpsilon0;
106 static const double InvTwoPiEpsilon0;
107
108 /// Default number elements per segment.
109 unsigned int m_nDivisions = 5;
110 unsigned int m_nCollocationPoints = 1;
111 bool m_autoSize = false;
112 bool m_randomCollocation = false;
113 unsigned int m_nMaxIterations = 3;
114
115 /// Background medium.
116 Medium* m_medium = nullptr;
117 /// Flag whether a z-range has been defined by the user.
118 bool m_useRangeZ = false;
119 /// Lower z limit.
120 double m_zmin = -1.;
121 /// Upper z limit.
122 double m_zmax = 1.;
123 /// Boundary condition type.
124 enum BC {
125 Voltage = 1, //< Fixed potential.
126 Charge, //< Fixed charge density (not implemented).
127 Floating, //< Floating conductor (not implemented).
128 Dielectric //< Dielectric-dielectric interface.
129 };
130 struct Region {
131 std::vector<double> xv; //< x-coordinates of the vertices.
132 std::vector<double> yv; //< y-coordinates of the vertices.
133 Medium* medium; //< Medium associated to the region.
134 std::pair<BC, double> bc; //< Applied boundary condition.
135 unsigned int depth; //< Level in the hierarchy.
136 int ndiv; //< Number of elements per edge segment.
137 };
138 /// Regions.
139 std::vector<Region> m_regions;
140
141 struct Segment {
142 std::array<double, 2> x0; //< Coordinates of the start point.
143 std::array<double, 2> x1; //< Coordinates of the end point.
144 int region1; //< Inner region.
145 int region2; //< Outer region.
146 std::pair<BC, double> bc; //< Applied boundary condition.
147 int ndiv; //< Number of elements.
148 };
149 /// User-specified conducting straight-line segments.
150 std::vector<Segment> m_segments;
151
152 struct Wire {
153 double x, y; //< Coordinates of the centre.
154 double r; //< Radius.
155 double v; //< Potential.
156 double q; //< Charge.
157 int ntrap; //< Trap radius (in units of the wire radius).
158 };
159 /// Wires.
160 std::vector<Wire> m_wires;
161
162 struct Element {
163 double x, y; //< Coordinates of the element centre (collocation point).
164 double a; //< Half-length.
165 double cphi; //< Rotation.
166 double sphi; //< Rotation.
167 double q; //< Charge density (solution).
168 std::pair<BC, double> bc; //< Boundary condition.
169 double lambda; //< Ratio of dielectric permittivities.
170 };
171 /// Straight-line boundary elements.
172 std::vector<Element> m_elements;
173
174 /// Split/merge overlapping segments.
175 void EliminateOverlaps(std::vector<Segment>& segments);
176 /// Create elements from a straight-line segment.
177 bool Discretise(const Segment& segment, std::vector<Element>& elements,
178 const double lambda, const unsigned int ndiv);
179
180 bool ComputeInfluenceMatrix(std::vector<std::vector<double> >& infmat) const;
181 bool InvertMatrix(std::vector<std::vector<double> >& influenceMatrix,
182 std::vector<std::vector<double> >& inverseMatrix) const;
183 bool LUDecomposition(std::vector<std::vector<double> >& mat,
184 std::vector<int>& index) const;
185 void LUSubstitution(const std::vector<std::vector<double> >& mat,
186 const std::vector<int>& index,
187 std::vector<double>& col) const;
188
189 bool Solve(const std::vector<std::vector<double> >& inverseMatrix,
190 const std::vector<double>& bc);
191 bool CheckConvergence(const double tol, std::vector<bool>& ok);
192 void SplitElement(Element& oldElement, std::vector<Element>& elements);
193
194 /// Potential of a line segment (half-width a) in local coordinates.
195 double LinePotential(const double a, const double x, const double y) const;
196 /// Potential of a thin wire with radius r0.
197 double WirePotential(const double r0, const double x, const double y) const;
198 /// Field of a line segment (half-width a) in local coordinates.
199 void LineFlux(const double a, const double x, const double y, double& ex,
200 double& ey) const;
201 /// Field of a thin wire with radius r0.
202 void WireFlux(const double r0, const double x, const double y, double& ex,
203 double& ey) const;
204
205 void Reset() override;
206 void UpdatePeriodicity() override;
207 /// Rotation from global to local coordinates.
208 void ToLocal(const double xIn, const double yIn,
209 const double cphi, const double sphi,
210 double& xOut, double& yOut) const;
211 /// Rotation from local to global coordinates.
212 void ToGlobal(const double xIn, const double yIn,
213 const double cphi, const double sphi,
214 double& xOut, double& yOut) const;
215
216};
217}
218
219#endif
Abstract base class for components.
Two-dimensional implementation of the nearly exact Boundary Element Method.
unsigned int GetNumberOfSegments() const
Return the number of conducting straight-line segments.
void SetNumberOfDivisions(const unsigned int ndiv)
Set the default number of elements per segment.
bool AddSegment(const double x0, const double y0, const double x1, const double y1, const double v, const int ndiv=-1)
void ElectricField(const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status) override
bool GetVoltageRange(double &vmin, double &vmax) override
Calculate the voltage range [V].
bool GetSegment(const unsigned int i, double &x0, double &y0, double &x1, double &x2, double &v) const
Return the coordinates and voltage of a given straight-line segment.
bool GetBoundingBox(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax) override
Get the bounding box coordinates.
void EnableRandomCollocation(const bool on=true)
void SetRangeZ(const double zmin, const double zmax)
Set the extent of the drift region along z.
void EnableAutoResizing(const bool on=true)
void SetNumberOfCollocationPoints(const unsigned int ncoll)
unsigned int GetNumberOfRegions() const
Return the number of regions.
bool IsWireCrossed(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, double &xc, double &yc, double &zc, const bool centre, double &rc) override
unsigned int GetNumberOfElements() const
Return the number of boundary elements.
bool Initialise()
Discretise the geometry and compute the solution.
bool GetRegion(const unsigned int i, std::vector< double > &xv, std::vector< double > &yv, Medium *&medium, unsigned int &bctype, double &v)
Return the properties of a given region.
void SetMedium(Medium *medium)
Set the "background" medium.
Medium * GetMedium(const double x, const double y, const double z) override
Get the medium at a given location (x, y, z).
bool GetWire(const unsigned int i, double &x, double &y, double &d, double &v, double &q) const
Return the coordinates, diameter, potential and charge of a given wire.
bool AddWire(const double x, const double y, const double d, const double v, const int ntrap=5)
void SetMaxNumberOfIterations(const unsigned int niter)
unsigned int GetNumberOfWires() const
Return the number of wires.
bool IsInTrapRadius(const double q0, const double x0, const double y0, const double z0, double &xw, double &yx, double &rw) override
bool AddRegion(const std::vector< double > &xp, const std::vector< double > &yp, Medium *medium, const unsigned int bctype=4, const double v=0., const int ndiv=-1)
bool GetElement(const unsigned int i, double &x0, double &y0, double &x1, double &y1, double &q) const
Return the coordinates and charge of a given boundary element.
Abstract base class for media.
Definition: Medium.hh:13