Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
ComponentElmer.hh
Go to the documentation of this file.
1// Copied and modified ComponentAnsys123.hh
2
3#ifndef G_COMPONENT_ELMER_H
4#define G_COMPONENT_ELMER_H
5
7
8namespace Garfield {
9
10/// Component for importing field maps computed by Elmer.
11
13 public:
14 /// Default constructor
16 /// Constructor with a set of field map files, see Initialise().
17 ComponentElmer(const std::string& header, const std::string& elist,
18 const std::string& nlist, const std::string& mplist,
19 const std::string& volt, const std::string& unit);
20 /// Destructor
22
23 void ElectricField(const double x, const double y, const double z, double& ex,
24 double& ey, double& ez, Medium*& m, int& status) override;
25 void ElectricField(const double x, const double y, const double z, double& ex,
26 double& ey, double& ez, double& v, Medium*& m,
27 int& status) override;
28
29 void WeightingField(const double x, const double y, const double z,
30 double& wx, double& wy, double& wz,
31 const std::string& label) override;
32 double WeightingPotential(const double x, const double y, const double z,
33 const std::string& label) override;
34
35 Medium* GetMedium(const double x, const double y, const double z) override;
36
37 /** Import a field map from a set of files.
38 * \param header name of the header file
39 (contains the number of elements and nodes).
40 * \param elist name of the file that contains the list of mesh elements
41 * \param nlist name of the file that contains the list of mesh nodes
42 * \param mplist name of the file that contains the material properties
43 * \param volt output of the field solver (list of voltages)
44 * \param unit length unit to be used
45 */
46 bool Initialise(const std::string& header = "mesh.header",
47 const std::string& elist = "mesh.elements",
48 const std::string& nlist = "mesh.nodes",
49 const std::string& mplist = "dielectrics.dat",
50 const std::string& volt = "out.result",
51 const std::string& unit = "cm");
52 /// Import a list of voltages to be used as weighting field.
53 bool SetWeightingField(std::string prnsol, std::string label);
54
55 protected:
56 // Verify periodicities
58
59 double GetElementVolume(const unsigned int i) override;
60 void GetAspectRatio(const unsigned int i, double& dmin,
61 double& dmax) override;
62};
63}
64#endif
Component for importing field maps computed by Elmer.
void WeightingField(const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label) override
bool Initialise(const std::string &header="mesh.header", const std::string &elist="mesh.elements", const std::string &nlist="mesh.nodes", const std::string &mplist="dielectrics.dat", const std::string &volt="out.result", const std::string &unit="cm")
void UpdatePeriodicity() override
Verify periodicities.
bool SetWeightingField(std::string prnsol, std::string label)
Import a list of voltages to be used as weighting field.
ComponentElmer()
Default constructor.
void ElectricField(const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status) override
double GetElementVolume(const unsigned int i) override
void GetAspectRatio(const unsigned int i, double &dmin, double &dmax) override
double WeightingPotential(const double x, const double y, const double z, const std::string &label) override
Medium * GetMedium(const double x, const double y, const double z) override
Get the medium at a given location (x, y, z).
Base class for components based on finite-element field maps.
Abstract base class for media.
Definition: Medium.hh:13