Reference documentation for deal.II version 8.4.2
petsc_matrix_free.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2012 - 2014 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 
17 #include <deal.II/lac/petsc_matrix_free.h>
18 
19 #ifdef DEAL_II_WITH_PETSC
20 
21 DEAL_II_NAMESPACE_OPEN
22 
23 namespace PETScWrappers
24 {
26  : communicator (PETSC_COMM_SELF)
27  {
28  const int m=0;
29  do_reinit (m, m, m, m);
30  }
31 
32 
33 
35  const unsigned int m,
36  const unsigned int n,
37  const unsigned int local_rows,
38  const unsigned int local_columns)
39  : communicator (communicator)
40  {
41  do_reinit (m, n, local_rows, local_columns);
42  }
43 
44 
45 
47  const unsigned int m,
48  const unsigned int n,
49  const std::vector<unsigned int> &local_rows_per_process,
50  const std::vector<unsigned int> &local_columns_per_process,
51  const unsigned int this_process)
52  : communicator (communicator)
53  {
54  Assert (local_rows_per_process.size() == local_columns_per_process.size(),
55  ExcDimensionMismatch (local_rows_per_process.size(),
56  local_columns_per_process.size()));
57  Assert (this_process < local_rows_per_process.size(),
58  ExcInternalError());
59 
60  do_reinit (m, n,
61  local_rows_per_process[this_process],
62  local_columns_per_process[this_process]);
63  }
64 
65 
66 
67  MatrixFree::MatrixFree (const unsigned int m,
68  const unsigned int n,
69  const unsigned int local_rows,
70  const unsigned int local_columns)
71  : communicator (MPI_COMM_WORLD)
72  {
73  do_reinit (m, n, local_rows, local_columns);
74  }
75 
76 
77 
78  MatrixFree::MatrixFree (const unsigned int m,
79  const unsigned int n,
80  const std::vector<unsigned int> &local_rows_per_process,
81  const std::vector<unsigned int> &local_columns_per_process,
82  const unsigned int this_process)
83  : communicator (MPI_COMM_WORLD)
84  {
85  Assert (local_rows_per_process.size() == local_columns_per_process.size(),
86  ExcDimensionMismatch (local_rows_per_process.size(),
87  local_columns_per_process.size()));
88  Assert (this_process < local_rows_per_process.size(),
89  ExcInternalError());
90 
91  do_reinit (m, n,
92  local_rows_per_process[this_process],
93  local_columns_per_process[this_process]);
94  }
95 
96 
97 
98  void MatrixFree::reinit (const MPI_Comm &communicator,
99  const unsigned int m,
100  const unsigned int n,
101  const unsigned int local_rows,
102  const unsigned int local_columns)
103  {
104  this->communicator = communicator;
105 
106  // destroy the matrix and
107  // generate a new one
108 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
109  int ierr = MatDestroy (matrix);
110 #else
111  int ierr = MatDestroy (&matrix);
112 #endif
113  AssertThrow (ierr == 0, ExcPETScError(ierr));
114 
115  do_reinit (m, n, local_rows, local_columns);
116  }
117 
118 
119 
120  void MatrixFree::reinit (const MPI_Comm &communicator,
121  const unsigned int m,
122  const unsigned int n,
123  const std::vector<unsigned int> &local_rows_per_process,
124  const std::vector<unsigned int> &local_columns_per_process,
125  const unsigned int this_process)
126  {
127  Assert (local_rows_per_process.size() == local_columns_per_process.size(),
128  ExcDimensionMismatch (local_rows_per_process.size(),
129  local_columns_per_process.size()));
130  Assert (this_process < local_rows_per_process.size(),
131  ExcInternalError());
132 
133  this->communicator = communicator;
134 
135 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
136  int ierr = MatDestroy (matrix);
137 #else
138  int ierr = MatDestroy (&matrix);
139 #endif
140  AssertThrow (ierr == 0, ExcPETScError(ierr));
141 
142  do_reinit (m, n,
143  local_rows_per_process[this_process],
144  local_columns_per_process[this_process]);
145  }
146 
147 
148 
149  void MatrixFree::reinit (const unsigned int m,
150  const unsigned int n,
151  const unsigned int local_rows,
152  const unsigned int local_columns)
153  {
154  reinit (MPI_COMM_WORLD, m, n, local_rows, local_columns);
155  }
156 
157 
158 
159  void MatrixFree::reinit (const unsigned int m,
160  const unsigned int n,
161  const std::vector<unsigned int> &local_rows_per_process,
162  const std::vector<unsigned int> &local_columns_per_process,
163  const unsigned int this_process)
164  {
165  reinit (MPI_COMM_WORLD, m, n, local_rows_per_process, local_columns_per_process, this_process);
166  }
167 
168 
169 
171  {
172 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
173  int ierr = MatDestroy (matrix);
174 #else
175  int ierr = MatDestroy (&matrix);
176 #endif
177  AssertThrow (ierr == 0, ExcPETScError(ierr));
178 
179  const int m=0;
180  do_reinit (m, m, m, m);
181  }
182 
183 
184 
185  void MatrixFree::vmult (Vec &dst, const Vec &src) const
186  {
187 
188 //TODO: Translate the given PETSc Vec* vector into a deal.II
189 // vector so we can call the vmult function with the usual
190 // interface; then convert back. This could be much more
191 // efficient, if the PETScWrappers::*::Vector classes
192 // had a way to simply generate such a vector object from
193 // a given PETSc Vec* object without allocating new memory
194 // and without taking ownership of the Vec*
195 
196  VectorBase *x = 0;
197  VectorBase *y = 0;
198  // because we do not know,
199  // if dst and src are sequential
200  // or distributed vectors,
201  // we ask for the vector-type
202  // and reinit x and y with
203  // ::PETScWrappers::*::Vector:
204  const char *vec_type;
205  int ierr = VecGetType (src, &vec_type);
206 
207  PetscInt local_size;
208  ierr = VecGetLocalSize (src, &local_size);
209  AssertThrow (ierr == 0, ExcPETScError(ierr));
210 
211  if (strcmp(vec_type,"mpi") == 0)
212  {
213  PetscInt size;
214  ierr = VecGetSize (src, &size);
215  AssertThrow (ierr == 0, ExcPETScError(ierr));
216 
217  x = new PETScWrappers::MPI::Vector (this->get_mpi_communicator (), size, local_size);
218  y = new PETScWrappers::MPI::Vector (this->get_mpi_communicator (), size, local_size);
219  }
220  else if (strcmp(vec_type,"seq") == 0)
221  {
222  x = new PETScWrappers::Vector (local_size);
223  y = new PETScWrappers::Vector (local_size);
224  }
225  else
226  AssertThrow (false, ExcMessage("PETScWrappers::MPI::MatrixFree::do_matrix_vector_action: "
227  "This only works for Petsc Vec Type = VECMPI | VECSEQ"));
228 
229  // copy src to x
230  x->equ(1., PETScWrappers::VectorBase(src));
231  // and call vmult(x,y) which must
232  // be reimplemented in derived classes
233  vmult (*y, *x);
234 
235  // copy the result back to dst
236  ierr = VecCopy (static_cast<const Vec &>(*y), dst);
237  AssertThrow (ierr == 0, ExcPETScError(ierr));
238 
239  delete (x);
240  delete (y);
241  }
242 
243 
244 
245  int MatrixFree::matrix_free_mult (Mat A, Vec src, Vec dst)
246  {
247  // create a pointer to this MatrixFree
248  // object and link the given matrix A
249  // to the matrix-vector multiplication
250  // of this MatrixFree object,
251  void *this_object;
252  int ierr = MatShellGetContext (A, &this_object);
253  AssertThrow (ierr == 0, ExcPETScError(ierr));
254 
255  // call vmult of this object:
256  reinterpret_cast<MatrixFree *>(this_object)->vmult (dst, src);
257 
258  return (0);
259  }
260 
261 
262 
263  void MatrixFree::do_reinit (const unsigned int m,
264  const unsigned int n,
265  const unsigned int local_rows,
266  const unsigned int local_columns)
267  {
268  Assert (local_rows <= m, ExcDimensionMismatch (local_rows, m));
269  Assert (local_columns <= n, ExcDimensionMismatch (local_columns, n));
270 
271  int ierr;
272  // create a PETSc MatShell matrix-type
273  // object of dimension m x n and local size
274  // local_rows x local_columns
275  ierr = MatCreateShell(communicator, local_rows, local_columns, m, n, (void *)this, &matrix);
276  AssertThrow (ierr == 0, ExcPETScError(ierr));
277  // register the MatrixFree::matrix_free_mult function
278  // as the matrix multiplication used by this matrix
279  ierr = MatShellSetOperation (matrix, MATOP_MULT,
281  AssertThrow (ierr == 0, ExcPETScError(ierr));
282 
283  ierr = MatSetFromOptions (matrix);
284  AssertThrow (ierr == 0, ExcPETScError(ierr));
285  }
286 }
287 
288 
289 DEAL_II_NAMESPACE_CLOSE
290 
291 #endif // DEAL_II_WITH_PETSC
::ExceptionBase & ExcMessage(std::string arg1)
virtual void vmult(VectorBase &dst, const VectorBase &src) const =0
#define AssertThrow(cond, exc)
Definition: exceptions.h:358
#define Assert(cond, exc)
Definition: exceptions.h:294
void reinit(const MPI_Comm &communicator, const unsigned int m, const unsigned int n, const unsigned int local_rows, const unsigned int local_columns)
const MPI_Comm & get_mpi_communicator() const
void equ(const PetscScalar a, const VectorBase &V)
static int matrix_free_mult(Mat A, Vec src, Vec dst)
void do_reinit(const unsigned int m, const unsigned int n, const unsigned int local_rows, const unsigned int local_columns)