Garfield++ 5.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
|
Component for importing and interpolating two-dimensional ANSYS field maps. More...
#include <ComponentAnsys121.hh>
Public Member Functions | |
ComponentAnsys121 () | |
Constructor. | |
~ComponentAnsys121 () | |
Destructor. | |
bool | Initialise (const std::string &elist="ELIST.lis", const std::string &nlist="NLIST.lis", const std::string &mplist="MPLIST.lis", const std::string &prnsol="PRNSOL.lis", const std::string &unit="cm") |
bool | SetWeightingPotential (const std::string &prnsol, const std::string &label) |
Import weighting potentials. | |
bool | SetWeightingField (const std::string &prnsol, const std::string &label) |
void | SetRangeZ (const double zmin, const double zmax) |
Set the limits of the active region along z. | |
![]() | |
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. | |
Medium * | GetMedium (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) |
Medium * | GetMedium (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). | |
![]() | |
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 ![]() | |
void | EnablePeriodicityY (const bool on=true) |
Enable simple periodicity in the ![]() | |
void | EnablePeriodicityZ (const bool on=true) |
Enable simple periodicity in the ![]() | |
void | IsPeriodic (bool &perx, bool &pery, bool &perz) |
Return periodicity flags. | |
void | EnableMirrorPeriodicityX (const bool on=true) |
Enable mirror periodicity in the ![]() | |
void | EnableMirrorPeriodicityY (const bool on=true) |
Enable mirror periodicity in the ![]() | |
void | EnableMirrorPeriodicityZ (const bool on=true) |
Enable mirror periodicity in the ![]() | |
void | IsMirrorPeriodic (bool &perx, bool &pery, bool &perz) |
Return mirror periodicity flags. | |
void | EnableAxialPeriodicityX (const bool on=true) |
Enable axial periodicity in the ![]() | |
void | EnableAxialPeriodicityY (const bool on=true) |
Enable axial periodicity in the ![]() | |
void | EnableAxialPeriodicityZ (const bool on=true) |
Enable axial periodicity in the ![]() | |
void | IsAxiallyPeriodic (bool &perx, bool &pery, bool &perz) |
Return axial periodicity flags. | |
void | EnableRotationSymmetryX (const bool on=true) |
Enable rotation symmetry around the ![]() | |
void | EnableRotationSymmetryY (const bool on=true) |
Enable rotation symmetry around the ![]() | |
void | EnableRotationSymmetryZ (const bool on=true) |
Enable rotation symmetry around the ![]() | |
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 | |
![]() | |
enum class | ElementType { Unknown = 0 , Serendipity = 5 , CurvedTetrahedron = 13 } |
![]() | |
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 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) |
![]() | |
bool | m_is3d = true |
ElementType | m_elementType = ElementType::CurvedTetrahedron |
std::vector< Element > | m_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< Node > | m_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< Material > | m_materials |
std::map< std::string, WeightingFieldCopy > | m_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 |
![]() | |
std::string | m_className = "Component" |
Class name. | |
Geometry * | m_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. | |
Component for importing and interpolating two-dimensional ANSYS field maps.
Definition at line 10 of file ComponentAnsys121.hh.
Garfield::ComponentAnsys121::ComponentAnsys121 | ( | ) |
Constructor.
Definition at line 10 of file ComponentAnsys121.cc.
|
inline |
bool Garfield::ComponentAnsys121::Initialise | ( | const std::string & | elist = "ELIST.lis", |
const std::string & | nlist = "NLIST.lis", | ||
const std::string & | mplist = "MPLIST.lis", | ||
const std::string & | prnsol = "PRNSOL.lis", | ||
const std::string & | unit = "cm" ) |
Import a field map.
elist | name of the file containing the list of elements |
nlist | name of the file containing the list of nodes |
mplist | name of the file containing the list of materials |
prnsol | name of the file containing the nodal solutions |
unit | length unit |
Definition at line 18 of file ComponentAnsys121.cc.
void Garfield::ComponentAnsys121::SetRangeZ | ( | const double | zmin, |
const double | zmax ) |
Set the limits of the active region along z.
Definition at line 541 of file ComponentAnsys121.cc.
bool Garfield::ComponentAnsys121::SetWeightingField | ( | const std::string & | prnsol, |
const std::string & | label ) |
Definition at line 519 of file ComponentAnsys121.cc.
Referenced by SetWeightingPotential().
|
inline |
Import weighting potentials.
Definition at line 30 of file ComponentAnsys121.hh.