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

Component for importing and interpolating field maps from CST. More...

#include <ComponentCST.hh>

Inheritance diagram for Garfield::ComponentCST:
Garfield::ComponentFieldMap

Public Member Functions

 ComponentCST ()
 Constructor.
 ~ComponentCST ()
 Destructor.
void ShiftComponent (const double xShift, const double yShift, const double zShift)
void GetNumberOfMeshLines (unsigned int &nx, unsigned int &ny, unsigned int &nz) const
size_t GetNumberOfElements () const override
bool GetElementNodes (const size_t i, std::vector< size_t > &nodes) const override
bool GetElementRegion (const size_t i, size_t &mat, bool &drift) const override
size_t GetNumberOfNodes () const override
bool GetNode (const size_t i, double &x, double &y, double &z) const override
void GetElementBoundaries (unsigned int element, double &xmin, double &xmax, double &ymin, double &ymax, double &zmin, double &zmax) const
Medium * GetMedium (const double x, const double y, const double z) override
void ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status) 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 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
bool Initialise (std::string elist, std::string nlist, std::string mplist, std::string prnsol, std::string unit="cm")
bool Initialise (std::string dataFile, std::string unit="cm")
 Import of field data based on binary files.
bool SetWeightingField (std::string prnsol, std::string label, bool isBinary=true)
 Initialise a weighting field.
void SetRangeZ (const double zmin, const double zmax)
void DisableXField ()
 Use these functions to disable a certain field component.
void DisableYField ()
void DisableZField ()
void EnableShaping ()
 If you calculate the electric field component in $x$ direction along a line in x direction this field component will be constant inside mesh elements (by construction).
void DisableShaping ()
int Index2Element (const unsigned int i, const unsigned int j, const unsigned int k) const
 Calculate the element index from the position in the x/y/z position vectors (m_xlines, m_ylines, m_zlines).
bool Coordinate2Index (const double x, const double y, const double z, unsigned int &i, unsigned int &j, unsigned int &k) const
 Find the positions in the x/y/z position vectors (m_xlines, m_ylines, m_zlines) for a given point.
Public Member Functions inherited from Garfield::ComponentFieldMap
 ComponentFieldMap ()=delete
 Default constructor.
 ComponentFieldMap (const std::string &name)
 Constructor.
virtual ~ComponentFieldMap ()
 Destructor.
bool Check ()
 Check element aspect ratio.
void PrintRange ()
 Show x, y, z, V and angular ranges.
void PrintMaterials ()
 List all currently defined materials.
void DriftMedium (const size_t imat)
 Flag a field map material as a drift medium.
void NotDriftMedium (const size_t imat)
 Flag a field map materials as a non-drift medium.
size_t GetNumberOfMaterials () const
 Return the number of materials in the field map.
double GetPermittivity (const size_t imat) const
 Return the relative permittivity of a field map material.
double GetConductivity (const size_t imat) const
 Return the conductivity of a field map material.
void SetMedium (const size_t imat, Medium *medium)
 Associate a field map material with a Medium object.
Medium * GetMedium (const size_t imat) const
 Return the Medium associated to a field map material.
void SetGas (Medium *medium)
 Associate all field map materials with a relative permittivity of unity to a given Medium class.
size_t GetNumberOfElements () const override
bool GetElementNodes (const size_t i, std::vector< size_t > &nodes) const override
bool GetElementRegion (const size_t i, size_t &mat, bool &drift) const override
bool GetElement (const size_t i, double &vol, double &dmin, double &dmax) const
 Return the volume and aspect ratio of a mesh element.
size_t GetNumberOfNodes () const override
bool GetNode (const size_t i, double &x, double &y, double &z) const override
double GetPotential (const size_t i) const
 Return the potential at a given node.
void EnableCheckMapIndices (const bool on=true)
void EnableDeleteBackgroundElements (const bool on=true)
 Option to eliminate mesh elements in conductors (default: on).
void EnableTetrahedralTreeForElementSearch (const bool on=true)
 Enable or disable the usage of the tetrahedral tree for searching the element in the mesh.
