17 #include <deal.II/lac/petsc_matrix_free.h> 19 #ifdef DEAL_II_WITH_PETSC 21 DEAL_II_NAMESPACE_OPEN
26 : communicator (PETSC_COMM_SELF)
37 const unsigned int local_rows,
38 const unsigned int local_columns)
39 : communicator (communicator)
41 do_reinit (m, n, local_rows, local_columns);
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)
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(),
61 local_rows_per_process[this_process],
62 local_columns_per_process[this_process]);
69 const unsigned int local_rows,
70 const unsigned int local_columns)
73 do_reinit (m, n, local_rows, local_columns);
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)
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(),
92 local_rows_per_process[this_process],
93 local_columns_per_process[this_process]);
100 const unsigned int n,
101 const unsigned int local_rows,
102 const unsigned int local_columns)
108 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 109 int ierr = MatDestroy (
matrix);
111 int ierr = MatDestroy (&
matrix);
115 do_reinit (m, n, local_rows, local_columns);
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)
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(),
135 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 136 int ierr = MatDestroy (
matrix);
138 int ierr = MatDestroy (&
matrix);
143 local_rows_per_process[this_process],
144 local_columns_per_process[this_process]);
150 const unsigned int n,
151 const unsigned int local_rows,
152 const unsigned int local_columns)
154 reinit (MPI_COMM_WORLD, m, n, local_rows, local_columns);
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)
165 reinit (MPI_COMM_WORLD, m, n, local_rows_per_process, local_columns_per_process, this_process);
172 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 173 int ierr = MatDestroy (
matrix);
175 int ierr = MatDestroy (&
matrix);
204 const char *vec_type;
205 int ierr = VecGetType (src, &vec_type);
208 ierr = VecGetLocalSize (src, &local_size);
211 if (strcmp(vec_type,
"mpi") == 0)
214 ierr = VecGetSize (src, &size);
220 else if (strcmp(vec_type,
"seq") == 0)
227 "This only works for Petsc Vec Type = VECMPI | VECSEQ"));
236 ierr = VecCopy (static_cast<const Vec &>(*y), dst);
252 int ierr = MatShellGetContext (A, &this_object);
264 const unsigned int n,
265 const unsigned int local_rows,
266 const unsigned int local_columns)
268 Assert (local_rows <= m, ExcDimensionMismatch (local_rows, m));
269 Assert (local_columns <= n, ExcDimensionMismatch (local_columns, n));
275 ierr = MatCreateShell(
communicator, local_rows, local_columns, m, n, (
void *)
this, &
matrix);
279 ierr = MatShellSetOperation (
matrix, MATOP_MULT,
283 ierr = MatSetFromOptions (
matrix);
289 DEAL_II_NAMESPACE_CLOSE
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)
size_type local_size() const
#define Assert(cond, exc)
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)