Garfield++ 5.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
Garfield::ComponentElmer2d Class Reference

Component for importing two-dimensional field maps computed by Elmer. More...

#include <ComponentElmer2d.hh>

+ Inheritance diagram for Garfield::ComponentElmer2d:

Public Member Functions

 ComponentElmer2d ()
 Default constructor.
 
 ComponentElmer2d (const std::string &header, const std::string &elist, const std::string &nlist, const std::string &mplist, const std::string &volt, const std::string &unit)
 Constructor with a set of field map files, see Initialise().
 
 ~ComponentElmer2d ()
 Destructor.
 
bool Initialise (const std::string &header="mesh.header", const std::string &elist="mesh.elements", const std::string &nlist="mesh.nodes", const std::string &mplist="dielectrics.dat", const std::string &volt="out.result", const std::string &unit="cm")
 
bool SetWeightingPotential (const std::string &prnsol, const std::string &label)
 Import a list of voltages to be used as weighting potentials.
 
bool SetWeightingField (const std::string &prnsol, const std::string &label)
 
void SetRangeZ (const double zmin, const double zmax)
 
- Public Member Functions inherited from Garfield::ComponentFieldMap
 ComponentFieldMap ()=delete
 Default constructor.
 
 ComponentFieldMap (const std::string &name)
 Constructor.
 
virtual ~ComponentFieldMap ()
 Destructor.
 
bool Check ()
 Check element aspect ratio.
 
void PrintRange ()
 Show x, y, z, V and angular ranges.
 
void PrintMaterials ()
 List all currently defined materials.
 
void DriftMedium (const size_t imat)
 Flag a field map material as a drift medium.
 
void NotDriftMedium (const size_t imat)
 Flag a field map materials as a non-drift medium.
 
size_t GetNumberOfMaterials () const
 Return the number of materials in the field map.
 
double GetPermittivity (const size_t imat) const
 Return the relative permittivity of a field map material.
 
double GetConductivity (const size_t imat) const
 Return the conductivity of a field map material.
 
void SetMedium (const size_t imat, Medium *medium)
 Associate a field map material with a Medium object.
 
MediumGetMedium (const size_t imat) const
 Return the Medium associated to a field map material.
 
void SetGas (Medium *medium)
 
virtual size_t GetNumberOfElements () const
 Return the number of mesh elements.
 
bool GetElement (const size_t i, double &vol, double &dmin, double &dmax) const
 Return the volume and aspect ratio of a mesh element.
 
virtual bool GetElement (const size_t i, size_t &mat, bool &drift, std::vector< size_t > &nodes) const
 Return the material and node indices of a mesh element.
 
virtual size_t GetNumberOfNodes () const
 
virtual bool GetNode (const size_t i, double &x, double &y, double &z) const
 
double GetPotential (const size_t i) const
 
void EnableCheckMapIndices (const bool on=true)
 
void EnableDeleteBackgroundElements (const bool on=true)
 Option to eliminate mesh elements in conductors (default: on).
 
void EnableTetrahedralTreeForElementSearch (const bool on=true)
 
void EnableConvergenceWarnings (const bool on=true)
 
MediumGetMedium (const double x, const double y, const double z) override
 Get the medium at a given location (x, y, z).
 
void ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status) override
 
void ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, double &v, Medium *&m, int &status) override
 Calculate the drift field [V/cm] and potential [V] at (x, y, z).
 
void WeightingField (const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label) override
 
double WeightingPotential (const double x, const double y, const double z, const std::string &label) override
 
double DelayedWeightingPotential (double x, double y, double z, const double t, const std::string &label) override
 
bool IsInBoundingBox (const double x, const double y, const double z) const
 
bool GetBoundingBox (double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax) override
 Get the bounding box coordinates.
 
bool GetElementaryCell (double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax) override
 Get the coordinates of the elementary cell.
 
bool GetVoltageRange (double &vmin, double &vmax) override
 Calculate the voltage range [V].
 