void EnableConvergenceWarnings (const bool on=true)
 Enable or disable warnings that the calculation of the local coordinates did not achieve the requested precision.
Medium * GetMedium (const double x, const double y, const double z) override
void ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status) 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
__DEVICE__ 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
double DelayedWeightingPotential (double x, double y, double z, const double t, const std::string &label) override
void DelayedWeightingPotentials (const double x, const double y, const double z, const std::string &label, std::vector< double > &dwp) override
bool IsInBoundingBox (const double x, const double y, const double z) const
bool Is3d () override
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
std::map< std::string, std::vector< double > > GetWeightingPotentials ()
bool GetVoltageRange (double &vmin, double &vmax) override
void CopyWeightingPotential (const std::string &label, const std::string &labelSource, const double x, const double y, const double z, const double alpha, const double beta, const double gamma)
 Makes a weighting potential copy of a imported map which can be translated and rotated.
double CreateGPUTransferObject (ComponentGPU *&comp_gpu) override
 Create and initialise GPU Transfer class.

Protected Member Functions

void SetRange () override
double GetElementVolume (const size_t i) const override
void GetAspectRatio (const size_t i, double &dmin, double &dmax) const override
bool Coordinate2Index (const double x, const double y, const double z, unsigned int &i, unsigned int &j, unsigned int &k, double *position_mapped, bool *mirrored) const
 Calculate the index in the vectors m_xlines, m_ylines, m_zlines, which is before the given coordinate.
Protected Member Functions inherited from Garfield::ComponentFieldMap
void Reset () override
void Prepare ()
void UpdatePeriodicity () override
void UpdatePeriodicity2d ()
void UpdatePeriodicityCommon ()
bool SetDefaultDriftMedium ()
 Find lowest epsilon, check for eps = 0, set default drift media flags.
__DEVICE__ int Field (const double x, const double y, const double z, double &fx, double &fy, double &fz, int &iel, const std::vector< double > &potentials) const
 Compute the electric/weighting field.
double Potential (const double x, const double y, const double z, const std::vector< double > &potentials) const
 Compute the electrostatic/weighting potential.
int FindElement5 (const double x, const double y, double &t1, double &t2, double &t3, double &t4, double jac[4][4], double &det) const
 Find the element for a point in curved quadratic quadrilaterals.
__DEVICE__ int FindElement13 (const double x, const double y, const double z, double &t1, double &t2, double &t3, double &t4, double jac[4][4], double &det) const
 Find the element for a point in curved quadratic tetrahedra.
int FindElementCube (const double x, const double y, const double z, double &t1, double &t2, double &t3, TMatrixD *&jac, std::vector< TMatrixD * > &dN) const
 Find the element for a point in a cube.
__DEVICE__ void MapCoordinates (double &xpos, double &ypos, double &zpos, bool &xmirrored, bool &ymirrored, bool &zmirrored, double &rcoordinate, double &rotation) const
 Move (xpos, ypos, zpos) to field map coordinates.
__DEVICE__ void UnmapFields (double &ex, double &ey, double &ez, const double xpos, const double ypos, const double zpos, const bool xmirrored, const bool ymirrored, const bool zmirrored, const double rcoordinate, const double rotation) const
 Move (ex, ey, ez) to global coordinates.
void PrintWarning (const std::string &header)
void PrintNotReady (const std::string &header) const
void PrintCouldNotOpen (const std::string &header, const std::string &filename) const
void PrintElement (const std::string &header, const double x, const double y, const double z, const double t1, const double t2, const double t3, const double t4, const size_t i, const std::vector< double > &potential) const
void TimeInterpolation (const double t, double &f0, double &f1, int &i0, int &i1)
 Interpolation of potential between two time slices.
int Coordinates3 (const double x, const double y, double &t1, double &t2, double &t3, double &t4, double jac[4][4], double &det, const std::array< double, 8 > &xn, const std::array< double, 8 > &yn) const
 Calculate local coordinates for curved quadratic triangles.
