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

Base class for components based on finite-element field maps. More...

#include <ComponentFieldMap.hh>

+ Inheritance diagram for Garfield::ComponentFieldMap:

Classes

struct  Element
 
struct  Material
 
struct  Node
 

Public Member Functions

 ComponentFieldMap ()=delete
 Default constructor.
 
 ComponentFieldMap (const std::string &name)
 Constructor.
 
virtual ~ComponentFieldMap ()
 Destructor.
 
virtual void SetRange ()
 Calculate x, y, z, V and angular ranges.
 
void PrintRange ()
 Show x, y, z, V and angular ranges.
 
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 PrintMaterials ()
 List all currently defined materials.
 
void DriftMedium (const unsigned int imat)
 Flag a field map material as a drift medium.
 
void NotDriftMedium (const unsigned int 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 unsigned int imat) const
 Return the permittivity of a field map material.
 
double GetConductivity (const unsigned int imat) const
 Return the conductivity of a field map material.
 
void SetMedium (const unsigned int imat, Medium *medium)
 Associate a field map material with a Medium class.
 
MediumGetMedium (const unsigned int i) const
 Return the Medium associated to a field map material.
 
virtual size_t GetNumberOfElements () const
 Return the number of mesh elements.
 
bool GetElement (const size_t i, double &vol, double &dmin, double &dmax)
 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
 
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)
 
virtual MediumGetMedium (const double x, const double y, const double z)
 Get the medium at a given location (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.
 
virtual MediumGetMedium (const double x, const double y, const double z)
 Get the medium at a given location (x, y, z).
 
virtual void ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status)=0
 
virtual void ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, double &v, Medium *&m, int &status)=0
 Calculate the drift field [V/cm] and potential [V] at (x, y, z).
 
virtual bool GetVoltageRange (double &vmin, double &vmax)=0
 Calculate the voltage range [V].
 
virtual void WeightingField (const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label)
 
virtual double WeightingPotential (const double x, const double y, const double z, const std::string &label)
 
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 double DelayedWeightingPotential (const double x, const double y, const double z, const double t, 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?
 
virtual bool GetBoundingBox (double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
 Get the bounding box coordinates.
 
virtual bool GetElementaryCell (double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
 Get the coordinates of the elementary cell.
 
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 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)
 
virtual bool IsInTrapRadius (const double q0, const double x0, const double y0, const double z0, double &xw, double &yw, double &rw)
 
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 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 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)
 

Protected Member Functions

void Reset () override
 Reset the component.
 
void Prepare ()
 
void UpdatePeriodicity2d ()
 
void UpdatePeriodicityCommon ()
 
bool SetDefaultDriftMedium ()
 Find lowest epsilon, check for eps = 0, set default drift media flags.
 
int FindElement5 (const double x, const double y, const double z, double &t1, double &t2, double &t3, double &t4, double jac[4][4], double &det)
 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)
 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)
 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, double &xpos, double &ypos, double &zpos, bool &xmirrored, bool &ymirrored, bool &zmirrored, double &rcoordinate, double &rotation) const
 Move (ex, ey, ez) to global coordinates.
 
int ReadInteger (char *token, int def, bool &error)
 
double ReadDouble (char *token, double def, bool &error)
 
virtual double GetElementVolume (const unsigned int i)=0
 
virtual void GetAspectRatio (const unsigned int i, double &dmin, double &dmax)=0
 
size_t GetWeightingFieldIndex (const std::string &label) const
 
size_t GetOrCreateWeightingFieldIndex (const std::string &label)
 
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 Element &element, const unsigned int n, const int iw=-1) const
 
- Protected Member Functions inherited from Garfield::Component
virtual void Reset ()=0
 Reset the component.
 
virtual void UpdatePeriodicity ()=0
 Verify periodicities.
 

Static Protected Member Functions

static double ScalingFactor (std::string unit)
 

Protected Attributes

bool m_is3d = true
 
std::vector< Elementm_elements
 
std::vector< Nodem_nodes
 
std::vector< Materialm_materials
 
std::vector< std::string > m_wfields
 
std::vector< bool > m_wfieldsOk
 
std::vector< bool > m_dwfieldsOk
 
std::vector< double > m_wdtimes
 
bool m_hasBoundingBox = false
 
std::array< double, 3 > m_minBoundingBox
 
std::array< double, 3 > m_maxBoundingBox
 
std::array< double, 3 > m_mapmin
 
std::array< double, 3 > m_mapmax
 
std::array< double, 3 > m_mapamin
 
std::array< double, 3 > m_mapamax
 
std::array< double, 3 > m_mapna
 
std::array< double, 3 > m_cells
 
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.
 

Friends

class ViewFEMesh
 

Detailed Description

Base class for components based on finite-element field maps.

Definition at line 16 of file ComponentFieldMap.hh.

Constructor & Destructor Documentation

◆ ComponentFieldMap() [1/2]

Garfield::ComponentFieldMap::ComponentFieldMap ( )
delete

Default constructor.

◆ ComponentFieldMap() [2/2]

Garfield::ComponentFieldMap::ComponentFieldMap ( const std::string &  name)

Constructor.

Definition at line 15 of file ComponentFieldMap.cc.

16 : Component(name) {}
Component()=delete
Default constructor.

◆ ~ComponentFieldMap()

Garfield::ComponentFieldMap::~ComponentFieldMap ( )
virtual

Destructor.

Definition at line 18 of file ComponentFieldMap.cc.

18{}

Member Function Documentation

◆ DriftMedium()

void Garfield::ComponentFieldMap::DriftMedium ( const unsigned int  imat)

Flag a field map material as a drift medium.

Definition at line 50 of file ComponentFieldMap.cc.

50 {
51 // Do not proceed if not properly initialised.
52 if (!m_ready) PrintNotReady("DriftMedium");
53
54 // Check value
55 if (imat >= m_materials.size()) {
56 std::cerr << m_className << "::DriftMedium: Index out of range.\n";
57 return;
58 }
59
60 // Make drift medium
61 m_materials[imat].driftmedium = true;
62}
std::vector< Material > m_materials
void PrintNotReady(const std::string &header) const
std::string m_className
Class name.
Definition: Component.hh:329
bool m_ready
Ready for use?
Definition: Component.hh:338

◆ EnableCheckMapIndices()

void Garfield::ComponentFieldMap::EnableCheckMapIndices ( const bool  on = true)
inline

Definition at line 76 of file ComponentFieldMap.hh.

76 {
77 m_checkMultipleElement = on;
78 if (on) m_lastElement = -1;
79 }

◆ EnableConvergenceWarnings()

void Garfield::ComponentFieldMap::EnableConvergenceWarnings ( const bool  on = true)
inline

Enable or disable warnings that the calculation of the local coordinates did not achieve the requested precision.

Definition at line 93 of file ComponentFieldMap.hh.

◆ EnableDeleteBackgroundElements()

void Garfield::ComponentFieldMap::EnableDeleteBackgroundElements ( const bool  on = true)
inline

Option to eliminate mesh elements in conductors (default: on).

Definition at line 81 of file ComponentFieldMap.hh.

◆ EnableTetrahedralTreeForElementSearch()

void Garfield::ComponentFieldMap::EnableTetrahedralTreeForElementSearch ( const bool  on = true)
inline

Enable or disable the usage of the tetrahedral tree for searching the element in the mesh.

Definition at line 87 of file ComponentFieldMap.hh.

87 {
88 m_useTetrahedralTree = on;
89 }

◆ FindElement13()

int Garfield::ComponentFieldMap::FindElement13 ( const double  x,
const double  y,
const double  z,
double &  t1,
double &  t2,
double &  t3,
double &  t4,
double  jac[4][4],
double &  det 
)
protected

Find the element for a point in curved quadratic tetrahedra.

Definition at line 337 of file ComponentFieldMap.cc.