void CopyWeightingPotential (const std::string &label, const std::string &labelSource, const double x, const double y, const double z, const double alpha, const double beta, const double gamma)
 
std::array< double, 3 > ElectricField (const double x, const double y, const double z)
 Calculate the drift field [V/cm] at (x, y, z).
 
- Public Member Functions inherited from Garfield::Component
 Component ()=delete
 Default constructor.
 
 Component (const std::string &name)
 Constructor.
 
virtual ~Component ()
 Destructor.
 
virtual void SetGeometry (Geometry *geo)
 Define the geometry.
 
virtual void Clear ()
 Reset.
 
std::array< double, 3 > ElectricField (const double x, const double y, const double z)
 Calculate the drift field [V/cm] at (x, y, z).
 
virtual double ElectricPotential (const double x, const double y, const double z)
 Calculate the (drift) electrostatic potential [V] at (x, y, z).
 
virtual void DelayedWeightingField (const double x, const double y, const double z, const double t, double &wx, double &wy, double &wz, const std::string &label)
 
virtual void MagneticField (const double x, const double y, const double z, double &bx, double &by, double &bz, int &status)
 
void SetMagneticField (const double bx, const double by, const double bz)
 Set a constant magnetic field.
 
virtual bool IsReady ()
 Ready for use?
 
double IntegrateFluxCircle (const double xc, const double yc, const double r, const unsigned int nI=50)
 
double IntegrateFluxSphere (const double xc, const double yc, const double zc, const double r, const unsigned int nI=20)
 
double IntegrateFluxParallelogram (const double x0, const double y0, const double z0, const double dx1, const double dy1, const double dz1, const double dx2, const double dy2, const double dz2, const unsigned int nU=20, const unsigned int nV=20)
 
double IntegrateWeightingFluxParallelogram (const std::string &label, const double x0, const double y0, const double z0, const double dx1, const double dy1, const double dz1, const double dx2, const double dy2, const double dz2, const unsigned int nU=20, const unsigned int nV=20)
 
double IntegrateFluxLine (const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, const double xp, const double yp, const double zp, const unsigned int nI, const int isign=0)
 
virtual bool CrossedWire (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)
 
virtual bool InTrapRadius (const double q0, const double x0, const double y0, const double z0, double &xw, double &yw, double &rw)
 
virtual bool CrossedPlane (const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, double &xc, double &yc, double &zc)
 
void EnablePeriodicityX (const bool on=true)
 Enable simple periodicity in the $x$ direction.
 
void EnablePeriodicityY (const bool on=true)
 Enable simple periodicity in the $y$ direction.
 
void EnablePeriodicityZ (const bool on=true)
 Enable simple periodicity in the $z$ direction.
 
void IsPeriodic (bool &perx, bool &pery, bool &perz)
 Return periodicity flags.
 
void EnableMirrorPeriodicityX (const bool on=true)
 Enable mirror periodicity in the $x$ direction.
 
void EnableMirrorPeriodicityY (const bool on=true)
 Enable mirror periodicity in the $y$ direction.
 
void EnableMirrorPeriodicityZ (const bool on=true)
 Enable mirror periodicity in the $y$ direction.
 
void IsMirrorPeriodic (bool &perx, bool &pery, bool &perz)
 Return mirror periodicity flags.
 
void EnableAxialPeriodicityX (const bool on=true)
 Enable axial periodicity in the $x$ direction.
 
void EnableAxialPeriodicityY (const bool on=true)
 Enable axial periodicity in the $y$ direction.
 
void EnableAxialPeriodicityZ (const bool on=true)
 Enable axial periodicity in the $z$ direction.
 
void IsAxiallyPeriodic (bool &perx, bool &pery, bool &perz)
 Return axial periodicity flags.
 
void EnableRotationSymmetryX (const bool on=true)
 Enable rotation symmetry around the $x$ axis.
 
void EnableRotationSymmetryY (const bool on=true)
 Enable rotation symmetry around the $y$ axis.
 
