Feellgood
electrostatSolver.h
Go to the documentation of this file.
1 
8 #include <iostream>
9 #include <map>
10 
11 #include "config.h"
12 #include "fem.h"
13 #include "mesh.h"
14 
15 #include "facette.h"
16 #include "tetra.h"
17 
18 #include "spinTransferTorque.h"
19 
21 inline void assembling_mat(Tetra::Tet const &tet, double Ke[Tetra::N][Tetra::N], std::vector<Eigen::Triplet<double>> &K)
22  {
23  for (int ie = 0; ie < Tetra::N; ie++)
24  {
25  for (int je = 0; je < Tetra::N; je++)
26  {
27  double val = Ke[ie][je];
28  if (val != 0)
29  { K.push_back( Eigen::Triplet(tet.ind[ie], tet.ind[je], val) ); }
30  }
31  }
32  }
33 
35 inline void assembling_vect(Facette::Fac const &fac, std::vector<double> const &Le, Eigen::Ref<Eigen::VectorXd> L)
36  {
37  for (int ie = 0; ie < Facette::N; ie++)
38  { L(fac.ind[ie]) += Le[ie]; }
39  }
40 
43 void integrales(Tetra::Tet const &tet, double sigma, double AE[Tetra::N][Tetra::N]);
44 
46 void integrales(Facette::Fac const &fac, double pot_val, std::vector<double> &BE);
47 
53  {
54 public:
57  Mesh::mesh const &_msh ,
58  STT const &_p_stt ,
59  const double _tol ,
60  const bool v ,
61  const int max_iter ,
62  const std::string
63  _fileName )
64  : msh(_msh), p_stt(_p_stt), verbose(v), MAXITER(max_iter), precision(PRECISION_STT),
65  fileName(_fileName)
66  {
68  D0 = 2.0 * p_stt.sigma / (Nodes::sq(CHARGE_ELECTRON) * p_stt.N0);
69  pf = Nodes::sq(p_stt.lJ) / (D0 * (1. + ksi * ksi)) * BOHRS_MUB * p_stt.beta / CHARGE_ELECTRON;
70 
71  if (verbose)
72  {
73  std::cout << "Dirichlet boundary conditions..." << std::endl;
74  infos();
75  }
76  V.resize(_msh.getNbNodes());
77  bool has_converged = solve(_tol);
78  if (has_converged)
79  {
80  if (p_stt.V_file)
81  {
82  if (verbose)
83  {
84  std::cout << "writing electrostatic potential solutions to file " << fileName
85  << std::endl;
86  }
87  bool iznogood = msh.savesol(precision, fileName, "## columns: index\tV\n", V);
88  if (verbose && iznogood)
89  {
90  std::cout << "file " << fileName << " status : " << iznogood << std::endl;
91  }
92  }
93  std::for_each(msh.tet.begin(), msh.tet.end(),
94  [this](Tetra::Tet const &tet)
95  {
96  Eigen::Matrix<double,Nodes::DIM,Tetra::NPI> _gradV;
97  calc_gradV(tet, _gradV);
98  gradV.push_back(_gradV);
99 
100  Eigen::Matrix<double,Nodes::DIM,Tetra::NPI> _Hm;
101  calc_Hm(tet, _gradV, _Hm);
102  Hm.push_back(_Hm);
103  });
104  prepareExtras();
105  }
106  else
107  {
108  std::cerr << "Solver (STT) has not converged" << std::endl;
109  exit(1);
110  }
111  }
112 
114  void calc_gradV(Tetra::Tet const &tet, Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,Tetra::NPI>> _gradV) const
115  {
116  for (int npi = 0; npi < Tetra::NPI; npi++)
117  {
118  Eigen::Vector3d v(0,0,0);
119  for (int i = 0; i < Tetra::N; i++)
120  {
121  v += V[tet.ind[i]] * tet.da.row(i);
122  }
123  _gradV.col(npi) = v;
124  }
125  }
126 
128  void calc_Hm(Tetra::Tet const &tet, Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,Tetra::NPI>> _gradV,
129  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,Tetra::NPI>> _Hm) const
130  {
131  Eigen::Matrix<double,Nodes::DIM,Tetra::NPI> p_g;
132  tet.getPtGauss(p_g);
133 
134  for (int npi = 0; npi < Tetra::NPI; npi++)
135  { _Hm.col(npi) = -p_stt.sigma * _gradV.col(npi).cross(p_g.col(npi)); }
136  }
137 
138 private:
140  void prepareExtras(void)
141  {
142  std::for_each(
143  msh.tet.begin(), msh.tet.end(),
144  [this](Tetra::Tet &tet)
145  {
146  const int _idx = tet.idx;
147  tet.extraField = [this, _idx](Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,Tetra::NPI>> H)
148  {
149  for(int npi = 0; npi<Tetra::NPI; npi++) { H.col(npi) += Hm[_idx].col(npi); }
150  };
151 
152  tet.extraCoeffs_BE = [this, &tet](double Js, Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,Tetra::NPI>> U,
153  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,Tetra::NPI>> dUdx,
154  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,Tetra::NPI>> dUdy,
155  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,Tetra::NPI>> dUdz,
156  Eigen::Ref<Eigen::Matrix<double,Nodes::DIM,Tetra::N>> BE)
157  {
158  const double prefactor = D0 / Nodes::sq(p_stt.lJ) / (gamma0*Js/mu0);
159  for (int npi = 0; npi < Tetra::NPI; npi++)
160  {
161  Eigen::Vector3d const &_gV = gradV[tet.idx].col(npi);
162  Eigen::Vector3d j_grad_u =
163  -p_stt.sigma
164  * Eigen::Vector3d(_gV.dot( Eigen::Vector3d(dUdx(Nodes::IDX_X,npi), dUdy(Nodes::IDX_X,npi), dUdz(Nodes::IDX_X,npi))),
165  _gV.dot( Eigen::Vector3d(dUdx(Nodes::IDX_Y,npi), dUdy(Nodes::IDX_Y,npi), dUdz(Nodes::IDX_Y,npi))),
166  _gV.dot( Eigen::Vector3d(dUdx(Nodes::IDX_Z,npi), dUdy(Nodes::IDX_Z,npi), dUdz(Nodes::IDX_Z,npi))));
167 
168  Eigen::Vector3d m = pf * (ksi * j_grad_u + U.col(npi).cross(j_grad_u));
169  for (int i = 0; i < Tetra::N; i++)
170  {
171  const double ai_w = tet.weight[npi] * Tetra::a[i][npi];
172  BE.col(i) += ai_w*(Hm[tet.idx].col(npi) + prefactor*m);
173  }
174  } // end loop on npi
175  }; //end lambda
176  });//end for_each
177  }
178 
180  double ksi;
181 
183  double D0;
184 
186  double pf;
187 
191 
194 
196  const bool verbose;
197 
199  const unsigned int MAXITER; // fixed to 5000 in ref code
200 
202  const int precision;
203 
205  const std::string fileName;
206 
209  std::vector<double> V;
210 
212  std::vector< Eigen::Matrix<double,Nodes::DIM,Tetra::NPI> > gradV;
213 
216  std::vector< Eigen::Matrix<double,Nodes::DIM,Tetra::NPI> > Hm;
217 
219  inline void infos(void)
220  {
221  std::cout << "sigma: " << p_stt.sigma << std::endl;
222 
223  std::for_each(p_stt.boundaryCond.begin(), p_stt.boundaryCond.end(),
224  [](std::pair<std::string, double> const &p)
225  { std::cout << "regName: " << p.first << "\tV :" << p.second << std::endl; });
226  }
227 
229  void prepareData(std::vector<Eigen::Triplet<double>> &Kw, Eigen::Ref<Eigen::VectorXd> Lw)
230  {
231  const double sigma = p_stt.sigma;
232 
233  std::for_each(msh.tet.begin(), msh.tet.end(),
234  [this, sigma, &Kw](Tetra::Tet const &tet)
235  {
236  double K[Tetra::N][Tetra::N] = {{0}};
237  integrales(tet, sigma, K);
238  assembling_mat(tet, K, Kw);
239  });
240 
241  double pot_val = 0; // we initialize pot_val to the average of the potentials set by the
242  // boundary conditions (does not seem to change convergence speed
243  // whatever value it is ..)
244  std::for_each(p_stt.boundaryCond.begin(), p_stt.boundaryCond.end(),
245  [&pot_val](auto const &it) { pot_val += it.second; });
246  pot_val /= p_stt.boundaryCond.size();
247 
248  std::for_each(msh.fac.begin(), msh.fac.end(),
249  [this, pot_val, &Lw](Facette::Fac const &fac)
250  {
251  std::vector<double> L(Facette::N);
252  integrales(fac, pot_val, L);
253  assembling_vect(fac, L, Lw);
254  });
255  }
256 
259  int solve(const double _tol)
260  {
261  const int NOD = msh.getNbNodes();
262 
263  std::vector<Eigen::Triplet<double>> Kw;
264  Eigen::VectorXd Lw = Eigen::VectorXd::Zero(NOD);
265  prepareData(Kw, Lw);
266 
267  Eigen::VectorXd Xw(NOD);
268  /* do something here with boundary conditions */
269 
270  if (verbose)
271  { std::cout << "line weighting..." << std::endl; }
272 
273  Eigen::SparseMatrix<double> Kr(NOD,NOD);
274  Kr.setFromTriplets(Kw.begin(),Kw.end());
275  std::vector<double> maxRow(NOD,0);
276 
277  //here we fill the vector maxRow with infinity norm : maxRow[i] = max{|Kr.row(i)|}
278  for(int i=0;i<NOD;i++)
279  {
280  for(int k=0;k<Kr.outerSize();++k)
281  for(Eigen::SparseMatrix<double>::InnerIterator it(Kr,k); it; ++it)
282  { if((it.row() == i)&&(fabs(it.value()) > maxRow[i])) { maxRow[i] = fabs(it.value()); } }
283  }
284 
285  for (int i = 0; i < NOD; i++)
286  {
287  double norme = maxRow[i];
288 
289  Lw[i] /= norme;
290  Kr.row(i) /= norme;
291  }
292 
293  if (verbose)
294  { std::cout << "solving ..." << std::endl; }
295 
296  Eigen::BiCGSTAB<Eigen::SparseMatrix<double>> _solver;
297 
298  _solver.setTolerance(_tol);
299  _solver.setMaxIterations(MAXITER);
300  _solver.compute(Kr);
301 
302  Eigen::VectorXd sol = _solver.solve(Lw);
303  for (int i=0;i<NOD;i++)
304  { V[i]= sol(i); }
305  return (_solver.iterations() < MAXITER);
306  }
307 
308  }; // end class electrostatSolver
Definition: facette.h:68
Definition: mesh.h:28
int getNbNodes(void) const
Definition: mesh.h:98
void savesol(const int precision, const std::string fileName, std::string const &metadata) const
Definition: save.cpp:153
std::vector< Tetra::Tet > tet
Definition: mesh.h:162
std::vector< Facette::Fac > fac
Definition: mesh.h:159
Definition: tetra.h:124
void getPtGauss(Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, NPI >> result) const
Definition: tetra.h:321
std::function< void(double Js, 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:318
Eigen::Matrix< double, N, Nodes::DIM > da
Definition: tetra.h:169
Definition: electrostatSolver.h:53
int solve(const double _tol)
Definition: electrostatSolver.h:259
std::vector< Eigen::Matrix< double, Nodes::DIM, Tetra::NPI > > Hm
Definition: electrostatSolver.h:216
const std::string fileName
Definition: electrostatSolver.h:205
const bool verbose
Definition: electrostatSolver.h:196
const int precision
Definition: electrostatSolver.h:202
const unsigned int MAXITER
Definition: electrostatSolver.h:199
Mesh::mesh msh
Definition: electrostatSolver.h:190
double D0
Definition: electrostatSolver.h:183
electrostatSolver(Mesh::mesh const &_msh, STT const &_p_stt, const double _tol, const bool v, const int max_iter, const std::string _fileName)
Definition: electrostatSolver.h:56
std::vector< double > V
Definition: electrostatSolver.h:209
void prepareExtras(void)
Definition: electrostatSolver.h:140
std::vector< Eigen::Matrix< double, Nodes::DIM, Tetra::NPI > > gradV
Definition: electrostatSolver.h:212
double ksi
Definition: electrostatSolver.h:180
STT p_stt
Definition: electrostatSolver.h:193
void infos(void)
Definition: electrostatSolver.h:219
double pf
Definition: electrostatSolver.h:186
void calc_gradV(Tetra::Tet const &tet, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, Tetra::NPI >> _gradV) const
Definition: electrostatSolver.h:114
void prepareData(std::vector< Eigen::Triplet< double >> &Kw, Eigen::Ref< Eigen::VectorXd > Lw)
Definition: electrostatSolver.h:229
void calc_Hm(Tetra::Tet const &tet, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, Tetra::NPI >> _gradV, Eigen::Ref< Eigen::Matrix< double, Nodes::DIM, Tetra::NPI >> _Hm) const
Definition: electrostatSolver.h:128
std::vector< int > ind
Definition: element.h:44
void integrales(Tetra::Tet const &tet, double sigma, double AE[Tetra::N][Tetra::N])
Definition: electrostatSolver.cpp:3
void assembling_vect(Facette::Fac const &fac, std::vector< double > const &Le, Eigen::Ref< Eigen::VectorXd > L)
Definition: electrostatSolver.h:35
void assembling_mat(Tetra::Tet const &tet, double Ke[Tetra::N][Tetra::N], std::vector< Eigen::Triplet< double >> &K)
Definition: electrostatSolver.h:21
contains namespace Facette header containing Fac class, and some constants and a less_than operator t...
principal header, contains the struct fem This file is called by many source files in feellgood....
class mesh, readMesh is expecting a mesh file in gmsh format either text or binary,...
const int N
Definition: facette.h:17
constexpr double v[NPI]
Definition: facette.h:22
double sq(const double x)
Definition: node.h:42
const int NPI
Definition: tetra.h:23
const int N
Definition: tetra.h:22
Definition: spinTransferTorque.h:13
double sigma
Definition: spinTransferTorque.h:21
double lsf
Definition: spinTransferTorque.h:27
double beta
Definition: spinTransferTorque.h:15
double lJ
Definition: spinTransferTorque.h:24
double N0
Definition: spinTransferTorque.h:18
std::vector< std::pair< std::string, double > > boundaryCond
Definition: spinTransferTorque.h:36
bool V_file
Definition: spinTransferTorque.h:31
namespace Tetra header containing Tet class, some constants, and integrales