340 {
341 // Backup
342 double jacbak[4][4];
343 double detbak = 1.;
344 double t1bak = 0., t2bak = 0., t3bak = 0., t4bak = 0.;
345 int imapbak = -1;
346
347 // Initial values.
348 t1 = t2 = t3 = t4 = 0.;
349
350 // Check previously used element
351 if (m_lastElement > -1 && !m_checkMultipleElement) {
352 const Element& element = m_elements[m_lastElement];
353 if (Coordinates13(x, y, z, t1, t2, t3, t4, jac, det, element) == 0) {
354 if (t1 >= 0 && t1 <= +1 && t2 >= 0 && t2 <= +1 && t3 >= 0 && t3 <= +1 &&
355 t4 >= 0 && t4 <= +1) {
356 return m_lastElement;
357 }
358 }
359 }
360
361 // Tetra list in the block that contains the input 3D point.
362 std::vector<int> tetList;
363 if (m_useTetrahedralTree && m_octree) {
364 tetList = m_octree->GetElementsInBlock(Vec3(x, y, z));
365 }
366 // Number of elements to scan.
367 // With tetra tree disabled, all elements are scanned.
368 const int numElemToSearch =
369 m_useTetrahedralTree ? tetList.size() : m_elements.size();
370 // Verify the count of volumes that contain the point.
371 int nfound = 0;
372 int imap = -1;
373
374 // Scan all elements
375 for (int i = 0; i < numElemToSearch; i++) {
376 const int idxToElemList = m_useTetrahedralTree ? tetList[i] : i;
377 const Element& element = m_elements[idxToElemList];
378 if (x < element.bbMin[0] || x > element.bbMax[0] ||
379 y < element.bbMin[1] || y > element.bbMax[1] ||
380 z < element.bbMin[2] || z > element.bbMax[2])
381 continue;
382 if (Coordinates13(x, y, z, t1, t2, t3, t4, jac, det, element) != 0) {
383 continue;
384 }
385 if (t1 < 0 || t1 > 1 || t2 < 0 || t2 > 1 || t3 < 0 || t3 > 1 || t4 < 0 ||
386 t4 > 1) {
387 continue;
388 }
389 ++nfound;
390 imap = idxToElemList;
391 m_lastElement = idxToElemList;
392 if (m_debug) {
393 std::cout << m_className << "::FindElement13:\n";
394 std::cout << " Found matching element " << i << ".\n";
395 }
396 if (!m_checkMultipleElement) return idxToElemList;
397 for (int j = 0; j < 4; ++j) {
398 for (int k = 0; k < 4; ++k) jacbak[j][k] = jac[j][k];
399 }
400 detbak = det;
401 t1bak = t1;
402 t2bak = t2;
403 t3bak = t3;
404 t4bak = t4;
405 imapbak = imap;
406 if (m_debug) {
407 PrintElement("FindElement13", x, y, z, t1, t2, t3, t4, element, 10);
408 }
409 }
410
411 // In checking mode, verify the tetrahedron/triangle count.
412 if (m_checkMultipleElement) {
413 if (nfound < 1) {
414 if (m_debug) {
415 std::cout << m_className << "::FindElement13:\n";
416 std::cout << " No element matching point (" << x << ", " << y << ", "
417 << z << ") found.\n";
418 }
419 m_lastElement = -1;
420 return -1;
421 }
422 if (nfound > 1) {
423 std::cerr << m_className << "::FindElement13:\n";
424 std::cerr << " Found << " << nfound << " elements matching point ("
425 << x << ", " << y << ", " << z << ").\n";
426 }
427 for (int j = 0; j < 4; ++j) {
428 for (int k = 0; k < 4; ++k) jac[j][k] = jacbak[j][k];
429 }
430 det = detbak;
431 t1 = t1bak;
432 t2 = t2bak;
433 t3 = t3bak;
434 t4 = t4bak;
435 imap = imapbak;
436 m_lastElement = imap;
437 return imap;
438 }
439
440 if (m_debug) {
441 std::cout << m_className << "::FindElement13:\n";
442 std::cout << " No element matching point (" << x << ", " << y << ", "
443 << z << ") found.\n";
444 }
445 return -1;
446}
std::vector< Element > m_elements
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 Element &element, const unsigned int n, const int iw=-1) const
bool m_debug
Switch on/off debugging messages.
Definition: Component.hh:341
Definition: neBEM.h:155

Referenced by Garfield::ComponentComsol::DelayedWeightingPotential(), Garfield::ComponentAnsys123::ElectricField(), Garfield::ComponentComsol::ElectricField(), Garfield::ComponentElmer::ElectricField(), Garfield::ComponentAnsys123::GetMedium(), Garfield::ComponentComsol::GetMedium(), Garfield::ComponentElmer::GetMedium(), Garfield::ComponentAnsys123::WeightingField(), Garfield::ComponentComsol::WeightingField(), Garfield::ComponentElmer::WeightingField(), Garfield::ComponentAnsys123::WeightingPotential(), Garfield::ComponentComsol::WeightingPotential(), and Garfield::ComponentElmer::WeightingPotential().

◆ FindElement5()

int Garfield::ComponentFieldMap::FindElement5 ( const double  x,
const double  y,
const double  z,
double &  t1,
double &  t2,
double &  t3,
double &  t4,
double  jac[4][4],
double &  det 
)
protected

Find the element for a point in curved quadratic quadrilaterals.

Definition at line 190 of file ComponentFieldMap.cc.

193 {
194 // Check if bounding boxes of elements have been computed
195 if (!m_cacheElemBoundingBoxes) {
196 m_cacheElemBoundingBoxes = true;
197 }
198
199 // Tetra list in the block that contains the input 3D point.
200 std::vector<int> tetList;
201 if (m_useTetrahedralTree && m_octree) {
202 tetList = m_octree->GetElementsInBlock(Vec3(x, y, z));
203 }
204 // Backup
205 double jacbak[4][4], detbak = 1.;
206 double t1bak = 0., t2bak = 0., t3bak = 0., t4bak = 0.;
207 int imapbak = -1;
208
209 // Initial values.
210 t1 = t2 = t3 = t4 = 0;
211
212 // Check previously used element
213 if (m_lastElement > -1 && !m_checkMultipleElement) {
214 const Element& element = m_elements[m_lastElement];
215 if (element.degenerate) {
216 if (Coordinates3(x, y, z, t1, t2, t3, t4, jac, det, element) == 0) {
217 if (t1 >= 0 && t1 <= +1 && t2 >= 0 && t2 <= +1 && t3 >= 0 && t3 <= +1) {
218 return m_lastElement;
219 }
220 }
221 } else {
222 if (Coordinates5(x, y, z, t1, t2, t3, t4, jac, det, element) == 0) {
223 if (t1 >= -1 && t1 <= +1 && t2 >= -1 && t2 <= +1) return m_lastElement;
224 }
225 }
226 }
227
228 // Verify the count of volumes that contain the point.
229 int nfound = 0;
230 int imap = -1;
231
232 // Number of elements to scan.
233 // With tetra tree disabled, all elements are scanned.
234 const int numElemToSearch =
235 m_useTetrahedralTree ? tetList.size() : m_elements.size();
236 for (int i = 0; i < numElemToSearch; ++i) {
237 const int idxToElemList = m_useTetrahedralTree ? tetList[i] : i;
238 const Element& element = m_elements[idxToElemList];
239 if (x < element.bbMin[0] || x > element.bbMax[0] ||
240 y < element.bbMin[1] || y > element.bbMax[1] ||
241 z < element.bbMin[2] || z > element.bbMax[2])
242 continue;
243 if (element.degenerate) {
244 // Degenerate element
245 if (Coordinates3(x, y, z, t1, t2, t3, t4, jac, det, element) != 0) {
246 continue;
247 }
248 if (t1 < 0 || t1 > 1 || t2 < 0 || t2 > 1 || t3 < 0 || t3 > 1) continue;
249 ++nfound;
250 imap = idxToElemList;
251 m_lastElement = idxToElemList;
252 if (m_debug) {
253 std::cout << m_className << "::FindElement5:\n";
254 std::cout << " Found matching degenerate element " << idxToElemList
255 << ".\n";
256 }
257 if (!m_checkMultipleElement) return idxToElemList;
258 for (int j = 0; j < 4; ++j) {
259 for (int k = 0; k < 4; ++k) jacbak[j][k] = jac[j][k];
260 }
261 detbak = det;
262 t1bak = t1;
263 t2bak = t2;
264 t3bak = t3;
265 t4bak = t4;
266 imapbak = imap;
267 if (m_debug) {
268 PrintElement("FindElement5", x, y, z, t1, t2, t3, t4, element, 6);
269 }
270 } else {
271 // Non-degenerate element
272 if (Coordinates5(x, y, z, t1, t2, t3, t4, jac, det, element) != 0) {
273 continue;
274 }
275 if (t1 < -1 || t1 > 1 || t2 < -1 || t2 > 1) continue;
276 ++nfound;
277 imap = idxToElemList;
278 m_lastElement = idxToElemList;
279 if (m_debug) {
280 std::cout << m_className << "::FindElement5:\n";
281 std::cout << " Found matching non-degenerate element "
282 << idxToElemList << ".\n";
283 }
284 if (!m_checkMultipleElement) return idxToElemList;
285 for (int j = 0; j < 4; ++j) {
286 for (int k = 0; k < 4; ++k) jacbak[j][k] = jac[j][k];
287 }
288 detbak = det;
289 t1bak = t1;
290 t2bak = t2;
291 t3bak = t3;
292 t4bak = t4;
293 imapbak = imap;
294 if (m_debug) {
295 PrintElement("FindElement5", x, y, z, t1, t2, t3, t4, element, 8);
296 }
297 }
298 }
299
300 // In checking mode, verify the tetrahedron/triangle count.
301 if (m_checkMultipleElement) {
302 if (nfound < 1) {
303 if (m_debug) {
304 std::cout << m_className << "::FindElement5:\n";
305 std::cout << " No element matching point (" << x << ", " << y
306 << ") found.\n";
307 }
308 m_lastElement = -1;
309 return -1;
310 }
311 if (nfound > 1) {
312 std::cout << m_className << "::FindElement5:\n";
313 std::cout << " Found " << nfound << " elements matching point (" << x
314 << ", " << y << ").\n";
315 }
316 for (int j = 0; j < 4; ++j) {
317 for (int k = 0; k < 4; ++k) jac[j][k] = jacbak[j][k];
318 }
319 det = detbak;
320 t1 = t1bak;
321 t2 = t2bak;
322 t3 = t3bak;
323 t4 = t4bak;
324 imap = imapbak;
325 m_lastElement = imap;
326 return imap;
327 }
328
329 if (m_debug) {
330 std::cout << m_className << "::FindElement5:\n";
331 std::cout << " No element matching point (" << x << ", " << y
332 << ") found.\n";
333 }
334 return -1;
335}