void EnableRotationSymmetryZ (const bool on=true)
 Enable rotation symmetry around the $z$ axis.
 
void IsRotationSymmetric (bool &rotx, bool &roty, bool &rotz)
 Return rotation symmetry flags.
 
void EnableDebugging ()
 Switch on debugging messages.
 
void DisableDebugging ()
 Switch off debugging messages.
 
virtual bool HasMagneticField () const
 Does the component have a non-zero magnetic field?
 
virtual bool HasTownsendMap () const
 Does the component have maps of the Townsend coefficient?
 
virtual bool HasAttachmentMap () const
 Does the component have attachment maps?
 
virtual bool HasVelocityMap () const
 Does the component have velocity maps?
 
virtual bool ElectronAttachment (const double, const double, const double, double &eta)
 Get the electron attachment coefficient.
 
virtual bool HoleAttachment (const double, const double, const double, double &eta)
 Get the hole attachment coefficient.
 
virtual bool ElectronTownsend (const double, const double, const double, double &alpha)
 Get the electron Townsend coefficient.
 
virtual bool HoleTownsend (const double, const double, const double, double &alpha)
 Get the hole Townsend coefficient.
 
virtual bool ElectronVelocity (const double, const double, const double, double &vx, double &vy, double &vz)
 Get the electron drift velocity.
 
virtual bool HoleVelocity (const double, const double, const double, double &vx, double &vy, double &vz)
 Get the hole drift velocity.
 
virtual bool GetElectronLifetime (const double, const double, const double, double &etau)
 
virtual bool GetHoleLifetime (const double, const double, const double, double &htau)
 
virtual double StepSizeHint ()
 

Additional Inherited Members

- Protected Types inherited from Garfield::ComponentFieldMap
enum class  ElementType { Unknown = 0 , Serendipity = 5 , CurvedTetrahedron = 13 }
 
- Protected Member Functions inherited from Garfield::ComponentFieldMap
void Reset () override
 Reset the component.
 
void Prepare ()
 
virtual void SetRange ()
 
void UpdatePeriodicity () override
 Verify periodicities.
 
void UpdatePeriodicity2d ()
 
void UpdatePeriodicityCommon ()
 
bool SetDefaultDriftMedium ()
 Find lowest epsilon, check for eps = 0, set default drift media flags.
 
int Field (const double x, const double y, const double z, double &fx, double &fy, double &fz, int &iel, const std::vector< double > &potentials) const
 Compute the electric/weighting field.
 
double Potential (const double x, const double y, const double z, const std::vector< double > &potentials) const
 Compute the electrostatic/weighting potential.
 
int FindElement5 (const double x, const double y, double &t1, double &t2, double &t3, double &t4, double jac[4][4], double &det) const
 Find the element for a point in curved quadratic quadrilaterals.
 
int FindElement13 (const double x, const double y, const double z, double &t1, double &t2, double &t3, double &t4, double jac[4][4], double &det) const
 Find the element for a point in curved quadratic tetrahedra.
 
int FindElementCube (const double x, const double y, const double z, double &t1, double &t2, double &t3, TMatrixD *&jac, std::vector< TMatrixD * > &dN) const
 Find the element for a point in a cube.
 
void MapCoordinates (double &xpos, double &ypos, double &zpos, bool &xmirrored, bool &ymirrored, bool &zmirrored, double &rcoordinate, double &rotation) const
 Move (xpos, ypos, zpos) to field map coordinates.
 
void UnmapFields (double &ex, double &ey, double &ez, const double xpos, const double ypos, const double zpos, const bool xmirrored, const bool ymirrored, const bool zmirrored, const double rcoordinate, const double rotation) const
 Move (ex, ey, ez) to global coordinates.
 
virtual double GetElementVolume (const size_t i) const
 
virtual void GetAspectRatio (const size_t i, double &dmin, double &dmax) const
 
void PrintWarning (const std::string &header)
 
