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  {
23 const double epsilon = EPSILON;
24 
25 const int N = 4;
27 #if ONE_GAUSS_POINT
28  const int NPI = 1;
31  constexpr double A = 1. / 4.;
32  constexpr double u[NPI] = {A};
33  constexpr double v[NPI] = {A};
34  constexpr double w[NPI] = {A};
35  constexpr double pds[NPI] = {1./6.};
38  constexpr double a[N][NPI] = {{1. - u[0] - v[0] - w[0]},
39  {u[0]},
40  {v[0]},
41  {w[0]}};
42 
44  const Eigen::Matrix<double,N,NPI> eigen_a = (Eigen::MatrixXd(N,NPI) << a[0][0],
45  a[1][0],
46  a[2][0],
47  a[3][0] ).finished();
48 #else
49  const int NPI = 5;
51  constexpr double A = 1. / 4.;
52  constexpr double B = 1. / 6.;
53  constexpr double C = 1. / 2.;
54  constexpr double D = -2. / 15.;
55  constexpr double E = 3. / 40.;
56  constexpr double u[NPI] = {A, B, B, B, C};
57  constexpr double v[NPI] = {A, B, B, C, B};
58  constexpr double w[NPI] = {A, B, C, B, B};
59  constexpr double pds[NPI] = {D, E, E, E, E};
67  constexpr double a[N][NPI] = {{1. - u[0] - v[0] - w[0], 1. - u[1] - v[1] - w[1],
68  1. - u[2] - v[2] - w[2], 1. - u[3] - v[3] - w[3],
69  1. - u[4] - v[4] - w[4]},
70  {u[0], u[1], u[2], u[3], u[4]},
71  {v[0], v[1], v[2], v[3], v[4]},
72  {w[0], w[1], w[2], w[3], w[4]}};
73 
75  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],
76  a[1][0], a[1][1], a[1][2], a[1][3], a[1][4],
77  a[2][0], a[2][1], a[2][2], a[2][3], a[2][4],
78  a[3][0], a[3][1], a[3][2], a[3][3], a[3][4] ).finished();
79 #endif
80 
84 struct prm
85  {
86  std::string regName;
87  double alpha_LLG;
88  double A;
89  double Ms;
90  double K;
91  Eigen::Vector3d uk;
93  double K3;
94  Eigen::Vector3d ex;
95  Eigen::Vector3d ey;
96  Eigen::Vector3d ez;
98  double P;
99  double N0;
100  double sigma;
101  double lsd;
102  double lsf;
103  double spinHall;
105  double volume = 0;
107  void infos(void);
108  };
109 
134 class Tet : public element<N,NPI>
135  {
136 public:
143  inline Tet(const std::vector<Nodes::Node> &_p_node ,
144  const int _idx ,
145  std::initializer_list<int> _i )
146  : element<N,NPI>(_p_node,_idx,_i), idx(0)
147  {
148  zeroBasing();
149  da.setZero();
150 
151  if (existNodes())
152  {
153  orientate();
154 
155  Eigen::Matrix3d J;
156  // we have to rebuild the jacobian in case of ill oriented tetrahedron
157  double detJ = Jacobian(J);
158 
159  if (fabs(detJ) < Tetra::epsilon)
160  {
161  std::cerr << "Singular jacobian in tetrahedron: |det(J)|= " << fabs(detJ) << std::endl;
162  element::infos();
163  SYSTEM_ERROR;
164  }
165  Eigen::Matrix<double,N,Nodes::DIM> dadu; //Shape function derivatives(constant for linear tetrahedron)
166  dadu << -1., -1., -1., 1., 0., 0., 0., 1., 0., 0., 0., 1.;
167  da = dadu * J.inverse();
168 
169  for (int j = 0; j < NPI; j++)
170  { weight[j] = detJ * Tetra::pds[j]; }// if NPI=1 weight[0] is the tetrahedron volume
171  }
172  // do nothing lambda's (usefull for spin transfer torque)
173  extraField = [] ( Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,Tetra::NPI>> ) {};
174  }
175 
178  Eigen::Matrix<double,N,Nodes::DIM> da;
179 
184  inline void interpolation(std::function<double(Nodes::Node)> getter,
185  Eigen::Ref<Eigen::Matrix<double,NPI,1>> result) const
186  {
187  Eigen::Matrix<double,N,1> scalar_nod;
188  for (int i = 0; i < N; i++) scalar_nod(i) = getter(getNode(i));
189  result = scalar_nod.transpose() * eigen_a;
190  }
191 
194  inline void interpolation(std::function<Eigen::Vector3d(Nodes::Node)> getter,
195  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> result,
196  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> Tx,
197  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> Ty,
198  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> Tz) const
199  {
200  Eigen::Matrix<double,Nodes::DIM,N> vec_nod;
201  for (int i = 0; i < N; i++) vec_nod.col(i) = getter(getNode(i));
202 
203  result = vec_nod * eigen_a;// interpolated value at Gauss point
204 
205  //Spatial derivatives: constant for linear elements, evaluated at any point including Gauss point
206  Tx = vec_nod * (da.col(Nodes::IDX_X)).replicate(1,NPI);
207  Ty = vec_nod * (da.col(Nodes::IDX_Y)).replicate(1,NPI);
208  Tz = vec_nod * (da.col(Nodes::IDX_Z)).replicate(1,NPI);
209  }
210 
214  inline void interpolation_field(std::function<double(Nodes::Node)> getter,
215  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> X) const
216  {
217  Eigen::Matrix<double,N,1> scalar_nod;
218  for (int i = 0; i < N; i++) scalar_nod(i) = getter(getNode(i));
219 
220  X.setZero();
221  for (int j = 0; j < NPI; j++)
222  {
223  for (int i = 0; i < N; i++)
224  {
225  X.col(j) -= (scalar_nod[i] * da.row(i));
226  }
227  }
228  }
229 
232  inline void interpolation(std::function<double(Nodes::Node, Nodes::index)> getter, Nodes::index idx,
233  Eigen::Ref<Eigen::Matrix<double,Tetra::NPI,1>> result) const
234  {
235  Eigen::Matrix<double,N,1> scalar_nod;
236 
237  for (int i = 0; i < N; i++)
238  {
239  scalar_nod[i] = getter(getNode(i),idx);
240  }
241  result = scalar_nod.transpose() * eigen_a;
242  }
243 
253  void lumping(Eigen::Ref<Eigen::Matrix<double,NPI,1>> alpha_eff, double prefactor,
254  Eigen::Ref<Eigen::Matrix<double,3*N,3*N>> AE ) const;
255 
260  Eigen::Matrix<double,N,N> calcDiagBlock(const double c, Eigen::Matrix<double,N,1> &x) const;
261 
266  Eigen::Matrix<double,N,1> calcOffDiagBlock(const Nodes::index idx) const;
267 
269  void add_drift_BE(double alpha, double s_dt, double Vdrift,
270  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> U,
271  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> V,
272  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dUd_,
273  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dVd_,
274  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,N>> BE) const;
275 
278  Eigen::Matrix<double,NPI,1> calc_aniso_uniax(Eigen::Ref<const Eigen::Vector3d> uk,
279  const double Kbis, const double s_dt,
280  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> U,
281  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> V,
282  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> H_aniso) const;
283 
286  Eigen::Matrix<double,NPI,1> calc_aniso_cub(Eigen::Ref<const Eigen::Vector3d> ex,
287  Eigen::Ref<const Eigen::Vector3d> ey,
288  Eigen::Ref<const Eigen::Vector3d> ez,
289  const double K3bis, const double s_dt,
290  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> U,
291  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> V,
292  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> H_aniso) const;
293 
297  void integrales( Tetra::prm const &param, timing const &prm_t,
298  std::function< Eigen::Matrix<double,Nodes::DIM,NPI> (void)> calc_Hext,
299  Nodes::index idx_dir, double Vdrift);
300 
302  double exchangeEnergy(Tetra::prm const &param,
303  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudx,
304  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudy,
305  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudz) const;
306 
308  double anisotropyEnergy(Tetra::prm const &param, Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> u) const;
309 
311  void charges(Tetra::prm const &param,
312  std::function<Eigen::Vector3d(Nodes::Node)> getter,
313  std::vector<double> &srcDen,
314  int &nsrc) const;
315 
317  double demagEnergy(Tetra::prm const &param,
318  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudx,
319  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudy,
320  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudz,
321  Eigen::Ref<Eigen::Matrix<double,NPI,1>> phi) const;
322 
324  double zeemanEnergy(Tetra::prm const &param ,
325  Eigen::Ref<Eigen::Vector3d> const Hext ,
326  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,Tetra::NPI>> const u ) const;
327 
329  double zeemanEnergy(Tetra::prm const &param ,
330  double fieldAmp ,
331  std::vector<Eigen::Matrix<double,Nodes::DIM,Tetra::NPI>> &spaceField,
333  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,Tetra::NPI>> const u) const;
335 
337  double Jacobian(Eigen::Ref<Eigen::Matrix3d> J);
338 
340  double calc_vol(void) const;
341 
343  int idx;
344 
346  std::function<void( Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> H)> extraField;
347 
349  void getPtGauss(Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> result) const
350  {
351  Eigen::Matrix<double,Nodes::DIM,N> vec_nod;
352  for (int i = 0; i < N; i++)
353  {
354  vec_nod.col(i) << getNode(i).p;
355  }
356  result = vec_nod * eigen_a;
357  }
358 
359 private:
361  void orientate(void)
362  {
363  if (calc_vol() < 0.0)
364  std::swap(ind[2], ind[3]);
365  }
366 
367  }; // end class Tetra
368 
371  Eigen::Matrix<double,NPI,1> calc_alpha_eff(const double dt, const double alpha,
372  Eigen::Ref<Eigen::Matrix<double,NPI,1>> uHeff);
373 
375  Eigen::Matrix<double,Nodes::DIM,NPI> calc_gradV(Tet const &tet,std::vector<double> &V);
376 
378  Eigen::Matrix<double,Nodes::DIM,Tetra::NPI> calc_Hst(Tetra::Tet const &tet,
379  const double prefactor,
380  std::vector<Eigen::Vector3d> &s);
381 } // end namespace Tetra
382 
383 #endif /* tetra_h */
Definition: tetra.h:135
void getPtGauss(Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> result) const
Definition: tetra.h:349
double calc_vol(void) const
Definition: tetra.cpp:400
void charges(Tetra::prm const &param, std::function< Eigen::Vector3d(Nodes::Node)> getter, std::vector< double > &srcDen, int &nsrc) const
Definition: tetra.cpp:333
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:109
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:168
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:179
Eigen::Matrix< double, N, Nodes::DIM > da
Definition: tetra.h:178
Tet(const std::vector< Nodes::Node > &_p_node, const int _idx, std::initializer_list< int > _i)
Definition: tetra.h:143
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:300
Eigen::Matrix< double, N, 1 > calcOffDiagBlock(const Nodes::index idx) const
Definition: tetra.cpp:140
double Jacobian(Eigen::Ref< Eigen::Matrix3d > J)
Definition: tetra.cpp:383
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:148
int idx
Definition: tetra.h:343
std::function< void(Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> H)> extraField
Definition: tetra.h:346
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:194
double zeemanEnergy(Tetra::prm const &param, Eigen::Ref< Eigen::Vector3d > const Hext, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, Tetra::NPI >> const u) const
Definition: tetra.cpp:365
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:232
void interpolation(std::function< double(Nodes::Node)> getter, Eigen::Ref< Eigen::Matrix< double, NPI, 1 >> result) const
Definition: tetra.h:184
void interpolation_field(std::function< double(Nodes::Node)> getter, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> X) const
Definition: tetra.h:214
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:352
void orientate(void)
Definition: tetra.h:361
Eigen::Matrix< double, N, N > calcDiagBlock(const double c, Eigen::Matrix< double, N, 1 > &x) const
Definition: tetra.cpp:132
double anisotropyEnergy(Tetra::prm const &param, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> u) const
Definition: tetra.cpp:311
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:203
Template abstract class, mother class for tetraedrons and facettes.
Definition: element.h:25
void zeroBasing(void)
Definition: element.h:108
Eigen::Matrix< double, NPI, 1 > weight
Definition: element.h:50
std::vector< int > ind
Definition: element.h:44
void infos(void)
Definition: element.h:92
bool existNodes(void) const
Definition: element.h:105
const Nodes::Node & getNode(const int i) const
Definition: element.h:102
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:23
constexpr double A
Definition: tetra.h:51
constexpr double a[N][NPI]
Definition: tetra.h:67
constexpr double B
Definition: tetra.h:52
const int NPI
Definition: tetra.h:49
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:25
constexpr double w[NPI]
Definition: tetra.h:58
constexpr double D
Definition: tetra.h:54
constexpr double C
Definition: tetra.h:53
constexpr double pds[NPI]
Definition: tetra.h:59
Eigen::Matrix< double, Nodes::DIM, NPI > calc_gradV(Tet const &tet, std::vector< double > &V)
Definition: tetra.cpp:78
constexpr double v[NPI]
Definition: tetra.h:57
const Eigen::Matrix< double, N, NPI > eigen_a
Definition: tetra.h:75
constexpr double E
Definition: tetra.h:55
constexpr double u[NPI]
Definition: tetra.h:56
Eigen::Matrix< double, Nodes::DIM, Tetra::NPI > calc_Hst(Tetra::Tet const &tet, const double prefactor, std::vector< Eigen::Vector3d > &s)
Definition: tetra.cpp:94
header to define struct Node
Definition: node.h:61
Eigen::Vector3d p
Definition: node.h:62
Definition: tetra.h:85
Eigen::Vector3d ex
Definition: tetra.h:94
double P
Definition: tetra.h:98
double volume
Definition: tetra.h:105
Eigen::Vector3d uk
Definition: tetra.h:91
std::string regName
Definition: tetra.h:86
double A
Definition: tetra.h:88
double alpha_LLG
Definition: tetra.h:87
double N0
Definition: tetra.h:99
double spinHall
Definition: tetra.h:103
double lsd
Definition: tetra.h:101
Eigen::Vector3d ez
Definition: tetra.h:96
void infos(void)
Definition: tetra.cpp:15
Eigen::Vector3d ey
Definition: tetra.h:95
double Ms
Definition: tetra.h:89
double lsf
Definition: tetra.h:102
double K
Definition: tetra.h:90
double K3
Definition: tetra.h:93
double sigma
Definition: tetra.h:100