Referenced by Garfield::ComponentAnsys121::ElectricField(), Garfield::ComponentElmer2D::ElectricField(), Garfield::ComponentAnsys121::GetMedium(), Garfield::ComponentElmer2D::GetMedium(), Garfield::ComponentAnsys121::WeightingField(), Garfield::ComponentElmer2D::WeightingField(), Garfield::ComponentAnsys121::WeightingPotential(), and Garfield::ComponentElmer2D::WeightingPotential().

◆ FindElementCube()

int Garfield::ComponentFieldMap::FindElementCube ( const double  x,
const double  y,
const double  z,
double &  t1,
double &  t2,
double &  t3,
TMatrixD *&  jac,
std::vector< TMatrixD * > &  dN 
)
protected

Find the element for a point in a cube.

Definition at line 448 of file ComponentFieldMap.cc.

451 {
452 int imap = -1;
453 if (m_lastElement >= 0) {
454 const Element& element = m_elements[m_lastElement];
455 const Node& n3 = m_nodes[element.emap[3]];
456 if (x >= n3.x && y >= n3.y && z >= n3.z) {
457 const Node& n0 = m_nodes[element.emap[0]];
458 const Node& n2 = m_nodes[element.emap[2]];
459 const Node& n7 = m_nodes[element.emap[7]];
460 if (x < n0.x && y < n2.y && z < n7.z) {
461 imap = m_lastElement;
462 }
463 }
464 }
465
466 // Default element loop
467 if (imap == -1) {
468 const size_t nElements = m_elements.size();
469 for (size_t i = 0; i < nElements; ++i) {
470 const Element& element = m_elements[i];
471 const Node& n3 = m_nodes[element.emap[3]];
472 if (x < n3.x || y < n3.y || z < n3.z) continue;
473 const Node& n0 = m_nodes[element.emap[0]];
474 const Node& n2 = m_nodes[element.emap[2]];
475 const Node& n7 = m_nodes[element.emap[7]];
476 if (x < n0.x && y < n2.y && z < n7.z) {
477 imap = i;
478 break;
479 }
480 }
481 }
482
483 if (imap < 0) {
484 if (m_debug) {
485 std::cout << m_className << "::FindElementCube:\n";
486 std::cout << " Point (" << x << "," << y << "," << z
487 << ") not in the mesh, it is background or PEC.\n";
488 const Node& first0 = m_nodes[m_elements.front().emap[0]];
489 const Node& first2 = m_nodes[m_elements.front().emap[2]];
490 const Node& first3 = m_nodes[m_elements.front().emap[3]];
491 const Node& first7 = m_nodes[m_elements.front().emap[7]];
492 std::cout << " First node (" << first3.x << "," << first3.y << ","
493 << first3.z << ") in the mesh.\n";
494 std::cout << " dx= " << (first0.x - first3.x)
495 << ", dy= " << (first2.y - first3.y)
496 << ", dz= " << (first7.z - first3.z) << "\n";
497 const Node& last0 = m_nodes[m_elements.back().emap[0]];
498 const Node& last2 = m_nodes[m_elements.back().emap[2]];
499 const Node& last3 = m_nodes[m_elements.back().emap[3]];
500 const Node& last5 = m_nodes[m_elements.back().emap[5]];
501 const Node& last7 = m_nodes[m_elements.back().emap[7]];
502 std::cout << " Last node (" << last5.x << "," << last5.y << ","
503 << last5.z << ") in the mesh.\n";
504 std::cout << " dx= " << (last0.x - last3.x)
505 << ", dy= " << (last2.y - last3.y)
506 << ", dz= " << (last7.z - last3.z) << "\n";
507 }
508 return -1;
509 }
510 CoordinatesCube(x, y, z, t1, t2, t3, jac, dN, m_elements[imap]);
511 if (m_debug) {
512 PrintElement("FindElementCube", x, y, z, t1, t2, t3, 0., m_elements[imap],
513 8);
514 }
515 return imap;
516}

◆ GetAspectRatio()

virtual void Garfield::ComponentFieldMap::GetAspectRatio ( const unsigned int  i,
double &  dmin,
double &  dmax 
)
protectedpure virtual

◆ GetBoundingBox()

bool Garfield::ComponentFieldMap::GetBoundingBox ( double &  xmin,
double &  ymin,
double &  zmin,
double &  xmax,
double &  ymax,
double &  zmax 
)
overridevirtual

Get the bounding box coordinates.

Reimplemented from Garfield::Component.

Definition at line 1993 of file ComponentFieldMap.cc.

1995 {
1996 if (!m_ready) return false;
1997
1998 xmin = m_minBoundingBox[0];
1999 xmax = m_maxBoundingBox[0];
2000 ymin = m_minBoundingBox[1];
2001 ymax = m_maxBoundingBox[1];
2002 zmin = m_minBoundingBox[2];
2003 zmax = m_maxBoundingBox[2];
2004 return true;
2005}
std::array< double, 3 > m_maxBoundingBox
std::array< double, 3 > m_minBoundingBox

◆ GetConductivity()

double Garfield::ComponentFieldMap::GetConductivity ( const unsigned int  imat) const

Return the conductivity of a field map material.

Definition at line 87 of file ComponentFieldMap.cc.

87 {
88 if (imat >= m_materials.size()) {
89 std::cerr << m_className << "::GetConductivity: Index out of range.\n";
90 return -1.;
91 }
92
93 return m_materials[imat].ohm;
94}

◆ GetElement() [1/2]

bool Garfield::ComponentFieldMap::GetElement ( const size_t  i,
double &  vol,
double &  dmin,
double &  dmax 
)

Return the volume and aspect ratio of a mesh element.

Definition at line 124 of file ComponentFieldMap.cc.

125 {
126 if (i >= m_elements.size()) {
127 std::cerr << m_className << "::GetElement: Index out of range.\n";
128 return false;
129 }
130
131 vol = GetElementVolume(i);
132 GetAspectRatio(i, dmin, dmax);
133 return true;
134}
virtual double GetElementVolume(const unsigned int i)=0
virtual void GetAspectRatio(const unsigned int i, double &dmin, double &dmax)=0

◆ GetElement() [2/2]

bool Garfield::ComponentFieldMap::GetElement ( const size_t  i,
size_t &  mat,
bool &  drift,
std::vector< size_t > &  nodes 
) const
virtual

Return the material and node indices of a mesh element.

Reimplemented in Garfield::ComponentCST.

Definition at line 136 of file ComponentFieldMap.cc.

137 {
138 if (i >= m_elements.size()) {
139 std::cerr << m_className << "::GetElement: Index out of range.\n";
140 return false;
141 }
142 const auto& element = m_elements[i];
143 mat = element.matmap;
144 drift = m_materials[mat].driftmedium;
145 nodes.resize(4);
146 for (size_t j = 0; j < 4; ++j) nodes[j] = element.emap[j];
147 return true;
148}

◆ GetElementaryCell()

bool Garfield::ComponentFieldMap::GetElementaryCell ( double &  xmin,
double &  ymin,
double &  zmin,
double &  xmax,
double &  ymax,
double &  zmax 
)
overridevirtual

Get the coordinates of the elementary cell.

Reimplemented from Garfield::Component.

Definition at line 2007 of file ComponentFieldMap.cc.

2009 {
2010
2011 if (!m_ready) return false;
2012 xmin = m_mapmin[0];
2013 xmax = m_mapmax[0];
2014 ymin = m_mapmin[1];
2015 ymax = m_mapmax[1];
2016 zmin = m_mapmin[2];
2017 zmax = m_mapmax[2];
2018 return true;
2019}
std::array< double, 3 > m_mapmin
std::array< double, 3 > m_mapmax

◆ GetElementVolume()

virtual double Garfield::ComponentFieldMap::GetElementVolume ( const unsigned int  i)
protectedpure virtual

◆ GetMedium() [1/2]

Medium * Garfield::Component::GetMedium ( const double  x,
const double  y,
const double  z 
)
virtual

Get the medium at a given location (x, y, z).

Reimplemented from Garfield::Component.

Definition at line 28 of file Component.cc.

22 {
23 if (!m_geometry) return nullptr;
24 return m_geometry->GetMedium(x, y, z);
25}
Geometry * m_geometry
Pointer to the geometry.
Definition: Component.hh:332
virtual Medium * GetMedium(const double x, const double y, const double z, const bool tesselated=false) const =0
Retrieve the medium at a given point.

◆ GetMedium() [2/2]

Medium * Garfield::ComponentFieldMap::GetMedium ( const unsigned int  i) const

Return the Medium associated to a field map material.

Definition at line 115 of file ComponentFieldMap.cc.

115 {
116 if (imat >= m_materials.size()) {
117 std::cerr << m_className << "::GetMedium: Index out of range.\n";
118 return nullptr;
119 }
120
121 return m_materials[imat].medium;
122}

◆ GetNode()

bool Garfield::ComponentFieldMap::GetNode ( const size_t  i,
double &  x,
double &  y,
double &  z 
) const
virtual

Reimplemented in Garfield::ComponentCST.

Definition at line 150 of file ComponentFieldMap.cc.

