16 #ifndef dealii_block_vector_h 17 #define dealii_block_vector_h 20 #include <deal.II/base/config.h> 22 #include <deal.II/base/exceptions.h> 24 #include <deal.II/lac/block_indices.h> 25 #include <deal.II/lac/block_vector_base.h> 26 #include <deal.II/lac/vector_operation.h> 27 #include <deal.II/lac/vector_type_traits.h> 32 DEAL_II_NAMESPACE_OPEN
35 #ifdef DEAL_II_WITH_TRILINOS 67 template <
typename Number>
84 using value_type =
typename BaseClass::value_type;
86 using pointer =
typename BaseClass::pointer;
87 using const_pointer =
typename BaseClass::const_pointer;
88 using reference =
typename BaseClass::reference;
89 using const_reference =
typename BaseClass::const_reference;
90 using size_type =
typename BaseClass::size_type;
105 const size_type block_size = 0);
132 template <typename OtherNumber>
135 #ifdef DEAL_II_WITH_TRILINOS 147 BlockVector(
const std::vector<size_type> &block_sizes);
164 template <
typename InputIterator>
165 BlockVector(
const std::vector<size_type> &block_sizes,
166 const InputIterator first,
167 const InputIterator
end);
222 template <
class Number2>
232 #ifdef DEAL_II_WITH_TRILINOS 255 const size_type block_size = 0,
256 const bool omit_zeroing_entries =
false);
275 reinit(
const std::vector<size_type> &block_sizes,
276 const bool omit_zeroing_entries =
false);
289 const bool omit_zeroing_entries =
false);
304 template <
typename Number2>
307 const bool omit_zeroing_entries =
false);
313 template <
class BlockVector2>
315 scale(
const BlockVector2 &v);
335 print(std::ostream & out,
336 const unsigned int precision = 3,
337 const bool scientific =
true,
338 const bool across =
true)
const;
381 template <
typename Number>
382 template <
typename InputIterator>
384 const InputIterator first,
385 const InputIterator end)
391 reinit(block_sizes,
true);
392 InputIterator start = first;
393 for (size_type b = 0; b < block_sizes.size(); ++b)
395 InputIterator end = start;
396 std::advance(end, static_cast<signed int>(block_sizes[b]));
397 std::copy(start, end, this->block(b).begin());
400 Assert(start == end, ExcIteratorRangeDoesNotMatchVectorSize());
405 template <
typename Number>
411 BaseClass::operator=(s);
417 template <
typename Number>
422 BaseClass::operator=(v);
428 template <
typename Number>
432 BaseClass::operator=(v);
438 template <
typename Number>
439 template <
typename Number2>
444 BaseClass::operator=(v);
448 template <
typename Number>
452 for (size_type i = 0; i < this->n_blocks(); ++i)
453 this->components[i].compress(operation);
458 template <
typename Number>
467 template <
typename Number>
468 template <
class BlockVector2>
486 template <
typename Number>
496 namespace LinearOperatorImplementation
505 template <
typename number>
509 template <
typename Matrix>
513 bool omit_zeroing_entries)
515 v.
reinit(matrix.get_row_indices(), omit_zeroing_entries);
518 template <
typename Matrix>
522 bool omit_zeroing_entries)
524 v.
reinit(matrix.get_column_indices(), omit_zeroing_entries);
537 template <
typename Number>
541 DEAL_II_NAMESPACE_CLOSE
BlockVector & operator=(const value_type s)
void block_write(std::ostream &out) const
bool has_ghost_elements() const
void scale(const BlockVector2 &v)
static void reinit_range_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
void swap(BlockVector< Number > &u, BlockVector< Number > &v)
__global__ void scale(Number *val, const Number *V_val, const size_type N)
static void reinit_domain_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
BlockIndices block_indices
static ::ExceptionBase & ExcIteratorRangeDoesNotMatchVectorSize()
#define Assert(cond, exc)
void compress(::VectorOperation::values operation=::VectorOperation::unknown)
#define DeclException0(Exception0)
void reinit(const unsigned int n_blocks, const size_type block_size=0, const bool omit_zeroing_entries=false)
void swap(BlockVector< Number > &v)
unsigned int n_blocks() const
BlockVector(const unsigned int n_blocks=0, const size_type block_size=0)
typename BlockType::real_type real_type
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
#define AssertIsFinite(number)
void block_read(std::istream &in)
~BlockVector() override=default