|
Feellgood
|
#include <electrostatSolver.h>
Public Member Functions | |
| electrostatSolver (Mesh::mesh &_msh, const std::vector< Tetra::prm > &_pTetra, const std::vector< Triangle::prm > &_pTri, const double _tol, const bool v, const int max_iter) | |
| void | infos (void) const |
| void | integrales (const Tetra::Tet &tet, Eigen::Ref< Eigen::Matrix< double, Tetra::N, Tetra::N > > AE) const |
| void | integrales (const Triangle::Tri &tri, std::vector< double > &BE) const |
| bool | save (const std::string &V_fileName, const std::string &metadata) const |
| double | getSigma (const Tetra::Tet &tet) const |
| double | getCurrentDensity (const Triangle::Tri &tri) const |
| void | compute (bool verbose, const std::string &V_fileName) |
| void | checkBoundaryConditions (void) const override |
Public Member Functions inherited from solver< DIM_PB_ELEC > | |
| solver (Mesh::mesh &_msh, const std::vector< Tetra::prm > &_pTetra, const std::vector< Triangle::prm > &_pTri, const std::string &name, const double _tol, const bool v, const int max_iter, const std::function< bool(Mesh::Edge)> &edge_filter=[](Mesh::Edge){ return true;}) | |
Public Attributes | |
| std::vector< double > | V |
Private Member Functions | |
| bool | solve (void) |
Private Attributes | |
| const int | precision = 8 |
Additional Inherited Members | |
Protected Member Functions inherited from solver< DIM_PB_ELEC > | |
| algebra::MatrixShape | build_shape (const std::function< bool(Mesh::Edge)> &edge_filter) const |
| void | buildMat (std::array< int, N > &ind, Eigen::Matrix< double, DIM_PROBLEM *N, DIM_PROBLEM *N > &Ke) |
| void | buildVect (std::array< int, N > &ind, std::vector< double > &Le) |
Protected Attributes inherited from solver< DIM_PB_ELEC > | |
| Mesh::mesh * | msh |
| const int | NOD |
| const std::vector< Tetra::prm > & | paramTet |
| const std::vector< Triangle::prm > & | paramTri |
| const bool | verbose |
| algebra::iteration< double > | iter |
| algebra::SparseMatrix | K |
| std::vector< double > | L_rhs |
Static Protected Attributes inherited from solver< DIM_PB_ELEC > | |
| static const int | DIM_PB |
this class is containing both data and a solver to compute potential from dirichlet boundary conditions problem for the current density flowing in the sample.
|
inline |
constructor
| [in] | _msh | reference to the mesh |
| [in] | _pTetra | ref to vector of param tetra (volume region parameters) |
| [in] | _pTri | ref to vector of param triangle (surface region parameters) |
| [in] | _tol | cg_dir tolerance |
| [in] | v | verbose mode for iteration monitor |
| [in] | max_iter | maximum number of iterations |
|
overridevirtual |
check boundary conditions: mesh and settings have to define a single surface with constant current density J (normal componant to the surface) and another single surface with constant potential V
Implements solver< DIM_PB_ELEC >.
| void electrostatSolver::compute | ( | bool | verbose, |
| const std::string & | V_fileName | ||
| ) |
solves the potential and stores result in V, save to text file if needed if verbose set to true, some printing are sent to terminal
| [in] | verbose | |
| [in] | V_fileName | output file name for V solution |
| double electrostatSolver::getCurrentDensity | ( | const Triangle::Tri & | tri | ) | const |
returns current density of the triangle if it is defined in the boundary conditions, else zero
| double electrostatSolver::getSigma | ( | const Tetra::Tet & | tet | ) | const |
returns sigma of the tetraedron, (conductivity in (Ohm.m)^-1
| void electrostatSolver::infos | ( | void | ) | const |
basic informations on boundary conditions
| void electrostatSolver::integrales | ( | const Tetra::Tet & | tet, |
| Eigen::Ref< Eigen::Matrix< double, Tetra::N, Tetra::N > > | AE | ||
| ) | const |
compute side problem (electrostatic potential on the nodes) integrales for matrix coefficients,inputs from tet
| void electrostatSolver::integrales | ( | const Triangle::Tri & | tri, |
| std::vector< double > & | BE | ||
| ) | const |
compute integrales for vector coefficients, input from triangle
| bool electrostatSolver::save | ( | const std::string & | V_fileName, |
| const std::string & | metadata | ||
| ) | const |
text file (tsv) writing function for the solution V over all volume regions of the mesh, node indices are zero based
| [in] | V_fileName | output file name |
| [in] | metadata |
|
private |
solver, using conjugate gradient with masking, with diagonal preconditionner and Dirichlet boundary conditions, returns true if has converged
|
private |
number of digits in the optional output file
| std::vector<double> electrostatSolver::V |
electrostatic potential values for boundary conditions, V.size() is the size of the vector of nodes