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

Component for interpolating field maps on a regular mesh. More...

#include <ComponentGrid.hh>

Inheritance diagram for Garfield::ComponentGrid:
Garfield::Component

Classes

struct  Node

Public Member Functions

 ComponentGrid ()
 Constructor.
 ~ComponentGrid ()
 Destructor.
bool SetMesh (const unsigned int nx, const unsigned int ny, const unsigned int nz, const double xmin, const double xmax, const double ymin, const double ymax, const double zmin, const double zmax)
 Define the grid.
bool GetMesh (unsigned int &nx, unsigned int &ny, unsigned int &nz, double &xmin, double &xmax, double &ymin, double &ymax, double &zmin, double &zmax) const
 Retrieve the parameters of the grid.
void SetCartesianCoordinates ()
 Use Cartesian coordinates (default).
void SetCylindricalCoordinates ()
 Use cylindrical coordinates.
bool LoadElectricField (const std::string &filename, const std::string &format, const bool withPotential, const bool withFlag, const double scaleX=1., const double scaleE=1., const double scaleP=1.)
 Import electric field and potential values from a file.
bool LoadWeightingField (const std::string &filename, const std::string &format, const bool withPotential, const double scaleX=1., const double scaleE=1., const double scaleP=1.)
 Import (prompt) weighting field from file.
bool LoadWeightingField (const std::string &filename, const std::string &format, const double time, const bool withPotential, const double scaleX=1., const double scaleE=1., const double scaleP=1.)
 Import delayed weighting field from file.
void SetWeightingFieldOffset (const double x, const double y, const double z)
 Offset coordinates in the weighting field, such that the same numerical weighting field map can be used for electrodes at different positions.
bool LoadMagneticField (const std::string &filename, const std::string &format, const double scaleX=1., const double scaleB=1.)
 Import magnetic field values from a file.
bool SaveElectricField (Component *cmp, const std::string &filename, const std::string &fmt)
 Export the electric field and potential of a component to a text file.
bool SaveElectricFieldROOT (Component *cmp, const std::string &filename, const std::string &fmt)
bool SaveElectricField (Component *cmp)
 Export the electric field and potential of a component.
bool SaveWeightingField (Component *cmp, const std::string &id, const std::string &filename, const std::string &fmt)
 Export the weighting field and potential of a component to a text file.
bool GetElectricField (const unsigned int i, const unsigned int j, const unsigned int k, double &v, double &ex, double &ey, double &ez) const
 Return the field at a given node.
void SetMedium (Medium *m)
 Set the medium.
MediumGetMedium () const
 Get the medium.
void Print ()
 Print information about the mesh and the available data.
bool LoadElectronAttachment (const std::string &fname, const std::string &fmt, const unsigned int col, const double scaleX=1.)
 Import electron attachment coefficients from a file.
bool LoadHoleAttachment (const std::string &fname, const std::string &fmt, const unsigned int col, const double scaleX=1.)
 Import hole attachment coefficients from a file.
bool LoadElectronVelocity (const std::string &fname, const std::string &fmt, const double scaleX=1., const double scaleV=1.e-9)
 Import a map of electron drift velocities from a file.
bool LoadHoleVelocity (const std::string &fname, const std::string &fmt, const double scaleX=1., const double scaleV=1.e-9)
 Import a map of hole drift velocities from a file.
void Clear () override
 Reset.
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.
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.
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.
void MagneticField (const double x, const double y, const double z, double &bx, double &by, double &bz, int &status) override
 Calculate the magnetic field at a given point.
MediumGetMedium (const double x, const double y, const double z) override
 Get the medium at a given location (x, y, z).
bool GetVoltageRange (double &vmin, double &vmax) override
 Calculate the voltage range [V].
bool GetElectricFieldRange (double &exmin, double &exmax, double &eymin, double &eymax, double &ezmin, double &ezmax)
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 HasMagneticField () const override
 Does the component have a non-zero magnetic field?
bool HasAttachmentMap () const override
 Does the component have maps of the attachment coefficient?
bool ElectronAttachment (const double x, const double y, const double z, double &att) override
 Get the electron attachment coefficient.
bool HoleAttachment (const double x, const double y, const double z, double &att) override
 Get the hole attachment coefficient.