int Coordinates4 (const double x, const double y, double &t1, double &t2, double &t3, double &t4, double &det, const std::array< double, 8 > &xn, const std::array< double, 8 > &yn) const
 Calculate local coordinates for linear quadrilaterals.
int Coordinates5 (const double x, const double y, double &t1, double &t2, double &t3, double &t4, double jac[4][4], double &det, const std::array< double, 8 > &xn, const std::array< double, 8 > &yn) const
 Calculate local coordinates for curved quadratic quadrilaterals.
__DEVICE__ void Coordinates12 (const double x, const double y, const double z, double &t1, double &t2, double &t3, double &t4, const std::array< double, 10 > &xn, const std::array< double, 10 > &yn, const std::array< double, 10 > &zn, const std::array< std::array< double, 3 >, 4 > &w) const
 Calculate local coordinates in linear tetrahedra.
__DEVICE__ int Coordinates13 (const double x, const double y, const double z, double &t1, double &t2, double &t3, double &t4, double jac[4][4], double &det, const std::array< double, 10 > &xn, const std::array< double, 10 > &yn, const std::array< double, 10 > &zn, const std::array< std::array< double, 3 >, 4 > &w) const
 Calculate local coordinates for curved quadratic tetrahedra.
int CoordinatesCube (const double x, const double y, const double z, double &t1, double &t2, double &t3, TMatrixD *&jac, std::vector< TMatrixD * > &dN, const Element &element) const
 Calculate local coordinates for a cube.
void JacobianCube (const Element &element, const double t1, const double t2, const double t3, TMatrixD *&jac, std::vector< TMatrixD * > &dN) const
 Calculate Jacobian for a cube.
void CalculateElementBoundingBoxes ()
 Calculate the bounding boxes of all elements after initialization.
bool InitializeTetrahedralTree ()
 Initialize the tetrahedral tree.

Private Member Functions

void ElectricFieldBinary (const double x, const double y, const double z, double &ex, double &ey, double &ez, double &v, Medium *&m, int &status, const bool calculatePotential=false) const
float GetFieldComponent (const unsigned int i, const unsigned int j, const unsigned int k, const double rx, const double ry, const double rz, const char component, const std::vector< float > &potentials) const
float GetPotential (const unsigned int i, const unsigned int j, const unsigned int k, const double rx, const double ry, const double rz, const std::vector< float > &potentials) const
void ShapeField (float &ex, float &ey, float &ez, const double rx, const double ry, const double rz, const unsigned int i, const unsigned int j, const unsigned int k, const std::vector< float > &potentials) const
void Element2Index (const size_t element, unsigned int &i, unsigned int &j, unsigned int &k) const
int Index2Node (const unsigned int i, const unsigned int j, const unsigned int k) const
void Node2Index (const size_t node, unsigned int &i, unsigned int &j, unsigned int &k) const

Private Attributes

std::vector< double > m_xlines
 x positions used in the CST mesh
std::vector< double > m_ylines
 y positions used in the CST mesh
std::vector< double > m_zlines
 z positions used in the CST mesh
std::vector< float > m_potential
 Potentials resulting from the CST simulation.
std::map< std::string, std::vector< float > > m_weightingFields
 Map of weighting field potentials.
std::vector< unsigned char > m_elementMaterial
 Material id for each element (unsigned char since it uses only 1 byte)
unsigned int m_nx = 0
 Number of mesh lines in x direction.
unsigned int m_ny = 0
 Number of mesh lines in y direction.
unsigned int m_nz = 0
 Number of mesh lines in z direction.
size_t m_nElements = 0
 Number of elements.
size_t m_nNodes = 0
 Number of nodes.
bool disableFieldComponent [3] = {false, false, false}
bool doShaping = false

Additional Inherited Members

Protected Types inherited from Garfield::ComponentFieldMap
enum class  ElementType { Unknown = 0 , Serendipity = 5 , CurvedTetrahedron = 13 }
Static Protected Member Functions inherited from Garfield::ComponentFieldMap
static double ScalingFactor (std::string unit)
static double Potential3 (const std::array< double, 6 > &v, const std::array< double, 3 > &t)
 Interpolate the potential in a triangle.
