![]() |
Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
|
Component for interpolating field maps on a regular mesh. More...
#include <ComponentGrid.hh>
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. | |
| 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 |
| 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. | |
| Medium * | GetMedium (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 | |
| void | EnablePeriodicityY (const bool on=true) |
| Enable simple periodicity in the | |
| void | EnablePeriodicityZ (const bool on=true) |
| Enable simple periodicity in the | |
| void | IsPeriodic (bool &perx, bool &pery, bool &perz) |
| Return periodicity flags. | |
| void | EnableMirrorPeriodicityX (const bool on=true) |
| Enable mirror periodicity in the | |
| void | EnableMirrorPeriodicityY (const bool on=true) |
| Enable mirror periodicity in the | |
| void | EnableMirrorPeriodicityZ (const bool on=true) |
| Enable mirror periodicity in the | |
| void | IsMirrorPeriodic (bool &perx, bool &pery, bool &perz) |
| Return mirror periodicity flags. | |
| void | EnableAxialPeriodicityX (const bool on=true) |
| Enable axial periodicity in the | |
| void | EnableAxialPeriodicityY (const bool on=true) |
| Enable axial periodicity in the | |
| void | EnableAxialPeriodicityZ (const bool on=true) |
| Enable axial periodicity in the | |
| void | IsAxiallyPeriodic (bool &perx, bool &pery, bool &perz) |
| Return axial periodicity flags. | |
| void | EnableRotationSymmetryX (const bool on=true) |
| Enable rotation symmetry around the | |
| void | EnableRotationSymmetryY (const bool on=true) |
| Enable rotation symmetry around the | |
| void | EnableRotationSymmetryZ (const bool on=true) |
| Enable rotation symmetry around the | |
| 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 | |
| void | EnableTriangleSymmetricXZ (const bool on=true, const bool oct=2) |
| Enable triangular periodicity in the | |
| void | EnableTriangleSymmetricYZ (const bool on=true, const bool oct=2) |
| Enable triangular periodicity in the | |
| 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 |
| 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. |
Additional Inherited Members | |
| Protected Attributes inherited from Garfield::Component | |
| std::string | m_className {"Component"} |
| Class name. | |
| Geometry * | m_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. | |
Component for interpolating field maps on a regular mesh.
Definition at line 15 of file ComponentGrid.hh.
|
strongprivate |
| Enumerator | |
|---|---|
| Cartesian | |
| Cylindrical | |
Definition at line 227 of file ComponentGrid.hh.
|
strongprivate |
| Enumerator | |
|---|---|
| Unknown | |
| XY | |
| XZ | |
| XYZ | |
| IJ | |
| IK | |
| IJK | |
| YXZ | |
Definition at line 226 of file ComponentGrid.hh.
| Garfield::ComponentGrid::ComponentGrid | ( | ) |
Constructor.
|
inline |
| bool Garfield::ComponentGrid::AddElectricField | ( | Component * | cmp, |
| const double | scale ) |
Add the field values of cmp to current grid.
| 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.
|
inlineoverridevirtual |
Reset.
Reimplemented from Garfield::Component.
Definition at line 149 of file ComponentGrid.hh.
|
overridevirtual |
Calculate the delayed weighting field at a given point and time and for a given electrode.
| x,y,z | coordinates [cm]. |
| t | time [ns]. |
| wx,wy,wz | components of the weighting field [1/cm]. |
| label | name of the electrode |
Reimplemented from Garfield::Component.
|
overridevirtual |
Calculate the delayed weighting potential at a given point and time and for a given electrode.
| x,y,z | coordinates [cm]. |
| t | time [ns]. |
| label | name of the electrode |
Reimplemented from Garfield::Component.
| 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).
|
overridevirtual |
Calculate the drift field [V/cm] and potential [V] at (x, y, z).
Implements Garfield::Component.
|
overridevirtual |
Calculate the drift field at given point.
| x,y,z | coordinates [cm]. |
| ex,ey,ez | components of the electric field [V/cm]. |
| m | pointer to the medium at this location. |
| status | status 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.
|
overridevirtual |
Get the electron attachment coefficient.
Reimplemented from Garfield::Component.
|
overridevirtual |
Reimplemented from Garfield::Component.
|
overridevirtual |
Get the electron drift velocity.
Reimplemented from Garfield::Component.
|
overridevirtual |
Get the bounding box coordinates.
Reimplemented from Garfield::Component.
|
private |
Interpolation in a table of scalars.
| 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.
| bool Garfield::ComponentGrid::GetElectricFieldRange | ( | double & | exmin, |
| double & | exmax, | ||
| double & | eymin, | ||
| double & | eymax, | ||
| double & | ezmin, | ||
| double & | ezmax ) |
|
overridevirtual |
Get the coordinates of the elementary cell.
Reimplemented from Garfield::Component.
|
private |
Interpolation of the field and potential at a given point.
|
inline |
Get field values on all nodes.
Definition at line 209 of file ComponentGrid.hh.
|
private |
Decode a format string.
|
inline |
|
overridevirtual |
Get the medium at a given location (x, y, z).
Reimplemented from Garfield::Component.
| 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.
| 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).
|
overridevirtual |
Calculate the voltage range [V].
Implements Garfield::Component.
|
inlineoverridevirtual |
Does the component have maps of the attachment coefficient?
Reimplemented from Garfield::Component.
Definition at line 182 of file ComponentGrid.hh.
|
overridevirtual |
Does the component have a non-zero magnetic field?
Reimplemented from Garfield::Component.
|
inlineoverridevirtual |
Does the component have maps of the low-field mobility?
Reimplemented from Garfield::Component.
Definition at line 189 of file ComponentGrid.hh.
|
inlineoverridevirtual |
Does the component have velocity maps?
Reimplemented from Garfield::Component.
Definition at line 197 of file ComponentGrid.hh.
|
overridevirtual |
Get the hole attachment coefficient.
Reimplemented from Garfield::Component.
|
overridevirtual |
Get the hole Mobility coefficient.
Reimplemented from Garfield::Component.
|
overridevirtual |
Get the hole drift velocity.
Reimplemented from Garfield::Component.
|
private |
Set the dimensions of a table according to the mesh.
|
private |
Read electric field and potential from file.
|
private |
Load scalar data (e. g. attachment coefficients) from file.
| 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
followed by
Format types are:
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).
| 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.
| fname | name of the text file. |
| fmt | format string, see LoadElectricField. |
| col | column in the file which has the attachment coefficient. |
| scaleX | scaling factor to be applied to the coordinates. |
| 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.
| fname | name of the text file. |
| fmt | format string, see LoadElectricField |
| scaleX | scaling factor to be applied to the coordinates. |
| scaleV | scaling factor to be applied to the velocity components. |
| 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.
| 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.
| 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.
|
private |
Read/determine mesh parameters from file.
| 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.
| 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.
|
overridevirtual |
Calculate the magnetic field at a given point.
| x,y,z | coordinates [cm]. |
| bx,by,bz | components of the magnetic field [Tesla]. |
| status | status flag. |
Reimplemented from Garfield::Component.
| void Garfield::ComponentGrid::Print | ( | ) |
Print information about the mesh and the available data.
|
private |
Reduce a coordinate to the basic cell (in case of periodicity).
|
overrideprivatevirtual |
Reset the component.
Implements Garfield::Component.
| bool Garfield::ComponentGrid::SaveElectricField | ( | Component * | cmp | ) |
Export the electric field and potential of a component.
| cmp | Component object for which to export the field/potential |
| 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.
| cmp | Component object for which to export the field/potential |
| filename | name of the text file |
| fmt | format string, see LoadElectricField |
| bool Garfield::ComponentGrid::SaveElectricFieldROOT | ( | Component * | cmp, |
| const std::string & | filename, | ||
| const std::string & | fmt ) |
| 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.
| cmp | Component object for which to export the field/potential |
| id | identifier of the weighting field |
| filename | name of the text file |
| fmt | format string, see LoadElectricField |
|
inline |
| void Garfield::ComponentGrid::SetCylindricalCoordinates | ( | ) |
Use cylindrical coordinates.
| void Garfield::ComponentGrid::SetMedium | ( | Medium * | m | ) |
Set the medium.
| 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.
| nx,ny,nz | number of nodes along |
| xmin,xmax | range along |
| ymin,ymax | range along |
| zmin,zmax | range along |
| 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.
|
overrideprivatevirtual |
Verify periodicities.
Implements Garfield::Component.
|
overridevirtual |
Calculate the weighting field at a given point and for a given electrode.
| x,y,z | coordinates [cm]. |
| wx,wy,wz | components of the weighting field [1/cm]. |
| label | name of the electrode |
Reimplemented from Garfield::Component.
|
overridevirtual |
Calculate the weighting potential at a given point.
| x,y,z | coordinates [cm]. |
| label | name of the electrode. |
Reimplemented from Garfield::Component.
|
private |
Active medium flag.
Definition at line 250 of file ComponentGrid.hh.
|
private |
Magnetic field values.
Definition at line 235 of file ComponentGrid.hh.
|
private |
Definition at line 228 of file ComponentGrid.hh.
|
private |
Attachment maps for electrons and holes.
Definition at line 241 of file ComponentGrid.hh.
|
private |
Electric field values and potentials.
Definition at line 233 of file ComponentGrid.hh.
|
private |
Mobility maps for electrons and holes.
Definition at line 244 of file ComponentGrid.hh.
|
private |
Velocity maps for electrons and holes.
Definition at line 247 of file ComponentGrid.hh.
|
private |
Definition at line 258 of file ComponentGrid.hh.
|
private |
Definition at line 259 of file ComponentGrid.hh.
|
private |
Definition at line 242 of file ComponentGrid.hh.
|
private |
Definition at line 245 of file ComponentGrid.hh.
|
private |
Definition at line 248 of file ComponentGrid.hh.
|
private |
Definition at line 230 of file ComponentGrid.hh.
|
private |
Definition at line 253 of file ComponentGrid.hh.
|
private |
Definition at line 265 of file ComponentGrid.hh.
|
private |
Definition at line 265 of file ComponentGrid.hh.
|
private |
Definition at line 256 of file ComponentGrid.hh.
|
private |
Delayed weighting field values and potentials.
Definition at line 239 of file ComponentGrid.hh.
|
private |
Definition at line 262 of file ComponentGrid.hh.
|
private |
Prompt weighting field values and potentials.
Definition at line 237 of file ComponentGrid.hh.
|
private |
Definition at line 255 of file ComponentGrid.hh.
|
private |
Definition at line 254 of file ComponentGrid.hh.