Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
Loading...
Searching...
No Matches
Garfield::ComponentTcadBase< N > Class Template Referenceabstract

Interpolation in a field map created by Sentaurus Device. More...

#include <ComponentTcadBase.hh>

Inheritance diagram for Garfield::ComponentTcadBase< N >:
Garfield::Component

Classes

struct  Defect
struct  Element
struct  Region

Public Member Functions

 ComponentTcadBase ()=delete
 Default constructor.
 ComponentTcadBase (const std::string &name)
 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.
void SetMedium (const std::string &material, Medium *m)
 Set the medium to be associated to all regions with a given material.
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.
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
 Calculate the drift field at given point.
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).
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 DelayedWeightingPotentials (const double x, const double y, const double z, const std::string &label, std::vector< double > &dwp)
 Calculate the delayed weighting potentials at a given point and for a given electrode, for a set of pre-defined times.
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)?
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 CellSizeX ()
double CellSizeY ()
double CellSizeZ ()
virtual bool GetNode (const size_t i, double &x, double &y, double &z) const
 Get the coordinates of a mesh node.
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.

Protected Member Functions

void UpdatePeriodicity () override
 Verify periodicities.
void Cleanup ()
virtual bool Interpolate (const double x, const double y, const double z, const std::vector< double > &field, double &f)=0
virtual bool Interpolate (const double x, const double y, const double z, const std::vector< std::array< double, N > > &field, double &fx, double &fy, double &fz)=0
virtual void FillTree ()=0
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)
Protected Member Functions inherited from Garfield::Component
virtual void Reset ()=0
 Reset the component.

Static Protected Member Functions

static unsigned int ElementVertices (const Element &element)

Protected Attributes

std::vector< Regionm_regions
std::vector< std::array< double, N > > m_vertices
std::vector< Elementm_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< Defectm_donors
std::vector< Defectm_acceptors
bool m_useVelocityMap = false
bool m_useTrapOccMap = false
bool m_useLifetimeMap = false
bool m_useAlphaMap = false
std::array< double, 3 > m_bbMin = {{0., 0., 0.}}
std::array< double, 3 > m_bbMax = {{0., 0., 0.}}
double m_pMin = 0.
double m_pMax = 0.
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

static constexpr size_t nMaxVertices = 4

Detailed Description

template<size_t N>
class Garfield::ComponentTcadBase< N >

Interpolation in a field map created by Sentaurus Device.

Definition at line 15 of file ComponentTcadBase.hh.

Constructor & Destructor Documentation

◆ ComponentTcadBase() [1/2]

template<size_t N>
Garfield::ComponentTcadBase< N >::ComponentTcadBase ( )
delete

Default constructor.

◆ ComponentTcadBase() [2/2]

template<size_t N>
Garfield::ComponentTcadBase< N >::ComponentTcadBase ( const std::string & name)
inline

Constructor.

Definition at line 20 of file ComponentTcadBase.hh.

20 : Component(name) {
21 m_regions.reserve(10);
22 m_vertices.reserve(10000);
23 m_elements.reserve(10000);
24 }
Interpolation in a field map created by Sentaurus Device.
std::vector< std::array< double, N > > m_vertices
std::vector< Element > m_elements
std::vector< Region > m_regions
Component()=delete
Default constructor.

◆ ~ComponentTcadBase()

template<size_t N>
virtual Garfield::ComponentTcadBase< N >::~ComponentTcadBase ( )
inlinevirtual

Destructor.

Definition at line 26 of file ComponentTcadBase.hh.

26{}

Member Function Documentation

◆ Cleanup()

template<size_t N>
void Garfield::ComponentTcadBase< N >::Cleanup ( )
protected

◆ ComputeEtaFromLifetime()

template<size_t N>
void Garfield::ComponentTcadBase< N >::ComputeEtaFromLifetime ( )
protected

◆ ComputeEtaFromTraps()

template<size_t N>
void Garfield::ComponentTcadBase< N >::ComputeEtaFromTraps ( )
protected

◆ DelayedSignalTimes()

template<size_t N>
const std::vector< double > & Garfield::ComponentTcadBase< N >::DelayedSignalTimes ( const std::string & )
inlineoverridevirtual

Return the time steps at which the delayed weighting potential/field are stored/evaluated.

Reimplemented from Garfield::Component.

Definition at line 135 of file ComponentTcadBase.hh.

