36#ifdef DEAL_II_WITH_UMFPACK
50 template <
typename SparseMatrixType>
52 parallel_grainsize(
const SparseMatrixType &matrix)
54 const unsigned int avg_entries_per_row =
56 return std::max(1000 / avg_entries_per_row, 1u);
72#ifdef DEAL_II_WITH_UMFPACK
77 , symbolic_decomposition(nullptr)
78 , numeric_decomposition(nullptr)
79 , control(UMFPACK_CONTROL)
81 umfpack_dl_defaults(
control.data());
103 std::vector<types::suitesparse_index> tmp;
108 std::vector<types::suitesparse_index> tmp;
113 std::vector<double> tmp;
118 std::vector<double> tmp;
122 umfpack_dl_defaults(
control.data());
127template <
typename number>
142 for (size_type row = row_begin; row < row_end; ++row)
155 long int cursor = Ap[row];
156 while ((cursor < Ap[row + 1] - 1) && (Ai[cursor] > Ai[cursor + 1]))
158 std::swap(Ai[cursor], Ai[cursor + 1]);
160 std::swap(Ax[cursor], Ax[cursor + 1]);
161 if (numbers::NumberTraits<number>::is_complex == true)
162 std::swap(Az[cursor], Az[cursor + 1]);
168 parallel_grainsize(matrix));
173template <
typename number>
182 for (size_type row = row_begin; row < row_end; ++row)
184 long int cursor = Ap[row];
185 while ((cursor < Ap[row + 1] - 1) && (Ai[cursor] > Ai[cursor + 1]))
187 std::swap(Ai[cursor], Ai[cursor + 1]);
189 std::swap(Ax[cursor], Ax[cursor + 1]);
190 if (numbers::NumberTraits<number>::is_complex == true)
191 std::swap(Az[cursor], Az[cursor + 1]);
197 parallel_grainsize(matrix));
202template <
typename number>
215 for (size_type row = row_begin; row < row_end; ++row)
217 long int cursor = Ap[row];
218 for (size_type block = 0; block < matrix.n_block_cols(); ++block)
221 while ((cursor < Ap[row + 1] - 1) &&
222 (Ai[cursor] < Ai[cursor + 1]))
226 if (cursor == Ap[row + 1] - 1)
231 long int element = cursor;
232 while ((element < Ap[row + 1] - 1) &&
233 (Ai[element] > Ai[element + 1]))
235 std::swap(Ai[element], Ai[element + 1]);
237 std::swap(Ax[element], Ax[element + 1]);
238 if (numbers::NumberTraits<number>::is_complex == true)
239 std::swap(Az[element], Az[element + 1]);
246 parallel_grainsize(matrix));
251template <
class Matrix>
259 using number =
typename Matrix::value_type;
285 Ai.resize(matrix.n_nonzero_elements());
286 Ax.resize(matrix.n_nonzero_elements());
288 Az.resize(matrix.n_nonzero_elements());
293 Ap[row] =
Ap[row - 1] + matrix.get_row_length(row - 1);
303 for (size_type row = row_begin; row < row_end; ++row)
305 long int index = Ap[row];
306 for (typename Matrix::const_iterator p = matrix.begin(row);
307 p != matrix.end(row);
311 Ai[index] = p->column();
312 Ax[index] = std::real(p->value());
313 if (numbers::NumberTraits<number>::is_complex == true)
314 Az[index] = std::imag(p->value());
319 Assert(index == Ap[row + 1], ExcInternalError());
322 parallel_grainsize(matrix));
331 status = umfpack_dl_symbolic(N,
336 &symbolic_decomposition,
340 status = umfpack_zl_symbolic(N,
346 &symbolic_decomposition,
350 ExcUMFPACKError(
"umfpack_dl_symbolic", status));
353 status = umfpack_dl_numeric(Ap.data(),
356 symbolic_decomposition,
357 &numeric_decomposition,
361 status = umfpack_zl_numeric(Ap.data(),
365 symbolic_decomposition,
366 &numeric_decomposition,
370 ExcUMFPACKError(
"umfpack_dl_numeric", status));
372 umfpack_dl_free_symbolic(&symbolic_decomposition);
387 ExcMessage(
"You have previously factored a matrix using this class "
388 "that had complex-valued entries. This then requires "
389 "applying the factored matrix to a complex-valued "
390 "vector, but you are only providing a real-valued vector "
394 rhs = rhs_and_solution;
402 const int status = umfpack_dl_solve(
transpose ? UMFPACK_A : UMFPACK_At,
406 rhs_and_solution.
begin(),
421# ifdef DEAL_II_WITH_COMPLEX_VALUES
451 for (
unsigned int i = 0; i < rhs_and_solution.size(); ++i)
453 rhs_re(i) = std::real(rhs_and_solution(i));
454 rhs_im(i) = std::imag(rhs_and_solution(i));
469 const int status = umfpack_zl_solve(
transpose ? UMFPACK_A : UMFPACK_Aat,
485 for (
unsigned int i = 0; i < rhs_and_solution.size(); ++i)
486 rhs_and_solution(i) = {solution_re(i), solution_im(i)};
501 for (
unsigned int i = 0; i < rhs.
size(); ++i)
502 rhs_real_or_imag(i) = std::real(rhs(i));
506 rhs_and_solution = rhs_real_or_imag;
511 for (
unsigned int i = 0; i < rhs.
size(); ++i)
512 rhs_real_or_imag(i) = std::imag(rhs(i));
516 for (
unsigned int i = 0; i < rhs.
size(); ++i)
517 rhs_and_solution(i).imag(rhs_real_or_imag(i));
522 (void)rhs_and_solution;
526 "This function can't be called if deal.II has been configured "
527 "with DEAL_II_WITH_COMPLEX_VALUES=FALSE."));
540 tmp = rhs_and_solution;
542 rhs_and_solution = tmp;
551# ifdef DEAL_II_WITH_COMPLEX_VALUES
556 tmp = rhs_and_solution;
558 rhs_and_solution = tmp;
561 (void)rhs_and_solution;
565 "This function can't be called if deal.II has been configured "
566 "with DEAL_II_WITH_COMPLEX_VALUES=FALSE."));
572template <
class Matrix>
584template <
class Matrix>
587 Vector<std::complex<double>> &rhs_and_solution,
590# ifdef DEAL_II_WITH_COMPLEX_VALUES
597 (void)rhs_and_solution;
601 "This function can't be called if deal.II has been configured "
602 "with DEAL_II_WITH_COMPLEX_VALUES=FALSE."));
608template <
class Matrix>
620template <
class Matrix>
623 BlockVector<std::complex<double>> &rhs_and_solution,
626# ifdef DEAL_II_WITH_COMPLEX_VALUES
633 (void)rhs_and_solution;
637 "This function can't be called if deal.II has been configured "
638 "with DEAL_II_WITH_COMPLEX_VALUES=FALSE."));
649 , symbolic_decomposition(nullptr)
650 , numeric_decomposition(nullptr)
660template <
class Matrix>
667 "To call this function you need UMFPACK, but you configured deal.II "
668 "without passing the necessary switch to 'cmake'. Please consult the "
669 "installation instructions at https://dealii.org/current/readme.html"));
679 "To call this function you need UMFPACK, but you configured deal.II "
680 "without passing the necessary switch to 'cmake'. Please consult the "
681 "installation instructions at https://dealii.org/current/readme.html"));
692 "To call this function you need UMFPACK, but you configured deal.II "
693 "without passing the necessary switch to 'cmake'. Please consult the "
694 "installation instructions at https://dealii.org/current/readme.html"));
705 "To call this function you need UMFPACK, but you configured deal.II "
706 "without passing the necessary switch to 'cmake'. Please consult the "
707 "installation instructions at https://dealii.org/current/readme.html"));
719 "To call this function you need UMFPACK, but you configured deal.II "
720 "without passing the necessary switch to 'cmake'. Please consult the "
721 "installation instructions at https://dealii.org/current/readme.html"));
726template <
class Matrix>
733 "To call this function you need UMFPACK, but you configured deal.II "
734 "without passing the necessary switch to 'cmake'. Please consult the "
735 "installation instructions at https://dealii.org/current/readme.html"));
740template <
class Matrix>
743 Vector<std::complex<double>> &,
749 "To call this function you need UMFPACK, but you configured deal.II "
750 "without passing the necessary switch to 'cmake'. Please consult the "
751 "installation instructions at https://dealii.org/current/readme.html"));
756template <
class Matrix>
763 "To call this function you need UMFPACK, but you configured deal.II "
764 "without passing the necessary switch to 'cmake'. Please consult the "
765 "installation instructions at https://dealii.org/current/readme.html"));
770template <
class Matrix>
779 "To call this function you need UMFPACK, but you configured deal.II "
780 "without passing the necessary switch to 'cmake'. Please consult the "
781 "installation instructions at https://dealii.org/current/readme.html"));
787template <
class Matrix>
818 this->
solve(dst,
true);
828 this->
solve(dst,
true);
847#define InstantiateUMFPACK(MatrixType) \
848 template void SparseDirectUMFPACK::factorize(const MatrixType &); \
849 template void SparseDirectUMFPACK::solve(const MatrixType &, \
852 template void SparseDirectUMFPACK::solve(const MatrixType &, \
853 Vector<std::complex<double>> &, \
855 template void SparseDirectUMFPACK::solve(const MatrixType &, \
856 BlockVector<double> &, \
858 template void SparseDirectUMFPACK::solve( \
859 const MatrixType &, BlockVector<std::complex<double>> &, const bool); \
860 template void SparseDirectUMFPACK::initialize(const MatrixType &, \
861 const AdditionalData)
872#ifdef DEAL_II_WITH_COMPLEX_VALUES
void * numeric_decomposition
~SparseDirectUMFPACK() override
void initialize(const SparsityPattern &sparsity_pattern)
void Tvmult(Vector< double > &dst, const Vector< double > &src) const
void * symbolic_decomposition
void solve(Vector< double > &rhs_and_solution, const bool transpose=false) const
void sort_arrays(const SparseMatrixEZ< number > &)
void factorize(const Matrix &matrix)
std::vector< double > control
void vmult(Vector< double > &dst, const Vector< double > &src) const
std::vector< types::suitesparse_index > Ap
std::vector< types::suitesparse_index > Ai
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcUMFPACKError(std::string arg1, int arg2)
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcNotQuadratic()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
@ matrix
Contents is actually a matrix.
void apply_to_subranges(const Iterator &begin, const std_cxx20::type_identity_t< Iterator > &end, const Function &f, const unsigned int grainsize)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
#define InstantiateUMFPACK(MatrixType)