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/trilinos_vector.h> 25 # include <deal.II/lac/block_indices.h> 26 # include <deal.II/lac/block_vector_base.h> 27 # include <deal.II/lac/exceptions.h> 29 # include <functional> 31 DEAL_II_NAMESPACE_OPEN
87 typedef BaseClass::pointer pointer;
88 typedef BaseClass::const_pointer const_pointer;
89 typedef BaseClass::reference reference;
90 typedef BaseClass::const_reference const_reference;
91 typedef BaseClass::size_type size_type;
106 explicit BlockVector (
const std::vector<IndexSet> ¶llel_partitioning,
107 const MPI_Comm &communicator = MPI_COMM_WORLD);
114 BlockVector (
const std::vector<IndexSet> ¶llel_partitioning,
115 const std::vector<IndexSet> &ghost_values,
116 const MPI_Comm &communicator,
117 const bool vector_writable =
false);
170 template <
typename Number>
181 void reinit (
const std::vector<IndexSet> ¶llel_partitioning,
182 const MPI_Comm &communicator = MPI_COMM_WORLD,
183 const bool omit_zeroing_entries =
false);
202 void reinit (
const std::vector<IndexSet> &partitioning,
203 const std::vector<IndexSet> &ghost_values,
204 const MPI_Comm &communicator = MPI_COMM_WORLD,
205 const bool vector_writable =
false);
223 const bool omit_zeroing_entries =
false);
231 void reinit (
const size_type num_blocks);
283 void print (std::ostream &out,
284 const unsigned int precision = 3,
285 const bool scientific =
true,
286 const bool across =
true)
const;
304 const MPI_Comm &communicator)
306 reinit (parallel_partitioning, communicator,
false);
313 const std::vector<IndexSet> &ghost_values,
314 const MPI_Comm &communicator,
315 const bool vector_writable)
317 reinit(parallel_partitioning, ghost_values, communicator,
339 for (size_type i=0; i<this->
n_blocks(); ++i)
355 template <
typename Number>
361 std::vector<size_type> block_sizes (v.n_blocks(), 0);
367 for (size_type i=0; i<this->
n_blocks(); ++i)
381 bool ghosted=
block(0).has_ghost_elements();
383 for (
unsigned int i=0; i<this->
n_blocks(); ++i)
426 namespace LinearOperatorImplementation
428 template <
typename>
class ReinitHelper;
438 template <
typename Matrix>
442 bool omit_zeroing_entries)
444 v.
reinit(matrix.locally_owned_range_indices(),
445 matrix.get_mpi_communicator(), omit_zeroing_entries);
448 template <
typename Matrix>
452 bool omit_zeroing_entries)
454 v.
reinit(matrix.locally_owned_domain_indices(),
455 matrix.get_mpi_communicator(), omit_zeroing_entries);
473 DEAL_II_NAMESPACE_CLOSE
475 #endif // DEAL_II_WITH_TRILINOS BaseClass::value_type value_type
static ::ExceptionBase & ExcIteratorRangeDoesNotMatchVectorSize()
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
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)
BaseClass::BlockType BlockType
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
::BlockVectorBase< MPI::Vector > BaseClass
static ::ExceptionBase & ExcInternalError()