Reference documentation for deal.II version 9.1.0-pre
petsc_matrix_free.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2012 - 2017 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.md at
12 // the top level directory of deal.II.
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 # include <deal.II/lac/exceptions.h>
24 # include <deal.II/lac/petsc_matrix_base.h>
25 # include <deal.II/lac/petsc_vector.h>
26 DEAL_II_NAMESPACE_OPEN
27 
28 
29 
30 namespace PETScWrappers
31 {
60  class MatrixFree : public MatrixBase
61  {
62  public:
66  MatrixFree();
67 
82  MatrixFree(const MPI_Comm & communicator,
83  const unsigned int m,
84  const unsigned int n,
85  const unsigned int local_rows,
86  const unsigned int local_columns);
87 
99  MatrixFree(const MPI_Comm & communicator,
100  const unsigned int m,
101  const unsigned int n,
102  const std::vector<unsigned int> &local_rows_per_process,
103  const std::vector<unsigned int> &local_columns_per_process,
104  const unsigned int this_process);
105 
111  MatrixFree(const unsigned int m,
112  const unsigned int n,
113  const unsigned int local_rows,
114  const unsigned int local_columns);
115 
121  MatrixFree(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 
132  void
133  reinit(const MPI_Comm & communicator,
134  const unsigned int m,
135  const unsigned int n,
136  const unsigned int local_rows,
137  const unsigned int local_columns);
138 
144  void
145  reinit(const MPI_Comm & communicator,
146  const unsigned int m,
147  const unsigned int n,
148  const std::vector<unsigned int> &local_rows_per_process,
149  const std::vector<unsigned int> &local_columns_per_process,
150  const unsigned int this_process);
151 
156  void
157  reinit(const unsigned int m,
158  const unsigned int n,
159  const unsigned int local_rows,
160  const unsigned int local_columns);
161 
166  void
167  reinit(const unsigned int m,
168  const unsigned int n,
169  const std::vector<unsigned int> &local_rows_per_process,
170  const std::vector<unsigned int> &local_columns_per_process,
171  const unsigned int this_process);
172 
177  void
178  clear();
179 
184  const MPI_Comm &
185  get_mpi_communicator() const override;
186 
198  virtual void
199  vmult(VectorBase &dst, const VectorBase &src) const = 0;
200 
213  virtual void
214  Tvmult(VectorBase &dst, const VectorBase &src) const = 0;
215 
227  virtual void
228  vmult_add(VectorBase &dst, const VectorBase &src) const = 0;
229 
242  virtual void
243  Tvmult_add(VectorBase &dst, const VectorBase &src) const = 0;
244 
253  virtual void
254  vmult(Vec &dst, const Vec &src) const;
255 
256  private:
261  MPI_Comm communicator;
262 
274  static int
275  matrix_free_mult(Mat A, Vec src, Vec dst);
276 
282  void
283  do_reinit(const unsigned int m,
284  const unsigned int n,
285  const unsigned int local_rows,
286  const unsigned int local_columns);
287  };
288 
289 
290 
291  // -------- template and inline functions ----------
292 
293  inline const MPI_Comm &
295  {
296  return communicator;
297  }
298 } // namespace PETScWrappers
299 
300 DEAL_II_NAMESPACE_CLOSE
301 
302 # endif // DEAL_II_WITH_PETSC
303 
304 #endif
305 /*---------------------------- petsc_matrix_free.h --------------------------*/
const MPI_Comm & get_mpi_communicator() const override
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)
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)