16 #include <deal.II/base/mpi.h> 17 #include <deal.II/lac/petsc_parallel_vector.h> 19 #ifdef DEAL_II_WITH_PETSC 21 # include <deal.II/lac/petsc_vector.h> 25 DEAL_II_NAMESPACE_OPEN
33 : communicator (MPI_COMM_SELF)
38 const PetscErrorCode ierr = VecCreateSeq (PETSC_COMM_SELF, n, &
vector);
49 communicator (communicator)
61 communicator (communicator)
79 const MPI_Comm &communicator)
81 communicator (communicator)
94 const MPI_Comm &communicator)
96 communicator (communicator)
114 const PetscErrorCode ierr = VecCreateSeq (PETSC_COMM_SELF, n, &
vector);
124 const bool omit_zeroing_entries)
133 const int ierr = MPI_Allreduce (&k, &k_global, 1,
149 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 150 const PetscErrorCode ierr = VecDestroy (
vector);
152 const PetscErrorCode ierr = VecDestroy (&
vector);
161 if (omit_zeroing_entries ==
false)
169 const bool omit_zeroing_entries)
174 if (!omit_zeroing_entries)
176 const PetscErrorCode ierr = VecSet(
vector, 0.0);
189 const MPI_Comm &comm)
191 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 192 const PetscErrorCode ierr = VecDestroy (
vector);
194 const PetscErrorCode ierr = VecDestroy (&
vector);
210 const MPI_Comm &comm)
212 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 213 const PetscErrorCode ierr = VecDestroy (
vector);
215 const PetscErrorCode ierr = VecDestroy (&
vector);
231 ExcMessage(
"Call to compress() required before calling operator=."));
238 PetscScalar *dest_array;
239 PetscErrorCode ierr = VecGetArray (
vector, &dest_array);
244 PetscScalar *src_array;
245 ierr = VecGetArray (static_cast<const Vec &>(v), &src_array);
249 const std::pair<size_type, size_type>
251 std::copy (src_array + local_elements.first,
252 src_array + local_elements.second,
256 ierr = VecRestoreArray (
vector, &dest_array);
259 ierr = VecRestoreArray (static_cast<const Vec &>(v), &src_array);
264 ierr = VecGhostUpdateBegin(
vector, INSERT_VALUES, SCATTER_FORWARD);
266 ierr = VecGhostUpdateEnd(
vector, INSERT_VALUES, SCATTER_FORWARD);
301 std::vector<size_type> ghostindices;
305 = (ghostindices.size() > 0
307 (
const PetscInt *)(&(ghostindices[0]))
327 ierr = VecGetOwnershipRange (
vector, &begin, &end);
333 ierr = VecGhostGetLocalForm(
vector, &l);
337 ierr = VecGetSize(l, &lsize);
340 ierr = VecGhostRestoreLocalForm(
vector, &l);
354 #if DEAL_II_PETSC_VERSION_LT(3,6,0) 368 #ifdef DEAL_II_WITH_MPI 372 return num_nonzero == 0;
374 return has_nonzero == 0;
381 const unsigned int precision,
382 const bool scientific,
383 const bool across)
const 390 PetscInt nlocal, istart, iend;
392 PetscErrorCode ierr = VecGetArray (
vector, &val);
395 ierr = VecGetLocalSize (
vector, &nlocal);
398 ierr = VecGetOwnershipRange (
vector, &istart, &iend);
402 std::ios::fmtflags old_flags = out.flags();
403 unsigned int old_precision = out.precision (precision);
405 out.precision (precision);
407 out.setf (std::ios::scientific, std::ios::floatfield);
409 out.setf (std::ios::fixed, std::ios::floatfield);
411 for (
unsigned int i = 0;
422 out <<
"[Proc" << i <<
" " << istart <<
"-" << iend-1 <<
"]" <<
' ';
423 for (PetscInt i=0; i<nlocal; ++i)
424 out << val[i] <<
' ';
428 out <<
"[Proc " << i <<
" " << istart <<
"-" << iend-1 <<
"]" << std::endl;
429 for (PetscInt i=0; i<nlocal; ++i)
430 out << val[i] << std::endl;
436 out.flags (old_flags);
437 out.precision(old_precision);
441 ierr = VecRestoreArray (
vector, &val);
451 DEAL_II_NAMESPACE_CLOSE
453 #endif // DEAL_II_WITH_PETSC
types::global_dof_index size_type
Vector & operator=(const Vector &v)
static ::ExceptionBase & ExcIO()
#define AssertThrow(cond, exc)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
std::pair< size_type, size_type > local_range() const
bool is_ascending_and_one_to_one(const MPI_Comm &communicator) const
static ::ExceptionBase & ExcMessage(std::string arg1)
size_type local_size() const
T sum(const T &t, const MPI_Comm &mpi_communicator)
void subtract_set(const IndexSet &other)
VectorOperation::values last_action
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
void fill_index_vector(std::vector< size_type > &indices) const
#define AssertThrowMPI(error_code)
bool has_ghost_elements() const
void reinit(const MPI_Comm &communicator, const size_type N, const size_type local_size, const bool omit_zeroing_entries=false)
IndexSet locally_owned_elements() const
virtual void create_vector(const size_type n, const size_type local_size)
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
static ::ExceptionBase & ExcNotImplemented()
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
size_type n_elements() const
static ::ExceptionBase & ExcInternalError()