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:

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 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.
Medium * GetMedium () 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
void ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, double &v, Medium *&m, int &status) override
void ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status) override
void WeightingField (const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label) override
double WeightingPotential (const double x, const double y, const double z, const std::string &label) override
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
double DelayedWeightingPotential (const double x, const double y, const double z, const double t, const std::string &label) override
void MagneticField (const double x, const double y, const double z, double &bx, double &by, double &bz, int &status) override
Medium * GetMedium (const double x, const double y, const double z) override
bool GetVoltageRange (double &vmin, double &vmax) override
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
bool GetElementaryCell (double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax) override
bool HasMagneticField () const override
bool HasAttachmentMap () const override
bool ElectronAttachment (const double x, const double y, const double z, double &att) override
bool HoleAttachment (const double x, const double y, const double z, double &att) override
bool HasMobilityMap () const override
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
bool HasVelocityMap () const override
bool ElectronVelocity (const double x, const double y, const double z, double &vx, double &vy, double &vz) override
bool HoleVelocity (const double x, const double y, const double z, double &vx, double &vy, double &vz) override

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
void UpdatePeriodicity () override
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
Medium * m_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.

Detailed Description

Component for interpolating field maps on a regular mesh.

Definition at line 14 of file ComponentGrid.hh.

Member Enumeration Documentation

◆ Coordinates

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

Definition at line 197 of file ComponentGrid.hh.

197{ Cartesian, Cylindrical };

◆ Format

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

Definition at line 196 of file ComponentGrid.hh.

196{ 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 19 of file ComponentGrid.hh.

19{}

Member Function Documentation

◆ Clear()

void Garfield::ComponentGrid::Clear ( )
inlineoverride

Definition at line 139 of file ComponentGrid.hh.

139{ Reset(); }
void Reset() override

◆ 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 )
override

◆ DelayedWeightingPotential()

double Garfield::ComponentGrid::DelayedWeightingPotential ( const double x,
const double y,
const double z,
const double t,
const std::string & label )
override

◆ ElectricField() [1/2]

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 )
override

◆ ElectricField() [2/2]

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

◆ ElectronAttachment()

bool Garfield::ComponentGrid::ElectronAttachment ( const double x,
const double y,
const double z,
double & att )
override

◆ ElectronMobility()

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

◆ ElectronVelocity()

bool Garfield::ComponentGrid::ElectronVelocity ( const double x,
const double y,
const double z,
double & vx,
double & vy,
double & vz )
override

◆ GetBoundingBox()

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

◆ 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 )
override

◆ 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.

◆ 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 109 of file ComponentGrid.hh.

109{ return m_medium; }

◆ GetMedium() [2/2]

Medium * Garfield::ComponentGrid::GetMedium ( const double x,
const double y,
const double z )
override

◆ 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.

◆ GetVoltageRange()

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

◆ HasAttachmentMap()

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

Definition at line 172 of file ComponentGrid.hh.

172 {
173 return !(m_eAttachment.empty() && m_hAttachment.empty());
174 }
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
override

◆ HasMobilityMap()

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

Definition at line 179 of file ComponentGrid.hh.

179 {
180 return !(m_eMobility.empty() && m_hMobility.empty());
181 }
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
inlineoverride

Definition at line 187 of file ComponentGrid.hh.

187 {
188 return !(m_eVelocity.empty() && m_hVelocity.empty());
189 }
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 x,
const double y,
const double z,
double & att )
override

◆ HoleMobility()

bool Garfield::ComponentGrid::HoleMobility ( const double x,
const double y,
const double z,
double & mu )
override

◆ HoleVelocity()

bool Garfield::ComponentGrid::HoleVelocity ( const double x,
const double y,
const double z,
double & vx,
double & vy,
double & vz )
override

◆ 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 )
override

◆ 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 ( )
overrideprivate

◆ SaveElectricField()

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

◆ 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 36 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 ( )
overrideprivate

◆ WeightingField()

void Garfield::ComponentGrid::WeightingField ( const double x,
const double y,
const double z,
double & wx,
double & wy,
double & wz,
const std::string & label )
override

◆ WeightingPotential()

double Garfield::ComponentGrid::WeightingPotential ( const double x,
const double y,
const double z,
const std::string & label )
override

Member Data Documentation

◆ m_active

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

Active medium flag.

Definition at line 224 of file ComponentGrid.hh.

◆ m_bfields

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

Magnetic field values.

Definition at line 209 of file ComponentGrid.hh.

◆ m_coordinates

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

Definition at line 198 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 215 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 207 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 218 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 221 of file ComponentGrid.hh.

◆ m_hasMesh

bool Garfield::ComponentGrid::m_hasMesh = false
private

Definition at line 232 of file ComponentGrid.hh.

◆ m_hasPotential

bool Garfield::ComponentGrid::m_hasPotential = false
private

Definition at line 233 of file ComponentGrid.hh.

◆ m_hAttachment

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

Definition at line 216 of file ComponentGrid.hh.

◆ m_hMobility

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

Definition at line 219 of file ComponentGrid.hh.

◆ m_hVelocity

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

Definition at line 222 of file ComponentGrid.hh.

◆ m_medium

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

Definition at line 200 of file ComponentGrid.hh.

◆ m_nX

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

Definition at line 227 of file ComponentGrid.hh.

227{{1, 1, 1}};

◆ m_pMax

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

Definition at line 239 of file ComponentGrid.hh.

◆ m_pMin

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

Definition at line 239 of file ComponentGrid.hh.

◆ m_sX

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

Definition at line 230 of file ComponentGrid.hh.

230{{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 213 of file ComponentGrid.hh.

◆ m_wFieldOffset

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

Definition at line 236 of file ComponentGrid.hh.

236{{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 211 of file ComponentGrid.hh.

◆ m_xMax

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

Definition at line 229 of file ComponentGrid.hh.

229{{0., 0., 0.}};

◆ m_xMin

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

Definition at line 228 of file ComponentGrid.hh.

228{{0., 0., 0.}};

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