void PrintNotReady (const std::string &header) const
 
void PrintCouldNotOpen (const std::string &header, const std::string &filename) const
 
void PrintElement (const std::string &header, const double x, const double y, const double z, const double t1, const double t2, const double t3, const double t4, const size_t i, const std::vector< double > &potential) const
 
void TimeInterpolation (const double t, double &f0, double &f1, int &i0, int &i1)
 Interpolation of potential between two time slices.
 
- Static Protected Member Functions inherited from Garfield::ComponentFieldMap
static double ScalingFactor (std::string unit)
 
static double Potential3 (const std::array< double, 6 > &v, const std::array< double, 3 > &t)
 Interpolate the potential in a triangle.
 
static void Field3 (const std::array< double, 6 > &v, const std::array< double, 3 > &t, double jac[4][4], const double det, double &ex, double &ey)
 Interpolate the field in a triangle.
 
static double Potential5 (const std::array< double, 8 > &v, const std::array< double, 2 > &t)
 Interpolate the potential in a curved quadrilateral.
 
static void Field5 (const std::array< double, 8 > &v, const std::array< double, 2 > &t, double jac[4][4], const double det, double &ex, double &ey)
 Interpolate the field in a curved quadrilateral.
 
static double Potential13 (const std::array< double, 10 > &v, const std::array< double, 4 > &t)
 Interpolate the potential in a curved quadratic tetrahedron.
 
static void Field13 (const std::array< double, 10 > &v, const std::array< double, 4 > &t, double jac[4][4], const double det, double &ex, double &ey, double &ez)
 Interpolate the field in a curved quadratic tetrahedron.
 
static int ReadInteger (char *token, int def, bool &error)
 
static double ReadDouble (char *token, double def, bool &error)
 
- Protected Attributes inherited from Garfield::ComponentFieldMap
bool m_is3d = true
 
ElementType m_elementType = ElementType::CurvedTetrahedron
 
std::vector< Elementm_elements
 
std::vector< int > m_elementIndices
 
std::vector< bool > m_degenerate
 
std::vector< std::array< double, 3 > > m_bbMin
 
std::vector< std::array< double, 3 > > m_bbMax
 
std::vector< std::array< std::array< double, 3 >, 4 > > m_w12
 
std::vector< Nodem_nodes
 
std::vector< double > m_pot
 
std::map< std::string, std::vector< double > > m_wpot
 
std::map< std::string, std::vector< std::vector< double > > > m_dwpot
 
std::vector< Materialm_materials
 
std::map< std::string, WeightingFieldCopym_wfieldCopies
 
std::vector< double > m_wdtimes
 
bool m_hasBoundingBox = false
 
std::array< double, 3 > m_minBoundingBox = {{0., 0., 0.}}
 
std::array< double, 3 > m_maxBoundingBox = {{0., 0., 0.}}
 
std::array< double, 3 > m_mapmin = {{0., 0., 0.}}
 
std::array< double, 3 > m_mapmax = {{0., 0., 0.}}
 
std::array< double, 3 > m_mapamin = {{0., 0., 0.}}
 
std::array< double, 3 > m_mapamax = {{0., 0., 0.}}
 
std::array< double, 3 > m_mapna = {{0., 0., 0.}}
 
std::array< double, 3 > m_cells = {{0., 0., 0.}}
 
double m_mapvmin = 0.
 
double m_mapvmax = 0.
 
std::array< bool, 3 > m_setang
 
bool m_deleteBackground = true
 
bool m_warning = false
 
unsigned int m_nWarnings = 0
 
bool m_printConvergenceWarnings = true
 
- Protected Attributes inherited from Garfield::Component
std::string m_className = "Component"
 Class name.
 
Geometrym_geometry = nullptr
 Pointer to the geometry.
 
std::array< double, 3 > m_b0 = {{0., 0., 0.}}
 Constant magnetic field.
 
bool m_ready = false
 Ready for use?
 