bool HasMobilityMap () const override
 Does the component have maps of the low-field mobility?
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 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.
void GetFieldOnGrid (std::vector< std::vector< std::vector< ComponentGrid::Node > > > &efields)
 Get field values on all nodes.
bool AddElectricField (ComponentGrid *cmp, const double scale=1., const double xShift=0., const double yShift=0., const double zShift=0.)
 Add the field values of cmp to current grid.
bool AddElectricField (Component *cmp, const double scale)
 Add the field values of cmp to current grid.
bool GetNodeIndex (double x, const double y, const double z, unsigned int &i, unsigned int &j, unsigned int &k)
 Gives the closest node index (i, j, k) to coordinate (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).
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.
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 const std::vector< double > & DelayedSignalTimes (const std::string &)
 Return the time steps at which the delayed weighting potential/field are stored/evaluated.
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.
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 ()
virtual size_t GetNumberOfElements () const
 Return the number of mesh elements.
virtual bool GetElementNodes (const size_t, std::vector< size_t > &) const
 Get the indices of the nodes constituting a given element.
virtual bool GetElementRegion (const size_t, size_t &, bool &) const
 Get the region/material of a mesh element and a flag whether it is associated to an active medium.
virtual size_t GetNumberOfNodes () const
 Return the number of mesh nodes.
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 HasTownsendMap () const
 Does the component have maps of the Townsend 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 double StepSizeHint ()
virtual double CreateGPUTransferObject (ComponentGPU *&comp_gpu)
 Create and initialise GPU Transfer class.

Private Types

enum class  Format {
  Unknown , XY , XZ , XYZ ,
  IJ , IK , IJK , YXZ
}
enum class  Coordinates { Cartesian , Cylindrical }

Private Member Functions

bool LoadMesh (const std::string &filename, std::string format, const double scaleX)
 Read/determine mesh parameters from file.
bool LoadData (const std::string &filename, std::string format, const bool withPotential, const bool withFlag, const double scaleX, const double scaleF, const double scaleP, std::vector< std::vector< std::vector< Node > > > &field)
 Read electric field and potential from file.
bool LoadData (const std::string &filename, std::string format, const double scaleX, std::vector< std::vector< std::vector< double > > > &tab, const unsigned int col)
 Load scalar data (e. g. attachment coefficients) from file.
void Reset () override
 Reset the component.
void UpdatePeriodicity () override
 Verify periodicities.
bool GetField (const double x, const double y, const double z, const std::vector< std::vector< std::vector< Node > > > &field, double &fx, double &fy, double &fz, double &p, bool &active)
 Interpolation of the field and potential at a given point.
bool GetData (const double x, const double y, const double z, const std::vector< std::vector< std::vector< double > > > &table, double &value)
 Interpolation in a table of scalars.
double Reduce (const double xin, const double xmin, const double xmax, const bool simplePeriodic, const bool mirrorPeriodic, bool &isMirrored) const
 Reduce a coordinate to the basic cell (in case of periodicity).
void Initialise (std::vector< std::vector< std::vector< Node > > > &fields)
 Set the dimensions of a table according to the mesh.
Format GetFormat (std::string fmt)
 Decode a format string.

Private Attributes

Coordinates m_coordinates = Coordinates::Cartesian
Mediumm_medium = nullptr
std::vector< std::vector< std::vector< Node > > > m_efields
 Electric field values and potentials.
std::vector< std::vector< std::vector< Node > > > m_bfields
 Magnetic field values.
std::vector< std::vector< std::vector< Node > > > m_wfields
 Prompt weighting field values and potentials.
std::vector< std::vector< std::vector< std::vector< Node > > > > m_wdfields
 Delayed weighting field values and potentials.
std::vector< std::vector< std::vector< double > > > m_eAttachment
 Attachment maps for electrons and holes.
std::vector< std::vector< std::vector< double > > > m_hAttachment
std::vector< std::vector< std::vector< double > > > m_eMobility
 Mobility maps for electrons and holes.
std::vector< std::vector< std::vector< double > > > m_hMobility
std::vector< std::vector< std::vector< Node > > > m_eVelocity
 Velocity maps for electrons and holes.
