Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
Loading...
Searching...
No Matches
Garfield::ComponentTcad2d Class Reference

Interpolation in a two-dimensional field map created by Sentaurus Device. More...

#include <ComponentTcad2d.hh>

Inheritance diagram for Garfield::ComponentTcad2d:
Garfield::ComponentTcadBase< 2 > Garfield::Component

Public Member Functions

 ComponentTcad2d ()
 Constructor.
 ~ComponentTcad2d ()
 Destructor.
void SetRangeZ (const double zmin, const double zmax)
 Set the z-extent of the bounding box (default: unlimited).
bool GetElement (const size_t i, double &vol, double &dmin, double &dmax, int &type) const
 Retrieve the properties of an element.
bool GetNode (const size_t i, double &x, double &y, double &z) const override
 Get the coordinates of a mesh node.
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 ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status) override
 Calculate the drift field at given point.
MediumGetMedium (const double x, const double y, const double z) override
 Get the medium at a given location (x, y, z).
void DelayedWeightingPotentials (const double x, const double y, const double z, const std::string &label, std::vector< double > &dwp) override
 Calculate the delayed weighting potentials at a given point and for a given electrode, for a set of pre-defined times.
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.
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::ComponentTcadBase< 2 >
 ComponentTcadBase ()=delete
 Default constructor.
virtual ~ComponentTcadBase ()
 Destructor.
bool Initialise (const std::string &gridfilename, const std::string &datafilename)
 Import mesh and field map from files.
bool SetWeightingField (const std::string &datfile1, const std::string &datfile2, const double dv, const std::string &label)
 Import field maps defining the prompt weighting field and potential.
bool SetWeightingPotential (const std::string &datfile1, const std::string &datfile2, const double dv, const std::string &label)
 Import field maps defining the prompt weighting field and potential.
bool SetWeightingFieldShift (const std::string &label, const double x, const double y, const double z)
 Shift the maps of weighting field/potential for a given electrode with respect to the original mesh.
bool SetDynamicWeightingPotential (const std::string &datfile1, const std::string &datfile2, const double dv, const double t, const std::string &label)
 Import time-dependent weighting potentials at t >= 0.
bool SetDynamicWeightingField (const std::string &datfile1, const std::string &datfile2, const double dv, const double t, const std::string &label)
 Import time-dependent weighting fields at t >= 0.
void PrintRegions () const
 List all currently defined regions.
size_t GetNumberOfRegions () const
 Get the number of regions in the device.
void GetRegion (const size_t ireg, std::string &name, bool &active) const
 Get the name and "active volume" flag of a region.
void SetDriftRegion (const size_t ireg)
 Make a region active ("driftable").
void UnsetDriftRegion (const size_t ireg)
 Make a region inactive.
void SetMedium (const size_t ireg, Medium *m)
 Set the medium to be associated to a given region.
size_t GetNumberOfElements () const override
 Return the number of mesh elements.
bool GetElementNodes (const size_t i, std::vector< size_t > &nodes) const override
 Get the indices of the nodes constituting a given element.
bool GetElementRegion (const size_t i, size_t &region, bool &active) const override
 Get the region/material of a mesh element and a flag whether it is associated to an active medium.
size_t GetNumberOfNodes () const override
 Return the number of mesh nodes.
void EnableVelocityMap (const bool on)
 Switch use of the imported velocity map on/off.
size_t GetNumberOfDonors ()
 Get the number of donor states found in the map.
size_t GetNumberOfAcceptors ()
 Get the number of acceptor states found in the map.
bool SetDonor (const size_t donorNumber, const double exsec, const double hxsec, const double concentration)
 Set the properties of a donor-type defect state.
bool SetAcceptor (const size_t acceptorNumber, const double exsec, const double hxsec, const double concentration)
 Set the properties of an acceptor-type defect state.
void EnableAlphaMap (const bool on=true)
 Use the imported impact ionisation map or not.
void EnableTrapOccupationMap (const bool on=true)
 Use the imported trapping map or not.
void EnableLifetimeMap (const bool on=true)
 Use the imported lifetime map or not.
bool GetElectronMobility (const double x, const double y, const double z, double &mob)
 Get the electron mobility at a given point in the mesh.
bool GetHoleMobility (const double x, const double y, const double z, double &mob)
 Get the hole mobility at a given point in the mesh.
double ElectronLifetime (const double x, const double y, const double z)
double HoleLifetime (const double x, const double y, const double z)
void WeightingField (const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label) override
 Calculate the weighting field at a given point and for a given electrode.
double WeightingPotential (const double x, const double y, const double z, const std::string &label) override
 Calculate the weighting potential at a given point.
