40 #ifndef GMM_PRECOND_ILUTP_H 41 #define GMM_PRECOND_ILUTP_H 56 template <
typename Matrix>
59 typedef typename linalg_traits<Matrix>::value_type value_type;
62 typedef row_matrix<_rsvector> LU_Matrix;
63 typedef col_matrix<_wsvector> CLU_Matrix;
67 gmm::unsorted_sub_index indperm;
68 gmm::unsorted_sub_index indperminv;
69 mutable std::vector<value_type> temporary;
75 template<
typename M>
void do_ilutp(
const M&, row_major);
76 void do_ilutp(
const Matrix&, col_major);
79 void build_with(
const Matrix& A,
int k_ = -1,
double eps_ =
double(-1)) {
81 if (eps_ >=
double(0)) eps = eps_;
85 do_ilutp(A,
typename principal_orientation_type<
typename 86 linalg_traits<Matrix>::sub_orientation>::potype());
89 : L(mat_nrows(A), mat_ncols(A)), U(mat_nrows(A), mat_ncols(A)),
90 K(k_), eps(eps_) { build_with(A); }
93 size_type memsize()
const {
94 return sizeof(*this) + (
nnz(U)+
nnz(L))*
sizeof(value_type);
99 template<
typename Matrix>
template<
typename M>
101 typedef value_type T;
102 typedef typename number_traits<T>::magnitude_type R;
104 size_type n = mat_nrows(A);
107 std::vector<T> indiag(n);
109 std::vector<size_type> ipvt(n), ipvtinv(n);
110 for (size_type i = 0; i < n; ++i) ipvt[i] = ipvtinv[i] = i;
111 indperm = unsorted_sub_index(ipvt);
112 indperminv = unsorted_sub_index(ipvtinv);
118 R prec = default_tol(R());
119 R max_pivot = gmm::abs(A(0,0)) * prec;
121 for (size_type i = 0; i < n; ++i) {
123 copy(sub_vector(mat_const_row(A, i), indperm), w);
126 typename _wsvector::iterator wkold = w.end();
127 for (
typename _wsvector::iterator wk = w.begin();
128 wk != w.end() && wk->first < i; ) {
129 size_type k = wk->first;
130 tmp = (wk->second) * indiag[k];
131 if (gmm::abs(tmp) < eps * norm_row) w.erase(k);
132 else { wk->second += tmp; gmm::add(scaled(mat_row(U, k), -tmp), w); }
133 if (wkold == w.end()) wk = w.begin();
else { wk = wkold; ++wk; }
134 if (wk != w.end() && wk->first == k)
135 {
if (wkold == w.end()) wkold = w.begin();
else ++wkold; ++wk; }
138 gmm::clean(w, eps * norm_row);
141 std::sort(ww.begin(), ww.end(), elt_rsvector_value_less_<T>());
142 typename _rsvector::const_iterator wit = ww.begin(), wite = ww.end();
145 for (; wit != wite; ++wit)
146 if (wit->c >= i) { ip = wit->c; tmp = wit->e;
break; }
147 if (ip ==
size_type(-1) || gmm::abs(tmp) <= max_pivot)
148 { GMM_WARNING2(
"pivot " << i <<
" too small"); ip=i; ww[i]=tmp=T(1); }
149 max_pivot = std::max(max_pivot, std::min(gmm::abs(tmp) * prec, R(1)));
150 indiag[i] = T(1) / tmp;
153 size_type nnl = 0, nnu = 0;
154 L[i].base_resize(K); U[i].base_resize(K+1);
155 typename _rsvector::iterator witL = L[i].begin(), witU = U[i].begin();
156 for (; wit != wite; ++wit) {
157 if (wit->c < i) {
if (nnl < K) { *witL++ = *wit; ++nnl; } }
158 else if (nnu < K || wit->c == i)
159 { CU(i, wit->c) = wit->e; *witU++ = *wit; ++nnu; }
161 L[i].base_resize(nnl); U[i].base_resize(nnu);
162 std::sort(L[i].begin(), L[i].end());
163 std::sort(U[i].begin(), U[i].end());
166 typename _wsvector::const_iterator iti = CU.col(i).begin();
167 typename _wsvector::const_iterator itie = CU.col(i).end();
168 typename _wsvector::const_iterator itp = CU.col(ip).begin();
169 typename _wsvector::const_iterator itpe = CU.col(ip).end();
171 while (iti != itie && itp != itpe) {
172 if (iti->first < itp->first)
173 { U.row(iti->first).swap_indices(i, ip); ++iti; }
174 else if (iti->first > itp->first)
175 { U.row(itp->first).swap_indices(i,ip);++itp; }
177 { U.row(iti->first).swap_indices(i, ip); ++iti; ++itp; }
180 for( ; iti != itie; ++iti) U.row(iti->first).swap_indices(i, ip);
181 for( ; itp != itpe; ++itp) U.row(itp->first).swap_indices(i, ip);
186 indperminv.swap(ipvt[i], ipvt[ip]);
187 std::swap(ipvtinv[ipvt[i]], ipvtinv[ipvt[ip]]);
188 std::swap(ipvt[i], ipvt[ip]);
193 template<
typename Matrix>
195 do_ilutp(gmm::transposed(A), row_major());
199 template <
typename Matrix,
typename V1,
typename V2>
inline 202 gmm::copy(gmm::sub_vector(v1, P.indperm), v2);
203 gmm::lower_tri_solve(gmm::transposed(P.U), v2,
false);
204 gmm::upper_tri_solve(gmm::transposed(P.L), v2,
true);
207 gmm::copy(v1, P.temporary);
208 gmm::lower_tri_solve(P.L, P.temporary,
true);
209 gmm::upper_tri_solve(P.U, P.temporary,
false);
210 gmm::copy(gmm::sub_vector(P.temporary, P.indperminv), v2);
214 template <
typename Matrix,
typename V1,
typename V2>
inline 217 gmm::copy(v1, P.temporary);
218 gmm::lower_tri_solve(P.L, P.temporary,
true);
219 gmm::upper_tri_solve(P.U, P.temporary,
false);
220 gmm::copy(gmm::sub_vector(P.temporary, P.indperminv), v2);
223 gmm::copy(gmm::sub_vector(v1, P.indperm), v2);
224 gmm::lower_tri_solve(gmm::transposed(P.U), v2,
false);
225 gmm::upper_tri_solve(gmm::transposed(P.L), v2,
true);
229 template <
typename Matrix,
typename V1,
typename V2>
inline 232 gmm::copy(gmm::sub_vector(v1, P.indperm), v2);
233 gmm::lower_tri_solve(gmm::transposed(P.U), v2,
false);
237 gmm::lower_tri_solve(P.L, v2,
true);
241 template <
typename Matrix,
typename V1,
typename V2>
inline 245 gmm::upper_tri_solve(gmm::transposed(P.L), v2,
true);
248 copy(v1, P.temporary);
249 gmm::upper_tri_solve(P.U, P.temporary,
false);
250 gmm::copy(gmm::sub_vector(P.temporary, P.indperminv), v2);
254 template <
typename Matrix,
typename V1,
typename V2>
inline 258 copy(v1, P.temporary);
259 gmm::upper_tri_solve(P.U, P.temporary,
false);
260 gmm::copy(gmm::sub_vector(P.temporary, P.indperminv), v2);
264 gmm::upper_tri_solve(gmm::transposed(P.L), v2,
true);
268 template <
typename Matrix,
typename V1,
typename V2>
inline 273 gmm::lower_tri_solve(P.L, v2,
true);
276 gmm::copy(gmm::sub_vector(v1, P.indperm), v2);
277 gmm::lower_tri_solve(gmm::transposed(P.U), v2,
false);
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norm2(const V &v)
Euclidean norm of a vector.
size_t size_type
used as the common size type in the library
void copy(const L1 &l1, L2 &l2)
*/
sparse vector built upon std::vector.
void clear(L &l)
clear (fill with zeros) a vector or matrix.
ILUTP: Incomplete LU with threshold and K fill-in Preconditioner and column pivoting.
ILUT: Incomplete LU with threshold and K fill-in Preconditioner.
void resize(M &v, size_type m, size_type n)
*/
sparse vector built upon std::map.
size_type nnz(const L &l)
count the number of non-zero entries of a vector or matrix.