136 {
137 if (m_dwtp.count(label) > 0) return m_dwtp[label];
138 static const std::vector<double> emptyVector;
139 return emptyVector;
140 }
std::map< std::string, std::vector< double > > m_dwtp

◆ DelayedWeightingField()

template<size_t N>
void Garfield::ComponentTcadBase< N >::DelayedWeightingField ( const double x,
const double y,
const double z,
const double t,
double & wx,
double & wy,
double & wz,
const std::string & label )
overridevirtual

Calculate the delayed weighting field at a given point and time and for a given electrode.

Parameters
x,y,zcoordinates [cm].
ttime [ns].
wx,wy,wzcomponents of the weighting field [1/cm].
labelname of the electrode

Reimplemented from Garfield::Component.

◆ DelayedWeightingPotential()

template<size_t N>
double Garfield::ComponentTcadBase< N >::DelayedWeightingPotential ( const double x,
const double y,
const double z,
const double t,
const std::string & label )
overridevirtual

Calculate the delayed weighting potential at a given point and time and for a given electrode.

Parameters
x,y,zcoordinates [cm].
ttime [ns].
labelname of the electrode

Reimplemented from Garfield::Component.

◆ ElectronAttachment()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::ElectronAttachment ( const double ,
const double ,
const double ,
double & eta )
overridevirtual

Get the electron attachment coefficient.

Reimplemented from Garfield::Component.

◆ ElectronLifetime()

template<size_t N>
double Garfield::ComponentTcadBase< N >::ElectronLifetime ( const double x,
const double y,
const double z )

◆ ElectronMobility()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::ElectronMobility ( const double x,
const double y,
const double z,
double & mu )
overridevirtual

Reimplemented from Garfield::Component.

◆ ElectronTownsend()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::ElectronTownsend ( const double ,
const double ,
const double ,
double & alpha )
overridevirtual

Get the electron Townsend coefficient.

Reimplemented from Garfield::Component.

◆ ElectronVelocity()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::ElectronVelocity ( const double ,
const double ,
const double ,
double & vx,
double & vy,
double & vz )
overridevirtual

Get the electron drift velocity.

Reimplemented from Garfield::Component.

◆ ElementVertices()

template<size_t N>
unsigned int Garfield::ComponentTcadBase< N >::ElementVertices ( const Element & element)
inlinestaticprotected

Definition at line 299 of file ComponentTcadBase.hh.

299 {
300 return std::min(element.type + 1, 4U);
301 }

◆ EnableAlphaMap()

template<size_t N>
void Garfield::ComponentTcadBase< N >::EnableAlphaMap ( const bool on = true)
inline

Use the imported impact ionisation map or not.

Definition at line 113 of file ComponentTcadBase.hh.

◆ EnableLifetimeMap()

template<size_t N>
void Garfield::ComponentTcadBase< N >::EnableLifetimeMap ( const bool on = true)

Use the imported lifetime map or not.

◆ EnableTrapOccupationMap()

template<size_t N>
void Garfield::ComponentTcadBase< N >::EnableTrapOccupationMap ( const bool on = true)

Use the imported trapping map or not.

◆ EnableVelocityMap()

template<size_t N>
void Garfield::ComponentTcadBase< N >::EnableVelocityMap ( const bool on)

Switch use of the imported velocity map on/off.

◆ FillTree()

template<size_t N>
virtual void Garfield::ComponentTcadBase< N >::FillTree ( )
protectedpure virtual

◆ FindRegion()

template<size_t N>
size_t Garfield::ComponentTcadBase< N >::FindRegion ( const std::string & name) const
protected

◆ GetElectronMobility()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::GetElectronMobility ( const double x,
const double y,
const double z,
double & mob )

Get the electron mobility at a given point in the mesh.

◆ GetElementNodes()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::GetElementNodes ( const size_t ,
std::vector< size_t > &  ) const
overridevirtual

Get the indices of the nodes constituting a given element.

Reimplemented from Garfield::Component.

◆ GetElementRegion()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::GetElementRegion ( const size_t ,
size_t & ,
bool &  ) const
overridevirtual

Get the region/material of a mesh element and a flag whether it is associated to an active medium.

Reimplemented from Garfield::Component.

◆ GetHoleMobility()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::GetHoleMobility ( const double x,
const double y,
const double z,
double & mob )

Get the hole mobility at a given point in the mesh.

◆ GetNumberOfAcceptors()

