18#ifdef DEAL_II_WITH_PETSC
38 const unsigned int local_rows,
39 const unsigned int local_columns)
41 do_reinit(communicator,
m,
n, local_rows, local_columns);
50 const std::vector<unsigned int> &local_rows_per_process,
51 const std::vector<unsigned int> &local_columns_per_process,
52 const unsigned int this_process)
54 Assert(local_rows_per_process.size() == local_columns_per_process.size(),
56 local_columns_per_process.size()));
62 local_rows_per_process[this_process],
63 local_columns_per_process[this_process]);
70 const unsigned int local_rows,
71 const unsigned int local_columns)
73 do_reinit(MPI_COMM_WORLD,
m,
n, local_rows, local_columns);
81 const std::vector<unsigned int> &local_rows_per_process,
82 const std::vector<unsigned int> &local_columns_per_process,
83 const unsigned int this_process)
85 Assert(local_rows_per_process.size() == local_columns_per_process.size(),
87 local_columns_per_process.size()));
93 local_rows_per_process[this_process],
94 local_columns_per_process[this_process]);
101 const unsigned int m,
102 const unsigned int n,
103 const unsigned int local_rows,
104 const unsigned int local_columns)
107 const PetscErrorCode ierr = MatDestroy(&
matrix);
110 do_reinit(communicator,
m,
n, local_rows, local_columns);
117 const unsigned int m,
118 const unsigned int n,
119 const std::vector<unsigned int> &local_rows_per_process,
120 const std::vector<unsigned int> &local_columns_per_process,
121 const unsigned int this_process)
123 Assert(local_rows_per_process.size() == local_columns_per_process.size(),
125 local_columns_per_process.size()));
128 const PetscErrorCode ierr = MatDestroy(&
matrix);
134 local_rows_per_process[this_process],
135 local_columns_per_process[this_process]);
142 const unsigned int n,
143 const unsigned int local_rows,
144 const unsigned int local_columns)
153 const unsigned int n,
154 const std::vector<unsigned int> &local_rows_per_process,
155 const std::vector<unsigned int> &local_columns_per_process,
156 const unsigned int this_process)
161 local_rows_per_process,
162 local_columns_per_process,
171 const PetscErrorCode ierr = MatDestroy(&
matrix);
201 const PetscErrorCode ierr = MatShellGetContext(A, &this_object);
214 const unsigned int m,
215 const unsigned int n,
216 const unsigned int local_rows,
217 const unsigned int local_columns)
225 PetscErrorCode ierr = MatCreateShell(communicator,
230 static_cast<void *
>(
this),
235 ierr = MatShellSetOperation(
238 reinterpret_cast<void (*)()
>(
242 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)