16 #ifndef dealii_trilinos_parallel_block_vector_h 17 #define dealii_trilinos_parallel_block_vector_h 20 #include <deal.II/base/config.h> 22 #ifdef DEAL_II_WITH_TRILINOS 24 # include <deal.II/lac/block_indices.h> 25 # include <deal.II/lac/block_vector_base.h> 26 # include <deal.II/lac/exceptions.h> 27 # include <deal.II/lac/trilinos_vector.h> 29 # include <functional> 31 DEAL_II_NAMESPACE_OPEN
34 template <
typename Number>
88 using pointer = BaseClass::pointer;
89 using const_pointer = BaseClass::const_pointer;
90 using reference = BaseClass::reference;
91 using const_reference = BaseClass::const_reference;
92 using size_type = BaseClass::size_type;
107 explicit BlockVector(
const std::vector<IndexSet> ¶llel_partitioning,
108 const MPI_Comm &communicator = MPI_COMM_WORLD);
115 BlockVector(
const std::vector<IndexSet> ¶llel_partitioning,
116 const std::vector<IndexSet> &ghost_values,
117 const MPI_Comm & communicator,
118 const bool vector_writable =
false);
174 template <
typename Number>
176 operator=(const ::BlockVector<Number> &v);
187 reinit(
const std::vector<IndexSet> ¶llel_partitioning,
188 const MPI_Comm & communicator = MPI_COMM_WORLD,
189 const bool omit_zeroing_entries =
false);
209 reinit(
const std::vector<IndexSet> &partitioning,
210 const std::vector<IndexSet> &ghost_values,
211 const MPI_Comm & communicator = MPI_COMM_WORLD,
212 const bool vector_writable =
false);
239 reinit(
const size_type num_blocks);
295 print(std::ostream & out,
296 const unsigned int precision = 3,
297 const bool scientific =
true,
298 const bool across =
true)
const;
315 const std::vector<IndexSet> ¶llel_partitioning,
316 const MPI_Comm & communicator)
318 reinit(parallel_partitioning, communicator,
false);
324 const std::vector<IndexSet> ¶llel_partitioning,
325 const std::vector<IndexSet> &ghost_values,
326 const MPI_Comm & communicator,
327 const bool vector_writable)
329 reinit(parallel_partitioning,
350 for (size_type i = 0; i < this->
n_blocks(); ++i)
365 template <
typename Number>
371 std::vector<size_type> block_sizes(v.n_blocks(), 0);
377 for (size_type i = 0; i < this->
n_blocks(); ++i)
390 bool ghosted =
block(0).has_ghost_elements();
392 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
433 namespace LinearOperatorImplementation
446 template <
typename Matrix>
450 bool omit_zeroing_entries)
452 v.
reinit(matrix.locally_owned_range_indices(),
453 matrix.get_mpi_communicator(),
454 omit_zeroing_entries);
457 template <
typename Matrix>
461 bool omit_zeroing_entries)
463 v.
reinit(matrix.locally_owned_domain_indices(),
464 matrix.get_mpi_communicator(),
465 omit_zeroing_entries);
482 DEAL_II_NAMESPACE_CLOSE
484 #endif // DEAL_II_WITH_TRILINOS
static ::ExceptionBase & ExcIteratorRangeDoesNotMatchVectorSize()
~BlockVector() override=default
static ::ExceptionBase & ExcNonMatchingBlockVectors()
void reinit(const std::vector< IndexSet > ¶llel_partitioning, const MPI_Comm &communicator=MPI_COMM_WORLD, const bool omit_zeroing_entries=false)
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)
BlockIndices block_indices
BaseClass::BlockType BlockType
void swap(BlockVector &u, BlockVector &v)
#define Assert(cond, exc)
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
#define DeclException0(Exception0)
std::vector< MPI::Vector > components
void swap(Vector< Number > &u, Vector< Number > &v)
BlockVector & operator=(const value_type s)
unsigned int n_blocks() const
bool has_ghost_elements() const
void import_nonlocal_data_for_fe(const TrilinosWrappers::BlockSparseMatrix &m, const BlockVector &v)
void swap(BlockVector &v)
LinearAlgebra::distributed::BlockVector< Number > BlockVector
BlockType & block(const unsigned int i)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
BaseClass::value_type value_type
static ::ExceptionBase & ExcInternalError()