151 {
152 if (i >= m_nodes.size()) {
153 std::cerr << m_className << "::GetNode: Index out of range.\n";
154 return false;
155 }
156 x = m_nodes[i].x;
157 y = m_nodes[i].y;
158 z = m_nodes[i].z;
159 return true;
160}

◆ GetNumberOfElements()

virtual size_t Garfield::ComponentFieldMap::GetNumberOfElements ( ) const
inlinevirtual

Return the number of mesh elements.

Reimplemented in Garfield::ComponentCST.

Definition at line 66 of file ComponentFieldMap.hh.

66{ return m_elements.size(); }

◆ GetNumberOfMaterials()

size_t Garfield::ComponentFieldMap::GetNumberOfMaterials ( ) const
inline

Return the number of materials in the field map.

Definition at line 54 of file ComponentFieldMap.hh.

54{ return m_materials.size(); }

◆ GetNumberOfNodes()

virtual size_t Garfield::ComponentFieldMap::GetNumberOfNodes ( ) const
inlinevirtual

Reimplemented in Garfield::ComponentCST.

Definition at line 72 of file ComponentFieldMap.hh.

72{ return m_nodes.size(); }

◆ GetOrCreateWeightingFieldIndex()

size_t Garfield::ComponentFieldMap::GetOrCreateWeightingFieldIndex ( const std::string &  label)
protected

Definition at line 2357 of file ComponentFieldMap.cc.

2358 {
2359 // Check if a weighting field with the same label already exists.
2360 size_t nWeightingFields = m_wfields.size();
2361 for (size_t i = 0; i < nWeightingFields; ++i) {
2362 if (m_wfields[i] == label) return i;
2363 }
2364 ++nWeightingFields;
2365 m_wfields.resize(nWeightingFields);
2366 m_wfieldsOk.resize(nWeightingFields);
2367 for (auto& node : m_nodes) {
2368 node.w.resize(nWeightingFields);
2369 node.dw.resize(nWeightingFields);
2370 }
2371 m_wfields.back() = label;
2372 return nWeightingFields - 1;
2373}
std::vector< bool > m_wfieldsOk
std::vector< std::string > m_wfields

Referenced by Garfield::ComponentComsol::SetDelayedWeightingPotential(), Garfield::ComponentComsol::SetWeightingField(), Garfield::ComponentAnsys121::SetWeightingField(), Garfield::ComponentAnsys123::SetWeightingField(), Garfield::ComponentElmer::SetWeightingField(), and Garfield::ComponentElmer2D::SetWeightingField().

◆ GetPermittivity()

double Garfield::ComponentFieldMap::GetPermittivity ( const unsigned int  imat) const

Return the permittivity of a field map material.

Definition at line 78 of file ComponentFieldMap.cc.

78 {
79 if (imat >= m_materials.size()) {
80 std::cerr << m_className << "::GetPermittivity: Index out of range.\n";
81 return -1.;
82 }
83
84 return m_materials[imat].eps;
85}

◆ GetVoltageRange()

bool Garfield::ComponentFieldMap::GetVoltageRange ( double &  vmin,
double &  vmax 
)
inlineoverridevirtual

Calculate the voltage range [V].

Implements Garfield::Component.

Definition at line 41 of file ComponentFieldMap.hh.

41 {
42 vmin = m_mapvmin;
43 vmax = m_mapvmax;
44 return true;
45 }

Referenced by Garfield::ViewFEMesh::Plot().

◆ GetWeightingFieldIndex()

◆ IsInBoundingBox()

bool Garfield::ComponentFieldMap::IsInBoundingBox ( const double  x,
const double  y,
const double  z 
) const
inline

Definition at line 30 of file ComponentFieldMap.hh.

30 {
31 return x >= m_minBoundingBox[0] && x <= m_maxBoundingBox[0] &&
32 y >= m_minBoundingBox[1] && y <= m_maxBoundingBox[1] &&
33 z >= m_minBoundingBox[2] && y <= m_maxBoundingBox[2];
34 }

◆ MapCoordinates()

void Garfield::ComponentFieldMap::MapCoordinates ( double &  xpos,
double &  ypos,
double &  zpos,
bool &  xmirrored,
bool &  ymirrored,
bool &  zmirrored,
double &  rcoordinate,
double &  rotation 
) const
protected

Move (xpos, ypos, zpos) to field map coordinates.

Definition at line 2021 of file ComponentFieldMap.cc.

2024 {
2025 // Initial values
2026 rotation = 0;
2027
2028 // If chamber is periodic, reduce to the cell volume.
2029 xmirrored = false;
2030 if (m_periodic[0]) {
2031 const double xrange = m_mapmax[0] - m_mapmin[0];
2032 xpos = m_mapmin[0] + fmod(xpos - m_mapmin[0], xrange);
2033 if (xpos < m_mapmin[0]) xpos += xrange;
2034 } else if (m_mirrorPeriodic[0]) {
2035 const double xrange = m_mapmax[0] - m_mapmin[0];
2036 double xnew = m_mapmin[0] + fmod(xpos - m_mapmin[0], xrange);
2037 if (xnew < m_mapmin[0]) xnew += xrange;
2038 int nx = int(floor(0.5 + (xnew - xpos) / xrange));
2039 if (nx != 2 * (nx / 2)) {
2040 xnew = m_mapmin[0] + m_mapmax[0] - xnew;
2041 xmirrored = true;
2042 }
2043 xpos = xnew;
2044 }
2045 if (m_axiallyPeriodic[0] && (zpos != 0 || ypos != 0)) {
2046 const double auxr = sqrt(zpos * zpos + ypos * ypos);
2047 double auxphi = atan2(zpos, ypos);
2048 const double phirange = m_mapamax[0] - m_mapamin[0];
2049 const double phim = 0.5 * (m_mapamin[0] + m_mapamax[0]);
2050 rotation = phirange * floor(0.5 + (auxphi - phim) / phirange);
2051 if (auxphi - rotation < m_mapamin[0]) rotation -= phirange;
2052 if (auxphi - rotation > m_mapamax[0]) rotation += phirange;
2053 auxphi = auxphi - rotation;
2054 ypos = auxr * cos(auxphi);
2055 zpos = auxr * sin(auxphi);
2056 }
2057
2058 ymirrored = false;
2059 if (m_periodic[1]) {
2060 const double yrange = m_mapmax[1] - m_mapmin[1];
2061 ypos = m_mapmin[1] + fmod(ypos - m_mapmin[1], yrange);
2062 if (ypos < m_mapmin[1]) ypos += yrange;
2063 } else if (m_mirrorPeriodic[1]) {
2064 const double yrange = m_mapmax[1] - m_mapmin[1];
2065 double ynew = m_mapmin[1] + fmod(ypos - m_mapmin[1], yrange);
2066 if (ynew < m_mapmin[1]) ynew += yrange;
2067 int ny = int(floor(0.5 + (ynew - ypos) / yrange));
2068 if (ny != 2 * (ny / 2)) {
2069 ynew = m_mapmin[1] + m_mapmax[1] - ynew;
2070 ymirrored = true;
2071 }
2072 ypos = ynew;
2073 }
2074 if (m_axiallyPeriodic[1] && (xpos != 0 || zpos != 0)) {
2075 const double auxr = sqrt(xpos * xpos + zpos * zpos);
2076 double auxphi = atan2(xpos, zpos);
2077 const double phirange = (m_mapamax[1] - m_mapamin[1]);
2078 const double phim = 0.5 * (m_mapamin[1] + m_mapamax[1]);
2079 rotation = phirange * floor(0.5 + (auxphi - phim) / phirange);
2080 if (auxphi - rotation < m_mapamin[1]) rotation -= phirange;
2081 if (auxphi - rotation > m_mapamax[1]) rotation += phirange;
2082 auxphi = auxphi - rotation;
2083 zpos = auxr * cos(auxphi);
2084 xpos = auxr * sin(auxphi);
2085 }
2086
2087 zmirrored = false;
2088 if (m_periodic[2]) {
2089 const double zrange = m_mapmax[2] - m_mapmin[2];
2090 zpos = m_mapmin[2] + fmod(zpos - m_mapmin[2], zrange);
2091 if (zpos < m_mapmin[2]) zpos += zrange;
2092 } else if (m_mirrorPeriodic[2]) {
2093 const double zrange = m_mapmax[2] - m_mapmin[2];
2094 double znew = m_mapmin[2] + fmod(zpos - m_mapmin[2], zrange);
2095 if (znew < m_mapmin[2]) znew += zrange;
2096 int nz = int(floor(0.5 + (znew - zpos) / zrange));
2097 if (nz != 2 * (nz / 2)) {
2098 znew = m_mapmin[2] + m_mapmax[2] - znew;
2099 zmirrored = true;
2100 }
2101 zpos = znew;
2102 }
2103 if (m_axiallyPeriodic[2] && (ypos != 0 || xpos != 0)) {
2104 const double auxr = sqrt(ypos * ypos + xpos * xpos);
2105 double auxphi = atan2(ypos, xpos);
2106 const double phirange = m_mapamax[2] - m_mapamin[2];
2107 const double phim = 0.5 * (m_mapamin[2] + m_mapamax[2]);
2108 rotation = phirange * floor(0.5 + (auxphi - phim) / phirange);
2109 if (auxphi - rotation < m_mapamin[2]) rotation -= phirange;
2110 if (auxphi - rotation > m_mapamax[2]) rotation += phirange;
2111 auxphi = auxphi - rotation;
2112 xpos = auxr * cos(auxphi);
2113 ypos = auxr * sin(auxphi);
2114 }
2115
2116 // If we have a rotationally symmetric field map, store coordinates.
2117 rcoordinate = 0;
2118 double zcoordinate = 0;
2119 if (m_rotationSymmetric[0]) {
2120 rcoordinate = sqrt(ypos * ypos + zpos * zpos);
2121 zcoordinate = xpos;
2122 } else if (m_rotationSymmetric[1]) {
2123 rcoordinate = sqrt(xpos * xpos + zpos * zpos);
2124 zcoordinate = ypos;
2125 } else if (m_rotationSymmetric[2]) {
2126 rcoordinate = sqrt(xpos * xpos + ypos * ypos);
2127 zcoordinate = zpos;
2128 }
2129
2132 xpos = rcoordinate;
2133 ypos = zcoordinate;
2134 zpos = 0;
2135 }
2136}
std::array< double, 3 > m_mapamin
std::array< double, 3 > m_mapamax
std::array< bool, 3 > m_rotationSymmetric
Rotation symmetry around x-axis, y-axis, z-axis.
Definition: Component.hh:350
std::array< bool, 3 > m_mirrorPeriodic
Mirror periodicity in x, y, z.
Definition: Component.hh:346
std::array< bool, 3 > m_periodic
Simple periodicity in x, y, z.
Definition: Component.hh:344
std::array< bool, 3 > m_axiallyPeriodic
Axial periodicity in x, y, z.
Definition: Component.hh:348
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:432
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:384
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314

