19#ifdef DEAL_II_WITH_PETSC
39 const unsigned int local_rows,
40 const unsigned int local_columns)
42 do_reinit(communicator,
m,
n, local_rows, local_columns);
51 const std::vector<unsigned int> &local_rows_per_process,
52 const std::vector<unsigned int> &local_columns_per_process,
53 const unsigned int this_process)
55 Assert(local_rows_per_process.size() == local_columns_per_process.size(),
57 local_columns_per_process.size()));
63 local_rows_per_process[this_process],
64 local_columns_per_process[this_process]);
71 const unsigned int local_rows,
72 const unsigned int local_columns)
74 do_reinit(MPI_COMM_WORLD,
m,
n, local_rows, local_columns);
82 const std::vector<unsigned int> &local_rows_per_process,
83 const std::vector<unsigned int> &local_columns_per_process,
84 const unsigned int this_process)
86 Assert(local_rows_per_process.size() == local_columns_per_process.size(),
88 local_columns_per_process.size()));
94 local_rows_per_process[this_process],
95 local_columns_per_process[this_process]);
102 const unsigned int m,
103 const unsigned int n,
104 const unsigned int local_rows,
105 const unsigned int local_columns)
108 const PetscErrorCode ierr = MatDestroy(&
matrix);
111 do_reinit(communicator,
m,
n, local_rows, local_columns);
118 const unsigned int m,
119 const unsigned int n,
120 const std::vector<unsigned int> &local_rows_per_process,
121 const std::vector<unsigned int> &local_columns_per_process,
122 const unsigned int this_process)
124 Assert(local_rows_per_process.size() == local_columns_per_process.size(),
126 local_columns_per_process.size()));
129 const PetscErrorCode ierr = MatDestroy(&
matrix);
135 local_rows_per_process[this_process],
136 local_columns_per_process[this_process]);
143 const unsigned int n,
144 const unsigned int local_rows,
145 const unsigned int local_columns)
154 const unsigned int n,
155 const std::vector<unsigned int> &local_rows_per_process,
156 const std::vector<unsigned int> &local_columns_per_process,
157 const unsigned int this_process)
162 local_rows_per_process,
163 local_columns_per_process,
172 const PetscErrorCode ierr = MatDestroy(&
matrix);
202 const PetscErrorCode ierr = MatShellGetContext(A, &this_object);
215 const unsigned int m,
216 const unsigned int n,
217 const unsigned int local_rows,
218 const unsigned int local_columns)
226 PetscErrorCode ierr = MatCreateShell(communicator,
231 static_cast<void *
>(
this),
236 ierr = MatShellSetOperation(
239 reinterpret_cast<void (*)()
>(
243 ierr = MatSetFromOptions(
matrix);
MPI_Comm get_mpi_communicator() const
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 MPI_Comm comm, const unsigned int m, const unsigned int n, const unsigned int local_rows, const unsigned int local_columns)
virtual void vmult(VectorBase &dst, const VectorBase &src) const =0
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define AssertThrow(cond, exc)