![]() |
Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
|
Component for interpolating field maps stored in a regular mesh. More...
#include <ComponentVoxel.hh>
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. | |
| Medium * | GetMedium (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? | |
| 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. | |
| 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 | 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 ®ion) |
| 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. | |
| 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 stored in a regular mesh.
Definition at line 10 of file ComponentVoxel.hh.
| Garfield::ComponentVoxel::ComponentVoxel | ( | ) |
Constructor.
|
inline |
|
inlineoverridevirtual |
Reset.
Reimplemented from Garfield::Component.
Definition at line 85 of file ComponentVoxel.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.
|
inline |
Interpolate between field values at the element centres.
Definition at line 18 of file ComponentVoxel.hh.
|
overridevirtual |
Get the bounding box coordinates.
Reimplemented from Garfield::Component.
| bool Garfield::ComponentVoxel::GetElectricFieldRange | ( | double & | exmin, |
| double & | exmax, | ||
| double & | eymin, | ||
| double & | eymax, | ||
| double & | ezmin, | ||
| double & | ezmax ) |
| 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.
| 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.
|
overridevirtual |
Get the coordinates of the elementary cell.
Reimplemented from Garfield::Component.
|
private |
Look up/interpolate the field at a given point.
|
overridevirtual |
Get the medium at a given location (x, y, z).
Reimplemented from Garfield::Component.
| Medium * Garfield::ComponentVoxel::GetMedium | ( | const unsigned int | i | ) | const |
Get the medium in region i.
|
overridevirtual |
Calculate the voltage range [V].
Implements Garfield::Component.
|
overridevirtual |
Does the component have a non-zero magnetic field?
Reimplemented from Garfield::Component.
|
private |
Set the dimensions of a table according to the mesh.
|
private |
|
private |
Read data from file.
| 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
followed by
Format types are:
| 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.
| 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.
| 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.
|
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::ComponentVoxel::PrintRegions | ( | ) | const |
Print all regions.
|
private |
Reduce a coordinate to the basic cell (in case of periodicity).
|
overrideprivatevirtual |
Reset the component.
Implements Garfield::Component.
| void Garfield::ComponentVoxel::SetMedium | ( | const unsigned int | i, |
| Medium * | m ) |
Set the medium in region i.
| 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.
| nx,ny,nz | number of bins along x, y, z. |
| xmin,xmax | range along |
| ymin,ymax | range along |
| zmin,zmax | range along |
| 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.
|
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 |
Magnetic field values at each mesh element.
Definition at line 129 of file ComponentVoxel.hh.
|
private |
Definition at line 139 of file ComponentVoxel.hh.
|
private |
Definition at line 139 of file ComponentVoxel.hh.
|
private |
Definition at line 139 of file ComponentVoxel.hh.
|
private |
Electric field values and potentials at each mesh element.
Definition at line 127 of file ComponentVoxel.hh.
|
private |
Definition at line 146 of file ComponentVoxel.hh.
|
private |
Definition at line 145 of file ComponentVoxel.hh.
|
private |
Definition at line 143 of file ComponentVoxel.hh.
|
private |
Definition at line 144 of file ComponentVoxel.hh.
|
private |
Definition at line 147 of file ComponentVoxel.hh.
|
private |
Definition at line 141 of file ComponentVoxel.hh.
|
private |
Definition at line 118 of file ComponentVoxel.hh.
|
private |
Definition at line 136 of file ComponentVoxel.hh.
|
private |
Definition at line 136 of file ComponentVoxel.hh.
|
private |
Definition at line 136 of file ComponentVoxel.hh.
|
private |
Definition at line 155 of file ComponentVoxel.hh.
|
private |
Definition at line 155 of file ComponentVoxel.hh.
|
private |
Region indices.
Definition at line 125 of file ComponentVoxel.hh.
|
private |
Delayed weighting field values and potentials at each mesh element.
Definition at line 133 of file ComponentVoxel.hh.
|
private |
Definition at line 150 of file ComponentVoxel.hh.
|
private |
Definition at line 151 of file ComponentVoxel.hh.
|
private |
Definition at line 152 of file ComponentVoxel.hh.
|
private |
Prompt weighting field values and potentials at each mesh element.
Definition at line 131 of file ComponentVoxel.hh.
|
private |
Definition at line 138 of file ComponentVoxel.hh.
|
private |
Definition at line 137 of file ComponentVoxel.hh.
|
private |
Definition at line 138 of file ComponentVoxel.hh.
|
private |
Definition at line 137 of file ComponentVoxel.hh.
|
private |
Definition at line 138 of file ComponentVoxel.hh.
|
private |
Definition at line 137 of file ComponentVoxel.hh.