Referenced by Garfield::ComponentCST::Coordinate2Index(), Garfield::ComponentComsol::DelayedWeightingPotential(), Garfield::ComponentAnsys121::ElectricField(), Garfield::ComponentAnsys123::ElectricField(), Garfield::ComponentComsol::ElectricField(), Garfield::ComponentElmer::ElectricField(), Garfield::ComponentElmer2D::ElectricField(), Garfield::ComponentAnsys121::GetMedium(), Garfield::ComponentAnsys123::GetMedium(), Garfield::ComponentComsol::GetMedium(), Garfield::ComponentElmer::GetMedium(), Garfield::ComponentElmer2D::GetMedium(), Garfield::ComponentAnsys121::WeightingField(), Garfield::ComponentAnsys123::WeightingField(), Garfield::ComponentComsol::WeightingField(), Garfield::ComponentElmer::WeightingField(), Garfield::ComponentElmer2D::WeightingField(), Garfield::ComponentAnsys121::WeightingPotential(), Garfield::ComponentAnsys123::WeightingPotential(), Garfield::ComponentComsol::WeightingPotential(), Garfield::ComponentElmer::WeightingPotential(), and Garfield::ComponentElmer2D::WeightingPotential().

◆ NotDriftMedium()

void Garfield::ComponentFieldMap::NotDriftMedium ( const unsigned int  imat)

Flag a field map materials as a non-drift medium.

Definition at line 64 of file ComponentFieldMap.cc.

64 {
65 // Do not proceed if not properly initialised.
66 if (!m_ready) PrintNotReady("NotDriftMedium");
67
68 // Check value
69 if (imat >= m_materials.size()) {
70 std::cerr << m_className << "::NotDriftMedium: Index out of range.\n";
71 return;
72 }
73
74 // Make drift medium
75 m_materials[imat].driftmedium = false;
76}

◆ Prepare()

void Garfield::ComponentFieldMap::Prepare ( )
protected

Definition at line 1713 of file ComponentFieldMap.cc.

1713 {
1714
1715 // Establish the ranges.
1716 SetRange();
1718 std::cout << m_className << "::Prepare:\n"
1719 << " Caching the bounding boxes of all elements...";
1720 CalculateElementBoundingBoxes();
1721 std::cout << " done.\n";
1722 // Initialize the tetrahedral tree.
1723 InitializeTetrahedralTree();
1724}
virtual void SetRange()
Calculate x, y, z, V and angular ranges.
virtual void UpdatePeriodicity()=0
Verify periodicities.

Referenced by Garfield::ComponentElmer::Initialise(), Garfield::ComponentElmer2D::Initialise(), Garfield::ComponentComsol::Initialise(), Garfield::ComponentCST::Initialise(), Garfield::ComponentAnsys121::Initialise(), and Garfield::ComponentAnsys123::Initialise().

◆ PrintCouldNotOpen()

void Garfield::ComponentFieldMap::PrintCouldNotOpen ( const std::string &  header,
const std::string &  filename 
) const
protected

◆ PrintElement()

void Garfield::ComponentFieldMap::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 Element element,
const unsigned int  n,
const int  iw = -1 
) const
protected

Definition at line 2395 of file ComponentFieldMap.cc.

2400 {
2401 std::cout << m_className << "::" << header << ":\n"
2402 << " Global = (" << x << ", " << y << ", " << z << ")\n"
2403 << " Local = (" << t1 << ", " << t2 << ", " << t3 << ", " << t4
2404 << ")\n";
2405 if (element.degenerate) std::cout << " Element is degenerate.\n";
2406 std::cout << " Node x y z V\n";
2407 for (unsigned int ii = 0; ii < n; ++ii) {
2408 const Node& node = m_nodes[element.emap[ii]];
2409 const double v = iw < 0 ? node.v : node.w[iw];
2410 printf(" %-5d %12g %12g %12g %12g\n", element.emap[ii], node.x, node.y,
2411 node.z, v);
2412 }
2413}

Referenced by Garfield::ComponentComsol::DelayedWeightingPotential(), Garfield::ComponentAnsys121::ElectricField(), Garfield::ComponentAnsys123::ElectricField(), Garfield::ComponentComsol::ElectricField(), Garfield::ComponentElmer::ElectricField(), Garfield::ComponentElmer2D::ElectricField(), FindElement13(), FindElement5(), FindElementCube(), Garfield::ComponentAnsys121::GetMedium(), Garfield::ComponentAnsys123::GetMedium(), Garfield::ComponentComsol::GetMedium(), Garfield::ComponentElmer::GetMedium(), Garfield::ComponentElmer2D::GetMedium(), Garfield::ComponentAnsys121::WeightingField(), Garfield::ComponentAnsys123::WeightingField(), Garfield::ComponentComsol::WeightingField(), Garfield::ComponentElmer::WeightingField(), Garfield::ComponentElmer2D::WeightingField(), Garfield::ComponentAnsys121::WeightingPotential(), Garfield::ComponentAnsys123::WeightingPotential(), Garfield::ComponentComsol::WeightingPotential(), Garfield::ComponentElmer::WeightingPotential(), and Garfield::ComponentElmer2D::WeightingPotential().

◆ PrintMaterials()

void Garfield::ComponentFieldMap::PrintMaterials ( )

List all currently defined materials.

Definition at line 20 of file ComponentFieldMap.cc.

20 {
21 // Do not proceed if not properly initialised.
22 if (!m_ready) PrintNotReady("PrintMaterials");
23
24 if (m_materials.empty()) {
25 std::cerr << m_className << "::PrintMaterials:\n"
26 << " No materials are currently defined.\n";
27 return;
28 }
29
30 const size_t nMaterials = m_materials.size();
31 std::cout << m_className << "::PrintMaterials:\n"
32 << " Currently " << nMaterials << " materials are defined.\n"
33 << " Index Permittivity Resistivity Notes\n";
34 for (size_t i = 0; i < nMaterials; ++i) {
35 printf(" %5zu %12g %12g", i, m_materials[i].eps, m_materials[i].ohm);
36 if (m_materials[i].medium) {
37 std::string name = m_materials[i].medium->GetName();
38 std::cout << " " << name;
39 if (m_materials[i].medium->IsDriftable()) std::cout << ", drift medium";
40 if (m_materials[i].medium->IsIonisable()) std::cout << ", ionisable";
41 }
42 if (m_materials[i].driftmedium) {
43 std::cout << " (drift medium)\n";
44 } else {
45 std::cout << "\n";
46 }
47 }
48}

Referenced by Garfield::ComponentCST::Initialise(), Garfield::ComponentAnsys121::Initialise(), and Garfield::ComponentAnsys123::Initialise().

◆ PrintNotReady()

◆ PrintRange()

void Garfield::ComponentFieldMap::PrintRange ( )

Show x, y, z, V and angular ranges.

Definition at line 1964 of file ComponentFieldMap.cc.