template<size_t N>
size_t Garfield::ComponentTcadBase< N >::GetNumberOfAcceptors ( )
inline

Get the number of acceptor states found in the map.

Definition at line 98 of file ComponentTcadBase.hh.

98{ return m_acceptors.size(); }
std::vector< Defect > m_acceptors

◆ GetNumberOfDonors()

template<size_t N>
size_t Garfield::ComponentTcadBase< N >::GetNumberOfDonors ( )
inline

Get the number of donor states found in the map.

Definition at line 96 of file ComponentTcadBase.hh.

96{ return m_donors.size(); }
std::vector< Defect > m_donors

◆ GetNumberOfElements()

template<size_t N>
size_t Garfield::ComponentTcadBase< N >::GetNumberOfElements ( ) const
inlineoverridevirtual

Return the number of mesh elements.

Reimplemented from Garfield::Component.

Definition at line 85 of file ComponentTcadBase.hh.

85{ return m_elements.size(); }

◆ GetNumberOfNodes()

template<size_t N>
size_t Garfield::ComponentTcadBase< N >::GetNumberOfNodes ( ) const
inlineoverridevirtual

Return the number of mesh nodes.

Reimplemented from Garfield::Component.

Definition at line 90 of file ComponentTcadBase.hh.

90{ return m_vertices.size(); }

◆ GetNumberOfRegions()

template<size_t N>
size_t Garfield::ComponentTcadBase< N >::GetNumberOfRegions ( ) const
inline

Get the number of regions in the device.

Definition at line 73 of file ComponentTcadBase.hh.

73{ return m_regions.size(); }

◆ GetOffset()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::GetOffset ( const std::string & label,
double & dx,
double & dy,
double & dz ) const
protected

◆ GetRegion()

template<size_t N>
void Garfield::ComponentTcadBase< N >::GetRegion ( const size_t ireg,
std::string & name,
bool & active ) const

Get the name and "active volume" flag of a region.

◆ GetVoltageRange()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::GetVoltageRange ( double & vmin,
double & vmax )
overridevirtual

Calculate the voltage range [V].

Implements Garfield::Component.

◆ HasAttachmentMap()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::HasAttachmentMap ( ) const
inlineoverridevirtual

Does the component have maps of the attachment coefficient?

Reimplemented from Garfield::Component.

Definition at line 161 of file ComponentTcadBase.hh.

161 {
162 return ((m_useTrapOccMap || m_useLifetimeMap) &&
163 !(m_eEta.empty() && m_hEta.empty()));
164 }

◆ HasMobilityMap()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::HasMobilityMap ( ) const
inlineoverridevirtual

Does the component have maps of the low-field mobility?

Reimplemented from Garfield::Component.

Definition at line 165 of file ComponentTcadBase.hh.

165 {
166 return !(m_eMobility.empty() && m_hMobility.empty());
167 }
std::vector< double > m_hMobility
std::vector< double > m_eMobility

◆ HasTownsendMap()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::HasTownsendMap ( ) const
inlineoverridevirtual

Does the component have maps of the Townsend coefficient?

Reimplemented from Garfield::Component.

Definition at line 158 of file ComponentTcadBase.hh.

158 {
159 return m_useAlphaMap && !(m_eAlpha.empty() && m_hAlpha.empty());
160 }
std::vector< double > m_hAlpha
std::vector< double > m_eAlpha

◆ HasVelocityMap()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::HasVelocityMap ( ) const
inlineoverridevirtual

Does the component have velocity maps?

Reimplemented from Garfield::Component.

Definition at line 150 of file ComponentTcadBase.hh.

150 {
151 return m_useVelocityMap && !(m_eVelocity.empty() && m_hVelocity.empty());
152 }
std::vector< std::array< double, N > > m_hVelocity
std::vector< std::array< double, N > > m_eVelocity

◆ HoleAttachment()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::HoleAttachment ( const double ,
const double ,
const double ,
double & eta )
overridevirtual

Get the hole attachment coefficient.

Reimplemented from Garfield::Component.

◆ HoleLifetime()

template<size_t N>
double Garfield::ComponentTcadBase< N >::HoleLifetime ( const double x,
const double y,
const double z )

◆ HoleMobility()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::HoleMobility ( const double ,
const double ,
const double ,
double & mu )
overridevirtual

Get the hole Mobility coefficient.

Reimplemented from Garfield::Component.

