Feellgood
|
#include <linear_algebra.h>
Public Member Functions | |
LinAlgebra (Settings &s, Mesh::mesh &my_msh) | |
algebra::MatrixShape | build_shape () |
void | prepareElements (Eigen::Vector3d const &Hext, timing const &t_prm) |
void | prepareElements (double const A_Hext, timing const &t_prm) |
void | buildInitGuess (std::vector< double > &G) const |
int | solver (timing const &t_prm) |
void | set_DW_vz (double vz) |
double | get_v_max (void) |
void | setExtSpaceField (Settings &s) |
void | base_projection () |
Private Attributes | |
Nodes::index | idx_dir |
const int | NOD |
Mesh::mesh * | refMsh |
algebra::iteration< double > | iter |
algebra::r_sparseMat | K |
std::vector< double > | L_rhs |
std::vector< double > | Xw |
const int | verbose |
const std::vector< Tetra::prm > & | prmTetra |
const std::vector< Facette::prm > & | prmFacette |
double | DW_vz |
double | v_max |
std::vector< Eigen::Matrix< double, Nodes::DIM, Tetra::NPI > > | extSpaceField |
convenient class to grab altogether some part of the calculations involved using algebra::bicg solver at each timestep. The solver is handled by solver method, and is using algebra sparse matrices(Row major). The write sparse matrix is prepared in 'batch mode', to add all non zero coefficients with a '+=' logic. Then it is turned into a read sparse matrix before being used by bicg algo. Be aware of time units: when entering solver method, division by gamma0 and multiplication by gamma0 when ending are mandatory. The bicg algorithm is monitored by iter object. When debugging it might be usefull to set iter verbosity differently from LinAlgebra Solver (see constructor list initialization)
|
inline |
constructor When debugging it might be usefull to set iter verbosity differently
[in] | s | |
[in] | my_msh |
void LinAlgebra::base_projection | ( | ) |
computes local vector basis {ep,eq} in the tangeant plane for projection on the elements
|
inline |
build a matrix shape suitable for our mesh
void LinAlgebra::buildInitGuess | ( | std::vector< double > & | G | ) | const |
build init guess for bicg solver
[out] | G |
|
inline |
getter for v_max
void LinAlgebra::prepareElements | ( | double const | A_Hext, |
timing const & | t_prm | ||
) |
computes inner data structures of tetraedrons and triangular facettes (K matrices and L vectors) this member function is overloaded to fit to two different situations, either if std::function passed to element = tetra is corresponding to the simple case of constant external field applied to the magnetic region or space dependant. Here is the variable space applied field.
[in] | A_Hext | amplitude applied field (might be time dependant) |
[in] | t_prm |
void LinAlgebra::prepareElements | ( | Eigen::Vector3d const & | Hext, |
timing const & | t_prm | ||
) |
computes inner data structures of tetraedrons and triangular facettes (K matrices and L vectors) this member function is overloaded to fit to two different situations, either if std::function passed to element = tetra is corresponding to the simple case of constant external field applied to the magnetic region or space dependant. Here is the constant space applied field.
[in] | Hext | applied field |
[in] | t_prm |
|
inline |
setter for DW_dz
[in] | vz |
void LinAlgebra::setExtSpaceField | ( | Settings & | s | ) |
when external applied field is of field_type R4toR3 values of field_space are stored in spaceField
[in] | s |
int LinAlgebra::solver | ( | timing const & | t_prm | ) |
solver, uses stabilized biconjugate gradient solver (bicg) with diagonal preconditionner, sparse matrix and vector are filled with multiThreading. Sparse matrix is row major.
[in] | t_prm |
|
private |
speed of the domain wall
|
private |
external applied space field, values on gauss points, size is number of tetraedrons
|
private |
recentering index direction if any
|
private |
monitor the solver
|
private |
matrix of the system to solve
|
private |
RHS vector of the system to solve
|
private |
number of nodes, also an offset for filling sparseMatrix, initialized by constructor
|
private |
material parameters of the facettes
|
private |
material parameters of the tetrahedrons
|
private |
direct access to the mesh
|
private |
maximum speed of the magnetization in the whole physical object
|
private |
verbosity
|
private |
solution of the system to solve