1964 {
1965 std::cout << m_className << "::PrintRange:\n";
1966 std::cout << " Dimensions of the elementary block\n";
1967 printf(" %15g < x < %-15g cm,\n", m_mapmin[0], m_mapmax[0]);
1968 printf(" %15g < y < %-15g cm,\n", m_mapmin[1], m_mapmax[1]);
1969 printf(" %15g < z < %-15g cm,\n", m_mapmin[2], m_mapmax[2]);
1970 printf(" %15g < V < %-15g V.\n", m_mapvmin, m_mapvmax);
1971
1972 std::cout << " Periodicities\n";
1973 const std::array<std::string, 3> axes = {{"x", "y", "z"}};
1974 for (unsigned int i = 0; i < 3; ++i) {
1975 std::cout << " " << axes[i] << ":";
1976 if (m_periodic[i]) {
1977 std::cout << " simple with length " << m_cells[i] << " cm";
1978 }
1979 if (m_mirrorPeriodic[i]) {
1980 std::cout << " mirror with length " << m_cells[i] << " cm";
1981 }
1982 if (m_axiallyPeriodic[i]) {
1983 std::cout << " axial " << int(0.5 + m_mapna[i]) << "-fold repetition";
1984 }
1985 if (m_rotationSymmetric[i]) std::cout << " rotational symmetry";
1986 if (!(m_periodic[i] || m_mirrorPeriodic[i] || m_axiallyPeriodic[i] ||
1988 std::cout << " none";
1989 std::cout << "\n";
1990 }
1991}
std::array< double, 3 > m_mapna
std::array< double, 3 > m_cells

Referenced by SetRange(), and UpdatePeriodicityCommon().

◆ PrintWarning()

void Garfield::ComponentFieldMap::PrintWarning ( const std::string &  header)
protected

◆ ReadDouble()

◆ ReadInteger()

int Garfield::ComponentFieldMap::ReadInteger ( char *  token,
int  def,
bool &  error 
)
protected

◆ Reset()

void Garfield::ComponentFieldMap::Reset ( )
overrideprotectedvirtual

Reset the component.

Implements Garfield::Component.

Definition at line 1695 of file ComponentFieldMap.cc.

1695 {
1696
1697 m_ready = false;
1698
1699 m_elements.clear();
1700 m_nodes.clear();
1701 m_materials.clear();
1702 m_wfields.clear();
1703 m_wfieldsOk.clear();
1704 m_hasBoundingBox = false;
1705 m_warning = false;
1706 m_nWarnings = 0;
1707
1708 m_octree.reset(nullptr);
1709 m_cacheElemBoundingBoxes = false;
1710 m_lastElement = -1;
1711}

Referenced by Garfield::ComponentElmer::Initialise(), Garfield::ComponentElmer2D::Initialise(), Garfield::ComponentComsol::Initialise(), Garfield::ComponentCST::Initialise(), Garfield::ComponentAnsys121::Initialise(), and Garfield::ComponentAnsys123::Initialise().

◆ ScalingFactor()

double Garfield::ComponentFieldMap::ScalingFactor ( std::string  unit)
staticprotected

Definition at line 2213 of file ComponentFieldMap.cc.

2213 {
2214 std::transform(unit.begin(), unit.end(), unit.begin(), toupper);
2215 if (unit == "MUM" || unit == "MICRON" || unit == "MICROMETER") {
2216 return 0.0001;
2217 } else if (unit == "MM" || unit == "MILLIMETER") {
2218 return 0.1;
2219 } else if (unit == "CM" || unit == "CENTIMETER") {
2220 return 1.0;
2221 } else if (unit == "M" || unit == "METER") {
2222 return 100.0;
2223 }
2224 return -1.;
2225}

Referenced by Garfield::ComponentElmer::Initialise(), Garfield::ComponentElmer2D::Initialise(), Garfield::ComponentComsol::Initialise(), Garfield::ComponentCST::Initialise(), Garfield::ComponentAnsys121::Initialise(), and Garfield::ComponentAnsys123::Initialise().

◆ SetDefaultDriftMedium()

bool Garfield::ComponentFieldMap::SetDefaultDriftMedium ( )
protected

Find lowest epsilon, check for eps = 0, set default drift media flags.

Definition at line 162 of file ComponentFieldMap.cc.

162 {
163
164 // Find lowest epsilon and set drift medium flags.
165 const size_t nMaterials = m_materials.size();
166 double epsmin = -1;
167 size_t iepsmin = 0;
168 for (size_t i = 0; i < nMaterials; ++i) {
169 m_materials[i].driftmedium = false;
170 if (m_materials[i].eps < 0) continue;
171 // Check for eps == 0.
172 if (m_materials[i].eps == 0) {
173 std::cerr << m_className << "::SetDefaultDriftMedium:\n"
174 << " Material " << i << " has zero permittivity.\n";
175 m_materials[i].eps = -1.;
176 } else if (epsmin < 0. || epsmin > m_materials[i].eps) {
177 epsmin = m_materials[i].eps;
178 iepsmin = i;
179 }
180 }
181 if (epsmin < 0.) {
182 std::cerr << m_className << "::SetDefaultDriftMedium:\n"
183 << " Found no material with positive permittivity.\n";
184 return false;
185 }
186 m_materials[iepsmin].driftmedium = true;
187 return true;
188}

Referenced by Garfield::ComponentElmer::Initialise(), Garfield::ComponentElmer2D::Initialise(), Garfield::ComponentComsol::Initialise(), Garfield::ComponentCST::Initialise(), Garfield::ComponentAnsys121::Initialise(), and Garfield::ComponentAnsys123::Initialise().

◆ SetMedium()

void Garfield::ComponentFieldMap::SetMedium ( const unsigned int  imat,
Medium medium 
)

Associate a field map material with a Medium class.

Definition at line 96 of file ComponentFieldMap.cc.

96 {
97 if (imat >= m_materials.size()) {
98 std::cerr << m_className << "::SetMedium: Index out of range.\n";
99 return;
100 }
101
102 if (!m) {
103 std::cerr << m_className << "::SetMedium: Null pointer.\n";
104 return;
105 }
106
107 if (m_debug) {
108 std::cout << m_className << "::SetMedium:\n Associated material " << imat
109 << " with medium " << m->GetName() << ".\n";
110 }
111
112 m_materials[imat].medium = m;
113}

Referenced by main().

◆ SetRange()

void Garfield::ComponentFieldMap::SetRange ( )
virtual

Calculate x, y, z, V and angular ranges.

Reimplemented in Garfield::ComponentCST.

Definition at line 1873 of file ComponentFieldMap.cc.

1873 {
1874 // Initial values
1875 m_mapmin.fill(0.);
1876 m_mapmax.fill(0.);
1877 m_mapamin.fill(0.);
1878 m_mapamax.fill(0.);
1879 m_mapvmin = m_mapvmax = 0.;
1880 m_setang.fill(false);
1881
1882 // Make sure the required data is available.
1883 if (!m_ready || m_nodes.empty()) {
1884 std::cerr << m_className << "::SetRange: Field map not yet set.\n";
1885 return;
1886 }
1887
1888 // Loop over the nodes.
1889 m_mapmin[0] = m_mapmax[0] = m_nodes[0].x;
1890 m_mapmin[1] = m_mapmax[1] = m_nodes[0].y;
1891 m_mapmin[2] = m_mapmax[2] = m_nodes[0].z;
1892 m_mapvmin = m_mapvmax = m_nodes[0].v;
1893
1894 for (const auto& node : m_nodes) {
1895 const std::array<double, 3> pos = {{node.x, node.y, node.z}};
1896 for (unsigned int i = 0; i < 3; ++i) {
1897 m_mapmin[i] = std::min(m_mapmin[i], pos[i]);
1898 m_mapmax[i] = std::max(m_mapmax[i], pos[i]);
1899 }
1900 m_mapvmin = std::min(m_mapvmin, node.v);
1901 m_mapvmax = std::max(m_mapvmax, node.v);
1902
1903 if (node.y != 0 || node.z != 0) {
1904 const double ang = atan2(node.z, node.y);
1905 if (m_setang[0]) {
1906 m_mapamin[0] = std::min(m_mapamin[0], ang);
1907 m_mapamax[0] = std::max(m_mapamax[0], ang);
1908 } else {
1909 m_mapamin[0] = m_mapamax[0] = ang;
1910 m_setang[0] = true;
1911 }
1912 }
1913
1914 if (node.z != 0 || node.x != 0) {
1915 const double ang = atan2(node.x, node.z);
1916 if (m_setang[1]) {
1917 m_mapamin[1] = std::min(m_mapamin[1], ang);
1918 m_mapamax[1] = std::max(m_mapamax[1], ang);
1919 } else {
1920 m_mapamin[1] = m_mapamax[1] = ang;
1921 m_setang[1] = true;
1922 }
1923 }
1924
1925 if (node.x != 0 || node.y != 0) {
1926 const double ang = atan2(node.y, node.x);
1927 if (m_setang[2]) {
1928 m_mapamin[2] = std::min(m_mapamin[2], ang);
1929 m_mapamax[2] = std::max(m_mapamax[2], ang);
1930 } else {
1931 m_mapamin[2] = m_mapamax[2] = ang;
1932 m_setang[2] = true;
1933 }
1934 }
1935 }
1936
1937 // Fix the angular ranges.
1938 for (unsigned int i = 0; i < 3; ++i) {
1939 if (m_mapamax[i] - m_mapamin[i] > Pi) {
1940 const double aux = m_mapamin[i];
1941 m_mapamin[i] = m_mapamax[i];
1942 m_mapamax[i] = aux + TwoPi;
1943 }
1944 }
1945
1946 // Set provisional cell dimensions.
1947 m_minBoundingBox[0] = m_mapmin[0];
1948 m_maxBoundingBox[0] = m_mapmax[0];
1949 m_minBoundingBox[1] = m_mapmin[1];
1950 m_maxBoundingBox[1] = m_mapmax[1];
1951 if (m_is3d) {
1952 m_minBoundingBox[2] = m_mapmin[2];
1953 m_maxBoundingBox[2] = m_mapmax[2];
1954 } else {
1955 m_mapmin[2] = m_minBoundingBox[2];
1956 m_mapmax[2] = m_maxBoundingBox[2];
1957 }
1958 m_hasBoundingBox = true;
1959
1960 // Display the range if requested.
1961 if (m_debug) PrintRange();
1962}
void PrintRange()
Show x, y, z, V and angular ranges.
std::array< bool, 3 > m_setang

Referenced by Prepare().

◆ UnmapFields()

void Garfield::ComponentFieldMap::UnmapFields ( double &  ex,
double &  ey,
double &  ez,
double &  xpos,
double &  ypos,
double &  zpos,
bool &  xmirrored,
bool &  ymirrored,
bool &  zmirrored,
double &  rcoordinate,
double &  rotation 
) const
protected

Move (ex, ey, ez) to global coordinates.

Definition at line 2138 of file ComponentFieldMap.cc.

2142 {
2143 // Apply mirror imaging.
2144 if (xmirrored) ex = -ex;
2145 if (ymirrored) ey = -ey;
2146 if (zmirrored) ez = -ez;
2147
2148 // Rotate the field.
2149 double er, theta;
2150 if (m_axiallyPeriodic[0]) {
2151 er = sqrt(ey * ey + ez * ez);
2152 theta = atan2(ez, ey);
2153 theta += rotation;
2154 ey = er * cos(theta);
2155 ez = er * sin(theta);
2156 }
2157 if (m_axiallyPeriodic[1]) {
2158 er = sqrt(ez * ez + ex * ex);
2159 theta = atan2(ex, ez);
2160 theta += rotation;
2161 ez = er * cos(theta);
2162 ex = er * sin(theta);
2163 }
2164 if (m_axiallyPeriodic[2]) {
2165 er = sqrt(ex * ex + ey * ey);
2166 theta = atan2(ey, ex);
2167 theta += rotation;
2168 ex = er * cos(theta);
2169 ey = er * sin(theta);
2170 }
2171
2172 // Take care of symmetry.
2173 double eaxis;
2174 er = ex;
2175 eaxis = ey;
2176
2177 // Rotational symmetry
2178 if (m_rotationSymmetric[0]) {
2179 if (rcoordinate <= 0) {
2180 ex = eaxis;
2181 ey = 0;
2182 ez = 0;
2183 } else {
2184 ex = eaxis;
2185 ey = er * ypos / rcoordinate;
2186 ez = er * zpos / rcoordinate;
2187 }
2188 }
2189 if (m_rotationSymmetric[1]) {
2190 if (rcoordinate <= 0) {
2191 ex = 0;
2192 ey = eaxis;
2193 ez = 0;
2194 } else {
2195 ex = er * xpos / rcoordinate;
2196 ey = eaxis;
2197 ez = er * zpos / rcoordinate;
2198 }
2199 }
2200 if (m_rotationSymmetric[2]) {
2201 if (rcoordinate <= 0) {
2202 ex = 0;
2203 ey = 0;
2204 ez = eaxis;
2205 } else {
2206 ex = er * xpos / rcoordinate;
2207 ey = er * ypos / rcoordinate;
2208 ez = eaxis;
2209 }
2210 }
2211}

Referenced by Garfield::ComponentAnsys121::ElectricField(), Garfield::ComponentAnsys123::ElectricField(), Garfield::ComponentComsol::ElectricField(), Garfield::ComponentElmer::ElectricField(), Garfield::ComponentElmer2D::ElectricField(), Garfield::ComponentAnsys121::WeightingField(), Garfield::ComponentAnsys123::WeightingField(), Garfield::ComponentComsol::WeightingField(), Garfield::ComponentElmer::WeightingField(), and Garfield::ComponentElmer2D::WeightingField().

◆ UpdatePeriodicity2d()

void Garfield::ComponentFieldMap::UpdatePeriodicity2d ( )
protected

Definition at line 1845 of file ComponentFieldMap.cc.

1845 {
1846 // Check the required data is available.
1847 if (!m_ready) {
1848 PrintNotReady("UpdatePeriodicity2d");
1849 return;
1850 }
1851
1852 // No z-periodicity in 2d
1853 if (m_periodic[2] || m_mirrorPeriodic[2]) {
1854 std::cerr << m_className << "::UpdatePeriodicity2d:\n"
1855 << " Simple or mirror periodicity along z\n"
1856 << " requested for a 2d map; reset.\n";
1857 m_periodic[2] = false;
1858 m_mirrorPeriodic[2] = false;
1859 m_warning = true;
1860 }
1861
1862 // Only z-axial periodicity in 2d maps
1863 if (m_axiallyPeriodic[0] || m_axiallyPeriodic[1]) {
1864 std::cerr << m_className << "::UpdatePeriodicity2d:\n"
1865 << " Axial symmetry has been requested \n"
1866 << " around x or y for a 2d map; reset.\n";
1867 m_axiallyPeriodic[0] = false;
1868 m_axiallyPeriodic[1] = false;
1869 m_warning = true;
1870 }
1871}

Referenced by Garfield::ComponentAnsys121::UpdatePeriodicity(), Garfield::ComponentCST::UpdatePeriodicity(), and Garfield::ComponentElmer2D::UpdatePeriodicity().

◆ UpdatePeriodicityCommon()

void Garfield::ComponentFieldMap::UpdatePeriodicityCommon ( )
protected

Definition at line 1726 of file ComponentFieldMap.cc.

1726 {
1727 // Check the required data is available.
1728 if (!m_ready) {
1729 PrintNotReady("UpdatePeriodicityCommon");
1730 return;
1731 }
1732
1733 for (size_t i = 0; i < 3; ++i) {
1734 // No regular and mirror periodicity at the same time.
1735 if (m_periodic[i] && m_mirrorPeriodic[i]) {
1736 std::cerr << m_className << "::UpdatePeriodicityCommon:\n"
1737 << " Both simple and mirror periodicity requested. Reset.\n";
1738 m_periodic[i] = false;
1739 m_mirrorPeriodic[i] = false;
1740 m_warning = true;
1741 }
1742 // In case of axial periodicity,
1743 // the range must be an integral part of two pi.
1744 if (m_axiallyPeriodic[i]) {
1745 if (m_mapamin[i] >= m_mapamax[i]) {
1746 m_mapna[i] = 0;
1747 } else {
1748 m_mapna[i] = TwoPi / (m_mapamax[i] - m_mapamin[i]);
1749 }
1750 if (fabs(m_mapna[i] - int(0.5 + m_mapna[i])) > 0.001 ||
1751 m_mapna[i] < 1.5) {
1752 std::cerr << m_className << "::UpdatePeriodicityCommon:\n"
1753 << " Axial symmetry has been requested but map does not\n"
1754 << " cover an integral fraction of 2 pi. Reset.\n";
1755 m_axiallyPeriodic[i] = false;
1756 m_warning = true;
1757 }
1758 }
1759 }
1760
1761 // Not more than 1 rotational symmetry
1765 std::cerr << m_className << "::UpdatePeriodicityCommon:\n"
1766 << " Only one rotational symmetry allowed; reset.\n";
1767 m_rotationSymmetric.fill(false);
1768 m_warning = true;
1769 }
1770
1771 // No rotational symmetry as well as axial periodicity
1773 m_rotationSymmetric[2]) &&
1775 std::cerr << m_className << "::UpdatePeriodicityCommon:\n"
1776 << " Not allowed to combine rotational symmetry\n"
1777 << " and axial periodicity; reset.\n";
1778 m_axiallyPeriodic.fill(false);
1779 m_rotationSymmetric.fill(false);
1780 m_warning = true;
1781 }
1782
1783 // In case of rotational symmetry, the x-range should not straddle 0.
1786 if (m_mapmin[0] * m_mapmax[0] < 0) {
1787 std::cerr << m_className << "::UpdatePeriodicityCommon:\n"
1788 << " Rotational symmetry requested, \n"
1789 << " but x-range straddles 0; reset.\n";
1790 m_rotationSymmetric.fill(false);
1791 m_warning = true;
1792 }
1793 }
1794
1795 // Recompute the cell ranges.
1796 for (size_t i = 0; i < 3; ++i) {
1797 m_minBoundingBox[i] = m_mapmin[i];
1798 m_maxBoundingBox[i] = m_mapmax[i];
1799 m_cells[i] = fabs(m_mapmax[i] - m_mapmin[i]);
1800 }
1801 for (size_t i = 0; i < 3; ++i) {
1802 if (!m_rotationSymmetric[i]) continue;
1803 const double r = std::max(fabs(m_mapmin[0]), fabs(m_mapmax[0]));
1804 m_minBoundingBox.fill(-r);
1805 m_maxBoundingBox.fill(+r);
1806 m_minBoundingBox[i] = m_mapmin[1];
1807 m_maxBoundingBox[i] = m_mapmax[1];
1808 break;
1809 }
1810
1811 if (m_axiallyPeriodic[0]) {
1812 const double yzmax = std::max({fabs(m_mapmin[1]), fabs(m_mapmax[1]),
1813 fabs(m_mapmin[2]), fabs(m_mapmax[2])});
1814 m_minBoundingBox[1] = -yzmax;
1815 m_maxBoundingBox[1] = +yzmax;
1816 m_minBoundingBox[2] = -yzmax;
1817 m_maxBoundingBox[2] = +yzmax;
1818 } else if (m_axiallyPeriodic[1]) {
1819 const double xzmax = std::max({fabs(m_mapmin[0]), fabs(m_mapmax[0]),
1820 fabs(m_mapmin[2]), fabs(m_mapmax[2])});
1821 m_minBoundingBox[0] = -xzmax;
1822 m_maxBoundingBox[0] = +xzmax;
1823 m_minBoundingBox[2] = -xzmax;
1824 m_maxBoundingBox[2] = +xzmax;
1825 } else if (m_axiallyPeriodic[2]) {
1826 const double xymax = std::max({fabs(m_mapmin[0]), fabs(m_mapmax[0]),
1827 fabs(m_mapmin[1]), fabs(m_mapmax[1])});
1828 m_minBoundingBox[0] = -xymax;
1829 m_maxBoundingBox[0] = +xymax;
1830 m_minBoundingBox[1] = -xymax;
1831 m_maxBoundingBox[1] = +xymax;
1832 }
1833
1834 for (size_t i = 0; i < 3; ++i) {
1835 if (m_periodic[i] || m_mirrorPeriodic[i]) {
1836 m_minBoundingBox[i] = -INFINITY;
1837 m_maxBoundingBox[i] = +INFINITY;
1838 }
1839 }
1840
1841 // Display the range if requested.
1842 if (m_debug) PrintRange();
1843}
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615

