1 #ifndef VIENNACL_LINALG_DETAIL_ILUT_HPP_
2 #define VIENNACL_LINALG_DETAIL_ILUT_HPP_
55 double drop_tolerance = 1e-4,
56 bool with_level_scheduling =
false)
57 : entries_per_row_(entries_per_row),
58 drop_tolerance_(drop_tolerance),
59 use_level_scheduling_(with_level_scheduling) {}
64 drop_tolerance_ = tol;
80 unsigned int entries_per_row_;
81 double drop_tolerance_;
82 bool use_level_scheduling_;
93 template<
typename NumericT>
119 template<
typename IndexT,
typename NumericT>
121 IndexT
const * u_coords,
NumericT const * u_elements, IndexT u_size,
NumericT alpha,
122 IndexT * z_coords,
NumericT * z_elements)
130 if (index_w < w_size && index_u < u_size)
132 if (w_coords[index_w] < u_coords[index_u])
134 z_coords[index_z] = w_coords[index_w];
135 z_elements[index_z++] = w_elements[index_w++];
137 else if (w_coords[index_w] == u_coords[index_u])
139 z_coords[index_z] = w_coords[index_w];
140 z_elements[index_z++] = w_elements[index_w++] - alpha * u_elements[index_u++];
144 z_coords[index_z] = u_coords[index_u];
145 z_elements[index_z++] = - alpha * u_elements[index_u++];
148 else if (index_w == w_size && index_u < u_size)
150 z_coords[index_z] = u_coords[index_u];
151 z_elements[index_z++] = - alpha * u_elements[index_u++];
153 else if (index_w < w_size && index_u == u_size)
155 z_coords[index_z] = w_coords[index_w];
156 z_elements[index_z++] = w_elements[index_w++];
163 template<
typename SizeT,
typename NumericT>
167 NumericT abs_value = std::fabs(value);
171 std::size_t first_smaller_index = 0;
172 while (first_smaller_index < map.size() && std::fabs(map[first_smaller_index].second) > abs_value)
173 ++first_smaller_index;
175 std::pair<SizeT, NumericT> tmp(index, value);
176 for (std::size_t j=first_smaller_index; j<map.size(); ++j)
192 template<
typename NumericT>
198 assert(A.
size1() == L.
size1() && bool(
"Output matrix size mismatch") );
199 assert(A.
size1() == U.
size1() && bool(
"Output matrix size mismatch") );
209 std::vector<NumericT> diagonal_U(A.
size1());
211 NumericT const * elements_A = viennacl::linalg::host_based::detail::extract_raw_pointer<NumericT>(A.
handle());
212 unsigned int const * row_buffer_A = viennacl::linalg::host_based::detail::extract_raw_pointer<unsigned int>(A.
handle1());
213 unsigned int const * col_buffer_A = viennacl::linalg::host_based::detail::extract_raw_pointer<unsigned int>(A.
handle2());
215 NumericT * elements_L = viennacl::linalg::host_based::detail::extract_raw_pointer<NumericT>(L.
handle());
216 unsigned int * row_buffer_L = viennacl::linalg::host_based::detail::extract_raw_pointer<unsigned int>(L.
handle1()); row_buffer_L[0] = 0;
217 unsigned int * col_buffer_L = viennacl::linalg::host_based::detail::extract_raw_pointer<unsigned int>(L.
handle2());
219 NumericT * elements_U = viennacl::linalg::host_based::detail::extract_raw_pointer<NumericT>(U.
handle());
220 unsigned int * row_buffer_U = viennacl::linalg::host_based::detail::extract_raw_pointer<unsigned int>(U.
handle1()); row_buffer_U[0] = 0;
221 unsigned int * col_buffer_U = viennacl::linalg::host_based::detail::extract_raw_pointer<unsigned int>(U.
handle2());
223 std::vector<std::pair<unsigned int, NumericT> > sorted_entries_L(tag.
get_entries_per_row());
224 std::vector<std::pair<unsigned int, NumericT> > sorted_entries_U(tag.
get_entries_per_row());
228 std::fill(sorted_entries_L.begin(), sorted_entries_L.end(), std::pair<unsigned int, NumericT>(0,
NumericT(0)));
229 std::fill(sorted_entries_U.begin(), sorted_entries_U.end(), std::pair<unsigned int, NumericT>(0,
NumericT(0)));
235 for (
unsigned int j = row_buffer_A[i]; j < row_buffer_A[i+1]; ++j, ++k)
240 row_norm += entry * entry;
242 row_norm = std::sqrt(row_norm);
247 unsigned int current_col = (row_buffer_A[i+1] > row_buffer_A[i]) ? w_in->
col_indices_[k] : static_cast<unsigned int>(i);
248 while (current_col < i)
251 NumericT a_kk = diagonal_U[current_col];
257 if ( std::fabs(w_k_entry) > tau_i)
260 unsigned int row_U_begin = row_buffer_U[current_col];
261 unsigned int row_U_end = row_buffer_U[current_col + 1];
263 if (row_U_end > row_U_begin)
267 col_buffer_U + row_U_begin + 1, elements_U + row_U_begin + 1, (row_U_end - row_U_begin) - 1, w_k_entry,
276 for (
unsigned int r = 0; r < k; ++r)
281 for (
unsigned int r = k+1; r < w_in->
size_; ++r)
294 current_col = (k < w_in->
size_) ? w_in->
col_indices_[k] : static_cast<unsigned int>(i);
299 for (
unsigned int r = 0; r < w_in->
size_; ++r)
308 diagonal_U[i] = value;
309 if (value <= 0 && value >= 0)
311 std::cerr <<
"ViennaCL: FATAL ERROR in ILUT(): Diagonal entry computed to zero (" << value <<
") in row " << i <<
"!" << std::endl;
320 unsigned int offset_L = row_buffer_L[i];
321 std::sort(sorted_entries_L.begin(), sorted_entries_L.end());
323 if (std::fabs(sorted_entries_L[j].second) > 0)
325 col_buffer_L[offset_L] = sorted_entries_L[j].first;
326 elements_L[offset_L] = sorted_entries_L[j].second;
329 row_buffer_L[i+1] = offset_L;
331 unsigned int offset_U = row_buffer_U[i];
332 col_buffer_U[offset_U] =
static_cast<unsigned int>(i);
333 elements_U[offset_U] = diagonal_U[i];
335 std::sort(sorted_entries_U.begin(), sorted_entries_U.end());
337 if (std::fabs(sorted_entries_U[j].second) > 0)
339 col_buffer_U[offset_U] = sorted_entries_U[j].first;
340 elements_U[offset_U] = sorted_entries_U[j].second;
343 row_buffer_U[i+1] = offset_U;
351 template<
typename MatrixT>
354 typedef typename MatrixT::value_type NumericType;
365 template<
typename VectorT>
370 unsigned int const * row_buffer = viennacl::linalg::host_based::detail::extract_raw_pointer<unsigned int>(L_.
handle1());
371 unsigned int const * col_buffer = viennacl::linalg::host_based::detail::extract_raw_pointer<unsigned int>(L_.
handle2());
372 NumericType
const * elements = viennacl::linalg::host_based::detail::extract_raw_pointer<NumericType>(L_.
handle());
374 viennacl::linalg::host_based::detail::csr_inplace_solve<NumericType>(row_buffer, col_buffer, elements, vec, L_.
size2(),
unit_lower_tag());
377 unsigned int const * row_buffer = viennacl::linalg::host_based::detail::extract_raw_pointer<unsigned int>(U_.
handle1());
378 unsigned int const * col_buffer = viennacl::linalg::host_based::detail::extract_raw_pointer<unsigned int>(U_.
handle2());
379 NumericType
const * elements = viennacl::linalg::host_based::detail::extract_raw_pointer<NumericType>(U_.
handle());
381 viennacl::linalg::host_based::detail::csr_inplace_solve<NumericType>(row_buffer, col_buffer, elements, vec, U_.
size2(),
upper_tag());
386 void init(MatrixT
const & mat)
409 template<
typename NumericT,
unsigned int AlignmentV>
434 multifrontal_L_row_index_arrays_,
435 multifrontal_L_row_buffers_,
436 multifrontal_L_col_buffers_,
437 multifrontal_L_element_buffers_,
438 multifrontal_L_row_elimination_num_list_);
443 multifrontal_U_row_index_arrays_,
444 multifrontal_U_row_buffers_,
445 multifrontal_U_col_buffers_,
446 multifrontal_U_element_buffers_,
447 multifrontal_U_row_elimination_num_list_);
468 void init(MatrixType
const & mat)
496 multifrontal_U_diagonal_.resize(U_.
size1(),
false);
500 multifrontal_U_diagonal_,
501 multifrontal_L_row_index_arrays_,
502 multifrontal_L_row_buffers_,
503 multifrontal_L_col_buffers_,
504 multifrontal_L_element_buffers_,
505 multifrontal_L_row_elimination_num_list_);
509 multifrontal_U_diagonal_,
510 multifrontal_U_row_index_arrays_,
511 multifrontal_U_row_buffers_,
512 multifrontal_U_col_buffers_,
513 multifrontal_U_element_buffers_,
514 multifrontal_U_row_elimination_num_list_);
522 for (
typename std::list< viennacl::backend::mem_handle >::iterator it = multifrontal_L_row_index_arrays_.begin();
523 it != multifrontal_L_row_index_arrays_.end();
527 for (
typename std::list< viennacl::backend::mem_handle >::iterator it = multifrontal_L_row_buffers_.begin();
528 it != multifrontal_L_row_buffers_.end();
532 for (
typename std::list< viennacl::backend::mem_handle >::iterator it = multifrontal_L_col_buffers_.begin();
533 it != multifrontal_L_col_buffers_.end();
537 for (
typename std::list< viennacl::backend::mem_handle >::iterator it = multifrontal_L_element_buffers_.begin();
538 it != multifrontal_L_element_buffers_.end();
547 for (
typename std::list< viennacl::backend::mem_handle >::iterator it = multifrontal_U_row_index_arrays_.begin();
548 it != multifrontal_U_row_index_arrays_.end();
552 for (
typename std::list< viennacl::backend::mem_handle >::iterator it = multifrontal_U_row_buffers_.begin();
553 it != multifrontal_U_row_buffers_.end();
557 for (
typename std::list< viennacl::backend::mem_handle >::iterator it = multifrontal_U_col_buffers_.begin();
558 it != multifrontal_U_col_buffers_.end();
562 for (
typename std::list< viennacl::backend::mem_handle >::iterator it = multifrontal_U_element_buffers_.begin();
563 it != multifrontal_U_element_buffers_.end();
574 std::list<viennacl::backend::mem_handle> multifrontal_L_row_index_arrays_;
575 std::list<viennacl::backend::mem_handle> multifrontal_L_row_buffers_;
576 std::list<viennacl::backend::mem_handle> multifrontal_L_col_buffers_;
577 std::list<viennacl::backend::mem_handle> multifrontal_L_element_buffers_;
578 std::list<vcl_size_t > multifrontal_L_row_elimination_num_list_;
581 std::list<viennacl::backend::mem_handle> multifrontal_U_row_index_arrays_;
582 std::list<viennacl::backend::mem_handle> multifrontal_U_row_buffers_;
583 std::list<viennacl::backend::mem_handle> multifrontal_U_col_buffers_;
584 std::list<viennacl::backend::mem_handle> multifrontal_U_element_buffers_;
585 std::list<vcl_size_t > multifrontal_U_row_elimination_num_list_;
const vcl_size_t & size2() const
Returns the number of columns.
void fill(MatrixType &matrix, vcl_size_t row_index, vcl_size_t col_index, NumericT value)
Generic filler routine for setting an entry of a matrix to a particular value.
viennacl::vector_expression< const vector_base< T >, const vector_base< T >, op_element_binary< op_div > > element_div(vector_base< T > const &v1, vector_base< T > const &v2)
void inplace_solve(const matrix_base< NumericT > &A, matrix_base< NumericT > &B, SolverTagT)
Direct inplace solver for triangular systems with multiple right hand sides, i.e. A \ B (MATLAB notat...
IndexT merge_subtract_sparse_rows(IndexT const *w_coords, NumericT const *w_elements, IndexT w_size, IndexT const *u_coords, NumericT const *u_elements, IndexT u_size, NumericT alpha, IndexT *z_coords, NumericT *z_elements)
Subtracts a scaled sparse vector u from a sparse vector w and writes the output to z: z = w - alpha *...
ilut_precond(MatrixT const &mat, ilut_tag const &tag)
const vcl_size_t & size1() const
Returns the number of rows.
vcl_size_t size1(MatrixType const &mat)
Generic routine for obtaining the number of rows of a matrix (ViennaCL, uBLAS, etc.)
std::vector< unsigned int > col_indices_
void set_entries_per_row(unsigned int e)
void precondition(viennacl::compressed_matrix< NumericT > &A, ilu0_tag const &)
Implementation of a ILU-preconditioner with static pattern. Optimized version for CSR matrices...
This file provides the forward declarations for the main types used within ViennaCL.
bool use_level_scheduling() const
void level_scheduling_setup_U(viennacl::compressed_matrix< NumericT, AlignmentV > const &LU, viennacl::vector< NumericT > const &diagonal_LU, std::list< viennacl::backend::mem_handle > &row_index_arrays, std::list< viennacl::backend::mem_handle > &row_buffers, std::list< viennacl::backend::mem_handle > &col_buffers, std::list< viennacl::backend::mem_handle > &element_buffers, std::list< vcl_size_t > &row_elimination_num_list)
result_of::size_type< MatrixType >::type size2(MatrixType const &mat)
Generic routine for obtaining the number of columns of a matrix (ViennaCL, uBLAS, etc...
void level_scheduling_substitute(viennacl::vector< NumericT > &vec, std::list< viennacl::backend::mem_handle > const &row_index_arrays, std::list< viennacl::backend::mem_handle > const &row_buffers, std::list< viennacl::backend::mem_handle > const &col_buffers, std::list< viennacl::backend::mem_handle > const &element_buffers, std::list< vcl_size_t > const &row_elimination_num_list)
const handle_type & handle() const
Returns the OpenCL handle to the matrix entry array.
const handle_type & handle1() const
Returns the OpenCL handle to the row index array.
const vcl_size_t & nnz() const
Returns the number of nonzero entries.
Represents a generic 'context' similar to an OpenCL context, but is backend-agnostic and thus also su...
void apply(VectorT &vec) const
void resize_if_bigger(vcl_size_t s)
A tag class representing an upper triangular matrix.
A tag for incomplete LU factorization with threshold (ILUT)
Implementation of the compressed_matrix class.
ilut_precond(MatrixType const &mat, ilut_tag const &tag)
unsigned int get_entries_per_row() const
const handle_type & handle2() const
Returns the OpenCL handle to the column index array.
ILUT preconditioner class, can be supplied to solve()-routines.
void insert_with_value_sort(std::vector< std::pair< SizeT, NumericT > > &map, SizeT index, NumericT value)
Common routines for single-threaded or OpenMP-enabled execution on CPU.
viennacl::enable_if< viennacl::is_scalar< ScalarT1 >::value &&viennacl::is_scalar< ScalarT2 >::value >::type swap(ScalarT1 &s1, ScalarT2 &s2)
Swaps the contents of two scalars, data is copied.
ilut_sparse_vector(vcl_size_t alloc_size=0)
void reserve(vcl_size_t new_nonzeros, bool preserve=true)
Allocate memory for the supplied number of nonzeros in the matrix. Old values are preserved...
void apply(viennacl::vector< NumericT > &vec) const
void set_drop_tolerance(double tol)
viennacl::context context(T const &t)
Returns an ID for the currently active memory domain of an object.
void copy(std::vector< NumericT > &cpu_vec, circulant_matrix< NumericT, AlignmentV > &gpu_mat)
Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...
void use_level_scheduling(bool b)
std::vector< NumericT > elements_
A tag class representing a lower triangular matrix with unit diagonal.
Helper struct for holding a sparse vector in linear memory. For internal use only.
void row_info(compressed_matrix< NumericT, AlignmentV > const &mat, vector_base< NumericT > &vec, viennacl::linalg::detail::row_info_types info_selector)
void level_scheduling_setup_L(viennacl::compressed_matrix< NumericT, AlignmentV > const &LU, viennacl::vector< NumericT > const &diagonal_LU, std::list< viennacl::backend::mem_handle > &row_index_arrays, std::list< viennacl::backend::mem_handle > &row_buffers, std::list< viennacl::backend::mem_handle > &col_buffers, std::list< viennacl::backend::mem_handle > &element_buffers, std::list< vcl_size_t > &row_elimination_num_list)
Common routines used within ILU-type preconditioners.
const handle_type & handle() const
Returns the memory handle.
ilut_tag(unsigned int entries_per_row=20, double drop_tolerance=1e-4, bool with_level_scheduling=false)
The constructor.
memory_types get_active_handle_id() const
Returns an ID for the currently active memory buffer. Other memory buffers might contain old or no da...
double get_drop_tolerance() const
void switch_memory_context(T &obj, viennacl::context new_ctx)
Generic convenience routine for migrating data of an object to a new memory domain.