bool m_debug = false
 Switch on/off debugging messages.
 
std::array< bool, 3 > m_periodic = {{false, false, false}}
 Simple periodicity in x, y, z.
 
std::array< bool, 3 > m_mirrorPeriodic = {{false, false, false}}
 Mirror periodicity in x, y, z.
 
std::array< bool, 3 > m_axiallyPeriodic = {{false, false, false}}
 Axial periodicity in x, y, z.
 
std::array< bool, 3 > m_rotationSymmetric = {{false, false, false}}
 Rotation symmetry around x-axis, y-axis, z-axis.
 

Detailed Description

Component for importing two-dimensional field maps computed by Elmer.

Definition at line 10 of file ComponentElmer2d.hh.

Constructor & Destructor Documentation

◆ ComponentElmer2d() [1/2]

Garfield::ComponentElmer2d::ComponentElmer2d ( )

Default constructor.

Definition at line 20 of file ComponentElmer2d.cc.

20 : ComponentFieldMap("Elmer2d") {
21 m_is3d = false;
23
24 // Default bounding box
25 m_minBoundingBox[2] = -50;
26 m_maxBoundingBox[2] = 50;
27}
std::array< double, 3 > m_maxBoundingBox
std::array< double, 3 > m_minBoundingBox
ComponentFieldMap()=delete
Default constructor.

Referenced by ComponentElmer2d().

◆ ComponentElmer2d() [2/2]

Garfield::ComponentElmer2d::ComponentElmer2d ( const std::string & header,
const std::string & elist,
const std::string & nlist,
const std::string & mplist,
const std::string & volt,
const std::string & unit )

Constructor with a set of field map files, see Initialise().

Definition at line 29 of file ComponentElmer2d.cc.

36
37 Initialise(header, elist, nlist, mplist, volt, unit);
38}
ComponentElmer2d()
Default constructor.
bool Initialise(const std::string &header="mesh.header", const std::string &elist="mesh.elements", const std::string &nlist="mesh.nodes", const std::string &mplist="dielectrics.dat", const std::string &volt="out.result", const std::string &unit="cm")

◆ ~ComponentElmer2d()

Garfield::ComponentElmer2d::~ComponentElmer2d ( )
inline

Destructor.

Definition at line 19 of file ComponentElmer2d.hh.

19{}

Member Function Documentation

◆ Initialise()

bool Garfield::ComponentElmer2d::Initialise ( const std::string & header = "mesh.header",
const std::string & elist = "mesh.elements",
const std::string & nlist = "mesh.nodes",
const std::string & mplist = "dielectrics.dat",
const std::string & volt = "out.result",
const std::string & unit = "cm" )

Import a field map from a set of files.

Parameters
headername of the header file (contains the number of elements and nodes).
elistname of the file that contains the list of mesh elements
nlistname of the file that contains the list of mesh nodes
mplistname of the file that contains the material properties
voltoutput of the field solver (list of voltages)
unitlength unit to be used

Definition at line 40 of file ComponentElmer2d.cc.