std::vector< std::vector< std::vector< Node > > > m_hVelocity
std::vector< std::vector< std::vector< bool > > > m_active
 Active medium flag.
std::array< unsigned int, 3 > m_nX = {{1, 1, 1}}
std::array< double, 3 > m_xMin = {{0., 0., 0.}}
std::array< double, 3 > m_xMax = {{0., 0., 0.}}
std::array< double, 3 > m_sX = {{0., 0., 0.}}
bool m_hasMesh = false
bool m_hasPotential = false
std::array< double, 3 > m_wFieldOffset = {{0., 0., 0.}}
double m_pMin = 0.
double m_pMax = 0.

Additional Inherited Members

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.

Detailed Description

Component for interpolating field maps on a regular mesh.

Definition at line 15 of file ComponentGrid.hh.

Member Enumeration Documentation

◆ Coordinates

enum class Garfield::ComponentGrid::Coordinates
strongprivate
Enumerator
Cartesian 
Cylindrical 

Definition at line 227 of file ComponentGrid.hh.

227{ Cartesian, Cylindrical };

◆ Format

enum class Garfield::ComponentGrid::Format
strongprivate
Enumerator
Unknown 
XY 
XZ 
XYZ 
IJ 
IK 
IJK 
YXZ 

Definition at line 226 of file ComponentGrid.hh.

226{ Unknown, XY, XZ, XYZ, IJ, IK, IJK, YXZ };

Constructor & Destructor Documentation

◆ ComponentGrid()

Garfield::ComponentGrid::ComponentGrid ( )

Constructor.

◆ ~ComponentGrid()

Garfield::ComponentGrid::~ComponentGrid ( )
inline

Destructor.

Definition at line 20 of file ComponentGrid.hh.

20{}

Member Function Documentation

◆ AddElectricField() [1/2]

bool Garfield::ComponentGrid::AddElectricField ( Component * cmp,
const double scale )

Add the field values of cmp to current grid.

◆ AddElectricField() [2/2]

bool Garfield::ComponentGrid::AddElectricField ( ComponentGrid * cmp,
const double scale = 1.,
const double xShift = 0.,
const double yShift = 0.,
const double zShift = 0. )

Add the field values of cmp to current grid.

◆ Clear()

void Garfield::ComponentGrid::Clear ( )
inlineoverridevirtual

Reset.

Reimplemented from Garfield::Component.

Definition at line 149 of file ComponentGrid.hh.

149{ Reset(); }
void Reset() override
Reset the component.

◆ DelayedWeightingField()

void Garfield::ComponentGrid::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()

double Garfield::ComponentGrid::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.

◆ 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::ComponentGrid::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::ComponentGrid::ElectricField ( const double x,
const double y,
const double z,
double & ex,
double & ey,
double & ez,
Medium *& m,
int & status )
overridevirtual

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.

◆ ElectronAttachment()

bool Garfield::ComponentGrid::ElectronAttachment ( const double ,
const double ,
const double ,
double & eta )
overridevirtual

Get the electron attachment coefficient.

Reimplemented from Garfield::Component.

◆ ElectronMobility()

bool Garfield::ComponentGrid::ElectronMobility ( const double x,
const double y,
const double z,
double & mu )
overridevirtual

Reimplemented from Garfield::Component.

◆ ElectronVelocity()

bool Garfield::ComponentGrid::ElectronVelocity ( const double ,
const double ,
const double ,
double & vx,
double & vy,
double & vz )
overridevirtual

Get the electron drift velocity.

Reimplemented from Garfield::Component.

◆ GetBoundingBox()

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

Get the bounding box coordinates.

Reimplemented from Garfield::Component.

◆ GetData()

bool Garfield::ComponentGrid::GetData ( const double x,
const double y,
const double z,
const std::vector< std::vector< std::vector< double > > > & table,
double & value )
private

Interpolation in a table of scalars.

◆ GetElectricField()

bool Garfield::ComponentGrid::GetElectricField ( const unsigned int i,
const unsigned int j,
const unsigned int k,
double & v,
double & ex,
double & ey,
double & ez ) const

Return the field at a given node.

◆ GetElectricFieldRange()

