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

Component for interpolating field maps stored in a regular mesh. More...

#include <ComponentVoxel.hh>

Inheritance diagram for Garfield::ComponentVoxel:
Garfield::Component

Classes

struct  Element

Public Member Functions

 ComponentVoxel ()
 Constructor.
 ~ComponentVoxel ()
 Destructor.
void EnableInterpolation (const bool on=true)
 Interpolate between field values at the element centres.
void 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 LoadElectricField (const std::string &filename, const std::string &format, const bool withPotential, const bool withRegion, 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.
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 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.
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.
void SetMedium (const unsigned int i, Medium *m)
 Set the medium in region i.
MediumGetMedium (const unsigned int i) const
 Get the medium in region i.
void PrintRegions () const
 Print all regions.
bool GetElement (const double xi, const double yi, const double zi, unsigned int &i, unsigned int &j, unsigned int &k, bool &xMirrored, bool &yMirrored, bool &zMirrored) const
 Return the indices of the element at a given point.
bool GetElement (const unsigned int i, const unsigned int j, const unsigned int k, double &v, double &ex, double &ey, double &ez) const
 Return the field for an element with given index.
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.
bool HasMagneticField () const override
 Does the component have a non-zero magnetic field?
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.
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 HasAttachmentMap () const
 Does the component have maps of the attachment coefficient?
virtual bool HasMobilityMap () const
 Does the component have maps of the low-field mobility?
virtual bool HasVelocityMap () const
 Does the component have velocity maps?
virtual bool ElectronAttachment (const double, const double, const double, double &eta)
 Get the electron attachment coefficient.
virtual bool HoleAttachment (const double, const double, const double, double &eta)
 Get the hole attachment coefficient.
virtual bool ElectronMobility (const double, const double, const double, double &mu)
virtual bool HoleMobility (const double, const double, const double, double &mu)
 Get the hole Mobility coefficient.
virtual bool ElectronTownsend (const double, const double, const double, double &alpha)
 Get the electron Townsend coefficient.
virtual bool HoleTownsend (const double, const double, const double, double &alpha)
 Get the hole Townsend coefficient.
virtual bool ElectronVelocity (const double, const double, const double, double &vx, double &vy, double &vz)
 Get the electron drift velocity.
virtual bool HoleVelocity (const double, const double, const double, double &vx, double &vy, double &vz)
 Get the hole drift velocity.
virtual double StepSizeHint ()
virtual double CreateGPUTransferObject (ComponentGPU *&comp_gpu)
 Create and initialise GPU Transfer class.

Private Member Functions

bool LoadData (const std::string &filename, std::string format, const bool withPotential, const bool withRegion, const double scaleX, const double scaleF, const double scaleP, std::vector< std::vector< std::vector< Element > > > &field)
 Read data 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< Element > > > &field, double &fx, double &fy, double &fz, double &p, int &region)
 Look up/interpolate the field at a given point.
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< Element > > > &fields)
 Set the dimensions of a table according to the mesh.
void InitialiseRegions ()

Private Attributes

std::vector< Medium * > m_media
std::vector< std::vector< std::vector< int > > > m_regions
 Region indices.
std::vector< std::vector< std::vector< Element > > > m_efields
 Electric field values and potentials at each mesh element.
std::vector< std::vector< std::vector< Element > > > m_bfields
 Magnetic field values at each mesh element.
std::vector< std::vector< std::vector< Element > > > m_wfields
 Prompt weighting field values and potentials at each mesh element.
std::vector< std::vector< std::vector< std::vector< Element > > > > m_wdfields
 Delayed weighting field values and potentials at each mesh element.
unsigned int m_nX = 0
unsigned int m_nY = 0
unsigned int m_nZ = 0
double m_xMin = 0.
double m_yMin = 0.
double m_zMin = 0.
double m_xMax = 0.
double m_yMax = 0.
double m_zMax = 0.
double m_dx = 0.
double m_dy = 0.
double m_dz = 0.
bool m_interpolate = false
bool m_hasMesh = false
bool m_hasPotential = false
bool m_hasEfield = false
bool m_hasBfield = false
bool m_hasWfield = false
double m_wField_xOffset = 0.
double m_wField_yOffset = 0.
double m_wField_zOffset = 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 stored in a regular mesh.

Definition at line 10 of file ComponentVoxel.hh.

Constructor & Destructor Documentation

◆ ComponentVoxel()

Garfield::ComponentVoxel::ComponentVoxel ( )

Constructor.

◆ ~ComponentVoxel()

