Reference documentation for deal.II version 9.0.0
petsc_sparse_matrix.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2018 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/petsc_sparse_matrix.h>
17 
18 #ifdef DEAL_II_WITH_PETSC
19 
20 # include <deal.II/lac/exceptions.h>
21 # include <deal.II/lac/petsc_compatibility.h>
22 # include <deal.II/lac/petsc_vector_base.h>
23 # include <deal.II/lac/sparsity_pattern.h>
24 # include <deal.II/lac/dynamic_sparsity_pattern.h>
25 
26 DEAL_II_NAMESPACE_OPEN
27 
28 namespace PETScWrappers
29 {
30 
32  {
33  const int m=0, n=0, n_nonzero_per_row=0;
34  const PetscErrorCode ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, m, n,
35  n_nonzero_per_row, nullptr, &matrix);
36  AssertThrow (ierr == 0, ExcPETScError(ierr));
37  }
38 
39 
40 
42  const size_type n,
43  const size_type n_nonzero_per_row,
44  const bool is_symmetric)
45  {
46  do_reinit (m, n, n_nonzero_per_row, is_symmetric);
47  }
48 
49 
50 
52  const size_type n,
53  const std::vector<size_type> &row_lengths,
54  const bool is_symmetric)
55  {
56  do_reinit (m, n, row_lengths, is_symmetric);
57  }
58 
59 
60 
61  template <typename SparsityPatternType>
63  SparseMatrix (const SparsityPatternType &sparsity_pattern,
64  const bool preset_nonzero_locations)
65  {
66  do_reinit (sparsity_pattern, preset_nonzero_locations);
67  }
68 
69 
70 
71  SparseMatrix &
72  SparseMatrix::operator = (const double d)
73  {
75  return *this;
76  }
77 
78 
79 
80  void
82  const size_type n,
83  const size_type n_nonzero_per_row,
84  const bool is_symmetric)
85  {
86  // get rid of old matrix and generate a
87  // new one
88  const PetscErrorCode ierr = destroy_matrix (matrix);
89  AssertThrow (ierr == 0, ExcPETScError (ierr));
90 
91  do_reinit (m, n, n_nonzero_per_row, is_symmetric);
92  }
93 
94 
95 
96  void
98  const size_type n,
99  const std::vector<size_type> &row_lengths,
100  const bool is_symmetric)
101  {
102  // get rid of old matrix and generate a
103  // new one
104  const PetscErrorCode ierr = destroy_matrix (matrix);
105  AssertThrow (ierr == 0, ExcPETScError (ierr));
106 
107  do_reinit (m, n, row_lengths, is_symmetric);
108  }
109 
110 
111 
112  template <typename SparsityPatternType>
113  void
115  reinit (const SparsityPatternType &sparsity_pattern,
116  const bool preset_nonzero_locations)
117  {
118  // get rid of old matrix and generate a
119  // new one
120  const PetscErrorCode ierr = destroy_matrix (matrix);
121  AssertThrow (ierr == 0, ExcPETScError (ierr));
122 
123  do_reinit (sparsity_pattern, preset_nonzero_locations);
124  }
125 
126 
127 
128  const MPI_Comm &
130  {
131  static MPI_Comm comm;
132  const PetscErrorCode ierr = PetscObjectGetComm((PetscObject)matrix, &comm);
133  AssertThrow (ierr == 0, ExcPETScError(ierr));
134  return comm;
135  }
136 
137 
138 
139  void
141  const size_type n,
142  const size_type n_nonzero_per_row,
143  const bool is_symmetric)
144  {
145  // use the call sequence indicating only
146  // a maximal number of elements per row
147  // for all rows globally
148  const PetscErrorCode ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, m, n,
149  n_nonzero_per_row,
150  nullptr, &matrix);
151  AssertThrow (ierr == 0, ExcPETScError(ierr));
152 
153  // set symmetric flag, if so requested
154  if (is_symmetric == true)
155  {
156  set_matrix_option (matrix, MAT_SYMMETRIC, PETSC_TRUE);
157  }
158  }
159 
160 
161 
162  void
164  const size_type n,
165  const std::vector<size_type> &row_lengths,
166  const bool is_symmetric)
167  {
168  Assert (row_lengths.size() == m,
169  ExcDimensionMismatch (row_lengths.size(), m));
170 
171  // use the call sequence indicating a
172  // maximal number of elements for each
173  // row individually. annoyingly, we
174  // always use unsigned ints for cases
175  // like this, while PETSc wants to see
176  // signed integers. so we have to
177  // convert, unless we want to play dirty
178  // tricks with conversions of pointers
179  const std::vector<PetscInt>
180  int_row_lengths (row_lengths.begin(), row_lengths.end());
181 
182  const PetscErrorCode ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, m, n, 0,
183  int_row_lengths.data(), &matrix);
184  AssertThrow (ierr == 0, ExcPETScError(ierr));
185 
186  // set symmetric flag, if so requested
187  if (is_symmetric == true)
188  {
189  set_matrix_option(matrix, MAT_SYMMETRIC, PETSC_TRUE);
190  }
191  }
192 
193 
194 
195  template <typename SparsityPatternType>
196  void
197  SparseMatrix::do_reinit (const SparsityPatternType &sparsity_pattern,
198  const bool preset_nonzero_locations)
199  {
200  std::vector<size_type> row_lengths (sparsity_pattern.n_rows());
201  for (size_type i=0; i<sparsity_pattern.n_rows(); ++i)
202  row_lengths[i] = sparsity_pattern.row_length (i);
203 
204  do_reinit (sparsity_pattern.n_rows(),
205  sparsity_pattern.n_cols(),
206  row_lengths, false);
207 
208  // next preset the exact given matrix
209  // entries with zeros, if the user
210  // requested so. this doesn't avoid any
211  // memory allocations, but it at least
212  // avoids some searches later on. the
213  // key here is that we can use the
214  // matrix set routines that set an
215  // entire row at once, not a single
216  // entry at a time
217  //
218  // for the usefulness of this option
219  // read the documentation of this
220  // class.
221  if (preset_nonzero_locations == true)
222  {
223  std::vector<PetscInt> row_entries;
224  std::vector<PetscScalar> row_values;
225  for (size_type i=0; i<sparsity_pattern.n_rows(); ++i)
226  {
227  row_entries.resize (row_lengths[i]);
228  row_values.resize (row_lengths[i], 0.0);
229  for (size_type j=0; j<row_lengths[i]; ++j)
230  row_entries[j] = sparsity_pattern.column_number (i,j);
231 
232  const PetscInt int_row = i;
233  const PetscErrorCode ierr = MatSetValues (matrix, 1, &int_row,
234  row_lengths[i], row_entries.data(),
235  row_values.data(), INSERT_VALUES);
236  AssertThrow (ierr == 0, ExcPETScError(ierr));
237  }
239 
242  }
243  }
244 
245  size_t
247  {
248  PetscInt m,n;
249  const PetscErrorCode ierr = MatGetSize(matrix, &m, &n);
250  AssertThrow (ierr == 0, ExcPETScError(ierr));
251 
252  return m;
253  }
254 
255  size_t
257  {
258  PetscInt m,n;
259  const PetscErrorCode ierr = MatGetSize(matrix, &m, &n);
260  AssertThrow (ierr == 0, ExcPETScError(ierr));
261 
262  return n;
263  }
264 
265  void
267  const SparseMatrix &B,
268  const MPI::Vector &V) const
269  {
270  // Simply forward to the protected member function of the base class
271  // that takes abstract matrix and vector arguments (to which the compiler
272  // automatically casts the arguments).
273  MatrixBase::mmult (C, B, V);
274  }
275 
276  void
278  const SparseMatrix &B,
279  const MPI::Vector &V) const
280  {
281  // Simply forward to the protected member function of the base class
282  // that takes abstract matrix and vector arguments (to which the compiler
283  // automatically casts the arguments).
284  MatrixBase::Tmmult (C, B, V);
285  }
286 
287  // Explicit instantiations
288  //
289  template
291  const bool);
292  template
294  const bool);
295 
296  template void
298  const bool);
299  template void
301  const bool);
302 
303  template void
305  const bool);
306  template void
308  const bool);
309 }
310 
311 
312 DEAL_II_NAMESPACE_CLOSE
313 
314 #endif // DEAL_II_WITH_PETSC
PetscErrorCode destroy_matrix(Mat &matrix)
void Tmmult(MatrixBase &C, const MatrixBase &B, const VectorBase &V) const
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
void set_matrix_option(Mat &matrix, const MatOption option_name, const PetscBool option_value=PETSC_FALSE)
PetscBool is_symmetric(const double tolerance=1.e-12)
void compress(const VectorOperation::values operation)
MatrixBase & operator=(const MatrixBase &)=delete
SparseMatrix & operator=(const double d)
void reinit(const size_type m, const size_type n, const size_type n_nonzero_per_row, const bool is_symmetric=false)
virtual const MPI_Comm & get_mpi_communicator() const
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void mmult(MatrixBase &C, const MatrixBase &B, const VectorBase &V) const
types::global_dof_index size_type
void set_keep_zero_rows(Mat &matrix)
void mmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
void do_reinit(const size_type m, const size_type n, const size_type n_nonzero_per_row, const bool is_symmetric=false)
void close_matrix(Mat &matrix)