const std::vector< double > & DelayedSignalTimes (const std::string &label) override
 Return the time steps at which the delayed weighting potential/field are stored/evaluated.
void DelayedWeightingField (const double x, const double y, const double z, const double t, double &wx, double &wy, double &wz, const std::string &label) override
 Calculate the delayed weighting field at a given point and time and for a given electrode.
double DelayedWeightingPotential (const double x, const double y, const double z, const double t, const std::string &label) override
 Calculate the delayed weighting potential at a given point and time and for a given electrode.
bool GetVoltageRange (double &vmin, double &vmax) override
 Calculate the voltage range [V].
bool HasVelocityMap () const override
 Does the component have velocity maps?
bool ElectronVelocity (const double x, const double y, const double z, double &vx, double &vy, double &vz) override
 Get the electron drift velocity.
bool HoleVelocity (const double x, const double y, const double z, double &vx, double &vy, double &vz) override
 Get the hole drift velocity.
bool HasTownsendMap () const override
 Does the component have maps of the Townsend coefficient?
bool HasAttachmentMap () const override
 Does the component have maps of the attachment coefficient?
bool HasMobilityMap () const override
 Does the component have maps of the low-field mobility?
bool ElectronAttachment (const double x, const double y, const double z, double &eta) override
 Get the electron attachment coefficient.
bool HoleAttachment (const double x, const double y, const double z, double &eta) override
 Get the hole attachment coefficient.
bool ElectronMobility (const double x, const double y, const double z, double &mu) override
bool HoleMobility (const double x, const double y, const double z, double &mu) override
 Get the hole Mobility coefficient.
bool ElectronTownsend (const double x, const double y, const double z, double &alpha) override
 Get the electron Townsend coefficient.
bool HoleTownsend (const double x, const double y, const double z, double &alpha) override
 Get the hole Townsend coefficient.
Public Member Functions inherited from Garfield::Component
 Component ()=delete
 Default constructor.
 Component (const std::string &name)
 Constructor.
virtual ~Component ()=default
 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 MagneticField (const double x, const double y, const double z, double &bx, double &by, double &bz, int &status)
 Calculate the magnetic field at a given point.
void SetMagneticField (const double bx, const double by, const double bz)
 Set a constant magnetic field.
virtual bool IsReady ()
 Ready for use?
virtual bool Is3d ()
 Does the component have a 3D field (map)?
double CellSizeX ()
double CellSizeY ()
double CellSizeZ ()
double IntegrateFluxCircle (const double xc, const double yc, const double r, const unsigned int nI=50)
 Integrate the normal component of the electric field over a circle.
double IntegrateFluxSphere (const double xc, const double yc, const double zc, const double r, const unsigned int nI=20)
 Integrate the normal component of the electric field over a sphere.
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)
 Integrate the normal component of the electric field over a parallelogram.
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)
 Integrate the normal component of the weighting field over a parallelogram.
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)
 Integrate the electric field flux through a line from (x0,y0,z0) to (x1,y1,z1) along a direction (xp,yp,zp).
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)
 Determine whether the line between two points crosses a wire.
virtual bool InTrapRadius (const double q0, const double x0, const double y0, const double z0, double &xw, double &yw, double &rw)
 Determine whether a particle is inside the trap radius of a wire.
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)
 Determine whether the line between two points crosses a plane.
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 EnableTriangleSymmetricXY (const bool on=true, const bool oct=2)
 Enable triangular periodicity in the $xy$ plane.
void EnableTriangleSymmetricXZ (const bool on=true, const bool oct=2)
 Enable triangular periodicity in the $xz$ plane.
void EnableTriangleSymmetricYZ (const bool on=true, const bool oct=2)
 Enable triangular periodicity in the $yz$ plane.
void EnableDebugging (const bool on=true)
 Switch on debugging messages.
void DisableDebugging ()
 Switch off debugging messages.
virtual bool HasMagneticField () const
 Does the component have a non-zero magnetic field?
virtual double StepSizeHint ()
virtual double CreateGPUTransferObject (ComponentGPU *&comp_gpu)
 Create and initialise GPU Transfer class.

Private Member Functions

void Reset () override
 Reset the component.
size_t FindElement (const double x, const double y, std::array< double, nMaxVertices > &w) const
bool InElement (const double x, const double y, const Element &element, std::array< double, nMaxVertices > &w) const
bool InRectangle (const double x, const double y, const Element &element, std::array< double, nMaxVertices > &w) const
bool InTriangle (const double x, const double y, const Element &element, std::array< double, nMaxVertices > &w) const
bool OnLine (const double x, const double y, const Element &element, std::array< double, nMaxVertices > &w) const
bool AtPoint (const double x, const double y, const Element &element, std::array< double, nMaxVertices > &w) const
bool Interpolate (const double x, const double y, const double z, const std::vector< double > &field, double &f) override
bool Interpolate (const double x, const double y, const double z, const std::vector< std::array< double, 2 > > &field, double &fx, double &fy, double &fz) override
void FillTree () override