Garfield::ComponentVoxel::~ComponentVoxel ( )
inline

Destructor.

Definition at line 15 of file ComponentVoxel.hh.

15{}

Member Function Documentation

◆ Clear()

void Garfield::ComponentVoxel::Clear ( )
inlineoverridevirtual

Reset.

Reimplemented from Garfield::Component.

Definition at line 85 of file ComponentVoxel.hh.

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

◆ DelayedWeightingField()

void Garfield::ComponentVoxel::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::ComponentVoxel::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::ComponentVoxel::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::ComponentVoxel::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.

◆ EnableInterpolation()

void Garfield::ComponentVoxel::EnableInterpolation ( const bool on = true)
inline

Interpolate between field values at the element centres.

Definition at line 18 of file ComponentVoxel.hh.

◆ GetBoundingBox()

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

Get the bounding box coordinates.

Reimplemented from Garfield::Component.

◆ GetElectricFieldRange()

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

◆ GetElement() [1/2]

bool Garfield::ComponentVoxel::GetElement ( const double xi,
const double yi,
const double zi,
unsigned int & i,
unsigned int & j,
unsigned int & k,
bool & xMirrored,
bool & yMirrored,
bool & zMirrored ) const

Return the indices of the element at a given point.

◆ GetElement() [2/2]

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

Return the field for an element with given index.

◆ GetElementaryCell()

bool Garfield::ComponentVoxel::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::ComponentVoxel::GetField ( const double x,
const double y,
const double z,
const std::vector< std::vector< std::vector< Element > > > & field,
double & fx,
double & fy,
double & fz,
double & p,
int & region )
private

Look up/interpolate the field at a given point.

◆ GetMedium() [1/2]

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

◆ GetMedium() [2/2]

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

Get the medium in region i.

◆ GetVoltageRange()

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

Calculate the voltage range [V].

Implements Garfield::Component.

◆ HasMagneticField()

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

Does the component have a non-zero magnetic field?

Reimplemented from Garfield::Component.

◆ Initialise()

void Garfield::ComponentVoxel::Initialise ( std::vector< std::vector< std::vector< Element > > > & fields)
private

Set the dimensions of a table according to the mesh.

◆ InitialiseRegions()

void Garfield::ComponentVoxel::InitialiseRegions ( )
private

◆ LoadData()

bool Garfield::ComponentVoxel::LoadData ( const std::string & filename,
std::string format,
const bool withPotential,
const bool withRegion,
const double scaleX,
const double scaleF,
const double scaleP,
std::vector< std::vector< std::vector< Element > > > & field )
private

Read data from file.

◆ LoadElectricField()