bool Garfield::ComponentGrid::GetElectricFieldRange ( double & exmin,
double & exmax,
double & eymin,
double & eymax,
double & ezmin,
double & ezmax )

◆ GetElementaryCell()

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

◆ GetField()

bool Garfield::ComponentGrid::GetField ( const double x,
const double y,
const double z,
const std::vector< std::vector< std::vector< Node > > > & field,
double & fx,
double & fy,
double & fz,
double & p,
bool & active )
private

Interpolation of the field and potential at a given point.

◆ GetFieldOnGrid()

void Garfield::ComponentGrid::GetFieldOnGrid ( std::vector< std::vector< std::vector< ComponentGrid::Node > > > & efields)
inline

Get field values on all nodes.

Definition at line 209 of file ComponentGrid.hh.

210 {
211 efields = m_efields;
212 };
std::vector< std::vector< std::vector< Node > > > m_efields
Electric field values and potentials.

◆ GetFormat()

Format Garfield::ComponentGrid::GetFormat ( std::string fmt)
private

Decode a format string.

◆ GetMedium() [1/2]

Medium * Garfield::ComponentGrid::GetMedium ( ) const
inline

Get the medium.

Definition at line 119 of file ComponentGrid.hh.

119{ return m_medium; }

◆ GetMedium() [2/2]

Medium * Garfield::ComponentGrid::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.

◆ GetMesh()

bool Garfield::ComponentGrid::GetMesh ( unsigned int & nx,
unsigned int & ny,
unsigned int & nz,
double & xmin,
double & xmax,
double & ymin,
double & ymax,
double & zmin,
double & zmax ) const

Retrieve the parameters of the grid.

◆ GetNodeIndex()

bool Garfield::ComponentGrid::GetNodeIndex ( double x,
const double y,
const double z,
unsigned int & i,
unsigned int & j,
unsigned int & k )

Gives the closest node index (i, j, k) to coordinate (x, y, z).

◆ GetVoltageRange()

bool Garfield::ComponentGrid::GetVoltageRange ( double & vmin,
double & vmax )
overridevirtual

Calculate the voltage range [V].

Implements Garfield::Component.

◆ HasAttachmentMap()

bool Garfield::ComponentGrid::HasAttachmentMap ( ) const
inlineoverridevirtual

Does the component have maps of the attachment coefficient?

Reimplemented from Garfield::Component.

Definition at line 182 of file ComponentGrid.hh.

182 {
183 return !(m_eAttachment.empty() && m_hAttachment.empty());
184 }
std::vector< std::vector< std::vector< double > > > m_hAttachment
std::vector< std::vector< std::vector< double > > > m_eAttachment
Attachment maps for electrons and holes.

◆ HasMagneticField()

bool Garfield::ComponentGrid::HasMagneticField ( ) const
overridevirtual

Does the component have a non-zero magnetic field?

Reimplemented from Garfield::Component.

◆ HasMobilityMap()

bool Garfield::ComponentGrid::HasMobilityMap ( ) const
inlineoverridevirtual

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

Reimplemented from Garfield::Component.

Definition at line 189 of file ComponentGrid.hh.

189 {
190 return !(m_eMobility.empty() && m_hMobility.empty());
191 }
std::vector< std::vector< std::vector< double > > > m_hMobility
std::vector< std::vector< std::vector< double > > > m_eMobility
Mobility maps for electrons and holes.

◆ HasVelocityMap()

bool Garfield::ComponentGrid::HasVelocityMap ( ) const
inlineoverridevirtual

Does the component have velocity maps?

Reimplemented from Garfield::Component.

Definition at line 197 of file ComponentGrid.hh.

197 {
198 return !(m_eVelocity.empty() && m_hVelocity.empty());
199 }
std::vector< std::vector< std::vector< Node > > > m_eVelocity
Velocity maps for electrons and holes.
std::vector< std::vector< std::vector< Node > > > m_hVelocity

◆ HoleAttachment()

bool Garfield::ComponentGrid::HoleAttachment ( const double ,
const double ,
const double ,
double & eta )
overridevirtual

Get the hole attachment coefficient.

Reimplemented from Garfield::Component.

◆ HoleMobility()