static void Field3 (const std::array< double, 6 > &v, const std::array< double, 3 > &t, double jac[4][4], const double det, double &ex, double &ey)
 Interpolate the field in a triangle.
static double Potential5 (const std::array< double, 8 > &v, const std::array< double, 2 > &t)
 Interpolate the potential in a curved quadrilateral.
static void Field5 (const std::array< double, 8 > &v, const std::array< double, 2 > &t, double jac[4][4], const double det, double &ex, double &ey)
 Interpolate the field in a curved quadrilateral.
static double Potential13 (const std::array< double, 10 > &v, const std::array< double, 4 > &t)
 Interpolate the potential in a curved quadratic tetrahedron.
static __DEVICE__ void Field13 (const std::array< double, 10 > &v, const std::array< double, 4 > &t, double jac[4][4], const double det, double &ex, double &ey, double &ez)
 Interpolate the field in a curved quadratic tetrahedron.
static int ReadInteger (char *token, int def, bool &error)
static double ReadDouble (char *token, double def, bool &error)
static void Jacobian3 (const std::array< double, 8 > &xn, const std::array< double, 8 > &yn, const double u, const double v, const double w, double &det, double jac[4][4])
 Calculate Jacobian for curved quadratic triangles.
static void Jacobian5 (const std::array< double, 8 > &xn, const std::array< double, 8 > &yn, const double u, const double v, double &det, double jac[4][4])
 Calculate Jacobian for curved quadratic quadrilaterals.
static __DEVICE__ void Jacobian13 (const std::array< double, 10 > &xn, const std::array< double, 10 > &yn, const std::array< double, 10 > &zn, const double fourt0, const double fourt1, const double fourt2, const double fourt3, double &det, double jac[4][4])
 Calculate Jacobian for curved quadratic tetrahedra.
static std::array< std::array< double, 3 >, 4 > Weights12 (const std::array< double, 10 > &xn, const std::array< double, 10 > &yn, const std::array< double, 10 > &zn)
Protected Attributes inherited from Garfield::ComponentFieldMap
bool m_is3d = true
ElementType m_elementType = ElementType::CurvedTetrahedron
std::vector< Elementm_elements
std::vector< int > m_elementIndices
std::vector< bool > m_degenerate
std::vector< std::array< double, 3 > > m_bbMin
std::vector< std::array< double, 3 > > m_bbMax
std::vector< std::array< std::array< double, 3 >, 4 > > m_w12
std::vector< Nodem_nodes
std::vector< double > m_pot
std::map< std::string, std::vector< double > > m_wpot
std::map< std::string, std::vector< std::vector< double > > > m_dwpot
std::vector< Materialm_materials
std::map< std::string, WeightingFieldCopym_wfieldCopies
bool m_hasBoundingBox = false
std::array< double, 3 > m_minBoundingBox = {{0., 0., 0.}}
std::array< double, 3 > m_maxBoundingBox = {{0., 0., 0.}}
bool m_cacheElemBoundingBoxes = false
 Flag to check if bounding boxes of elements are cached.
std::array< double, 3 > m_mapmin = {{0., 0., 0.}}
std::array< double, 3 > m_mapmax = {{0., 0., 0.}}
std::array< double, 3 > m_mapamin = {{0., 0., 0.}}
std::array< double, 3 > m_mapamax = {{0., 0., 0.}}
std::array< double, 3 > m_mapna = {{0., 0., 0.}}
std::array< double, 3 > m_cells = {{0., 0., 0.}}
double m_mapvmin = 0.
double m_mapvmax = 0.
std::array< bool, 3 > m_setang
bool m_deleteBackground = true
bool m_warning = false
unsigned int m_nWarnings = 0
bool m_printConvergenceWarnings = true
bool m_checkMultipleElement = false
 Scan for multiple elements that contain a point.
bool m_useTetrahedralTree = true
std::unique_ptr< GARFIELD_CLASS_NAME(TetrahedralTree)> m_octree

Detailed Description

Component for importing and interpolating field maps from CST.

