35#ifdef DEAL_II_WITH_UMFPACK
49 template <
typename SparseMatrixType>
51 parallel_grainsize(
const SparseMatrixType &matrix)
53 const unsigned int avg_entries_per_row =
55 return std::max(1000 / avg_entries_per_row, 1u);
71#ifdef DEAL_II_WITH_UMFPACK
76 , symbolic_decomposition(nullptr)
77 , numeric_decomposition(nullptr)
78 , control(UMFPACK_CONTROL)
80 umfpack_dl_defaults(
control.data());
102 std::vector<types::suitesparse_index> tmp;
107 std::vector<types::suitesparse_index> tmp;
112 std::vector<double> tmp;
117 std::vector<double> tmp;
121 umfpack_dl_defaults(
control.data());
126template <
typename number>
141 for (size_type row = row_begin; row < row_end; ++row)
154 long int cursor = Ap[row];
155 while ((cursor < Ap[row + 1] - 1) && (Ai[cursor] > Ai[cursor + 1]))
157 std::swap(Ai[cursor], Ai[cursor + 1]);
159 std::swap(Ax[cursor], Ax[cursor + 1]);
160 if (numbers::NumberTraits<number>::is_complex == true)
161 std::swap(Az[cursor], Az[cursor + 1]);
167 parallel_grainsize(matrix));
172template <
typename number>
181 for (size_type row = row_begin; row < row_end; ++row)
183 long int cursor = Ap[row];
184 while ((cursor < Ap[row + 1] - 1) && (Ai[cursor] > Ai[cursor + 1]))
186 std::swap(Ai[cursor], Ai[cursor + 1]);
188 std::swap(Ax[cursor], Ax[cursor + 1]);
189 if (numbers::NumberTraits<number>::is_complex == true)
190 std::swap(Az[cursor], Az[cursor + 1]);
196 parallel_grainsize(matrix));
201template <
typename number>
214 for (size_type row = row_begin; row < row_end; ++row)
216 long int cursor = Ap[row];
217 for (size_type block = 0; block < matrix.n_block_cols(); ++block)
220 while ((cursor < Ap[row + 1] - 1) &&
221 (Ai[cursor] < Ai[cursor + 1]))
225 if (cursor == Ap[row + 1] - 1)
230 long int element = cursor;
231 while ((element < Ap[row + 1] - 1) &&
232 (Ai[element] > Ai[element + 1]))
234 std::swap(Ai[element], Ai[element + 1]);
236 std::swap(Ax[element], Ax[element + 1]);
237 if (numbers::NumberTraits<number>::is_complex == true)
238 std::swap(Az[element], Az[element + 1]);
245 parallel_grainsize(matrix));
250template <
class Matrix>
258 using number =
typename Matrix::value_type;
284 Ai.resize(matrix.n_nonzero_elements());
285 Ax.resize(matrix.n_nonzero_elements());
287 Az.resize(matrix.n_nonzero_elements());
292 Ap[row] =
Ap[row - 1] + matrix.get_row_length(row - 1);
302 for (size_type row = row_begin; row < row_end; ++row)
304 long int index = Ap[row];
305 for (typename Matrix::const_iterator p = matrix.begin(row);
306 p != matrix.end(row);
310 Ai[index] = p->column();
311 Ax[index] = std::real(p->value());
312 if (numbers::NumberTraits<number>::is_complex == true)
313 Az[index] = std::imag(p->value());
318 Assert(index == Ap[row + 1], ExcInternalError());
321 parallel_grainsize(matrix));
330 status = umfpack_dl_symbolic(N,
335 &symbolic_decomposition,
339 status = umfpack_zl_symbolic(N,
345 &symbolic_decomposition,
349 ExcUMFPACKError(
"umfpack_dl_symbolic", status));
352 status = umfpack_dl_numeric(Ap.data(),
355 symbolic_decomposition,
356 &numeric_decomposition,
360 status = umfpack_zl_numeric(Ap.data(),
364 symbolic_decomposition,
365 &numeric_decomposition,
369 ExcUMFPACKError(
"umfpack_dl_numeric", status));
371 umfpack_dl_free_symbolic(&symbolic_decomposition);
386 ExcMessage(
"You have previously factored a matrix using this class "
387 "that had complex-valued entries. This then requires "
388 "applying the factored matrix to a complex-valued "
389 "vector, but you are only providing a real-valued vector "
393 rhs = rhs_and_solution;
401 const int status = umfpack_dl_solve(
transpose ? UMFPACK_A : UMFPACK_At,
405 rhs_and_solution.
begin(),
420# ifdef DEAL_II_WITH_COMPLEX_VALUES
450 for (
unsigned int i = 0; i < rhs_and_solution.size(); ++i)
452 rhs_re(i) = std::real(rhs_and_solution(i));
453 rhs_im(i) = std::imag(rhs_and_solution(i));
468 const int status = umfpack_zl_solve(
transpose ? UMFPACK_A : UMFPACK_Aat,
484 for (
unsigned int i = 0; i < rhs_and_solution.size(); ++i)
485 rhs_and_solution(i) = {solution_re(i), solution_im(i)};
500 for (
unsigned int i = 0; i < rhs.
size(); ++i)
501 rhs_real_or_imag(i) = std::real(rhs(i));
505 rhs_and_solution = rhs_real_or_imag;
510 for (
unsigned int i = 0; i < rhs.
size(); ++i)
511 rhs_real_or_imag(i) = std::imag(rhs(i));
515 for (
unsigned int i = 0; i < rhs.
size(); ++i)
516 rhs_and_solution(i).imag(rhs_real_or_imag(i));
521 (void)rhs_and_solution;
525 "This function can't be called if deal.II has been configured "
526 "with DEAL_II_WITH_COMPLEX_VALUES=FALSE."));
539 tmp = rhs_and_solution;
541 rhs_and_solution = tmp;
550# ifdef DEAL_II_WITH_COMPLEX_VALUES
555 tmp = rhs_and_solution;
557 rhs_and_solution = tmp;
560 (void)rhs_and_solution;
564 "This function can't be called if deal.II has been configured "
565 "with DEAL_II_WITH_COMPLEX_VALUES=FALSE."));
571template <
class Matrix>
583template <
class Matrix>
586 Vector<std::complex<double>> &rhs_and_solution,
589# ifdef DEAL_II_WITH_COMPLEX_VALUES
596 (void)rhs_and_solution;
600 "This function can't be called if deal.II has been configured "
601 "with DEAL_II_WITH_COMPLEX_VALUES=FALSE."));
607template <
class Matrix>
619template <
class Matrix>
622 BlockVector<std::complex<double>> &rhs_and_solution,
625# ifdef DEAL_II_WITH_COMPLEX_VALUES
632 (void)rhs_and_solution;
636 "This function can't be called if deal.II has been configured "
637 "with DEAL_II_WITH_COMPLEX_VALUES=FALSE."));
648 , symbolic_decomposition(nullptr)
649 , numeric_decomposition(nullptr)
659template <
class Matrix>
666 "To call this function you need UMFPACK, but you configured deal.II "
667 "without passing the necessary switch to 'cmake'. Please consult the "
668 "installation instructions at https://dealii.org/current/readme.html"));
678 "To call this function you need UMFPACK, but you configured deal.II "
679 "without passing the necessary switch to 'cmake'. Please consult the "
680 "installation instructions at https://dealii.org/current/readme.html"));
691 "To call this function you need UMFPACK, but you configured deal.II "
692 "without passing the necessary switch to 'cmake'. Please consult the "
693 "installation instructions at https://dealii.org/current/readme.html"));
704 "To call this function you need UMFPACK, but you configured deal.II "
705 "without passing the necessary switch to 'cmake'. Please consult the "
706 "installation instructions at https://dealii.org/current/readme.html"));
718 "To call this function you need UMFPACK, but you configured deal.II "
719 "without passing the necessary switch to 'cmake'. Please consult the "
720 "installation instructions at https://dealii.org/current/readme.html"));
725template <
class Matrix>
732 "To call this function you need UMFPACK, but you configured deal.II "
733 "without passing the necessary switch to 'cmake'. Please consult the "
734 "installation instructions at https://dealii.org/current/readme.html"));
739template <
class Matrix>
742 Vector<std::complex<double>> &,
748 "To call this function you need UMFPACK, but you configured deal.II "
749 "without passing the necessary switch to 'cmake'. Please consult the "
750 "installation instructions at https://dealii.org/current/readme.html"));
755template <
class Matrix>
762 "To call this function you need UMFPACK, but you configured deal.II "
763 "without passing the necessary switch to 'cmake'. Please consult the "
764 "installation instructions at https://dealii.org/current/readme.html"));
769template <
class Matrix>
778 "To call this function you need UMFPACK, but you configured deal.II "
779 "without passing the necessary switch to 'cmake'. Please consult the "
780 "installation instructions at https://dealii.org/current/readme.html"));
786template <
class Matrix>
817 this->
solve(dst,
true);
827 this->
solve(dst,
true);
846#define InstantiateUMFPACK(MatrixType) \
847 template void SparseDirectUMFPACK::factorize(const MatrixType &); \
848 template void SparseDirectUMFPACK::solve(const MatrixType &, \
851 template void SparseDirectUMFPACK::solve(const MatrixType &, \
852 Vector<std::complex<double>> &, \
854 template void SparseDirectUMFPACK::solve(const MatrixType &, \
855 BlockVector<double> &, \
857 template void SparseDirectUMFPACK::solve( \
858 const MatrixType &, BlockVector<std::complex<double>> &, const bool); \
859 template void SparseDirectUMFPACK::initialize(const MatrixType &, \
860 const AdditionalData)
871#ifdef DEAL_II_WITH_COMPLEX_VALUES
virtual size_type size() const override
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
virtual size_type size() const override
#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)