bool Garfield::ComponentGrid::HoleMobility ( const double ,
const double ,
const double ,
double & mu )
overridevirtual

Get the hole Mobility coefficient.

Reimplemented from Garfield::Component.

◆ HoleVelocity()

bool Garfield::ComponentGrid::HoleVelocity ( const double ,
const double ,
const double ,
double & vx,
double & vy,
double & vz )
overridevirtual

Get the hole drift velocity.

Reimplemented from Garfield::Component.

◆ Initialise()

void Garfield::ComponentGrid::Initialise ( std::vector< std::vector< std::vector< Node > > > & fields)
private

Set the dimensions of a table according to the mesh.

◆ LoadData() [1/2]

bool Garfield::ComponentGrid::LoadData ( const std::string & filename,
std::string format,
const bool withPotential,
const bool withFlag,
const double scaleX,
const double scaleF,
const double scaleP,
std::vector< std::vector< std::vector< Node > > > & field )
private

Read electric field and potential from file.

◆ LoadData() [2/2]

bool Garfield::ComponentGrid::LoadData ( const std::string & filename,
std::string format,
const double scaleX,
std::vector< std::vector< std::vector< double > > > & tab,
const unsigned int col )
private

Load scalar data (e. g. attachment coefficients) from file.

◆ LoadElectricField()

bool Garfield::ComponentGrid::LoadElectricField ( const std::string & filename,
const std::string & format,
const bool withPotential,
const bool withFlag,
const double scaleX = 1.,
const double scaleE = 1.,
const double scaleP = 1. )

Import electric field and potential values from a file.

The file is supposed to contain one line for each grid point starting with

  • either two or three floating point numbers, specifying the coordinates (in cm) of the grid node or
  • two or three integers specifying the index of the node,

followed by

  • two or three floating point numbers for the electric field (in V/cm), and (depending on the value of withPotential and withFlag),
  • a floating point number specifying the potential (in V), and
  • an integer flag indicating whether the point is in an active region (1) or not (0).

Format types are:

  • "xy", "xz", "xyz": nodes are specified by their coordinates
  • "ij", "ik", "ijk": nodes are specified by their indices

If cylindrical coordinates are used, the first coordinate (x) corresponds to the radial distance and the second coordinate (y) corresponds to the azimuth (in radian).

◆ LoadElectronAttachment()

bool Garfield::ComponentGrid::LoadElectronAttachment ( const std::string & fname,
const std::string & fmt,
const unsigned int col,
const double scaleX = 1. )

Import electron attachment coefficients from a file.

Parameters
fnamename of the text file.
fmtformat string, see LoadElectricField.
colcolumn in the file which has the attachment coefficient.
scaleXscaling factor to be applied to the coordinates.

◆ LoadElectronVelocity()

bool Garfield::ComponentGrid::LoadElectronVelocity ( const std::string & fname,
const std::string & fmt,
const double scaleX = 1.,
const double scaleV = 1.e-9 )

Import a map of electron drift velocities from a file.

Parameters
fnamename of the text file.
fmtformat string, see LoadElectricField
scaleXscaling factor to be applied to the coordinates.
scaleVscaling factor to be applied to the velocity components.

◆ LoadHoleAttachment()

bool Garfield::ComponentGrid::LoadHoleAttachment ( const std::string & fname,
const std::string & fmt,
const unsigned int col,
const double scaleX = 1. )

Import hole attachment coefficients from a file.

◆ LoadHoleVelocity()

bool Garfield::ComponentGrid::LoadHoleVelocity ( const std::string & fname,
const std::string & fmt,
const double scaleX = 1.,
const double scaleV = 1.e-9 )

Import a map of hole drift velocities from a file.

◆ LoadMagneticField()

bool Garfield::ComponentGrid::LoadMagneticField ( const std::string & filename,
const std::string & format,
const double scaleX = 1.,
const double scaleB = 1. )

Import magnetic field values from a file.

◆ LoadMesh()

bool Garfield::ComponentGrid::LoadMesh ( const std::string & filename,
std::string format,
const double scaleX )
private

Read/determine mesh parameters from file.

◆ LoadWeightingField() [1/2]