This interface assumes a certain format of the ascii files Please find the tools to extract the field data from CST in the correct way here: http://www.desy.de/~zenker/garfieldpp.html

Definition at line 16 of file ComponentCST.hh.

Constructor & Destructor Documentation

◆ ComponentCST()

Garfield::ComponentCST::ComponentCST ( )

Constructor.

◆ ~ComponentCST()

Garfield::ComponentCST::~ComponentCST ( )
inline

Destructor.

Definition at line 21 of file ComponentCST.hh.

21{}

Member Function Documentation

◆ Coordinate2Index() [1/2]

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

Find the positions in the x/y/z position vectors (m_xlines, m_ylines, m_zlines) for a given point.

The operator used for the comparison is <=. Therefore, the last entry in the vector will never be returned for a point inside the mesh.

◆ Coordinate2Index() [2/2]

bool Garfield::ComponentCST::Coordinate2Index ( const double x,
const double y,
const double z,
unsigned int & i,
unsigned int & j,
unsigned int & k,
double * position_mapped,
bool * mirrored ) const
protected

Calculate the index in the vectors m_xlines, m_ylines, m_zlines, which is before the given coordinate.

Remarks
x, y, z need to be mapped before using this function!
Parameters
xThe x coordinate mapped to the basic cell.
yThe y coordinate mapped to the basic cell.
zThe z coordinate mapped to the basic cell.
iThe index of the m_xlines vector, where m_xlines.at(i) < x < m_xlines.at(i+1).
jThe index of the m_ylines vector, where m_ylines.at(j) < y < m_ylines.at(j+1).
kThe index of the m_zlines vector, where m_zlines.at(k) < z < m_zlines.at(k+1).
position_mappedThe calculated mapped position (x,y,z) -> (x_mapped, y_mapped, z_mapped)
mirroredInformation if x, y, or z direction is mirrored.

◆ DisableShaping()

void Garfield::ComponentCST::DisableShaping ( )
inline

Definition at line 186 of file ComponentCST.hh.

186{ doShaping = false; }

◆ DisableXField()

void Garfield::ComponentCST::DisableXField ( )
inline

Use these functions to disable a certain field component.

Is a field component is disabled ElectricField and WeightingField will return 0 for this component. This is useful if you want to have calculated global field distortions and you want to add the field of a GEM. If you would simply add both components the field component in drift direction would be added twice!

Definition at line 139 of file ComponentCST.hh.

139{ disableFieldComponent[0] = true; }

◆ DisableYField()

void Garfield::ComponentCST::DisableYField ( )
inline

Definition at line 140 of file ComponentCST.hh.

140{ disableFieldComponent[1] = true; }

◆ DisableZField()

void Garfield::ComponentCST::DisableZField ( )
inline

Definition at line 141 of file ComponentCST.hh.

141{ disableFieldComponent[2] = true; }

◆ ElectricField() [1/2]

void Garfield::ComponentCST::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::ComponentCST::ElectricField ( const double x,
const double y,
const double z,
double & ex,
double & ey,
double & ez,
Medium *& m,
int & status )
override

◆ ElectricFieldBinary()

void Garfield::ComponentCST::ElectricFieldBinary ( const double x,
const double y,
const double z,
double & ex,
double & ey,
double & ez,
double & v,
Medium *& m,
int & status,
const bool calculatePotential = false ) const
private

◆ Element2Index()

void Garfield::ComponentCST::Element2Index ( const size_t element,
unsigned int & i,
unsigned int & j,
unsigned int & k ) const
private

◆ EnableShaping()

void Garfield::ComponentCST::EnableShaping ( )
inline

If you calculate the electric field component in $x$ direction along a line in x direction this field component will be constant inside mesh elements (by construction).

This can be observed by plotting $E_x$ in $x$ direction. If you plot $E_x$ in y direction the field will be smooth (also by construction). Yuri Piadyk proposed also to shape the electric field. This is done as follows. The field component calculated as described above is assumed to appear in the center of the mesh element.


| M1 P | M2 | x direction |-—x—x-|--—x---—|-->

__________ ____________

element 1 element 2