Referenced by Garfield::ComponentAnsys121::UpdatePeriodicity(), Garfield::ComponentAnsys123::UpdatePeriodicity(), Garfield::ComponentComsol::UpdatePeriodicity(), Garfield::ComponentCST::UpdatePeriodicity(), Garfield::ComponentElmer::UpdatePeriodicity(), and Garfield::ComponentElmer2D::UpdatePeriodicity().

Friends And Related Function Documentation

◆ ViewFEMesh

friend class ViewFEMesh
friend

Definition at line 97 of file ComponentFieldMap.hh.

Member Data Documentation

◆ m_cells

std::array<double, 3> Garfield::ComponentFieldMap::m_cells
protected

Definition at line 157 of file ComponentFieldMap.hh.

Referenced by PrintRange(), and UpdatePeriodicityCommon().

◆ m_deleteBackground

bool Garfield::ComponentFieldMap::m_deleteBackground = true
protected

◆ m_dwfieldsOk

std::vector<bool> Garfield::ComponentFieldMap::m_dwfieldsOk
protected

Definition at line 142 of file ComponentFieldMap.hh.

◆ m_elements

std::vector<Element> Garfield::ComponentFieldMap::m_elements
protected

Definition at line 113 of file ComponentFieldMap.hh.

Referenced by Garfield::ComponentComsol::DelayedWeightingPotential(), Garfield::ComponentAnsys121::ElectricField(), Garfield::ComponentAnsys123::ElectricField(), Garfield::ComponentComsol::ElectricField(), Garfield::ComponentElmer::ElectricField(), Garfield::ComponentElmer2D::ElectricField(), FindElement13(), FindElement5(), FindElementCube(), Garfield::ComponentAnsys121::GetAspectRatio(), Garfield::ComponentAnsys123::GetAspectRatio(), Garfield::ComponentComsol::GetAspectRatio(), Garfield::ComponentElmer::GetAspectRatio(), Garfield::ComponentElmer2D::GetAspectRatio(), GetElement(), Garfield::ComponentAnsys121::GetElementVolume(), Garfield::ComponentAnsys123::GetElementVolume(), Garfield::ComponentComsol::GetElementVolume(), Garfield::ComponentElmer::GetElementVolume(), Garfield::ComponentElmer2D::GetElementVolume(), Garfield::ComponentAnsys121::GetMedium(), Garfield::ComponentAnsys123::GetMedium(), Garfield::ComponentComsol::GetMedium(), Garfield::ComponentElmer::GetMedium(), Garfield::ComponentElmer2D::GetMedium(), GetNumberOfElements(), Garfield::ComponentElmer::Initialise(), Garfield::ComponentElmer2D::Initialise(), Garfield::ComponentComsol::Initialise(), Garfield::ComponentAnsys121::Initialise(), Garfield::ComponentAnsys123::Initialise(), Reset(), Garfield::ComponentAnsys121::WeightingField(), Garfield::ComponentAnsys123::WeightingField(), Garfield::ComponentComsol::WeightingField(), Garfield::ComponentElmer::WeightingField(), Garfield::ComponentElmer2D::WeightingField(), Garfield::ComponentAnsys121::WeightingPotential(), Garfield::ComponentAnsys123::WeightingPotential(), Garfield::ComponentComsol::WeightingPotential(), Garfield::ComponentElmer::WeightingPotential(), and Garfield::ComponentElmer2D::WeightingPotential().