45 {
46 const std::string hdr = m_className + "::Initialise:";
48
49 // Keep track of the success.
50 bool ok = true;
51
52 // Buffer for reading
53 constexpr int size = 100;
54 char line[size];
55
56 // Open the header.
57 std::ifstream fheader(header);
58 if (!fheader) {
59 PrintCouldNotOpen("Initialise", header);
60 return false;
61 }
62
63 // Read the header to get the number of nodes and elements.
64 fheader.getline(line, size, '\n');
65 char* token = strtok(line, " ");
66 bool readerror = false;
67 const int nNodes = ReadInteger(token, 0, readerror);
68 token = strtok(nullptr, " ");
69 const int nElements = ReadInteger(token, 0, readerror);
70 std::cout << hdr << "\n Read " << nNodes << " nodes and " << nElements
71 << " elements from file " << header << ".\n";
72 if (readerror) {
73 PrintErrorReadingFile(hdr, header, 0);
74 fheader.close();
75 return false;
76 }
77
78 // Close the header file.
79 fheader.close();
80
81 // Open the nodes list.
82 std::ifstream fnodes(nlist);
83 if (!fnodes) {
84 PrintCouldNotOpen("Initialise", nlist);
85 return false;
86 }
87
88 // Check the value of the unit.
89 double funit = ScalingFactor(unit);
90 if (funit <= 0.) {
91 std::cerr << hdr << " Unknown length unit " << unit << " Will use cm.\n";
92 funit = 1.0;
93 }
94 if (m_debug) std::cout << " Unit scaling factor = " << funit << ".\n";
95
96 // Read the nodes from the file.
97 for (int il = 0; il < nNodes; il++) {
98 // Get a line from the nodes file.
99 fnodes.getline(line, size, '\n');
100
101 // Ignore the first two characters.
102 token = strtok(line, " ");
103 token = strtok(nullptr, " ");
104
105 // Get the node coordinates.
106 token = strtok(nullptr, " ");
107 double xnode = ReadDouble(token, -1, readerror);
108 token = strtok(nullptr, " ");
109 double ynode = ReadDouble(token, -1, readerror);
110 token = strtok(nullptr, " ");
111 double znode = ReadDouble(token, -1, readerror);
112 if (readerror) {
113 PrintErrorReadingFile(hdr, nlist, il);
114 fnodes.close();
115 return false;
116 }
117
118 // Set up and create a new node.
119 Node node;
120 node.x = xnode * funit;
121 node.y = ynode * funit;
122 node.z = znode * funit;
123 m_nodes.push_back(std::move(node));
124 }
125
126 // Close the nodes file.
127 fnodes.close();
128
129 // Read the potentials.
130 if (!LoadPotentials(volt, m_pot)) return false;
131
132 // Open the materials file.
133 std::ifstream fmplist(mplist);
134 if (!fmplist) {
135 PrintCouldNotOpen("Initialise", mplist);
136 return false;
137 }
138
139 // Read the dielectric constants from the materials file.
140 fmplist.getline(line, size, '\n');
141 token = strtok(line, " ");
142 if (readerror) {
143 std::cerr << hdr << "\n Error reading number of materials from "
144 << mplist << ".\n";
145 fmplist.close();
146 return false;
147 }
148 const unsigned int nMaterials = ReadInteger(token, 0, readerror);
149 m_materials.resize(nMaterials);
150 for (auto& material : m_materials) {
151 material.ohm = -1;
152 material.eps = -1;
153 material.medium = nullptr;
154 }
155 for (int il = 2; il < ((int)nMaterials + 2); il++) {
156 fmplist.getline(line, size, '\n');
157 token = strtok(line, " ");
158 ReadInteger(token, -1, readerror);
159 token = strtok(nullptr, " ");
160 double dc = ReadDouble(token, -1.0, readerror);
161 if (readerror) {
162 PrintErrorReadingFile(hdr, mplist, il);
163 fmplist.close();
164 return false;
165 }
166 m_materials[il - 2].eps = dc;
167 std::cout << " Set material " << il - 2 << " of "
168 << nMaterials << " to eps " << dc << ".\n";
169 }
170
171 // Close the materials file.
172 fmplist.close();
173
174 // Find lowest epsilon, check for eps = 0, set default drift medium.
175 if (!SetDefaultDriftMedium()) return false;
176
177 // Open the elements file.
178 std::ifstream felems(elist);
179 if (!felems) {
180 PrintCouldNotOpen("Initialise", elist);
181 return false;
182 }
183
184 // Read the elements and their material indices.
185 for (int il = 0; il < nElements; il++) {
186 // Get a line
187 felems.getline(line, size, '\n');
188
189 // Split into tokens.
190 token = strtok(line, " ");
191 // Read the 2nd-order element.
192 token = strtok(nullptr, " ");
193 int imat = ReadInteger(token, -1, readerror) - 1;
194 token = strtok(nullptr, " ");
195 std::vector<int> inode;
196 for (size_t k = 0; k < 8; ++k) {
197 token = strtok(nullptr, " ");
198 const int in = ReadInteger(token, -1, readerror);
199 if (!readerror) inode.push_back(in);
200 }
201 token = strtok(nullptr, " ");
202 if (inode.size() != 8) {
203 PrintErrorReadingFile(hdr, elist, il);
204 std::cerr << " Read " << inode.size() << " node indices for element"
205 << il << " (expected 8).\n";
206 felems.close();
207 return false;
208 }
209
210 if (m_debug && il < 10) {
211 std::cout << " Read nodes ";
212 for (size_t k = 0; k < 8; ++k) {
213 std::cout << inode[k];
214 if (k < 7) std::cout << ", ";
215 }
216 std::cout << "\n from element " << il + 1 << " of "
217 << nElements << " with material " << imat << ".\n";
218 }
219
220 // Check the material number and ensure that epsilon is non-negative.
221 if (imat < 0 || imat > (int)nMaterials) {
222 std::cerr << hdr << "\n Out-of-range material number on file " << elist
223 << " (line " << il << ").\n";
224 std::cerr << " Element: " << il << ", material: " << imat << "\n";
225 ok = false;
226 break;
227 }
228 if (m_materials[imat].eps < 0) {
229 std::cerr << hdr << "\n Element " << il << " in " << elist << "\n"
230 << " uses material " << imat << " which does not have\n"
231 << " a positive permittivity in " << mplist << ".\n";
232 ok = false;
233 break;
234 }
235
236 // Check the node numbers.
237 bool degenerate = false;
238 for (size_t k = 0; k < 8; ++k) {
239 if (inode[k] < 1) {
240 std::cerr << hdr << "\n Found a node number < 1 on file " << elist
241 << " (line " << il << ").\n Element: " << il
242 << ", material: " << imat << ".\n";
243 ok = false;
244 }
245 // These elements must not be degenerate.
246 for (size_t kk = k + 1; kk < 8; ++kk) {
247 if (inode[k] == inode[kk]) degenerate = true;
248 }
249 }
250 if (degenerate) {
251 std::cerr << hdr << "\n Element " << il << " of file " << elist
252 << " is degenerate,\n"
253 << " no such elements are allowed in this type of map.\n";
254 ok = false;
255 }
256 if (!ok) break;
257 Element element;
258 // Store the material reference.
259 element.matmap = imat;
260 // Store the node references.
261 // 3 -- 6 -- 2
262 // | |
263 // 7 5
264 // | |
265 // 0 -- 4 -- 1
266 // Note: Ordering of Elmer elements is described in the
267 // ElmerSolver manual (appendix D. at the time of this comment).
268 // In order to work properly with Coordinates4 and Coordinates5, the
269 // Elmer ordering 0,1,2,3,4,5,6,7 (counter-clockwise about the element)
270 // will need to be changed to 3,2,1,0,6,5,4,7 (clockwise about the element)
271 // if the normal of the defined element points in the -Z direction.
272
273 // Check the direction of the element normal, +Z or -Z.
274 double x01 = m_nodes[inode[1] - 1].x - m_nodes[inode[0] - 1].x;
275 double y01 = m_nodes[inode[1] - 1].y - m_nodes[inode[0] - 1].y;
276 double x12 = m_nodes[inode[2] - 1].x - m_nodes[inode[1] - 1].x;
277 double y12 = m_nodes[inode[2] - 1].y - m_nodes[inode[1] - 1].y;
278 const double crossprod = x01 * y12 - y01 * x12;
279 if (crossprod < 0) {
280 element.emap[3] = inode[0] - 1;
281 element.emap[2] = inode[1] - 1;
282 element.emap[1] = inode[2] - 1;
283 element.emap[0] = inode[3] - 1;
284 element.emap[6] = inode[4] - 1;
285 element.emap[5] = inode[5] - 1;
286 element.emap[4] = inode[6] - 1;
287 element.emap[7] = inode[7] - 1;
288 } else {
289 element.emap[0] = inode[0] - 1;
290 element.emap[1] = inode[1] - 1;
291 element.emap[2] = inode[2] - 1;
292 element.emap[3] = inode[3] - 1;
293 element.emap[4] = inode[4] - 1;
294 element.emap[5] = inode[5] - 1;
295 element.emap[6] = inode[6] - 1;
296 element.emap[7] = inode[7] - 1;
297 }
298 m_elements.push_back(std::move(element));
299 }
300
301 // Close the elements file.
302 felems.close();
303 if (!ok) return false;
304
305 m_degenerate.assign(m_elements.size(), false);
306 // Set the ready flag.
307 m_ready = true;
308 std::cout << " Finished.\n";
309
310 Prepare();
311 return true;
312}
bool SetDefaultDriftMedium()
Find lowest epsilon, check for eps = 0, set default drift media flags.
static double ReadDouble(char *token, double def, bool &error)
static double ScalingFactor(std::string unit)
static int ReadInteger(char *token, int def, bool &error)
std::vector< Material > m_materials
void PrintCouldNotOpen(const std::string &header, const std::string &filename) const
std::vector< Element > m_elements
void Reset() override
Reset the component.
bool m_debug
Switch on/off debugging messages.
Definition Component.hh:371
std::string m_className
Class name.
Definition Component.hh:359
bool m_ready
Ready for use?
Definition Component.hh:368