Lets consider only the $x$ direction and we want to calculate $E_x(P)$. The field in the center of the element containing $P$ is $E_x(M_1) = E_1$. Without shaping it is $E_1$ along the $x$ direction in everywhere in element 1. The idea of the shaping is to do a linear interpolation of the $E_x$ between the field $E_1$ and $E_x(M_2)=E_2$. This results in a smooth electric field $E_x$ also in $x$ direction. If P would be left from $M_1$ the field in the left neighboring element would be considered. In addition it is also checked if the material in both elements used for the interpolation is the same. Else no interpolation is done.

Remarks
This shaping gives you a nice and smooth field, but you introduce additional information. This information is not coming from the CST simulation, but from the assumption that the field between elements changes in a linear way, which might be wrong! So you might consider to increase the number of mesh cells used in the simulation rather than using this smoothing.

Definition at line 185 of file ComponentCST.hh.

185{ doShaping = true; }

◆ GetAspectRatio()

void Garfield::ComponentCST::GetAspectRatio ( const size_t i,
double & dmin,
double & dmax ) const
overrideprotectedvirtual

Reimplemented from Garfield::ComponentFieldMap.

◆ GetElementBoundaries()

void Garfield::ComponentCST::GetElementBoundaries ( unsigned int element,
double & xmin,
double & xmax,
double & ymin,
double & ymax,
double & zmin,
double & zmax ) const

◆ GetElementNodes()

bool Garfield::ComponentCST::GetElementNodes ( const size_t i,
std::vector< size_t > & nodes ) const
override

◆ GetElementRegion()

bool Garfield::ComponentCST::GetElementRegion ( const size_t i,
size_t & mat,
bool & drift ) const
override

◆ GetElementVolume()

double Garfield::ComponentCST::GetElementVolume ( const size_t i) const
overrideprotectedvirtual

Reimplemented from Garfield::ComponentFieldMap.

◆ GetFieldComponent()

float Garfield::ComponentCST::GetFieldComponent ( const unsigned int i,
const unsigned int j,
const unsigned int k,
const double rx,
const double ry,
const double rz,
const char component,
const std::vector< float > & potentials ) const
private

◆ GetMedium()

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

◆ GetNode()

bool Garfield::ComponentCST::GetNode ( const size_t i,
double & x,
double & y,
double & z ) const
override

◆ GetNumberOfElements()

size_t Garfield::ComponentCST::GetNumberOfElements ( ) const
inlineoverride

Definition at line 28 of file ComponentCST.hh.

28{ return m_nElements; }
size_t m_nElements
Number of elements.

◆ GetNumberOfMeshLines()

void Garfield::ComponentCST::GetNumberOfMeshLines ( unsigned int & nx,
unsigned int & ny,
unsigned int & nz ) const

◆ GetNumberOfNodes()

size_t Garfield::ComponentCST::GetNumberOfNodes ( ) const
inlineoverride

Definition at line 33 of file ComponentCST.hh.

33{ return m_nNodes; }
size_t m_nNodes
Number of nodes.

◆ GetPotential()

float Garfield::ComponentCST::GetPotential ( const unsigned int i,
const unsigned int j,
const unsigned int k,
const double rx,
const double ry,
const double rz,
const std::vector< float > & potentials ) const
private

◆ Index2Element()

int Garfield::ComponentCST::Index2Element ( const unsigned int i,
const unsigned int j,
const unsigned int k ) const

Calculate the element index from the position in the x/y/z position vectors (m_xlines, m_ylines, m_zlines).

This is public since it is used in ViewFEMesh::DrawCST.

◆ Index2Node()

int Garfield::ComponentCST::Index2Node ( const unsigned int i,
const unsigned int j,
const unsigned int k ) const
private

◆ Initialise() [1/2]

bool Garfield::ComponentCST::Initialise ( std::string dataFile,
std::string unit = "cm" )

Import of field data based on binary files.

See http://www.desy.de/~zenker/garfieldpp.html to get information about the binary files export from CST.

Parameters
dataFileThe binary file containing the field data exported from CST.
unitThe units used in the binary file. They are not necessarily equal to CST units.