◆ m_hasBoundingBox

bool Garfield::ComponentFieldMap::m_hasBoundingBox = false
protected

Definition at line 147 of file ComponentFieldMap.hh.

Referenced by Reset(), SetRange(), and Garfield::ComponentCST::SetRange().

◆ m_is3d

bool Garfield::ComponentFieldMap::m_is3d = true
protected

◆ m_mapamax

std::array<double, 3> Garfield::ComponentFieldMap::m_mapamax
protected

Definition at line 155 of file ComponentFieldMap.hh.

Referenced by MapCoordinates(), SetRange(), and UpdatePeriodicityCommon().

◆ m_mapamin

std::array<double, 3> Garfield::ComponentFieldMap::m_mapamin
protected

Definition at line 154 of file ComponentFieldMap.hh.

Referenced by MapCoordinates(), SetRange(), and UpdatePeriodicityCommon().

◆ m_mapmax

◆ m_mapmin

◆ m_mapna

std::array<double, 3> Garfield::ComponentFieldMap::m_mapna
protected

Definition at line 156 of file ComponentFieldMap.hh.

Referenced by PrintRange(), and UpdatePeriodicityCommon().

◆ m_mapvmax

double Garfield::ComponentFieldMap::m_mapvmax = 0.
protected

◆ m_mapvmin

double Garfield::ComponentFieldMap::m_mapvmin = 0.
protected

◆ m_materials

◆ m_maxBoundingBox

◆ m_minBoundingBox

◆ m_nodes

std::vector<Node> Garfield::ComponentFieldMap::m_nodes
protected

Definition at line 126 of file ComponentFieldMap.hh.

Referenced by Garfield::ComponentComsol::DelayedWeightingPotential(), Garfield::ComponentAnsys121::ElectricField(), Garfield::ComponentAnsys123::ElectricField(), Garfield::ComponentComsol::ElectricField(), Garfield::ComponentElmer::ElectricField(), Garfield::ComponentElmer2D::ElectricField(), FindElementCube(), Garfield::ComponentAnsys121::GetAspectRatio(), Garfield::ComponentAnsys123::GetAspectRatio(), Garfield::ComponentComsol::GetAspectRatio(), Garfield::ComponentElmer::GetAspectRatio(), Garfield::ComponentElmer2D::GetAspectRatio(), Garfield::ComponentAnsys121::GetElementVolume(), Garfield::ComponentAnsys123::GetElementVolume(), Garfield::ComponentComsol::GetElementVolume(), Garfield::ComponentElmer::GetElementVolume(), Garfield::ComponentElmer2D::GetElementVolume(), GetNode(), GetNumberOfNodes(), GetOrCreateWeightingFieldIndex(), Garfield::ComponentElmer::Initialise(), Garfield::ComponentElmer2D::Initialise(), Garfield::ComponentComsol::Initialise(), Garfield::ComponentAnsys121::Initialise(), Garfield::ComponentAnsys123::Initialise(), PrintElement(), Reset(), Garfield::ComponentComsol::SetDelayedWeightingPotential(), SetRange(), Garfield::ComponentComsol::SetWeightingField(), Garfield::ComponentAnsys121::SetWeightingField(), Garfield::ComponentAnsys123::SetWeightingField(), Garfield::ComponentElmer::SetWeightingField(), Garfield::ComponentElmer2D::SetWeightingField(), Garfield::ComponentAnsys121::WeightingField(), Garfield::ComponentAnsys123::WeightingField(), Garfield::ComponentComsol::WeightingField(), Garfield::ComponentElmer::WeightingField(), Garfield::ComponentElmer2D::WeightingField(), Garfield::ComponentAnsys121::WeightingPotential(), Garfield::ComponentAnsys123::WeightingPotential(), Garfield::ComponentComsol::WeightingPotential(), Garfield::ComponentElmer::WeightingPotential(), and Garfield::ComponentElmer2D::WeightingPotential().

◆ m_nWarnings

unsigned int Garfield::ComponentFieldMap::m_nWarnings = 0
protected

Definition at line 169 of file ComponentFieldMap.hh.

Referenced by PrintWarning(), and Reset().

◆ m_printConvergenceWarnings

bool Garfield::ComponentFieldMap::m_printConvergenceWarnings = true
protected

Definition at line 173 of file ComponentFieldMap.hh.

Referenced by EnableConvergenceWarnings().

◆ m_setang

std::array<bool, 3> Garfield::ComponentFieldMap::m_setang
protected

Definition at line 162 of file ComponentFieldMap.hh.

Referenced by SetRange().

◆ m_warning

◆ m_wdtimes

std::vector<double> Garfield::ComponentFieldMap::m_wdtimes
protected

◆ m_wfields

◆ m_wfieldsOk


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