36#ifdef DEAL_II_WITH_UMFPACK
53#ifdef DEAL_II_WITH_UMFPACK
58 , symbolic_decomposition(nullptr)
59 , numeric_decomposition(nullptr)
60 , control(UMFPACK_CONTROL)
62 umfpack_dl_defaults(
control.data());
84 std::vector<long int> tmp;
89 std::vector<long int> tmp;
94 std::vector<double> tmp;
99 std::vector<double> tmp;
103 umfpack_dl_defaults(
control.data());
108template <
typename number>
132 long int cursor =
Ap[row];
133 while ((cursor <
Ap[row + 1] - 1) && (
Ai[cursor] >
Ai[cursor + 1]))
148template <
typename number>
155 long int cursor =
Ap[row];
156 while ((cursor <
Ap[row + 1] - 1) && (
Ai[cursor] >
Ai[cursor + 1]))
171template <
typename number>
182 long int cursor =
Ap[row];
186 while ((cursor <
Ap[row + 1] - 1) && (
Ai[cursor] <
Ai[cursor + 1]))
190 if (cursor ==
Ap[row + 1] - 1)
195 long int element = cursor;
196 while ((element <
Ap[row + 1] - 1) && (
Ai[element] >
Ai[element + 1]))
212template <
class Matrix>
220 using number =
typename Matrix::value_type;
246 Ai.resize(
matrix.n_nonzero_elements());
247 Ax.resize(
matrix.n_nonzero_elements());
249 Az.resize(
matrix.n_nonzero_elements());
254 Ap[row] =
Ap[row - 1] +
matrix.get_row_length(row - 1);
264 std::vector<long int> row_pointers =
Ap;
270 for (
typename Matrix::const_iterator p =
matrix.begin(row);
275 Ai[row_pointers[row]] = p->column();
276 Ax[row_pointers[row]] = std::real(p->value());
278 Az[row_pointers[row]] = std::imag(p->value());
297 status = umfpack_dl_symbolic(
N,
306 status = umfpack_zl_symbolic(
N,
319 status = umfpack_dl_numeric(
Ap.data(),
327 status = umfpack_zl_numeric(
Ap.data(),
353 ExcMessage(
"You have previously factored a matrix using this class "
354 "that had complex-valued entries. This then requires "
355 "applying the factored matrix to a complex-valued "
356 "vector, but you are only providing a real-valued vector "
360 rhs = rhs_and_solution;
368 const int status = umfpack_dl_solve(
transpose ? UMFPACK_A : UMFPACK_At,
372 rhs_and_solution.
begin(),
387# ifdef DEAL_II_WITH_COMPLEX_VALUES
417 for (
unsigned int i = 0; i < rhs_and_solution.size(); ++i)
419 rhs_re(i) = std::real(rhs_and_solution(i));
420 rhs_im(i) = std::imag(rhs_and_solution(i));
435 const int status = umfpack_zl_solve(
transpose ? UMFPACK_A : UMFPACK_Aat,
451 for (
unsigned int i = 0; i < rhs_and_solution.size(); ++i)
452 rhs_and_solution(i) = {solution_re(i), solution_im(i)};
467 for (
unsigned int i = 0; i < rhs.
size(); ++i)
468 rhs_real_or_imag(i) = std::real(rhs(i));
472 rhs_and_solution = rhs_real_or_imag;
477 for (
unsigned int i = 0; i < rhs.
size(); ++i)
478 rhs_real_or_imag(i) = std::imag(rhs(i));
482 for (
unsigned int i = 0; i < rhs.
size(); ++i)
483 rhs_and_solution(i).imag(rhs_real_or_imag(i));
488 (void)rhs_and_solution;
492 "This function can't be called if deal.II has been configured "
493 "with DEAL_II_WITH_COMPLEX_VALUES=FALSE."));
506 tmp = rhs_and_solution;
508 rhs_and_solution = tmp;
517# ifdef DEAL_II_WITH_COMPLEX_VALUES
522 tmp = rhs_and_solution;
524 rhs_and_solution = tmp;
527 (void)rhs_and_solution;
531 "This function can't be called if deal.II has been configured "
532 "with DEAL_II_WITH_COMPLEX_VALUES=FALSE."));
538template <
class Matrix>
550template <
class Matrix>
553 Vector<std::complex<double>> &rhs_and_solution,
556# ifdef DEAL_II_WITH_COMPLEX_VALUES
563 (void)rhs_and_solution;
567 "This function can't be called if deal.II has been configured "
568 "with DEAL_II_WITH_COMPLEX_VALUES=FALSE."));
574template <
class Matrix>
586template <
class Matrix>
589 BlockVector<std::complex<double>> &rhs_and_solution,
592# ifdef DEAL_II_WITH_COMPLEX_VALUES
599 (void)rhs_and_solution;
603 "This function can't be called if deal.II has been configured "
604 "with DEAL_II_WITH_COMPLEX_VALUES=FALSE."));
615 , symbolic_decomposition(nullptr)
616 , numeric_decomposition(nullptr)
626template <
class Matrix>
633 "To call this function you need UMFPACK, but you configured deal.II "
634 "without passing the necessary switch to 'cmake'. Please consult the "
635 "installation instructions in doc/readme.html."));
645 "To call this function you need UMFPACK, but you configured deal.II "
646 "without passing the necessary switch to 'cmake'. Please consult the "
647 "installation instructions in doc/readme.html."));
658 "To call this function you need UMFPACK, but you configured deal.II "
659 "without passing the necessary switch to 'cmake'. Please consult the "
660 "installation instructions in doc/readme.html."));
671 "To call this function you need UMFPACK, but you configured deal.II "
672 "without passing the necessary switch to 'cmake'. Please consult the "
673 "installation instructions in doc/readme.html."));
685 "To call this function you need UMFPACK, but you configured deal.II "
686 "without passing the necessary switch to 'cmake'. Please consult the "
687 "installation instructions in doc/readme.html."));
692template <
class Matrix>
699 "To call this function you need UMFPACK, but you configured deal.II "
700 "without passing the necessary switch to 'cmake'. Please consult the "
701 "installation instructions in doc/readme.html."));
706template <
class Matrix>
709 Vector<std::complex<double>> &,
715 "To call this function you need UMFPACK, but you configured deal.II "
716 "without passing the necessary switch to 'cmake'. Please consult the "
717 "installation instructions in doc/readme.html."));
722template <
class Matrix>
729 "To call this function you need UMFPACK, but you configured deal.II "
730 "without passing the necessary switch to 'cmake'. Please consult the "
731 "installation instructions in doc/readme.html."));
736template <
class Matrix>
745 "To call this function you need UMFPACK, but you configured deal.II "
746 "without passing the necessary switch to 'cmake'. Please consult the "
747 "installation instructions in doc/readme.html."));
753template <
class Matrix>
784 this->
solve(dst,
true);
794 this->
solve(dst,
true);
813#define InstantiateUMFPACK(MatrixType) \
814 template void SparseDirectUMFPACK::factorize(const MatrixType &); \
815 template void SparseDirectUMFPACK::solve(const MatrixType &, \
818 template void SparseDirectUMFPACK::solve(const MatrixType &, \
819 Vector<std::complex<double>> &, \
821 template void SparseDirectUMFPACK::solve(const MatrixType &, \
822 BlockVector<double> &, \
824 template void SparseDirectUMFPACK::solve( \
825 const MatrixType &, BlockVector<std::complex<double>> &, const bool); \
826 template void SparseDirectUMFPACK::initialize(const MatrixType &, \
827 const AdditionalData)
838#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 swap(MemorySpaceData< Number, MemorySpace > &, MemorySpaceData< Number, MemorySpace > &)
#define InstantiateUMFPACK(MatrixType)