◆ HoleTownsend()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::HoleTownsend ( const double ,
const double ,
const double ,
double & alpha )
overridevirtual

Get the hole Townsend coefficient.

Reimplemented from Garfield::Component.

◆ HoleVelocity()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::HoleVelocity ( const double ,
const double ,
const double ,
double & vx,
double & vy,
double & vz )
overridevirtual

Get the hole drift velocity.

Reimplemented from Garfield::Component.

◆ InBoundingBox()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::InBoundingBox ( const std::array< double, N > & x) const
inlineprotected

Definition at line 312 of file ComponentTcadBase.hh.

312 {
313 for (size_t i = 0; i < N; ++i) {
314 if (x[i] < m_bbMin[i] || x[i] > m_bbMax[i]) return false;
315 }
316 return true;
317 }
std::array< double, 3 > m_bbMax
std::array< double, 3 > m_bbMin

◆ Initialise()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::Initialise ( const std::string & gridfilename,
const std::string & datafilename )

Import mesh and field map from files.

Parameters
gridfilenamename of the .grd file containing the mesh
datafilenamename of the .dat file containing the nodal solution

◆ Interpolate() [1/2]

template<size_t N>
virtual bool Garfield::ComponentTcadBase< N >::Interpolate ( const double x,
const double y,
const double z,
const std::vector< double > & field,
double & f )
protectedpure virtual

◆ Interpolate() [2/2]

template<size_t N>
virtual bool Garfield::ComponentTcadBase< N >::Interpolate ( const double x,
const double y,
const double z,
const std::vector< std::array< double, N > > & field,
double & fx,
double & fy,
double & fz )
protectedpure virtual

◆ LoadData()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::LoadData ( const std::string & datafilename)
protected

◆ LoadGrid()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::LoadGrid ( const std::string & gridfilename)
protected

◆ LoadWeightingField()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::LoadWeightingField ( const std::string & datafilename,
std::vector< std::array< double, N > > & wf,
std::vector< double > & wp )
protected

◆ MapCoordinates()

template<size_t N>
void Garfield::ComponentTcadBase< N >::MapCoordinates ( std::array< double, N > & x,
std::array< bool, N > & mirr ) const
protected

◆ PrintRegions()

template<size_t N>
void Garfield::ComponentTcadBase< N >::PrintRegions ( ) const

List all currently defined regions.

◆ ReadDataset()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::ReadDataset ( std::ifstream & datafile,
const std::string & dataset )
protected

◆ SetAcceptor()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::SetAcceptor ( const size_t acceptorNumber,
const double exsec,
const double hxsec,
const double concentration )

Set the properties of an acceptor-type defect state.

◆ SetDonor()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::SetDonor ( const size_t donorNumber,
const double exsec,
const double hxsec,
const double concentration )

Set the properties of a donor-type defect state.

Parameters
donorNumberindex of the donor
exseccross-section [cm2] for electrons
hxseccross-section [cm2] for holes
concentrationdefect density [cm-3]

◆ SetDriftRegion()

template<size_t N>
void Garfield::ComponentTcadBase< N >::SetDriftRegion ( const size_t ireg)

Make a region active ("driftable").

◆ SetDynamicWeightingField()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::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.

◆ SetDynamicWeightingPotential()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::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.

◆ SetMedium() [1/2]

template<size_t N>
void Garfield::ComponentTcadBase< N >::SetMedium ( const size_t ireg,
Medium * m )

Set the medium to be associated to a given region.

◆ SetMedium() [2/2]

template<size_t N>
void Garfield::ComponentTcadBase< N >::SetMedium ( const std::string & material,
Medium * m )

Set the medium to be associated to all regions with a given material.

◆ SetWeightingField()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::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.

Parameters
datfile1.dat file containing the field map at nominal bias.
datfile2.dat file containing the field map for a configuration with the potential at the electrode to be read out increased by a small voltage dv.
dvincrease in electrode potential between the two field maps.
labelname of the electrode

The field maps must use the same mesh as the drift field.

◆ SetWeightingFieldShift()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::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.

If the electrode does not exist yet, a new one will be added to the list.

◆ SetWeightingPotential()

template<size_t N>
bool Garfield::ComponentTcadBase< N >::SetWeightingPotential ( const std::string & datfile1,
const std::string & datfile2,
const double dv,
const std::string & label )
inline

Import field maps defining the prompt weighting field and potential.

