16 #ifndef dealii_petsc_parallel_vector_h 17 #define dealii_petsc_parallel_vector_h 20 #include <deal.II/base/config.h> 22 #ifdef DEAL_II_WITH_PETSC 24 # include <deal.II/base/subscriptor.h> 25 # include <deal.II/lac/exceptions.h> 26 # include <deal.II/lac/vector.h> 27 # include <deal.II/lac/vector_operation.h> 28 # include <deal.II/lac/petsc_vector_base.h> 29 # include <deal.II/base/index_set.h> 30 # include <deal.II/lac/vector_type_traits.h> 32 DEAL_II_NAMESPACE_OPEN
202 template <
typename Number>
204 const ::Vector<Number> &v,
295 template <
typename number>
317 const bool omit_zeroing_entries =
false);
329 const bool omit_zeroing_entries =
false);
369 void print (std::ostream &out,
370 const unsigned int precision = 3,
371 const bool scientific =
true,
372 const bool across =
true)
const;
432 template <
typename number>
434 const ::Vector<number> &v,
435 const size_type local_size)
437 communicator (communicator)
457 template <
typename number>
462 Assert (size() == v.size(),
491 for (size_type i=0; i<v.size(); ++i)
503 Vector::get_mpi_communicator ()
const 514 namespace LinearOperatorImplementation
516 template <
typename>
class ReinitHelper;
526 template <
typename Matrix>
532 v.
reinit(matrix.locally_owned_range_indices(), matrix.get_mpi_communicator());
535 template <
typename Matrix>
541 v.
reinit(matrix.locally_owned_domain_indices(), matrix.get_mpi_communicator());
562 DEAL_II_NAMESPACE_CLOSE
564 #endif // DEAL_II_WITH_PETSC
types::global_dof_index size_type
VectorBase & operator=(const PetscScalar s)
Vector< Number > & operator=(const Number s)
static void reinit_range_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
static void reinit_domain_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
void swap(BlockVector &u, BlockVector &v)
size_type local_size() const
unsigned int global_dof_index
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
const MPI_Comm & get_mpi_communicator() const
Vector & operator=(const Vector &v)
void reinit(const MPI_Comm &communicator, const size_type N, const size_type local_size, const bool omit_zeroing_entries=false)
virtual void create_vector(const size_type n, const size_type local_size)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
void swap(Vector &u, Vector &v)