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::ComponentElmer Class Reference

Component for importing field maps computed by Elmer. More...

#include <ComponentElmer.hh>

+ Inheritance diagram for Garfield::ComponentElmer:

Public Member Functions

 ComponentElmer ()
 Default constructor.
 
 ComponentElmer (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().
 
 ~ComponentElmer ()
 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)
 
- 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 field maps computed by Elmer.

Definition at line 10 of file ComponentElmer.hh.

Constructor & Destructor Documentation

◆ ComponentElmer() [1/2]

Garfield::ComponentElmer::ComponentElmer ( )

Default constructor.

Definition at line 20 of file ComponentElmer.cc.

20: ComponentFieldMap("Elmer") {}
ComponentFieldMap()=delete
Default constructor.

◆ ComponentElmer() [2/2]

Garfield::ComponentElmer::ComponentElmer ( 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 22 of file ComponentElmer.cc.

27 : ComponentFieldMap("Elmer") {
28 Initialise(header, elist, nlist, mplist, volt, unit);
29}
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")

◆ ~ComponentElmer()

Garfield::ComponentElmer::~ComponentElmer ( )
inline

Destructor.

Definition at line 19 of file ComponentElmer.hh.

19{}

Member Function Documentation

◆ Initialise()

bool Garfield::ComponentElmer::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 31 of file ComponentElmer.cc.

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

◆ SetWeightingField()

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

Definition at line 281 of file ComponentElmer.cc.

282 {
283 const std::string hdr = m_className + "::SetWeightingField:";
284 if (!m_ready) {
285 PrintNotReady("SetWeightingField");
286 std::cerr << " Weighting field cannot be added.\n";
287 return false;
288 }
289
290 std::cout << m_className << "::SetWeightingField:\n"
291 << " Loading field map for electrode " << label << ".\n";
292 if (m_wpot.count(label) > 0) {
293 std::cout << " Replacing existing weighting field.\n";
294 m_wpot[label].clear();
295 }
296 std::vector<double> pot(m_nodes.size(), 0.);
297 if (!LoadPotentials(wvolt, pot)) return false;
298 m_wpot[label] = std::move(pot);
299 return true;
300}
std::map< std::string, std::vector< double > > m_wpot
void PrintNotReady(const std::string &header) const

Referenced by SetWeightingPotential().

◆ SetWeightingPotential()

bool Garfield::ComponentElmer::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 ComponentElmer.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: