16 #include <deal.II/lac/sparse_direct.h> 17 #include <deal.II/base/memory_consumption.h> 18 #include <deal.II/base/thread_management.h> 19 #include <deal.II/lac/sparse_matrix.h> 20 #include <deal.II/lac/block_sparse_matrix.h> 21 #include <deal.II/lac/vector.h> 30 DEAL_II_NAMESPACE_OPEN
34 #ifdef DEAL_II_WITH_UMFPACK 52 #ifdef DEAL_II_WITH_UMFPACK 58 symbolic_decomposition (0),
59 numeric_decomposition (0),
60 control (UMFPACK_CONTROL)
62 umfpack_dl_defaults (&
control[0]);
77 if (numeric_decomposition != 0)
79 umfpack_dl_free_numeric (&numeric_decomposition);
80 numeric_decomposition = 0;
84 std::vector<long int> tmp;
89 std::vector<long int> tmp;
94 std::vector<double> tmp;
98 umfpack_dl_defaults (&
control[0]);
103 template <
typename number>
115 for (
size_type row=0; row<matrix.m(); ++row)
128 long int cursor =
Ap[row];
129 while ((cursor <
Ap[row+1]-1) &&
130 (Ai[cursor] > Ai[cursor+1]))
132 std::swap (Ai[cursor], Ai[cursor+1]);
133 std::swap (Ax[cursor], Ax[cursor+1]);
141 template <
typename number>
147 for (
size_type row=0; row<matrix.m(); ++row)
149 long int cursor =
Ap[row];
150 while ((cursor <
Ap[row+1]-1) &&
151 (Ai[cursor] > Ai[cursor+1]))
153 std::swap (Ai[cursor], Ai[cursor+1]);
154 std::swap (Ax[cursor], Ax[cursor+1]);
162 template <
typename number>
172 for (
size_type row=0; row<matrix.m(); ++row)
174 long int cursor =
Ap[row];
175 for (
size_type block=0; block<matrix.n_block_cols(); ++block)
178 while ((cursor <
Ap[row+1]-1) &&
179 (Ai[cursor] < Ai[cursor+1]))
183 if (cursor ==
Ap[row+1]-1)
188 long int element = cursor;
189 while ((element <
Ap[row+1]-1) &&
190 (Ai[element] > Ai[element+1]))
192 std::swap (Ai[element], Ai[element+1]);
193 std::swap (Ax[element], Ax[element+1]);
202 template <
class Matrix>
235 Ai.resize (matrix.n_nonzero_elements());
236 Ax.resize (matrix.n_nonzero_elements());
241 Ap[row] =
Ap[row-1] + matrix.get_row_length(row-1);
242 Assert (static_cast<size_type>(
Ap.back()) == Ai.size(),
252 std::vector<long int> row_pointers =
Ap;
256 for (
size_type row = 0; row < matrix.m(); ++row)
258 for (
typename Matrix::const_iterator p=matrix.begin(row);
259 p!=matrix.end(row); ++p)
262 Ai[row_pointers[row]] = p->column();
263 Ax[row_pointers[row]] = p->value();
281 status = umfpack_dl_symbolic (N, N,
282 &
Ap[0], &Ai[0], &Ax[0],
288 status = umfpack_dl_numeric (&
Ap[0], &Ai[0], &Ax[0],
290 &numeric_decomposition,
310 rhs = rhs_and_solution;
319 = umfpack_dl_solve (transpose ? UMFPACK_A : UMFPACK_At,
320 &
Ap[0], &Ai[0], &Ax[0],
321 rhs_and_solution.
begin(), rhs.begin(),
322 numeric_decomposition,
336 tmp = rhs_and_solution;
337 solve (tmp, transpose);
338 rhs_and_solution = tmp;
343 template <
class Matrix>
350 solve (rhs_and_solution, transpose);
354 template <
class Matrix>
361 solve (rhs_and_solution, transpose);
372 symbolic_decomposition (0),
373 numeric_decomposition (0),
383 template <
class Matrix>
386 AssertThrow(
false,
ExcMessage(
"To call this function you need UMFPACK, but you configured deal.II without passing the necessary switch to 'cmake'. Please consult the installation instructions in doc/readme.html."));
393 AssertThrow(
false,
ExcMessage(
"To call this function you need UMFPACK, but you configured deal.II without passing the necessary switch to 'cmake'. Please consult the installation instructions in doc/readme.html."));
401 AssertThrow(
false,
ExcMessage(
"To call this function you need UMFPACK, but you configured deal.II without passing the necessary switch to 'cmake'. Please consult the installation instructions in doc/readme.html."));
405 template <
class Matrix>
411 AssertThrow(
false,
ExcMessage(
"To call this function you need UMFPACK, but you configured deal.II without passing the necessary switch to 'cmake'. Please consult the installation instructions in doc/readme.html."));
416 template <
class Matrix>
422 AssertThrow(
false,
ExcMessage(
"To call this function you need UMFPACK, but you configured deal.II without passing the necessary switch to 'cmake'. Please consult the installation instructions in doc/readme.html."));
428 template <
class Matrix>
464 this->
solve(dst,
true);
475 this->
solve(dst,
true);
494 #define InstantiateUMFPACK(MatrixType) \ 496 void SparseDirectUMFPACK::factorize (const MatrixType &); \ 498 void SparseDirectUMFPACK::solve (const MatrixType &, \ 502 void SparseDirectUMFPACK::solve (const MatrixType &, \ 503 BlockVector<double> &, \ 506 void SparseDirectUMFPACK::initialize (const MatrixType &, \ 507 const AdditionalData); 516 DEAL_II_NAMESPACE_CLOSE
void * symbolic_decomposition
types::global_dof_index size_type
static ::ExceptionBase & ExcNotInitialized()
#define AssertThrow(cond, exc)
void factorize(const Matrix &matrix)
void vmult(Vector< double > &dst, const Vector< double > &src) const
static ::ExceptionBase & ExcMessage(std::string arg1)
void solve(Vector< double > &rhs_and_solution, const bool transpose=false) const
#define Assert(cond, exc)
std::vector< SuiteSparse_long > Ap
static ::ExceptionBase & ExcUMFPACKError(char *arg1, int arg2)
std::vector< double > control
static ::ExceptionBase & ExcNotQuadratic()
void initialize(const SparsityPattern &sparsity_pattern)
void Tvmult(Vector< double > &dst, const Vector< double > &src) const
void sort_arrays(const SparseMatrixEZ< number > &)
static ::ExceptionBase & ExcInternalError()