bool Garfield::ComponentGrid::LoadWeightingField ( const std::string & filename,
const std::string & format,
const bool withPotential,
const double scaleX = 1.,
const double scaleE = 1.,
const double scaleP = 1. )

Import (prompt) weighting field from file.

◆ LoadWeightingField() [2/2]

bool Garfield::ComponentGrid::LoadWeightingField ( const std::string & filename,
const std::string & format,
const double time,
const bool withPotential,
const double scaleX = 1.,
const double scaleE = 1.,
const double scaleP = 1. )

Import delayed weighting field from file.

◆ MagneticField()

void Garfield::ComponentGrid::MagneticField ( const double x,
const double y,
const double z,
double & bx,
double & by,
double & bz,
int & status )
overridevirtual

Calculate the magnetic field at a given point.

Parameters
x,y,zcoordinates [cm].
bx,by,bzcomponents of the magnetic field [Tesla].
statusstatus flag.

Reimplemented from Garfield::Component.

◆ Print()

void Garfield::ComponentGrid::Print ( )

Print information about the mesh and the available data.

◆ Reduce()

double Garfield::ComponentGrid::Reduce ( const double xin,
const double xmin,
const double xmax,
const bool simplePeriodic,
const bool mirrorPeriodic,
bool & isMirrored ) const
private

Reduce a coordinate to the basic cell (in case of periodicity).

◆ Reset()

void Garfield::ComponentGrid::Reset ( )
overrideprivatevirtual

Reset the component.

Implements Garfield::Component.

◆ SaveElectricField() [1/2]

bool Garfield::ComponentGrid::SaveElectricField ( Component * cmp)

Export the electric field and potential of a component.

Parameters
cmpComponent object for which to export the field/potential

◆ SaveElectricField() [2/2]

bool Garfield::ComponentGrid::SaveElectricField ( Component * cmp,
const std::string & filename,
const std::string & fmt )

Export the electric field and potential of a component to a text file.

Parameters
cmpComponent object for which to export the field/potential
filenamename of the text file
fmtformat string, see LoadElectricField

◆ SaveElectricFieldROOT()

bool Garfield::ComponentGrid::SaveElectricFieldROOT ( Component * cmp,
const std::string & filename,
const std::string & fmt )

◆ SaveWeightingField()

bool Garfield::ComponentGrid::SaveWeightingField ( Component * cmp,
const std::string & id,
const std::string & filename,
const std::string & fmt )

Export the weighting field and potential of a component to a text file.

Parameters
cmpComponent object for which to export the field/potential
ididentifier of the weighting field
filenamename of the text file
fmtformat string, see LoadElectricField

◆ SetCartesianCoordinates()

void Garfield::ComponentGrid::SetCartesianCoordinates ( )
inline

Use Cartesian coordinates (default).

Definition at line 37 of file ComponentGrid.hh.

◆ SetCylindricalCoordinates()

void Garfield::ComponentGrid::SetCylindricalCoordinates ( )

Use cylindrical coordinates.

◆ SetMedium()

void Garfield::ComponentGrid::SetMedium ( Medium * m)

Set the medium.

◆ SetMesh()

bool Garfield::ComponentGrid::SetMesh ( const unsigned int nx,
const unsigned int ny,
const unsigned int nz,
const double xmin,
const double xmax,
const double ymin,
const double ymax,
const double zmin,
const double zmax )

Define the grid.

Parameters
nx,ny,nznumber of nodes along $x, y, z$.
xmin,xmaxrange along $x$.
ymin,ymaxrange along $y$.
zmin,zmaxrange along $z$.

◆ SetWeightingFieldOffset()

void Garfield::ComponentGrid::SetWeightingFieldOffset ( const double x,
const double y,
const double z )

Offset coordinates in the weighting field, such that the same numerical weighting field map can be used for electrodes at different positions.

◆ UpdatePeriodicity()

void Garfield::ComponentGrid::UpdatePeriodicity ( )
overrideprivatevirtual

Verify periodicities.

Implements Garfield::Component.

◆ WeightingField()

void Garfield::ComponentGrid::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()

double Garfield::ComponentGrid::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_active

std::vector<std::vector<std::vector<bool> > > Garfield::ComponentGrid::m_active
private

Active medium flag.

