Feellgood
tetra.h
Go to the documentation of this file.
1 #ifndef tetra_h
2 #define tetra_h
3 
9 #include <set>
10 
11 #include "facette.h"
12 #include "node.h"
13 #include "time_integration.h"
14 #include "element.h"
15 
19 namespace Tetra
20  {
21 const int N = 4;
22 const int NPI = 5;
24 const double epsilon =
25  EPSILON;
28 constexpr double A = 1. / 4.;
29 constexpr double B = 1. / 6.;
30 constexpr double C = 1. / 2.;
31 constexpr double D = -2. / 15.;
32 constexpr double E = 3. / 40.;
33 constexpr double u[NPI] = {A, B, B, B, C};
34 constexpr double v[NPI] = {A, B, B, C, B};
35 constexpr double w[NPI] = {A, B, C, B, B};
36 constexpr double pds[NPI] = {D, E, E, E, E};
44 constexpr double a[N][NPI] = {{1. - u[0] - v[0] - w[0], 1. - u[1] - v[1] - w[1],
45  1. - u[2] - v[2] - w[2], 1. - u[3] - v[3] - w[3],
46  1. - u[4] - v[4] - w[4]},
47  {u[0], u[1], u[2], u[3], u[4]},
48  {v[0], v[1], v[2], v[3], v[4]},
49  {w[0], w[1], w[2], w[3], w[4]}};
50 
52 const Eigen::Matrix<double,N,NPI> eigen_a = (Eigen::MatrixXd(N,NPI) << a[0][0], a[0][1], a[0][2], a[0][3], a[0][4],
53  a[1][0], a[1][1], a[1][2], a[1][3], a[1][4],
54  a[2][0], a[2][1], a[2][2], a[2][3], a[2][4],
55  a[3][0], a[3][1], a[3][2], a[3][3], a[3][4] ).finished();
56 
60 struct prm
61  {
62  std::string regName;
63  double alpha_LLG;
64  double A;
65  double J = 0;
66  double K;
67  Eigen::Vector3d uk;
69  double K3;
70  Eigen::Vector3d ex;
71  Eigen::Vector3d ey;
72  Eigen::Vector3d ez;
74  double P;
75  double N0;
76  double sigma;
77  double lsd;
78  double lsf;
79  double spinHall;
81  double volume = 0;
83  void infos(void);
84  };
85 
110 class Tet : public element<N,NPI>
111  {
112 public:
117  inline Tet(const std::vector<Nodes::Node> &_p_node ,
118  const int _idx ,
119  std::initializer_list<int> _i )
120  : element<N,NPI>(_p_node,_idx,_i), idx(0)
121  {
122  zeroBasing();
123  da.setZero();
124 
125  if (existNodes())
126  {
127  orientate();
128 
129  Eigen::Matrix3d J;
130  // we have to rebuild the jacobian in case of ill oriented tetrahedron
131  double detJ = Jacobian(J);
132 
133  if (fabs(detJ) < Tetra::epsilon)
134  {
135  std::cerr << "Singular jacobian in tetrahedron: |det(J)|= " << fabs(detJ) << std::endl;
136  element::infos();
137  SYSTEM_ERROR;
138  }
139  Eigen::Matrix<double,N,Nodes::DIM> dadu;
140  dadu << -1., -1., -1., 1., 0., 0., 0., 1., 0., 0., 0., 1.;
141  da = dadu * J.inverse();
142 
143  for (int j = 0; j < NPI; j++)
144  { weight[j] = detJ * Tetra::pds[j]; }
145  }
146  // do nothing lambda's (usefull for spin transfer torque)
147  extraField = [] ( Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,Tetra::NPI>> ) {};
148  extraCoeffs_BE = [](Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>>,
149  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>>,
150  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>>,
151  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>>,
152  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,N>>) {};
153  }
154 
156  Eigen::Matrix<double,N,Nodes::DIM> da;
157 
160  inline void interpolation(std::function<double(Nodes::Node)> getter,
161  Eigen::Ref<Eigen::Matrix<double,NPI,1>> result) const
162  {
163  Eigen::Matrix<double,N,1> scalar_nod;
164  for (int i = 0; i < N; i++) scalar_nod(i) = getter(getNode(i));
165  result = scalar_nod.transpose() * eigen_a;
166  }
167 
170  inline void interpolation(std::function<Eigen::Vector3d(Nodes::Node)> getter,
171  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> result,
172  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> Tx,
173  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> Ty,
174  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> Tz) const
175  {
176  Eigen::Matrix<double,Nodes::DIM,N> vec_nod;
177  for (int i = 0; i < N; i++) vec_nod.col(i) = getter(getNode(i));
178 
179  result = vec_nod * eigen_a;
180  Tx = vec_nod * (da.col(Nodes::IDX_X)).replicate(1,NPI);
181  Ty = vec_nod * (da.col(Nodes::IDX_Y)).replicate(1,NPI);
182  Tz = vec_nod * (da.col(Nodes::IDX_Z)).replicate(1,NPI);
183  }
184 
187  inline void interpolation_field(std::function<double(Nodes::Node)> getter,
188  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> X) const
189  {
190  Eigen::Matrix<double,N,1> scalar_nod;
191  for (int i = 0; i < N; i++) scalar_nod(i) = getter(getNode(i));
192 
193  X.setZero();
194  for (int j = 0; j < NPI; j++)
195  {
196  for (int i = 0; i < N; i++)
197  {
198  X.col(j) -= (scalar_nod[i] * da.row(i));
199  }
200  }
201  }
202 
205  inline void interpolation(std::function<double(Nodes::Node, Nodes::index)> getter, Nodes::index idx,
206  Eigen::Ref<Eigen::Matrix<double,Tetra::NPI,1>> result) const
207  {
208  Eigen::Matrix<double,N,1> scalar_nod;
209 
210  for (int i = 0; i < N; i++)
211  {
212  scalar_nod[i] = getter(getNode(i),idx);
213  }
214  result = scalar_nod.transpose() * eigen_a;
215  }
216 
226  void lumping(Eigen::Ref<Eigen::Matrix<double,NPI,1>> alpha_eff, double prefactor,
227  Eigen::Ref<Eigen::Matrix<double,3*N,3*N>> AE ) const;
228 
230  void add_drift_BE(double alpha, double s_dt, double Vdrift,
231  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> U,
232  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> V,
233  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dUd_,
234  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dVd_,
235  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,N>> BE) const;
236 
239  Eigen::Matrix<double,NPI,1> calc_aniso_uniax(Eigen::Ref<const Eigen::Vector3d> uk,
240  const double Kbis, const double s_dt,
241  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> U,
242  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> V,
243  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> H_aniso) const;
244 
247  Eigen::Matrix<double,NPI,1> calc_aniso_cub(Eigen::Ref<const Eigen::Vector3d> ex,
248  Eigen::Ref<const Eigen::Vector3d> ey,
249  Eigen::Ref<const Eigen::Vector3d> ez,
250  const double K3bis, const double s_dt,
251  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> U,
252  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> V,
253  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> H_aniso) const;
254 
258  void integrales( Tetra::prm const &param, timing const &prm_t,
259  std::function< Eigen::Matrix<double,Nodes::DIM,NPI> (void)> calc_Hext,
260  Nodes::index idx_dir, double Vdrift);
261 
263  double exchangeEnergy(Tetra::prm const &param,
264  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudx,
265  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudy,
266  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudz) const;
267 
269  double anisotropyEnergy(Tetra::prm const &param, Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> u) const;
270 
272  void charges(Tetra::prm const &param,
273  std::function<Eigen::Vector3d(Nodes::Node)> getter,
274  std::vector<double> &srcDen,
275  int &nsrc) const;
276 
278  double demagEnergy(Tetra::prm const &param,
279  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudx,
280  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudy,
281  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudz,
282  Eigen::Ref<Eigen::Matrix<double,NPI,1>> phi) const;
283 
285  double zeemanEnergy(Tetra::prm const &param, Eigen::Ref<Eigen::Vector3d> const Hext,
286  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> const u) const;
287 
289  double Jacobian(Eigen::Ref<Eigen::Matrix3d> J);
290 
292  double calc_vol(void) const;
293 
295  int idx;
296 
298  std::function<void( Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> H)> extraField;
299 
301  std::function<void(Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> U,
302  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dUdx,
303  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dUdy,
304  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dUdz,
305  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,N>> BE)> extraCoeffs_BE;
306 
308  void getPtGauss(Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> result) const
309  {
310  Eigen::Matrix<double,Nodes::DIM,N> vec_nod;
311  for (int i = 0; i < N; i++)
312  {
313  vec_nod.col(i) << getNode(i).p;
314  }
315  result = vec_nod * eigen_a;
316  }
317 
318 private:
320  void orientate(void)
321  {
322  if (calc_vol() < 0.0)
323  std::swap(ind[2], ind[3]);
324  }
325 
326  }; // end class Tetra
327 
330  Eigen::Matrix<double,NPI,1> calc_alpha_eff(const double dt, const double alpha,
331  Eigen::Ref<Eigen::Matrix<double,NPI,1>> uHeff);
332  } // end namespace Tetra
333 
334 #endif /* tetra_h */
Definition: tetra.h:111
void getPtGauss(Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> result) const
Definition: tetra.h:308
double calc_vol(void) const
Definition: tetra.cpp:335
void charges(Tetra::prm const &param, std::function< Eigen::Vector3d(Nodes::Node)> getter, std::vector< double > &srcDen, int &nsrc) const
Definition: tetra.cpp:278
void lumping(Eigen::Ref< Eigen::Matrix< double, NPI, 1 >> alpha_eff, double prefactor, Eigen::Ref< Eigen::Matrix< double, 3 *N, 3 *N >> AE) const
Definition: tetra.cpp:78
Eigen::Matrix< double, NPI, 1 > calc_aniso_uniax(Eigen::Ref< const Eigen::Vector3d > uk, const double Kbis, const double s_dt, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> U, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> V, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> H_aniso) const
Definition: tetra.cpp:126
Eigen::Matrix< double, NPI, 1 > calc_aniso_cub(Eigen::Ref< const Eigen::Vector3d > ex, Eigen::Ref< const Eigen::Vector3d > ey, Eigen::Ref< const Eigen::Vector3d > ez, const double K3bis, const double s_dt, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> U, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> V, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> H_aniso) const
Definition: tetra.cpp:137
Eigen::Matrix< double, N, Nodes::DIM > da
Definition: tetra.h:156
Tet(const std::vector< Nodes::Node > &_p_node, const int _idx, std::initializer_list< int > _i)
Definition: tetra.h:117
double zeemanEnergy(Tetra::prm const &param, Eigen::Ref< Eigen::Vector3d > const Hext, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> const u) const
Definition: tetra.cpp:310
double exchangeEnergy(Tetra::prm const &param, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> dudx, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> dudy, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> dudz) const
Definition: tetra.cpp:245
double Jacobian(Eigen::Ref< Eigen::Matrix3d > J)
Definition: tetra.cpp:318
void add_drift_BE(double alpha, double s_dt, double Vdrift, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> U, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> V, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> dUd_, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> dVd_, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, N >> BE) const
Definition: tetra.cpp:106
int idx
Definition: tetra.h:295
std::function< void(Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> U, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> dUdx, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> dUdy, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> dUdz, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, N >> BE)> extraCoeffs_BE
Definition: tetra.h:305
std::function< void(Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> H)> extraField
Definition: tetra.h:298
void interpolation(std::function< Eigen::Vector3d(Nodes::Node)> getter, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> result, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> Tx, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> Ty, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> Tz) const
Definition: tetra.h:170
void interpolation(std::function< double(Nodes::Node, Nodes::index)> getter, Nodes::index idx, Eigen::Ref< Eigen::Matrix< double, Tetra::NPI, 1 >> result) const
Definition: tetra.h:205
void interpolation(std::function< double(Nodes::Node)> getter, Eigen::Ref< Eigen::Matrix< double, NPI, 1 >> result) const
Definition: tetra.h:160
void interpolation_field(std::function< double(Nodes::Node)> getter, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> X) const
Definition: tetra.h:187
double demagEnergy(Tetra::prm const &param, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> dudx, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> dudy, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> dudz, Eigen::Ref< Eigen::Matrix< double, NPI, 1 >> phi) const
Definition: tetra.cpp:297
void orientate(void)
Definition: tetra.h:320
double anisotropyEnergy(Tetra::prm const &param, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> u) const
Definition: tetra.cpp:256
void integrales(Tetra::prm const &param, timing const &prm_t, std::function< Eigen::Matrix< double, Nodes::DIM, NPI >(void)> calc_Hext, Nodes::index idx_dir, double Vdrift)
Definition: tetra.cpp:161
Template abstract class, mother class for tetraedrons and facettes.
Definition: element.h:25
void zeroBasing(void)
Definition: element.h:140
bool existNodes(void)
Definition: element.h:136
void infos() const
Definition: element.h:120
Eigen::Matrix< double, NPI, 1 > weight
Definition: element.h:50
std::vector< int > ind
Definition: element.h:44
const Nodes::Node & getNode(const int i) const
Definition: element.h:133
Definition: time_integration.h:6
contains namespace Facette header containing Fac class, and some constants and a less_than operator t...
index
Definition: node.h:34
const double epsilon
Definition: tetra.h:24
constexpr double A
Definition: tetra.h:28
constexpr double a[N][NPI]
Definition: tetra.h:44
constexpr double B
Definition: tetra.h:29
const int NPI
Definition: tetra.h:22
Eigen::Matrix< double, NPI, 1 > calc_alpha_eff(const double dt, const double alpha, Eigen::Ref< Eigen::Matrix< double, NPI, 1 >> uHeff)
Definition: tetra.cpp:48
const int N
Definition: tetra.h:21
constexpr double w[NPI]
Definition: tetra.h:35
constexpr double D
Definition: tetra.h:31
constexpr double C
Definition: tetra.h:30
constexpr double pds[NPI]
Definition: tetra.h:36
constexpr double v[NPI]
Definition: tetra.h:34
const Eigen::Matrix< double, N, NPI > eigen_a
Definition: tetra.h:52
constexpr double E
Definition: tetra.h:32
constexpr double u[NPI]
Definition: tetra.h:33
header to define struct Node
Definition: node.h:61
Eigen::Vector3d p
Definition: node.h:62
Definition: tetra.h:61
Eigen::Vector3d ex
Definition: tetra.h:70
double P
Definition: tetra.h:74
double volume
Definition: tetra.h:81
Eigen::Vector3d uk
Definition: tetra.h:67
std::string regName
Definition: tetra.h:62
double A
Definition: tetra.h:64
double alpha_LLG
Definition: tetra.h:63
double N0
Definition: tetra.h:75
double spinHall
Definition: tetra.h:79
double lsd
Definition: tetra.h:77
Eigen::Vector3d ez
Definition: tetra.h:72
void infos(void)
Definition: tetra.cpp:15
Eigen::Vector3d ey
Definition: tetra.h:71
double J
Definition: tetra.h:65
double lsf
Definition: tetra.h:78
double K
Definition: tetra.h:66
double K3
Definition: tetra.h:69
double sigma
Definition: tetra.h:76