Reference documentation for deal.II version 8.5.1
petsc_matrix_free.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2012 - 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 #ifndef dealii__petsc_matrix_free_h
17 #define dealii__petsc_matrix_free_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #ifdef DEAL_II_WITH_PETSC
23 
24 # include <deal.II/lac/exceptions.h>
25 # include <deal.II/lac/petsc_matrix_base.h>
26 # include <deal.II/lac/petsc_vector.h>
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 
31 
32 namespace PETScWrappers
33 {
64  class MatrixFree : public MatrixBase
65  {
66  public:
67 
71  MatrixFree ();
72 
87  MatrixFree (const MPI_Comm &communicator,
88  const unsigned int m,
89  const unsigned int n,
90  const unsigned int local_rows,
91  const unsigned int local_columns);
92 
104  MatrixFree (const MPI_Comm &communicator,
105  const unsigned int m,
106  const unsigned int n,
107  const std::vector<unsigned int> &local_rows_per_process,
108  const std::vector<unsigned int> &local_columns_per_process,
109  const unsigned int this_process);
110 
116  MatrixFree (const unsigned int m,
117  const unsigned int n,
118  const unsigned int local_rows,
119  const unsigned int local_columns);
120 
126  MatrixFree (const unsigned int m,
127  const unsigned int n,
128  const std::vector<unsigned int> &local_rows_per_process,
129  const std::vector<unsigned int> &local_columns_per_process,
130  const unsigned int this_process);
131 
137  void reinit (const MPI_Comm &communicator,
138  const unsigned int m,
139  const unsigned int n,
140  const unsigned int local_rows,
141  const unsigned int local_columns);
142 
148  void reinit (const MPI_Comm &communicator,
149  const unsigned int m,
150  const unsigned int n,
151  const std::vector<unsigned int> &local_rows_per_process,
152  const std::vector<unsigned int> &local_columns_per_process,
153  const unsigned int this_process);
154 
159  void reinit (const unsigned int m,
160  const unsigned int n,
161  const unsigned int local_rows,
162  const unsigned int local_columns);
163 
168  void reinit (const unsigned int m,
169  const unsigned int n,
170  const std::vector<unsigned int> &local_rows_per_process,
171  const std::vector<unsigned int> &local_columns_per_process,
172  const unsigned int this_process);
173 
178  void clear ();
179 
184  const MPI_Comm &get_mpi_communicator () const;
185 
197  virtual
198  void vmult (VectorBase &dst,
199  const VectorBase &src) const = 0;
200 
213  virtual
214  void Tvmult (VectorBase &dst,
215  const VectorBase &src) const = 0;
216 
228  virtual
229  void vmult_add (VectorBase &dst,
230  const VectorBase &src) const = 0;
231 
244  virtual
245  void Tvmult_add (VectorBase &dst,
246  const VectorBase &src) const = 0;
247 
256  virtual
257  void vmult (Vec &dst, const Vec &src) const;
258 
259  private:
260 
265  MPI_Comm communicator;
266 
278  static int matrix_free_mult (Mat A, Vec src, Vec dst);
279 
285  void do_reinit (const unsigned int m,
286  const unsigned int n,
287  const unsigned int local_rows,
288  const unsigned int local_columns);
289  };
290 
291 
292 
293 // -------- template and inline functions ----------
294 
295  inline
296  const MPI_Comm &
298  {
299  return communicator;
300  }
301 }
302 
303 DEAL_II_NAMESPACE_CLOSE
304 
305 #endif // DEAL_II_WITH_PETSC
306 
307 
308 /*---------------------------- petsc_matrix_free.h ---------------------------*/
309 
310 #endif
311 /*---------------------------- petsc_matrix_free.h ---------------------------*/
virtual void vmult_add(VectorBase &dst, const VectorBase &src) const =0
virtual void Tvmult_add(VectorBase &dst, const VectorBase &src) const =0
virtual void vmult(VectorBase &dst, const VectorBase &src) const =0
virtual void Tvmult(VectorBase &dst, const VectorBase &src) const =0
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
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)