Definition at line 50 of file ComponentTcadBase.hh.

52 {
54 }
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.

◆ UnsetDriftRegion()

template<size_t N>
void Garfield::ComponentTcadBase< N >::UnsetDriftRegion ( const size_t ireg)

Make a region inactive.

◆ UpdateAttachment()

template<size_t N>
void Garfield::ComponentTcadBase< N >::UpdateAttachment ( )
protected

◆ UpdatePeriodicity()

template<size_t N>
void Garfield::ComponentTcadBase< N >::UpdatePeriodicity ( )
overrideprotectedvirtual

Verify periodicities.

Implements Garfield::Component.

◆ WeightingField()

template<size_t N>
void Garfield::ComponentTcadBase< N >::WeightingField ( const double x,
const double y,
const double z,
double & wx,
double & wy,
double & wz,
const std::string & label )
overridevirtual

Calculate the weighting field at a given point and for a given electrode.

Parameters
x,y,zcoordinates [cm].
wx,wy,wzcomponents of the weighting field [1/cm].
labelname of the electrode

Reimplemented from Garfield::Component.

◆ WeightingPotential()

template<size_t N>
double Garfield::ComponentTcadBase< N >::WeightingPotential ( const double x,
const double y,
const double z,
const std::string & label )
overridevirtual

Calculate the weighting potential at a given point.

Parameters
x,y,zcoordinates [cm].
labelname of the electrode.
Returns
weighting potential [dimensionless].

Reimplemented from Garfield::Component.

Member Data Documentation

◆ m_acceptorOcc

template<size_t N>
std::vector<std::vector<float> > Garfield::ComponentTcadBase< N >::m_acceptorOcc
protected

Definition at line 262 of file ComponentTcadBase.hh.

◆ m_acceptors

template<size_t N>
std::vector<Defect> Garfield::ComponentTcadBase< N >::m_acceptors
protected

Definition at line 276 of file ComponentTcadBase.hh.

◆ m_bbMax

template<size_t N>
std::array<double, 3> Garfield::ComponentTcadBase< N >::m_bbMax = {{0., 0., 0.}}
protected

Definition at line 289 of file ComponentTcadBase.hh.

289{{0., 0., 0.}};

◆ m_bbMin

template<size_t N>
std::array<double, 3> Garfield::ComponentTcadBase< N >::m_bbMin = {{0., 0., 0.}}
protected

Definition at line 288 of file ComponentTcadBase.hh.

288{{0., 0., 0.}};

◆ m_donorOcc

template<size_t N>
std::vector<std::vector<float> > Garfield::ComponentTcadBase< N >::m_donorOcc
protected

Definition at line 261 of file ComponentTcadBase.hh.

◆ m_donors

template<size_t N>
std::vector<Defect> Garfield::ComponentTcadBase< N >::m_donors
protected

Definition at line 275 of file ComponentTcadBase.hh.

◆ m_dwf

template<size_t N>
std::map<std::string, std::vector<std::vector<std::array<double, N> > > > Garfield::ComponentTcadBase< N >::m_dwf
protected

Definition at line 242 of file ComponentTcadBase.hh.

◆ m_dwp

template<size_t N>
std::map<std::string, std::vector<std::vector<double> > > Garfield::ComponentTcadBase< N >::m_dwp
protected

Definition at line 243 of file ComponentTcadBase.hh.

◆ m_dwtf

template<size_t N>
std::map<std::string, std::vector<double> > Garfield::ComponentTcadBase< N >::m_dwtf
protected

Definition at line 245 of file ComponentTcadBase.hh.

◆ m_dwtp

template<size_t N>
std::map<std::string, std::vector<double> > Garfield::ComponentTcadBase< N >::m_dwtp
protected

Definition at line 246 of file ComponentTcadBase.hh.

◆ m_eAlpha

template<size_t N>
std::vector<double> Garfield::ComponentTcadBase< N >::m_eAlpha
protected

Definition at line 255 of file ComponentTcadBase.hh.

◆ m_eEta

template<size_t N>
std::vector<double> Garfield::ComponentTcadBase< N >::m_eEta
protected

Definition at line 264 of file ComponentTcadBase.hh.

◆ m_efield

template<size_t N>
std::vector<std::array<double, N> > Garfield::ComponentTcadBase< N >::m_efield
protected

Definition at line 232 of file ComponentTcadBase.hh.

◆ m_elements

