Reference documentation for deal.II version 8.5.1
sparse_direct.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 2016 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
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>
22 
23 #include <cerrno>
24 #include <iostream>
25 #include <list>
26 #include <typeinfo>
27 #include <vector>
28 
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 
33 // include UMFPACK file.
34 #ifdef DEAL_II_WITH_UMFPACK
35 # include <umfpack.h>
36 #endif
37 
38 
39 
41 {
42  clear ();
43 }
44 
45 
46 void
49 {}
50 
51 
52 #ifdef DEAL_II_WITH_UMFPACK
53 
55  :
56  _m (0),
57  _n (0),
58  symbolic_decomposition (0),
59  numeric_decomposition (0),
60  control (UMFPACK_CONTROL)
61 {
62  umfpack_dl_defaults (&control[0]);
63 }
64 
65 
66 
67 void
69 {
70  // delete objects that haven't been deleted yet
71  if (symbolic_decomposition != 0)
72  {
73  umfpack_dl_free_symbolic (&symbolic_decomposition);
75  }
76 
77  if (numeric_decomposition != 0)
78  {
79  umfpack_dl_free_numeric (&numeric_decomposition);
80  numeric_decomposition = 0;
81  }
82 
83  {
84  std::vector<long int> tmp;
85  tmp.swap (Ap);
86  }
87 
88  {
89  std::vector<long int> tmp;
90  tmp.swap (Ai);
91  }
92 
93  {
94  std::vector<double> tmp;
95  tmp.swap (Ax);
96  }
97 
98  umfpack_dl_defaults (&control[0]);
99 }
100 
101 
102 
103 template <typename number>
104 void
106 sort_arrays (const SparseMatrix<number> &matrix)
107 {
108  // do the copying around of entries so that the diagonal entry is in the
109  // right place. note that this is easy to detect: since all entries apart
110  // from the diagonal entry are sorted, we know that the diagonal entry is
111  // in the wrong place if and only if its column index is larger than the
112  // column index of the second entry in a row
113  //
114  // ignore rows with only one or no entry
115  for (size_type row=0; row<matrix.m(); ++row)
116  {
117  // we may have to move some elements that are left of the diagonal
118  // but presently after the diagonal entry to the left, whereas the
119  // diagonal entry has to move to the right. we could first figure out
120  // where to move everything to, but for simplicity we just make a
121  // series of swaps instead (this is kind of a single run of
122  // bubble-sort, which gives us the desired result since the array is
123  // already "almost" sorted)
124  //
125  // in the first loop, the condition in the while-header also checks
126  // that the row has at least two entries and that the diagonal entry
127  // is really in the wrong place
128  long int cursor = Ap[row];
129  while ((cursor < Ap[row+1]-1) &&
130  (Ai[cursor] > Ai[cursor+1]))
131  {
132  std::swap (Ai[cursor], Ai[cursor+1]);
133  std::swap (Ax[cursor], Ax[cursor+1]);
134  ++cursor;
135  }
136  }
137 }
138 
139 
140 
141 template <typename number>
142 void
145 {
146  //same thing for SparseMatrixEZ
147  for (size_type row=0; row<matrix.m(); ++row)
148  {
149  long int cursor = Ap[row];
150  while ((cursor < Ap[row+1]-1) &&
151  (Ai[cursor] > Ai[cursor+1]))
152  {
153  std::swap (Ai[cursor], Ai[cursor+1]);
154  std::swap (Ax[cursor], Ax[cursor+1]);
155  ++cursor;
156  }
157  }
158 }
159 
160 
161 
162 template <typename number>
163 void
166 {
167  // the case for block matrices is a bit more difficult, since all we know
168  // is that *within each block*, the diagonal of that block may come
169  // first. however, that means that there may be as many entries per row
170  // in the wrong place as there are block columns. we can do the same
171  // thing as above, but we have to do it multiple times
172  for (size_type row=0; row<matrix.m(); ++row)
173  {
174  long int cursor = Ap[row];
175  for (size_type block=0; block<matrix.n_block_cols(); ++block)
176  {
177  // find the next out-of-order element
178  while ((cursor < Ap[row+1]-1) &&
179  (Ai[cursor] < Ai[cursor+1]))
180  ++cursor;
181 
182  // if there is none, then just go on
183  if (cursor == Ap[row+1]-1)
184  break;
185 
186  // otherwise swap this entry with successive ones as long as
187  // necessary
188  long int element = cursor;
189  while ((element < Ap[row+1]-1) &&
190  (Ai[element] > Ai[element+1]))
191  {
192  std::swap (Ai[element], Ai[element+1]);
193  std::swap (Ax[element], Ax[element+1]);
194  ++element;
195  }
196  }
197  }
198 }
199 
200 
201 
202 template <class Matrix>
203 void
205 factorize (const Matrix &matrix)
206 {
207  Assert (matrix.m() == matrix.n(), ExcNotQuadratic())
208 
209  clear ();
210 
211  _m = matrix.m();
212  _n = matrix.n();
213 
214  const size_type N = matrix.m();
215 
216  // copy over the data from the matrix to the data structures UMFPACK
217  // wants. note two things: first, UMFPACK wants compressed column storage
218  // whereas we always do compressed row storage; we work around this by,
219  // rather than shuffling things around, copy over the data we have, but
220  // then call the umfpack_dl_solve function with the UMFPACK_At argument,
221  // meaning that we want to solve for the transpose system
222  //
223  // second: the data we have in the sparse matrices is "almost" right
224  // already; UMFPACK wants the entries in each row (i.e. really: column)
225  // to be sorted in ascending order. we almost have that, except that we
226  // usually store the diagonal first in each row to allow for some
227  // optimizations. thus, we have to resort things a little bit, but only
228  // within each row
229  //
230  // final note: if the matrix has entries in the sparsity pattern that are
231  // actually occupied by entries that have a zero numerical value, then we
232  // keep them anyway. people are supposed to provide accurate sparsity
233  // patterns.
234  Ap.resize (N+1);
235  Ai.resize (matrix.n_nonzero_elements());
236  Ax.resize (matrix.n_nonzero_elements());
237 
238  // first fill row lengths array
239  Ap[0] = 0;
240  for (size_type row=1; row<=N; ++row)
241  Ap[row] = Ap[row-1] + matrix.get_row_length(row-1);
242  Assert (static_cast<size_type>(Ap.back()) == Ai.size(),
243  ExcInternalError());
244 
245  // then copy over matrix elements. note that for sparse matrices,
246  // iterators are sorted so that they traverse each row from start to end
247  // before moving on to the next row. however, this isn't true for block
248  // matrices, so we have to do a bit of book keeping
249  {
250  // have an array that for each row points to the first entry not yet
251  // written to
252  std::vector<long int> row_pointers = Ap;
253 
254  // loop over the elements of the matrix row by row, as suggested in the
255  // documentation of the sparse matrix iterator class
256  for (size_type row = 0; row < matrix.m(); ++row)
257  {
258  for (typename Matrix::const_iterator p=matrix.begin(row);
259  p!=matrix.end(row); ++p)
260  {
261  // write entry into the first free one for this row
262  Ai[row_pointers[row]] = p->column();
263  Ax[row_pointers[row]] = p->value();
264 
265  // then move pointer ahead
266  ++row_pointers[row];
267  }
268  }
269 
270  // at the end, we should have written all rows completely
271  for (size_type i=0; i<Ap.size()-1; ++i)
272  Assert (row_pointers[i] == Ap[i+1], ExcInternalError());
273  }
274 
275  // make sure that the elements in each row are sorted. we have to be more
276  // careful for block sparse matrices, so ship this task out to a
277  // different function
278  sort_arrays (matrix);
279 
280  int status;
281  status = umfpack_dl_symbolic (N, N,
282  &Ap[0], &Ai[0], &Ax[0],
284  &control[0], 0);
285  AssertThrow (status == UMFPACK_OK,
286  ExcUMFPACKError("umfpack_dl_symbolic", status));
287 
288  status = umfpack_dl_numeric (&Ap[0], &Ai[0], &Ax[0],
290  &numeric_decomposition,
291  &control[0], 0);
292  AssertThrow (status == UMFPACK_OK,
293  ExcUMFPACKError("umfpack_dl_numeric", status));
294 
295  umfpack_dl_free_symbolic (&symbolic_decomposition) ;
296 }
297 
298 
299 
300 void
302  bool transpose /*=false*/) const
303 {
304  // make sure that some kind of factorize() call has happened before
305  Assert (Ap.size() != 0, ExcNotInitialized());
306  Assert (Ai.size() != 0, ExcNotInitialized());
307  Assert (Ai.size() == Ax.size(), ExcNotInitialized());
308 
309  Vector<double> rhs (rhs_and_solution.size());
310  rhs = rhs_and_solution;
311 
312  // solve the system. note that since UMFPACK wants compressed column
313  // storage instead of the compressed row storage format we use in
314  // deal.II's SparsityPattern classes, we solve for UMFPACK's A^T instead
315 
316  // Conversely, if we solve for the transpose, we have to use UMFPACK_A
317  // instead.
318  const int status
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,
323  &control[0], 0);
324  AssertThrow (status == UMFPACK_OK, ExcUMFPACKError("umfpack_dl_solve", status));
325 }
326 
327 
328 void
330  bool transpose /*=false*/) const
331 {
332  // the UMFPACK functions want a contiguous array of elements, so
333  // there is no way around copying data around. thus, just copy the
334  // data into a regular vector and back
335  Vector<double> tmp (rhs_and_solution.size());
336  tmp = rhs_and_solution;
337  solve (tmp, transpose);
338  rhs_and_solution = tmp;
339 }
340 
341 
342 
343 template <class Matrix>
344 void
345 SparseDirectUMFPACK::solve (const Matrix &matrix,
346  Vector<double> &rhs_and_solution,
347  bool transpose /*=false*/)
348 {
349  factorize (matrix);
350  solve (rhs_and_solution, transpose);
351 }
352 
353 
354 template <class Matrix>
355 void
356 SparseDirectUMFPACK::solve (const Matrix &matrix,
357  BlockVector<double> &rhs_and_solution,
358  bool transpose /*=false*/)
359 {
360  factorize (matrix);
361  solve (rhs_and_solution, transpose);
362 }
363 
364 
365 #else
366 
367 
369  :
370  _m(0),
371  _n(0),
372  symbolic_decomposition (0),
373  numeric_decomposition (0),
374  control (0)
375 {}
376 
377 
378 void
380 {}
381 
382 
383 template <class Matrix>
384 void SparseDirectUMFPACK::factorize (const Matrix &)
385 {
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."));
387 }
388 
389 
390 void
392 {
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."));
394 }
395 
396 
397 
398 void
400 {
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."));
402 }
403 
404 
405 template <class Matrix>
406 void
407 SparseDirectUMFPACK::solve (const Matrix &,
408  Vector<double> &,
409  bool)
410 {
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."));
412 }
413 
414 
415 
416 template <class Matrix>
417 void
418 SparseDirectUMFPACK::solve (const Matrix &,
420  bool)
421 {
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."));
423 }
424 
425 #endif
426 
427 
428 template <class Matrix>
429 void
431  const AdditionalData)
432 {
433  this->factorize(M);
434 }
435 
436 
437 void
439  Vector<double> &dst,
440  const Vector<double> &src) const
441 {
442  dst = src;
443  this->solve(dst);
444 }
445 
446 
447 
448 void
450  BlockVector<double> &dst,
451  const BlockVector<double> &src) const
452 {
453  dst = src;
454  this->solve(dst);
455 }
456 
457 
458 void
460  Vector<double> &dst,
461  const Vector<double> &src) const
462 {
463  dst = src;
464  this->solve(dst, /*transpose=*/ true);
465 }
466 
467 
468 
469 void
471  BlockVector<double> &dst,
472  const BlockVector<double> &src) const
473 {
474  dst = src;
475  this->solve(dst, /*transpose=*/ true);
476 }
477 
480 {
481  Assert (_m!=0, ExcNotInitialized());
482  return _m;
483 }
484 
487 {
488  Assert (_n!=0, ExcNotInitialized());
489  return _n;
490 }
491 
492 
493 // explicit instantiations for SparseMatrixUMFPACK
494 #define InstantiateUMFPACK(MatrixType) \
495  template \
496  void SparseDirectUMFPACK::factorize (const MatrixType &); \
497  template \
498  void SparseDirectUMFPACK::solve (const MatrixType &, \
499  Vector<double> &, \
500  bool); \
501  template \
502  void SparseDirectUMFPACK::solve (const MatrixType &, \
503  BlockVector<double> &, \
504  bool); \
505  template \
506  void SparseDirectUMFPACK::initialize (const MatrixType &, \
507  const AdditionalData);
508 
509 InstantiateUMFPACK(SparseMatrix<double>)
510 InstantiateUMFPACK(SparseMatrix<float>)
511 InstantiateUMFPACK(SparseMatrixEZ<double>)
512 InstantiateUMFPACK(SparseMatrixEZ<float>)
513 InstantiateUMFPACK(BlockSparseMatrix<double>)
514 InstantiateUMFPACK(BlockSparseMatrix<float>)
515 
516 DEAL_II_NAMESPACE_CLOSE
types::global_dof_index size_type
Definition: sparse_direct.h:88
size_type n() const
static ::ExceptionBase & ExcNotInitialized()
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
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)
Definition: exceptions.h:313
std::size_t size() const
std::vector< SuiteSparse_long > Ap
static ::ExceptionBase & ExcUMFPACKError(char *arg1, int arg2)
iterator begin()
std::vector< double > control
static ::ExceptionBase & ExcNotQuadratic()
void initialize(const SparsityPattern &sparsity_pattern)
size_type m() const
void Tvmult(Vector< double > &dst, const Vector< double > &src) const
void sort_arrays(const SparseMatrixEZ< number > &)
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
static ::ExceptionBase & ExcInternalError()