Private Attributes

bool m_hasRangeZ = false
std::unique_ptr< QuadTreem_tree

Additional Inherited Members

Protected Member Functions inherited from Garfield::ComponentTcadBase< 2 >
void UpdatePeriodicity () override
 Verify periodicities.
void Cleanup ()
size_t FindRegion (const std::string &name) const
void MapCoordinates (std::array< double, N > &x, std::array< bool, N > &mirr) const
bool InBoundingBox (const std::array< double, N > &x) const
void UpdateAttachment ()
void ComputeEtaFromLifetime ()
void ComputeEtaFromTraps ()
bool GetOffset (const std::string &label, double &dx, double &dy, double &dz) const
bool LoadGrid (const std::string &gridfilename)
bool LoadData (const std::string &datafilename)
bool ReadDataset (std::ifstream &datafile, const std::string &dataset)
bool LoadWeightingField (const std::string &datafilename, std::vector< std::array< double, N > > &wf, std::vector< double > &wp)
Static Protected Member Functions inherited from Garfield::ComponentTcadBase< 2 >
static unsigned int ElementVertices (const Element &element)
Protected Attributes inherited from Garfield::ComponentTcadBase< 2 >
std::vector< Region > m_regions
std::vector< std::array< double, N > > m_vertices
std::vector< Element > m_elements
std::vector< double > m_epot
std::vector< std::array< double, N > > m_efield
std::map< std::string, std::vector< std::array< double, N > > > m_wfield
std::map< std::string, std::vector< double > > m_wpot
std::map< std::string, std::vector< double > > m_wshift
std::map< std::string, std::vector< std::vector< std::array< double, N > > > > m_dwf
std::map< std::string, std::vector< std::vector< double > > > m_dwp
std::map< std::string, std::vector< double > > m_dwtf
std::map< std::string, std::vector< double > > m_dwtp
std::vector< std::array< double, N > > m_eVelocity
std::vector< std::array< double, N > > m_hVelocity
std::vector< double > m_eMobility
std::vector< double > m_hMobility
std::vector< double > m_eAlpha
std::vector< double > m_hAlpha
std::vector< double > m_eLifetime
std::vector< double > m_hLifetime
std::vector< std::vector< float > > m_donorOcc
std::vector< std::vector< float > > m_acceptorOcc
std::vector< double > m_eEta
std::vector< double > m_hEta
std::vector< Defect > m_donors
std::vector< Defect > m_acceptors
bool m_useVelocityMap
bool m_useTrapOccMap
bool m_useLifetimeMap
bool m_useAlphaMap
std::array< double, 3 > m_bbMin
std::array< double, 3 > m_bbMax
double m_pMin
double m_pMax
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.
std::array< bool, 3 > m_triangleSymmetric = {{false, false, false}}
 Triangle symmetric in the xy, xz, and yz plane.
int m_triangleSymmetricOct = 0
 Triangle symmetric octant of imported map (0 < phi < Pi/4 --> octant 1).
const std::array< int, 4 > m_triangleOctRules = {1, 4, 5, 8}
 Octants where |x| >= |y|.
bool m_outsideCone = false
std::vector< double > m_wdtimes
 Time steps at which the delayed weighting potentials/fields are stored.
Static Protected Attributes inherited from Garfield::ComponentTcadBase< 2 >
static constexpr size_t nMaxVertices

Detailed Description

Interpolation in a two-dimensional field map created by Sentaurus Device.

Definition at line 14 of file ComponentTcad2d.hh.

Constructor & Destructor Documentation

◆ ComponentTcad2d()

Garfield::ComponentTcad2d::ComponentTcad2d ( )

Constructor.

◆ ~ComponentTcad2d()

Garfield::ComponentTcad2d::~ComponentTcad2d ( )
inline

Destructor.

Definition at line 19 of file ComponentTcad2d.hh.

19{}

Member Function Documentation

◆ AtPoint()

bool Garfield::ComponentTcad2d::AtPoint ( const double x,
const double y,
const Element & element,
std::array< double, nMaxVertices > & w ) const
private

◆ DelayedWeightingPotentials()

void Garfield::ComponentTcad2d::DelayedWeightingPotentials ( const double x,
const double y,
const double z,
const std::string & label,
std::vector< double > & dwp )
overridevirtual

Calculate the delayed weighting potentials at a given point and for a given electrode, for a set of pre-defined times.

Reimplemented from Garfield::Component.

