Rheolef  7.1
an efficient C++ finite element environment
mic.cc
Go to the documentation of this file.
1 // preconditioner eigen, seq implementations
22 //
23 #include "rheolef/mic.h"
24 
25 #ifdef _RHEOLEF_HAVE_EIGEN
26 namespace rheolef {
27 using namespace std;
28 
29 // =========================================================================
30 // the class interface
31 // =========================================================================
32 template<class T, class M>
33 void
35 {
36  using namespace Eigen;
37  Matrix<int,Dynamic,1> nnz_row (a.nrow());
38  typename csr<T,M>::const_iterator ia = a.begin();
39  for (size_type i = 0, q = 0, n = a.nrow(); i < n; ++i) {
40  nnz_row[i] = ia[i+1] - ia[i];
41  }
42  SparseMatrix<T> a_tmp (a.nrow(),a.ncol());
43  if (a.nrow() != 0) {
44  a_tmp.reserve (nnz_row);
45  }
46  for (size_type i = 0, n = a.nrow(); i < n; ++i) {
47  for (typename csr<T,M>::const_data_iterator p = ia[i]; p < ia[i+1]; ++p) {
48  a_tmp.insert (i, (*p).first) = (*p).second;
49  }
50  }
51  a_tmp.makeCompressed();
52  if (a.nrow() != 0) {
53  _mic_a.setInitialShift (_shift);
54  _mic_a.compute (a_tmp);
55  check_macro (_mic_a.info() == Success, "eigen MIC incomplete factorization failed");
56  }
57 }
58 template<class T, class M>
61 {
62  if (b.dis_size() == 0) return b; // empty matrix
63  vec<T,M> x(b.ownership());
64  using namespace Eigen;
65  Map<Matrix<T,Dynamic,1> > b_map ((T*)(&(*b.begin())), b.size()),
66  x_map ( &(*x.begin()), x.size());
67  x_map = _mic_a.solve (b_map);
68  return x;
69 }
70 template<class T, class M>
73 {
74  return solve(b);
75 }
76 // ----------------------------------------------------------------------------
77 // instanciation in library
78 // ----------------------------------------------------------------------------
79 
81 
82 #ifdef _RHEOLEF_HAVE_MPI
84 #endif // _RHEOLEF_HAVE_MPI
85 
86 } // namespace rheolef
87 #endif // HAVE_EIGEN
void update_values(const csr< T, M > &a)
Definition: mic.cc:34
base::size_type size_type
Definition: mic.h:98
vec< T, M > trans_solve(const vec< T, M > &rhs) const
Definition: mic.cc:72
vec< T, M > solve(const vec< T, M > &rhs) const
Definition: mic.cc:60
Expr1::float_type T
Definition: field_expr.h:261
check_macro(expr1.have_homogeneous_space(Xh1), "dual(expr1,expr2); expr1 should have homogeneous space. HINT: use dual(interpolate(Xh, expr1),expr2)")
This file is part of Rheolef.
void solve(tiny_matrix< T > &a, tiny_vector< size_t > &piv, const tiny_vector< T > &b, tiny_vector< T > &x)
Definition: tiny_lu.h:92
Definition: sphere.icc:25