template<size_t N>
std::vector<Element> Garfield::ComponentTcadBase< N >::m_elements
protected

Definition at line 227 of file ComponentTcadBase.hh.

◆ m_eLifetime

template<size_t N>
std::vector<double> Garfield::ComponentTcadBase< N >::m_eLifetime
protected

Definition at line 258 of file ComponentTcadBase.hh.

◆ m_eMobility

template<size_t N>
std::vector<double> Garfield::ComponentTcadBase< N >::m_eMobility
protected

Definition at line 252 of file ComponentTcadBase.hh.

◆ m_epot

template<size_t N>
std::vector<double> Garfield::ComponentTcadBase< N >::m_epot
protected

Definition at line 230 of file ComponentTcadBase.hh.

◆ m_eVelocity

template<size_t N>
std::vector<std::array<double, N> > Garfield::ComponentTcadBase< N >::m_eVelocity
protected

Definition at line 249 of file ComponentTcadBase.hh.

◆ m_hAlpha

template<size_t N>
std::vector<double> Garfield::ComponentTcadBase< N >::m_hAlpha
protected

Definition at line 256 of file ComponentTcadBase.hh.

◆ m_hEta

template<size_t N>
std::vector<double> Garfield::ComponentTcadBase< N >::m_hEta
protected

Definition at line 265 of file ComponentTcadBase.hh.

◆ m_hLifetime

template<size_t N>
std::vector<double> Garfield::ComponentTcadBase< N >::m_hLifetime
protected

Definition at line 259 of file ComponentTcadBase.hh.

◆ m_hMobility

template<size_t N>
std::vector<double> Garfield::ComponentTcadBase< N >::m_hMobility
protected

Definition at line 253 of file ComponentTcadBase.hh.

◆ m_hVelocity

template<size_t N>
std::vector<std::array<double, N> > Garfield::ComponentTcadBase< N >::m_hVelocity
protected

Definition at line 250 of file ComponentTcadBase.hh.

◆ m_pMax

template<size_t N>
double Garfield::ComponentTcadBase< N >::m_pMax = 0.
protected

Definition at line 293 of file ComponentTcadBase.hh.

◆ m_pMin

template<size_t N>
double Garfield::ComponentTcadBase< N >::m_pMin = 0.
protected

Definition at line 292 of file ComponentTcadBase.hh.

◆ m_regions

template<size_t N>
std::vector<Region> Garfield::ComponentTcadBase< N >::m_regions
protected

Definition at line 197 of file ComponentTcadBase.hh.

◆ m_useAlphaMap

template<size_t N>
bool Garfield::ComponentTcadBase< N >::m_useAlphaMap = false
protected

Definition at line 285 of file ComponentTcadBase.hh.

◆ m_useLifetimeMap

template<size_t N>
bool Garfield::ComponentTcadBase< N >::m_useLifetimeMap = false
protected

Definition at line 283 of file ComponentTcadBase.hh.

◆ m_useTrapOccMap

template<size_t N>
bool Garfield::ComponentTcadBase< N >::m_useTrapOccMap = false
protected

Definition at line 281 of file ComponentTcadBase.hh.

◆ m_useVelocityMap

template<size_t N>
bool Garfield::ComponentTcadBase< N >::m_useVelocityMap = false
protected

Definition at line 279 of file ComponentTcadBase.hh.

◆ m_vertices

template<size_t N>
std::vector<std::array<double, N> > Garfield::ComponentTcadBase< N >::m_vertices
protected

Definition at line 200 of file ComponentTcadBase.hh.

◆ m_wfield

template<size_t N>
std::map<std::string, std::vector<std::array<double, N> > > Garfield::ComponentTcadBase< N >::m_wfield
protected

Definition at line 235 of file ComponentTcadBase.hh.

◆ m_wpot

template<size_t N>
std::map<std::string, std::vector<double> > Garfield::ComponentTcadBase< N >::m_wpot
protected

Definition at line 236 of file ComponentTcadBase.hh.

◆ m_wshift

template<size_t N>
std::map<std::string, std::vector<double> > Garfield::ComponentTcadBase< N >::m_wshift
protected

Definition at line 238 of file ComponentTcadBase.hh.

◆ nMaxVertices

template<size_t N>
size_t Garfield::ComponentTcadBase< N >::nMaxVertices = 4
staticconstexprprotected

Definition at line 184 of file ComponentTcadBase.hh.


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