Feellgood
spinAccumulationSolver.h
1 #ifndef spinAccumulationSolver_h
2 #define spinAccumulationSolver_h
3 
4 #include <vector>
5 #include "config.h"
6 #include "node.h"
7 #include "tetra.h"
8 #include "electrostatSolver.h"
9 #include "solver.h"
10 #include "meshUtils.h"
11 
13 const int DIM_PB_SPIN_ACC = 3;
14 
22 class spinAcc : public solver<DIM_PB_SPIN_ACC>
23  {
24  public:
27  std::vector<Tetra::prm> & _pTetra ,
29  std::vector<Facette::prm> & _pFac ,
31  const double _tol , // _tol could be 1e-6
32  const bool v ,
33  const int max_iter ):
34  solver<DIM_PB_SPIN_ACC>(_msh,_pTetra,_pFac,"bicg_dir",_tol,v,max_iter)
35  {
36  valDirichlet.resize(DIM_PB*NOD);
38  }
39 
44  void boundaryConditions(void); // should be private
45 
47  void preCompute(std::vector<double> &V);
48 
50  bool compute(void);
51 
53  std::vector<Eigen::Vector3d> s;
54 
58  void checkBoundaryConditions(void) const;
59 
60  private:
63  std::vector<double> valDirichlet;
64 
67  std::vector<int> idxDirichlet;
68 
71  void fillDirichletData(const int k, Eigen::Vector3d &s_value);
72 
76  std::vector< Eigen::Matrix<double,Nodes::DIM,Tetra::NPI> > gradV;
77 
79  const int precision = 8;
80 
82  double getMs(Tetra::Tet const &tet) const;
83 
85  double getSigma(Tetra::Tet const &tet) const;
86 
88  double getPolarizationRate(Tetra::Tet const &tet) const;
89 
91  double getDiffusionCst(Tetra::Tet const &tet) const;
92 
94  double getLsd(Tetra::Tet const &tet) const;
95 
97  double getLsf(Tetra::Tet const &tet) const;
98 
100  double getSpinHall(Tetra::Tet const &tet) const;
101 
103  void prepareExtras(void);
104 
107  bool solve(void);
108 
111  void integrales(Tetra::Tet &tet, Eigen::Matrix<double,DIM_PB*Tetra::N,DIM_PB*Tetra::N> &AE);
112 
115  void integrales(Tetra::Tet &tet,std::vector<double> &BE);
116  };
117 
118 #endif
Definition: mesh.h:32
Definition: tetra.h:141
template class for the different solvers template parameter DIM_PROBLEM: dimensionnality of the probl...
Definition: solver.h:25
const int NOD
Definition: solver.h:54
static const int DIM_PB
Definition: solver.h:48
Definition: spinAccumulationSolver.h:23
void prepareExtras(void)
Definition: spinAccumulationSolver.cpp:106
void boundaryConditions(void)
Definition: spinAccumulationSolver.cpp:56
spinAcc(Mesh::mesh &_msh, std::vector< Tetra::prm > &_pTetra, std::vector< Facette::prm > &_pFac, const double _tol, const bool v, const int max_iter)
Definition: spinAccumulationSolver.h:26
bool compute(void)
Definition: spinAccumulationSolver.cpp:137
void preCompute(std::vector< double > &V)
Definition: spinAccumulationSolver.cpp:126
void checkBoundaryConditions(void) const
Definition: spinAccumulationSolver.cpp:9
bool solve(void)
Definition: spinAccumulationSolver.cpp:151
std::vector< int > idxDirichlet
Definition: spinAccumulationSolver.h:67
double getSigma(Tetra::Tet const &tet) const
Definition: spinAccumulationSolver.cpp:85
const int precision
Definition: spinAccumulationSolver.h:79
void integrales(Tetra::Tet &tet, Eigen::Matrix< double, DIM_PB *Tetra::N, DIM_PB *Tetra::N > &AE)
Definition: spinAccumulationSolver.cpp:202
std::vector< double > valDirichlet
Definition: spinAccumulationSolver.h:63
double getDiffusionCst(Tetra::Tet const &tet) const
Definition: spinAccumulationSolver.cpp:88
double getMs(Tetra::Tet const &tet) const
Definition: spinAccumulationSolver.cpp:82
std::vector< Eigen::Vector3d > s
Definition: spinAccumulationSolver.h:53
double getPolarizationRate(Tetra::Tet const &tet) const
Definition: spinAccumulationSolver.cpp:94
std::vector< Eigen::Matrix< double, Nodes::DIM, Tetra::NPI > > gradV
Definition: spinAccumulationSolver.h:76
void fillDirichletData(const int k, Eigen::Vector3d &s_value)
Definition: spinAccumulationSolver.cpp:46
double getLsd(Tetra::Tet const &tet) const
Definition: spinAccumulationSolver.cpp:97
double getLsf(Tetra::Tet const &tet) const
Definition: spinAccumulationSolver.cpp:100
double getSpinHall(Tetra::Tet const &tet) const
Definition: spinAccumulationSolver.cpp:103
solver for electrostatic problem when spin accumulation is required header containing electrostatSolv...
constexpr double v[NPI]
Definition: facette.h:49
header to define struct Node
two templates to fill matrix and vectors in various dimensionnality situations. DIM_PROBLEM = 1 is us...
namespace Tetra header containing Tet class, some constants, and integrales