Definition at line 250 of file ComponentGrid.hh.

◆ m_bfields

std::vector<std::vector<std::vector<Node> > > Garfield::ComponentGrid::m_bfields
private

Magnetic field values.

Definition at line 235 of file ComponentGrid.hh.

◆ m_coordinates

Coordinates Garfield::ComponentGrid::m_coordinates = Coordinates::Cartesian
private

Definition at line 228 of file ComponentGrid.hh.

◆ m_eAttachment

std::vector<std::vector<std::vector<double> > > Garfield::ComponentGrid::m_eAttachment
private

Attachment maps for electrons and holes.

Definition at line 241 of file ComponentGrid.hh.

◆ m_efields

std::vector<std::vector<std::vector<Node> > > Garfield::ComponentGrid::m_efields
private

Electric field values and potentials.

Definition at line 233 of file ComponentGrid.hh.

◆ m_eMobility

std::vector<std::vector<std::vector<double> > > Garfield::ComponentGrid::m_eMobility
private

Mobility maps for electrons and holes.

Definition at line 244 of file ComponentGrid.hh.

◆ m_eVelocity

std::vector<std::vector<std::vector<Node> > > Garfield::ComponentGrid::m_eVelocity
private

Velocity maps for electrons and holes.

Definition at line 247 of file ComponentGrid.hh.

◆ m_hasMesh

bool Garfield::ComponentGrid::m_hasMesh = false
private

Definition at line 258 of file ComponentGrid.hh.

◆ m_hasPotential

bool Garfield::ComponentGrid::m_hasPotential = false
private

Definition at line 259 of file ComponentGrid.hh.

◆ m_hAttachment

std::vector<std::vector<std::vector<double> > > Garfield::ComponentGrid::m_hAttachment
private

Definition at line 242 of file ComponentGrid.hh.

◆ m_hMobility

std::vector<std::vector<std::vector<double> > > Garfield::ComponentGrid::m_hMobility
private

Definition at line 245 of file ComponentGrid.hh.

◆ m_hVelocity

std::vector<std::vector<std::vector<Node> > > Garfield::ComponentGrid::m_hVelocity
private

Definition at line 248 of file ComponentGrid.hh.

◆ m_medium

Medium* Garfield::ComponentGrid::m_medium = nullptr
private

Definition at line 230 of file ComponentGrid.hh.

◆ m_nX

std::array<unsigned int, 3> Garfield::ComponentGrid::m_nX = {{1, 1, 1}}
private

Definition at line 253 of file ComponentGrid.hh.

253{{1, 1, 1}};

◆ m_pMax

double Garfield::ComponentGrid::m_pMax = 0.
private

Definition at line 265 of file ComponentGrid.hh.

◆ m_pMin

double Garfield::ComponentGrid::m_pMin = 0.
private

Definition at line 265 of file ComponentGrid.hh.

◆ m_sX

std::array<double, 3> Garfield::ComponentGrid::m_sX = {{0., 0., 0.}}
private

Definition at line 256 of file ComponentGrid.hh.

256{{0., 0., 0.}};

◆ m_wdfields

std::vector<std::vector<std::vector<std::vector<Node> > > > Garfield::ComponentGrid::m_wdfields
private

Delayed weighting field values and potentials.

Definition at line 239 of file ComponentGrid.hh.

◆ m_wFieldOffset

std::array<double, 3> Garfield::ComponentGrid::m_wFieldOffset = {{0., 0., 0.}}
private

Definition at line 262 of file ComponentGrid.hh.

262{{0., 0., 0.}};

◆ m_wfields

std::vector<std::vector<std::vector<Node> > > Garfield::ComponentGrid::m_wfields
private

Prompt weighting field values and potentials.

Definition at line 237 of file ComponentGrid.hh.

◆ m_xMax

std::array<double, 3> Garfield::ComponentGrid::m_xMax = {{0., 0., 0.}}
private

Definition at line 255 of file ComponentGrid.hh.

255{{0., 0., 0.}};

◆ m_xMin

std::array<double, 3> Garfield::ComponentGrid::m_xMin = {{0., 0., 0.}}
private

Definition at line 254 of file ComponentGrid.hh.

254{{0., 0., 0.}};

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