Referenced by ComponentElmer2d().

◆ SetRangeZ()

void Garfield::ComponentElmer2d::SetRangeZ ( const double zmin,
const double zmax )

Definition at line 396 of file ComponentElmer2d.cc.

396 {
397 if (fabs(zmax - zmin) <= 0.) {
398 std::cerr << m_className << "::SetRangeZ: Zero range is not permitted.\n";
399 return;
400 }
401 m_minBoundingBox[2] = m_mapmin[2] = std::min(zmin, zmax);
402 m_maxBoundingBox[2] = m_mapmax[2] = std::max(zmin, zmax);
403}
std::array< double, 3 > m_mapmin
std::array< double, 3 > m_mapmax
DoubleAc fabs(const DoubleAc &f)
Definition DoubleAc.h:615

Referenced by main().

◆ SetWeightingField()

bool Garfield::ComponentElmer2d::SetWeightingField ( const std::string & prnsol,
const std::string & label )

Definition at line 314 of file ComponentElmer2d.cc.

315 {
316 const std::string hdr = m_className + "::SetWeightingField:";
317 if (!m_ready) {
318 PrintNotReady("SetWeightingField");
319 std::cerr << " Weighting field cannot be added.\n";
320 return false;
321 }
322 std::cout << m_className << "::SetWeightingField:\n"
323 << " Loading field map for electrode " << label << ".\n";
324 // Check if a weighting field with the same label already exists.
325 if (m_wpot.count(label) > 0) {
326 std::cout << " Replacing existing weighting field.\n";
327 m_wpot[label].clear();
328 }
329 std::vector<double> pot(m_nodes.size(), 0.);
330 if (!LoadPotentials(wvolt, pot)) return false;
331 m_wpot[label] = std::move(pot);
332 return true;
333}
std::map< std::string, std::vector< double > > m_wpot
void PrintNotReady(const std::string &header) const

Referenced by SetWeightingPotential().

◆ SetWeightingPotential()

bool Garfield::ComponentElmer2d::SetWeightingPotential ( const std::string & prnsol,
const std::string & label )
inline

Import a list of voltages to be used as weighting potentials.

Definition at line 37 of file ComponentElmer2d.hh.

38 {
39 return SetWeightingField(prnsol, label);
40 }
bool SetWeightingField(const std::string &prnsol, const std::string &label)

The documentation for this class was generated from the following files: