16 #include <deal.II/base/mpi.h> 17 #include <deal.II/lac/petsc_parallel_vector.h> 19 #ifdef DEAL_II_WITH_PETSC 24 DEAL_II_NAMESPACE_OPEN
32 : communicator (MPI_COMM_SELF)
46 communicator (communicator)
58 communicator (communicator)
76 const MPI_Comm &communicator)
78 communicator (communicator)
91 const MPI_Comm &communicator)
93 communicator (communicator)
127 ierr = VecGhostUpdateBegin(
vector, INSERT_VALUES, SCATTER_FORWARD);
129 ierr = VecGhostUpdateEnd(
vector, INSERT_VALUES, SCATTER_FORWARD);
152 const bool omit_zeroing_entries)
161 const int ierr = MPI_Allreduce (&k, &k_global, 1,
177 const PetscErrorCode ierr = VecDestroy (&
vector);
185 if (omit_zeroing_entries ==
false)
193 const bool omit_zeroing_entries)
198 if (!omit_zeroing_entries)
200 const PetscErrorCode ierr = VecSet(
vector, 0.0);
213 const MPI_Comm &comm)
215 const PetscErrorCode ierr = VecDestroy (&
vector);
230 const MPI_Comm &comm)
232 const PetscErrorCode ierr = VecDestroy (&
vector);
271 std::vector<size_type> ghostindices;
275 = (ghostindices.size() > 0
277 (
const PetscInt *)(&(ghostindices[0]))
297 ierr = VecGetOwnershipRange (
vector, &begin, &end);
303 ierr = VecGhostGetLocalForm(
vector, &l);
307 ierr = VecGetSize(l, &lsize);
310 ierr = VecGhostRestoreLocalForm(
vector, &l);
324 #if DEAL_II_PETSC_VERSION_LT(3,6,0) 338 #ifdef DEAL_II_WITH_MPI 342 return num_nonzero == 0;
344 return has_nonzero == 0;
351 const unsigned int precision,
352 const bool scientific,
353 const bool across)
const 360 PetscInt nlocal, istart, iend;
362 PetscErrorCode ierr = VecGetArray (
vector, &val);
365 ierr = VecGetLocalSize (
vector, &nlocal);
368 ierr = VecGetOwnershipRange (
vector, &istart, &iend);
372 std::ios::fmtflags old_flags = out.flags();
373 unsigned int old_precision = out.precision (precision);
375 out.precision (precision);
377 out.setf (std::ios::scientific, std::ios::floatfield);
379 out.setf (std::ios::fixed, std::ios::floatfield);
386 for (
unsigned int i = 0;
397 out <<
"[Proc" << i <<
" " << istart <<
"-" << iend-1 <<
"]" <<
' ';
398 for (PetscInt i=0; i<nlocal; ++i)
399 out << val[i] <<
' ';
403 out <<
"[Proc " << i <<
" " << istart <<
"-" << iend-1 <<
"]" << std::endl;
404 for (PetscInt i=0; i<nlocal; ++i)
405 out << val[i] << std::endl;
411 out.flags (old_flags);
412 out.precision(old_precision);
416 ierr = VecRestoreArray (
vector, &val);
426 DEAL_II_NAMESPACE_CLOSE
428 #endif // DEAL_II_WITH_PETSC
types::global_dof_index size_type
static ::ExceptionBase & ExcIO()
#define AssertThrow(cond, exc)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
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
Vector & operator=(const Vector &v)
#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()