◆ Initialise() [2/2]

bool Garfield::ComponentCST::Initialise ( std::string elist,
std::string nlist,
std::string mplist,
std::string prnsol,
std::string unit = "cm" )

◆ Node2Index()

void Garfield::ComponentCST::Node2Index ( const size_t node,
unsigned int & i,
unsigned int & j,
unsigned int & k ) const
private

◆ SetRange()

void Garfield::ComponentCST::SetRange ( )
overrideprotectedvirtual

Reimplemented from Garfield::ComponentFieldMap.

◆ SetRangeZ()

void Garfield::ComponentCST::SetRangeZ ( const double zmin,
const double zmax )

◆ SetWeightingField()

bool Garfield::ComponentCST::SetWeightingField ( std::string prnsol,
std::string label,
bool isBinary = true )

Initialise a weighting field.

This function can handle the deprecated text based file format (see Initialise( std::string elist...) for the expected file format, which is similar to prnsol. It also also handles binary files including the weighting field.

Parameters
prnsolThe input file (binary/text file)
labelThe name of the weighting field to be added. If a weighting field with same name already exist it is replaced by the new one.
isBinaryDepending on the file type you use, adapt this switch.

◆ ShapeField()

void Garfield::ComponentCST::ShapeField ( float & ex,
float & ey,
float & ez,
const double rx,
const double ry,
const double rz,
const unsigned int i,
const unsigned int j,
const unsigned int k,
const std::vector< float > & potentials ) const
private

◆ ShiftComponent()

void Garfield::ComponentCST::ShiftComponent ( const double xShift,
const double yShift,
const double zShift )

◆ WeightingField()

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

◆ WeightingPotential()

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

Member Data Documentation

◆ disableFieldComponent

bool Garfield::ComponentCST::disableFieldComponent[3] = {false, false, false}
private

Definition at line 250 of file ComponentCST.hh.

250{false, false, false};

◆ doShaping

bool Garfield::ComponentCST::doShaping = false
private

Definition at line 251 of file ComponentCST.hh.

◆ m_elementMaterial

std::vector<unsigned char> Garfield::ComponentCST::m_elementMaterial
private

Material id for each element (unsigned char since it uses only 1 byte)

Definition at line 242 of file ComponentCST.hh.

◆ m_nElements

size_t Garfield::ComponentCST::m_nElements = 0
private

Number of elements.

Definition at line 247 of file ComponentCST.hh.

◆ m_nNodes

size_t Garfield::ComponentCST::m_nNodes = 0
private

Number of nodes.

Definition at line 248 of file ComponentCST.hh.

◆ m_nx

unsigned int Garfield::ComponentCST::m_nx = 0
private

Number of mesh lines in x direction.

Definition at line 244 of file ComponentCST.hh.

◆ m_ny

unsigned int Garfield::ComponentCST::m_ny = 0
private

Number of mesh lines in y direction.

Definition at line 245 of file ComponentCST.hh.

◆ m_nz

unsigned int Garfield::ComponentCST::m_nz = 0
private

Number of mesh lines in z direction.

Definition at line 246 of file ComponentCST.hh.

◆ m_potential

std::vector<float> Garfield::ComponentCST::m_potential
private

Potentials resulting from the CST simulation.

Definition at line 238 of file ComponentCST.hh.

◆ m_weightingFields

std::map<std::string, std::vector<float> > Garfield::ComponentCST::m_weightingFields
private

Map of weighting field potentials.

Definition at line 240 of file ComponentCST.hh.

◆ m_xlines

std::vector<double> Garfield::ComponentCST::m_xlines
private

x positions used in the CST mesh

Definition at line 234 of file ComponentCST.hh.

◆ m_ylines

std::vector<double> Garfield::ComponentCST::m_ylines
private

y positions used in the CST mesh

Definition at line 235 of file ComponentCST.hh.

◆ m_zlines

std::vector<double> Garfield::ComponentCST::m_zlines
private

z positions used in the CST mesh

Definition at line 236 of file ComponentCST.hh.


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