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 "triangle.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.};
39  constexpr double a[N][NPI] = {{1. - u[0] - v[0] - w[0]},
40  {u[0]},
41  {v[0]},
42  {w[0]}};
43 
45  const Eigen::Matrix<double,N,NPI> eigen_a = (Eigen::MatrixXd(N,NPI) << a[0][0],
46  a[1][0],
47  a[2][0],
48  a[3][0] ).finished();
49 #else
50  const int NPI = 5;
52  constexpr double A = 1. / 4.;
53  constexpr double B = 1. / 6.;
54  constexpr double C = 1. / 2.;
55  constexpr double D = -2. / 15.;
56  constexpr double E = 3. / 40.;
57  constexpr double u[NPI] = {A, B, B, B, C};
58  constexpr double v[NPI] = {A, B, B, C, B};
59  constexpr double w[NPI] = {A, B, C, B, B};
60  constexpr double pds[NPI] = {D, E, E, E, E};
68  constexpr double a[N][NPI] = {{1. - u[0] - v[0] - w[0], 1. - u[1] - v[1] - w[1],
69  1. - u[2] - v[2] - w[2], 1. - u[3] - v[3] - w[3],
70  1. - u[4] - v[4] - w[4]},
71  {u[0], u[1], u[2], u[3], u[4]},
72  {v[0], v[1], v[2], v[3], v[4]},
73  {w[0], w[1], w[2], w[3], w[4]}};
74 
76  const Eigen::Matrix<double,N,NPI> eigen_a =
77  (Eigen::MatrixXd(N,NPI) << a[0][0], a[0][1], a[0][2], a[0][3], a[0][4],
78  a[1][0], a[1][1], a[1][2], a[1][3], a[1][4],
79  a[2][0], a[2][1], a[2][2], a[2][3], a[2][4],
80  a[3][0], a[3][1], a[3][2], a[3][3], a[3][4] ).finished();
81 #endif
82 
86 struct prm
87  {
88  std::string regName;
89  double alpha_LLG;
90  double A;
91  double Ms;
93  double K;
94  Eigen::Vector3d uk;
96  double K3;
97  Eigen::Vector3d ex;
98  Eigen::Vector3d ey;
99  Eigen::Vector3d ez;
101  double P;
102  double N0;
104  double sigma;
106  double lsd;
107  double lsf;
109  double volume = 0;
111  void infos(void);
112  };
113 
131 class Tet : public element<N,NPI>
132  {
133 public:
140  inline Tet(const std::vector<Nodes::Node> &_p_node ,
141  const int _idx ,
142  const std::initializer_list<int> _i )
143  : element<N,NPI>(_p_node,_idx,_i), idx(0)
144  {
145  if (existNodes())
146  {
147  orientate(); // enforce the correct orientation
148 
149  Eigen::Matrix3d J;
150  double detJ = Jacobian(J);
151  Eigen::Matrix<double,N,Nodes::DIM> dadu; // Shape function derivatives
152  // (constant for linear tetrahedron)
153  dadu << -1., -1., -1., 1., 0., 0., 0., 1., 0., 0., 0., 1.;
154  da = dadu * J.inverse();
155 
156  for (int j = 0; j < NPI; j++)
157  { weight[j] = detJ * Tetra::pds[j]; }// if NPI=1 weight[0] is the tetrahedron volume
158  // do nothing lambda's (usefull for spin transfer torque)
159  extraField = [] ( const Eigen::Ref<const Eigen::Matrix<double,Nodes::DIM,Tetra::NPI>> ) {};
160  }
161  else
162  { std::cerr<<"Error: Tet constructor has out of bound index\n"; }
163  }
164 
167  Eigen::Matrix<double,N,Nodes::DIM> da;
168 
173  inline void interpolation(const std::function<double(Nodes::Node)>& getter,
174  Eigen::Ref<Eigen::Matrix<double,NPI,1>> result) const
175  {
176  Eigen::Matrix<double,N,1> scalar_nod;
177  for (int i = 0; i < N; i++) scalar_nod(i) = getter(getNode(i));
178  result = scalar_nod.transpose() * eigen_a;
179  }
180 
183  inline void interpolation(const std::function<Eigen::Vector3d(Nodes::Node)>& getter,
184  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> result,
185  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> Tx,
186  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> Ty,
187  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> Tz) const
188  {
189  Eigen::Matrix<double,Nodes::DIM,N> vec_nod;
190  for (int i = 0; i < N; i++) vec_nod.col(i) = getter(getNode(i));
191 
192  result = vec_nod * eigen_a;// interpolated value at Gauss point
193 
194  //Spatial derivatives: constant for linear elements, evaluated at any point including Gauss
195  //point
196  Tx = vec_nod * (da.col(Nodes::IDX_X)).replicate(1,NPI);
197  Ty = vec_nod * (da.col(Nodes::IDX_Y)).replicate(1,NPI);
198  Tz = vec_nod * (da.col(Nodes::IDX_Z)).replicate(1,NPI);
199  }
200 
204  inline void interpolation_field(const std::function<double(Nodes::Node)>& getter,
205  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> X) const
206  {
207  Eigen::Matrix<double,N,1> scalar_nod;
208  for (int i = 0; i < N; i++) scalar_nod(i) = getter(getNode(i));
209 
210  X.setZero();
211  for (int j = 0; j < NPI; j++)
212  {
213  for (int i = 0; i < N; i++)
214  {
215  X.col(j) -= (scalar_nod[i] * da.row(i));
216  }
217  }
218  }
219 
222  inline void interpolation(const std::function<double(Nodes::Node, Nodes::index)>& getter,
223  const Nodes::index idx, Eigen::Ref<Eigen::Matrix<double,Tetra::NPI,1>> result) const
224  {
225  Eigen::Matrix<double,N,1> scalar_nod;
226 
227  for (int i = 0; i < N; i++)
228  {
229  scalar_nod[i] = getter(getNode(i),idx);
230  }
231  result = scalar_nod.transpose() * eigen_a;
232  }
233 
244  void lumping(const Eigen::Ref<const Eigen::Matrix<double,NPI,1>> alpha_eff, double prefactor,
245  Eigen::Ref<Eigen::Matrix<double,3*N,3*N>> AE ) const;
246 
251  Eigen::Matrix<double,N,N> calcDiagBlock(const double c,
252  const Eigen::Matrix<double,N,1> &x) const;
253 
258  Eigen::Matrix<double,N,1> calcOffDiagBlock(const Nodes::index idx) const;
259 
261  void add_drift_BE(double alpha, double s_dt, double Vdrift,
262  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> U,
263  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> V,
264  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dUd_,
265  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dVd_,
266  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,N>> BE) const;
267 
270  Eigen::Matrix<double,NPI,1> calc_aniso_uniax(const Eigen::Ref<const Eigen::Vector3d> uk,
271  const double Kbis, const double s_dt,
272  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> U,
273  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> V,
274  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> H_aniso) const;
275 
278  Eigen::Matrix<double,NPI,1> calc_aniso_cub(const Eigen::Ref<const Eigen::Vector3d> ex,
279  const Eigen::Ref<const Eigen::Vector3d> ey,
280  const Eigen::Ref<const Eigen::Vector3d> ez,
281  const double K3bis, const double s_dt,
282  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> U,
283  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> V,
284  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> H_aniso) const;
285 
289  void integrales( const Tetra::prm &param, const timing &prm_t,
290  const std::function< Eigen::Matrix<double,Nodes::DIM,NPI> (void)>& calc_Hext,
291  Nodes::index idx_dir, double Vdrift);
292 
294  double exchangeEnergy(const Tetra::prm &param,
295  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudx,
296  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudy,
297  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudz) const;
298 
300  double uniaxialAnisotropyEnergy(const Tetra::prm &param,
301  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> u) const;
302 
304  double cubicAnisotropyEnergy(const Tetra::prm &param,
305  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> u) const;
306 
308  Eigen::Matrix<double,Tetra::NPI,1> charges(const double &Ms,
309  const std::function<Eigen::Vector3d(const Nodes::Node&)> &getter) const override;
310 
312  double demagEnergy(const Tetra::prm &param,
313  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudx,
314  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudy,
315  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> dudz,
316  Eigen::Ref<Eigen::Matrix<double,NPI,1>> phi) const;
317 
319  double zeemanEnergy(const Tetra::prm &param ,
320  const Eigen::Ref<const Eigen::Vector3d> Hext ,
321  const Eigen::Ref<const Eigen::Matrix<double,Nodes::DIM,Tetra::NPI>> u ) const;
322 
324  double zeemanEnergy(const Tetra::prm &param ,
325  double fieldAmp ,
326  const std::vector<Eigen::Matrix<double,Nodes::DIM,Tetra::NPI>> &spaceField,
328  const Eigen::Ref<const Eigen::Matrix<double,Nodes::DIM,Tetra::NPI>> u) const;
330 
332  double Jacobian(Eigen::Ref<Eigen::Matrix3d> J) const;
333 
335  double calc_vol(void) const;
336 
338  int idx;
339 
341  std::function<void( Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,NPI>> H)> extraField;
342 
344  Eigen::Matrix<double,Nodes::DIM,NPI> getPtGauss(void) const override
345  {
346  Eigen::Matrix<double,Nodes::DIM,N> vec_nod;
347  for (int i = 0; i < N; i++)
348  {
349  vec_nod.col(i) << getNode(i).p;
350  }
351  return vec_nod*eigen_a;
352  }
353 
354 private:
356  void orientate(void) override;
357  }; // end class Tetra
358 
361  Eigen::Matrix<double,NPI,1> calc_alpha_eff(const double dt, const double alpha,
362  Eigen::Ref<Eigen::Matrix<double,NPI,1>> uHeff);
363 
365  Eigen::Matrix<double,Nodes::DIM,NPI> calc_gradV(Tet const &tet, const std::vector<double> &V);
366 
368  Eigen::Matrix<double,Nodes::DIM,Tetra::NPI> calc_Hst(const Tetra::Tet &tet,
369  const double prefactor,
370  const std::vector<Eigen::Vector3d> &s);
371 } // end namespace Tetra
372 
373 #endif /* tetra_h */
Definition: tetra.h:132
Eigen::Matrix< double, NPI, 1 > calc_aniso_cub(const Eigen::Ref< const Eigen::Vector3d > ex, const Eigen::Ref< const Eigen::Vector3d > ey, const 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:183
double calc_vol(void) const
Definition: tetra.cpp:426
Eigen::Matrix< double, N, N > calcDiagBlock(const double c, const Eigen::Matrix< double, N, 1 > &x) const
Definition: tetra.cpp:133
double zeemanEnergy(const Tetra::prm &param, const Eigen::Ref< const Eigen::Vector3d > Hext, const Eigen::Ref< const Eigen::Matrix< double, Nodes::DIM, Tetra::NPI >> u) const
Definition: tetra.cpp:374
double cubicAnisotropyEnergy(const Tetra::prm &param, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> u) const
Definition: tetra.cpp:330
double Jacobian(Eigen::Ref< Eigen::Matrix3d > J) const
Definition: tetra.cpp:393
Eigen::Matrix< double, N, Nodes::DIM > da
Definition: tetra.h:167
double uniaxialAnisotropyEnergy(const Tetra::prm &param, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> u) const
Definition: tetra.cpp:320
void integrales(const Tetra::prm &param, const timing &prm_t, const std::function< Eigen::Matrix< double, Nodes::DIM, NPI >(void)> &calc_Hext, Nodes::index idx_dir, double Vdrift)
Definition: tetra.cpp:210
double exchangeEnergy(const Tetra::prm &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:309
Eigen::Matrix< double, Tetra::NPI, 1 > charges(const double &Ms, const std::function< Eigen::Vector3d(const Nodes::Node &)> &getter) const override
Definition: tetra.cpp:347
Eigen::Matrix< double, N, 1 > calcOffDiagBlock(const Nodes::index idx) const
Definition: tetra.cpp:142
Eigen::Matrix< double, Nodes::DIM, NPI > getPtGauss(void) const override
Definition: tetra.h:344
void interpolation(const std::function< double(Nodes::Node)> &getter, Eigen::Ref< Eigen::Matrix< double, NPI, 1 >> result) const
Definition: tetra.h:173
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:150
void interpolation(const 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:183
int idx
Definition: tetra.h:338
void lumping(const Eigen::Ref< const Eigen::Matrix< double, NPI, 1 >> alpha_eff, double prefactor, Eigen::Ref< Eigen::Matrix< double, 3 *N, 3 *N >> AE) const
Definition: tetra.cpp:108
std::function< void(Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> H)> extraField
Definition: tetra.h:341
void orientate(void) override
Definition: tetra.cpp:410
Tet(const std::vector< Nodes::Node > &_p_node, const int _idx, const std::initializer_list< int > _i)
Definition: tetra.h:140
void interpolation_field(const std::function< double(Nodes::Node)> &getter, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> X) const
Definition: tetra.h:204
void interpolation(const std::function< double(Nodes::Node, Nodes::index)> &getter, const Nodes::index idx, Eigen::Ref< Eigen::Matrix< double, Tetra::NPI, 1 >> result) const
Definition: tetra.h:222
Eigen::Matrix< double, NPI, 1 > calc_aniso_uniax(const 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:171
double demagEnergy(const Tetra::prm &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:361
Template abstract class, mother class for tetraedrons and triangles.
Definition: element.h:30
Eigen::Matrix< double, NPI, 1 > weight
Definition: element.h:59
bool existNodes(void) const
Definition: element.h:119
const Nodes::Node & getNode(const int i) const
Definition: element.h:131
Definition: time_integration.h:7
index
Definition: node.h:33
const double epsilon
Definition: tetra.h:23
constexpr double A
Definition: tetra.h:52
constexpr double a[N][NPI]
Definition: tetra.h:68
constexpr double B
Definition: tetra.h:53
const int NPI
Definition: tetra.h:50
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:47
const int N
Definition: tetra.h:25
constexpr double w[NPI]
Definition: tetra.h:59
constexpr double D
Definition: tetra.h:55
constexpr double C
Definition: tetra.h:54
constexpr double pds[NPI]
Definition: tetra.h:60
constexpr double v[NPI]
Definition: tetra.h:58
const Eigen::Matrix< double, N, NPI > eigen_a
Definition: tetra.h:76
constexpr double E
Definition: tetra.h:56
constexpr double u[NPI]
Definition: tetra.h:57
Eigen::Matrix< double, Nodes::DIM, NPI > calc_gradV(Tet const &tet, const std::vector< double > &V)
Definition: tetra.cpp:77
Eigen::Matrix< double, Nodes::DIM, Tetra::NPI > calc_Hst(const Tetra::Tet &tet, const double prefactor, const std::vector< Eigen::Vector3d > &s)
Definition: tetra.cpp:93
header to define struct Node
Definition: node.h:60
Eigen::Vector3d p
Definition: node.h:61
Definition: tetra.h:87
Eigen::Vector3d ex
Definition: tetra.h:97
double P
Definition: tetra.h:101
double volume
Definition: tetra.h:109
Eigen::Vector3d uk
Definition: tetra.h:94
std::string regName
Definition: tetra.h:88
double A
Definition: tetra.h:90
double alpha_LLG
Definition: tetra.h:89
double N0
Definition: tetra.h:102
double lsd
Definition: tetra.h:106
Eigen::Vector3d ez
Definition: tetra.h:99
void infos(void)
Definition: tetra.cpp:14
Eigen::Vector3d ey
Definition: tetra.h:98
double Ms
Definition: tetra.h:91
double lsf
Definition: tetra.h:107
double K
Definition: tetra.h:93
double K3
Definition: tetra.h:96
double sigma
Definition: tetra.h:104
contains namespace Triangle header containing Tri class, and some constants and a less_than operator ...