◆ ElectricField() [1/3]

std::array< double, 3 > Garfield::Component::ElectricField ( const double x,
const double y,
const double z )

Calculate the drift field [V/cm] at (x, y, z).

◆ ElectricField() [2/3]

void Garfield::ComponentTcad2d::ElectricField ( const double x,
const double y,
const double z,
double & ex,
double & ey,
double & ez,
double & v,
Medium *& m,
int & status )
overridevirtual

Calculate the drift field [V/cm] and potential [V] at (x, y, z).

Implements Garfield::Component.

◆ ElectricField() [3/3]

void Garfield::ComponentTcad2d::ElectricField ( const double x,
const double y,
const double z,
double & ex,
double & ey,
double & ez,
Medium *& m,
int & status )
inlineoverridevirtual

Calculate the drift field at given point.

Parameters
x,y,zcoordinates [cm].
ex,ey,ezcomponents of the electric field [V/cm].
mpointer to the medium at this location.
statusstatus flag

Status flags:

        0: Inside an active medium
      > 0: Inside a wire of type X
-4 ... -1: On the side of a plane where no wires are
       -5: Inside the mesh but not in an active medium
       -6: Outside the mesh
      -10: Unknown potential type (should not occur)
    other: Other cases (should not occur)

Implements Garfield::Component.

Definition at line 40 of file ComponentTcad2d.hh.

41 {
42 double v = 0.;
43 ElectricField(x, y, z, ex, ey, ez, v, m, status);
44 }
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).

◆ FillTree()

void Garfield::ComponentTcad2d::FillTree ( )
overrideprivatevirtual

◆ FindElement()

size_t Garfield::ComponentTcad2d::FindElement ( const double x,
const double y,
std::array< double, nMaxVertices > & w ) const
private

◆ GetBoundingBox()

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

Get the bounding box coordinates.

Reimplemented from Garfield::Component.

◆ GetElement()

bool Garfield::ComponentTcad2d::GetElement ( const size_t i,
double & vol,
double & dmin,
double & dmax,
int & type ) const

Retrieve the properties of an element.

Parameters
iindex of the element
volvolume
dminsmallest length in the element
dmaxlargest length in the element
typeelement type

◆ GetElementaryCell()

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

◆ GetMedium()

Medium * Garfield::ComponentTcad2d::GetMedium ( const double x,
const double y,
const double z )
overridevirtual

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

Reimplemented from Garfield::Component.

◆ GetNode()

bool Garfield::ComponentTcad2d::GetNode ( const size_t i,
double & x,
double & y,
double & z ) const
overridevirtual

Get the coordinates of a mesh node.

Reimplemented from Garfield::Component.

◆ InElement()

bool Garfield::ComponentTcad2d::InElement ( const double x,
const double y,
const Element & element,
std::array< double, nMaxVertices > & w ) const
private

◆ InRectangle()

bool Garfield::ComponentTcad2d::InRectangle ( const double x,
const double y,
const Element & element,
std::array< double, nMaxVertices > & w ) const
private

◆ Interpolate() [1/2]

bool Garfield::ComponentTcad2d::Interpolate ( const double x,
const double y,
const double z,
const std::vector< double > & field,
double & f )
overrideprivatevirtual

◆ Interpolate() [2/2]

bool Garfield::ComponentTcad2d::Interpolate ( const double x,
const double y,
const double z,
const std::vector< std::array< double, 2 > > & field,
double & fx,
double & fy,
double & fz )
overrideprivate

◆ InTriangle()

bool Garfield::ComponentTcad2d::InTriangle ( const double x,
const double y,
const Element & element,
std::array< double, nMaxVertices > & w ) const
private

◆ OnLine()

bool Garfield::ComponentTcad2d::OnLine ( const double x,
const double y,
const Element & element,
std::array< double, nMaxVertices > & w ) const
private

◆ Reset()

void Garfield::ComponentTcad2d::Reset ( )
inlineoverrideprivatevirtual

Reset the component.

Implements Garfield::Component.

Definition at line 63 of file ComponentTcad2d.hh.

63 {
64 Cleanup();
65 m_hasRangeZ = false;
66 m_ready = false;
67 }
bool m_ready
Ready for use?
Definition Component.hh:441

◆ SetRangeZ()

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

Set the z-extent of the bounding box (default: unlimited).

Member Data Documentation

◆ m_hasRangeZ

bool Garfield::ComponentTcad2d::m_hasRangeZ = false
private

Definition at line 58 of file ComponentTcad2d.hh.

◆ m_tree

std::unique_ptr<QuadTree> Garfield::ComponentTcad2d::m_tree
private

Definition at line 61 of file ComponentTcad2d.hh.


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