16 #ifndef dealii_la_parallel_block_vector_h 17 #define dealii_la_parallel_block_vector_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/exceptions.h> 22 #include <deal.II/lac/block_indices.h> 23 #include <deal.II/lac/block_vector_base.h> 24 #include <deal.II/lac/la_parallel_vector.h> 25 #include <deal.II/lac/vector_operation.h> 26 #include <deal.II/lac/vector_type_traits.h> 31 DEAL_II_NAMESPACE_OPEN
34 #ifdef DEAL_II_WITH_PETSC 44 #ifdef DEAL_II_WITH_TRILINOS 80 template <
typename Number>
113 typedef typename BaseClass::pointer pointer;
114 typedef typename BaseClass::const_pointer const_pointer;
115 typedef typename BaseClass::reference reference;
116 typedef typename BaseClass::const_reference const_reference;
117 typedef typename BaseClass::size_type size_type;
136 explicit BlockVector (
const size_type num_blocks = 0,
137 const size_type block_size = 0);
157 template <
typename OtherNumber>
165 BlockVector (
const std::vector<size_type> &block_sizes);
171 BlockVector (
const std::vector<IndexSet> &local_ranges,
172 const std::vector<IndexSet> &ghost_indices,
173 const MPI_Comm communicator);
178 BlockVector (
const std::vector<IndexSet> &local_ranges,
179 const MPI_Comm communicator);
198 template <
class Number2>
208 #ifdef DEAL_II_WITH_PETSC 220 #ifdef DEAL_II_WITH_TRILINOS 246 void reinit (
const size_type num_blocks,
247 const size_type block_size = 0,
248 const bool omit_zeroing_entries =
false);
269 void reinit (
const std::vector<size_type> &N,
270 const bool omit_zeroing_entries=
false);
286 template <
typename Number2>
288 const bool omit_zeroing_entries=
false);
345 template <
typename OtherNumber>
346 void add (
const std::vector<size_type> &indices,
347 const ::Vector<OtherNumber> &values);
353 void sadd (
const Number s,
371 void sadd (
const Number s,
382 virtual bool all_zero ()
const override;
393 real_type
lp_norm (
const real_type p)
const;
425 const bool omit_zeroing_entries =
false)
override;
457 std::shared_ptr<const CommunicationPatternBase> communication_pattern =
458 std::shared_ptr<const CommunicationPatternBase> ())
override;
480 template <
typename FullMatrixType>
483 const bool symmetric =
false)
const;
500 template <
typename FullMatrixType>
503 const bool symmetric =
false)
const;
514 template <
typename FullMatrixType>
516 const FullMatrixType &matrix,
517 const Number s = Number(0.),
518 const Number b = Number(1.))
const;
523 virtual void add(
const Number a)
override;
540 virtual void add (
const std::vector<size_type> &indices,
541 const std::vector<Number> &values);
547 virtual void sadd(
const Number s,
const Number a,
566 virtual real_type
l1_norm()
const override;
572 virtual real_type
l2_norm()
const override;
612 virtual size_type
size()
const override;
630 virtual void print(std::ostream &out,
631 const unsigned int precision=3,
632 const bool scientific=
true,
633 const bool across=
true)
const override;
675 template <
typename Number>
689 template <
typename Number>
695 DEAL_II_NAMESPACE_CLOSE
698 #include <deal.II/lac/la_parallel_block_vector.templates.h> static ::ExceptionBase & ExcIteratorRangeDoesNotMatchVectorSize()
virtual std::size_t memory_consumption() const override
virtual BlockVector< Number > & operator-=(const VectorSpaceVector< Number > &V) override
virtual void compress(::VectorOperation::values operation) override
virtual BlockVector & operator=(const value_type s) override
real_type lp_norm(const real_type p) const
virtual BlockVector< Number > & operator/=(const Number factor) override
virtual real_type l2_norm() const override
virtual Number mean_value() const override
virtual real_type linfty_norm() const override
virtual bool all_zero() const override
void reinit(const size_type num_blocks, const size_type block_size=0, const bool omit_zeroing_entries=false)
void equ(const Number a, const BlockVector< Number > &u, const Number b, const BlockVector< Number > &v)
BaseClass::BlockType BlockType
BlockVectorBase< Vector< Number > > BaseClass
static ::ExceptionBase & ExcVectorTypeNotCompatible()
virtual Number operator*(const VectorSpaceVector< Number > &V) const override
BaseClass::value_type value_type
BlockType::real_type real_type
Vector< Number > BlockType
void swap(LinearAlgebra::distributed::BlockVector< Number > &u, LinearAlgebra::distributed::BlockVector< Number > &v)
real_type norm_sqr() const
#define DeclException0(Exception0)
virtual real_type l1_norm() const override
void swap(BlockVector< Number > &v)
void update_ghost_values() const
void mmult(BlockVector< Number > &V, const FullMatrixType &matrix, const Number s=Number(0.), const Number b=Number(1.)) const
virtual void scale(const VectorSpaceVector< Number > &scaling_factors) override
BlockVector(const size_type num_blocks=0, const size_type block_size=0)
Number multivector_inner_product_with_metric(const FullMatrixType &matrix, const BlockVector< Number > &V, const bool symmetric=false) const
virtual size_type size() const override
void zero_out_ghosts() const
virtual void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const override
void add(const std::vector< size_type > &indices, const ::Vector< OtherNumber > &values)
virtual BlockVector< Number > & operator+=(const VectorSpaceVector< Number > &V) override
virtual ::IndexSet locally_owned_elements() const override
void multivector_inner_product(FullMatrixType &matrix, const BlockVector< Number > &V, const bool symmetric=false) const
void sadd(const Number s, const BlockVector< Number > &V)
static constexpr unsigned int communication_block_size
virtual BlockVector< Number > & operator*=(const Number factor) override
virtual Number add_and_dot(const Number a, const VectorSpaceVector< Number > &V, const VectorSpaceVector< Number > &W) override
bool has_ghost_elements() const