bool Garfield::ComponentVoxel::LoadElectricField ( const std::string & filename,
const std::string & format,
const bool withPotential,
const bool withRegion,
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 mesh point starting with

  • either two or three floating point numbers, specifying the coordinates (in cm) of the element centre or
  • two or three integers specifying the index of the element in the mesh,

followed by

  • two or three floating point numbers for the electric field (in V/cm), and (depending on the values of withPotential and withRegion),
  • a floating point number specifying the potential (in V), and
  • an integer specifying the "region" of the element.

Format types are:

  • "xy", "xyz": elements are specified by the coordinates of their centres
  • "ij", "ijk": elements are specified by their indices

◆ LoadMagneticField()

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

Import magnetic field values from a file.

◆ LoadWeightingField() [1/2]

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

◆ PrintRegions()

void Garfield::ComponentVoxel::PrintRegions ( ) const

Print all regions.

◆ Reduce()

double Garfield::ComponentVoxel::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::ComponentVoxel::Reset ( )
overrideprivatevirtual

Reset the component.

Implements Garfield::Component.

◆ SetMedium()

void Garfield::ComponentVoxel::SetMedium ( const unsigned int i,
Medium * m )

Set the medium in region i.

◆ SetMesh()

void Garfield::ComponentVoxel::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 bins along x, y, z.
xmin,xmaxrange along $x$.
ymin,ymaxrange along $y$.
zmin,zmaxrange along $z$.

◆ SetWeightingFieldOffset()

void Garfield::ComponentVoxel::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::ComponentVoxel::UpdatePeriodicity ( )
overrideprivatevirtual

Verify periodicities.

Implements Garfield::Component.

◆ WeightingField()

void Garfield::ComponentVoxel::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::ComponentVoxel::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_bfields

std::vector<std::vector<std::vector<Element> > > Garfield::ComponentVoxel::m_bfields
private

Magnetic field values at each mesh element.

Definition at line 129 of file ComponentVoxel.hh.

◆ m_dx

double Garfield::ComponentVoxel::m_dx = 0.
private

Definition at line 139 of file ComponentVoxel.hh.

◆ m_dy

double Garfield::ComponentVoxel::m_dy = 0.
private

Definition at line 139 of file ComponentVoxel.hh.

◆ m_dz

double Garfield::ComponentVoxel::m_dz = 0.
private

Definition at line 139 of file ComponentVoxel.hh.

◆ m_efields

std::vector<std::vector<std::vector<Element> > > Garfield::ComponentVoxel::m_efields
private

Electric field values and potentials at each mesh element.

Definition at line 127 of file ComponentVoxel.hh.

◆ m_hasBfield

bool Garfield::ComponentVoxel::m_hasBfield = false
private

Definition at line 146 of file ComponentVoxel.hh.

◆ m_hasEfield

bool Garfield::ComponentVoxel::m_hasEfield = false
private

Definition at line 145 of file ComponentVoxel.hh.

◆ m_hasMesh

bool Garfield::ComponentVoxel::m_hasMesh = false
private

Definition at line 143 of file ComponentVoxel.hh.

◆ m_hasPotential

bool Garfield::ComponentVoxel::m_hasPotential = false
private

Definition at line 144 of file ComponentVoxel.hh.

◆ m_hasWfield

bool Garfield::ComponentVoxel::m_hasWfield = false
private

Definition at line 147 of file ComponentVoxel.hh.

◆ m_interpolate

bool Garfield::ComponentVoxel::m_interpolate = false
private

Definition at line 141 of file ComponentVoxel.hh.

◆ m_media

std::vector<Medium*> Garfield::ComponentVoxel::m_media
private

Definition at line 118 of file ComponentVoxel.hh.

◆ m_nX

unsigned int Garfield::ComponentVoxel::m_nX = 0
private

Definition at line 136 of file ComponentVoxel.hh.

◆ m_nY

unsigned int Garfield::ComponentVoxel::m_nY = 0
private

Definition at line 136 of file ComponentVoxel.hh.

◆ m_nZ

unsigned int Garfield::ComponentVoxel::m_nZ = 0
private

Definition at line 136 of file ComponentVoxel.hh.

◆ m_pMax

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

Definition at line 155 of file ComponentVoxel.hh.

◆ m_pMin

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

Definition at line 155 of file ComponentVoxel.hh.

◆ m_regions

std::vector<std::vector<std::vector<int> > > Garfield::ComponentVoxel::m_regions
private

Region indices.

Definition at line 125 of file ComponentVoxel.hh.

◆ m_wdfields

std::vector<std::vector<std::vector<std::vector<Element> > > > Garfield::ComponentVoxel::m_wdfields
private

Delayed weighting field values and potentials at each mesh element.

Definition at line 133 of file ComponentVoxel.hh.

◆ m_wField_xOffset

double Garfield::ComponentVoxel::m_wField_xOffset = 0.
private

Definition at line 150 of file ComponentVoxel.hh.

◆ m_wField_yOffset

double Garfield::ComponentVoxel::m_wField_yOffset = 0.
private

Definition at line 151 of file ComponentVoxel.hh.

◆ m_wField_zOffset

double Garfield::ComponentVoxel::m_wField_zOffset = 0.
private

Definition at line 152 of file ComponentVoxel.hh.

◆ m_wfields

std::vector<std::vector<std::vector<Element> > > Garfield::ComponentVoxel::m_wfields
private

Prompt weighting field values and potentials at each mesh element.

Definition at line 131 of file ComponentVoxel.hh.

◆ m_xMax

double Garfield::ComponentVoxel::m_xMax = 0.
private

Definition at line 138 of file ComponentVoxel.hh.

◆ m_xMin

double Garfield::ComponentVoxel::m_xMin = 0.
private

Definition at line 137 of file ComponentVoxel.hh.

◆ m_yMax

double Garfield::ComponentVoxel::m_yMax = 0.
private

Definition at line 138 of file ComponentVoxel.hh.

◆ m_yMin

double Garfield::ComponentVoxel::m_yMin = 0.
private

Definition at line 137 of file ComponentVoxel.hh.

◆ m_zMax

double Garfield::ComponentVoxel::m_zMax = 0.
private

Definition at line 138 of file ComponentVoxel.hh.

◆ m_zMin

double Garfield::ComponentVoxel::m_zMin = 0.
private

Definition at line 137 of file ComponentVoxel.hh.


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