Reference documentation for deal.II version 8.5.1
trilinos_block_vector.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__trilinos_block_vector_h
17 #define dealii__trilinos_block_vector_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #ifdef DEAL_II_WITH_TRILINOS
23 
24 # include <deal.II/lac/trilinos_vector.h>
25 # include <deal.II/lac/trilinos_parallel_block_vector.h>
26 # include <deal.II/lac/block_indices.h>
27 # include <deal.II/lac/block_vector_base.h>
28 # include <deal.II/lac/exceptions.h>
29 # include <deal.II/lac/vector_type_traits.h>
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 // forward declaration
34 template <typename Number> class BlockVector;
35 
40 namespace TrilinosWrappers
41 {
42  // forward declaration
43  namespace MPI
44  {
45  class BlockVector;
46  }
47  class BlockVector;
48  class BlockSparseMatrix;
49 
50 
70  class BlockVector : public BlockVectorBase<Vector>
71  {
72  public:
77 
82 
86  typedef BaseClass::value_type value_type;
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;
94 
98  BlockVector ();
99 
105  explicit BlockVector (const std::vector<Epetra_Map> &partitioner) DEAL_II_DEPRECATED;
106 
112  explicit BlockVector (const std::vector<IndexSet> &partitioner,
113  const MPI_Comm &communicator = MPI_COMM_WORLD) DEAL_II_DEPRECATED;
114 
119  BlockVector (const MPI::BlockVector &V) DEAL_II_DEPRECATED;
120 
125  BlockVector (const BlockVector &V) DEAL_II_DEPRECATED;
126 
132  explicit BlockVector (const size_type num_blocks) DEAL_II_DEPRECATED;
133 
140  explicit BlockVector (const std::vector<size_type> &N) DEAL_II_DEPRECATED;
141 
150  template <typename InputIterator>
151  BlockVector (const std::vector<size_type> &n,
152  const InputIterator first,
153  const InputIterator end) DEAL_II_DEPRECATED;
154 
158  ~BlockVector ();
159 
164  BlockVector &
165  operator = (const value_type s);
166 
170  BlockVector &
171  operator = (const MPI::BlockVector &V);
172 
176  BlockVector &
177  operator = (const BlockVector &V);
178 
189  template <typename Number>
190  BlockVector &
191  operator = (const ::BlockVector<Number> &V);
192 
205  void reinit (const std::vector<Epetra_Map> &partitioning,
206  const bool omit_zeroing_entries = false);
207 
219  void reinit (const std::vector<IndexSet> &partitioning,
220  const MPI_Comm &communicator = MPI_COMM_WORLD,
221  const bool omit_zeroing_entries = false);
222 
231  void reinit (const std::vector<size_type> &N,
232  const bool omit_zeroing_entries=false);
233 
238  void reinit (const MPI::BlockVector &V);
239 
254  void reinit (const BlockVector &V,
255  const bool omit_zeroing_entries = false);
256 
263  void reinit (const size_type num_blocks);
264 
281  void swap (BlockVector &v);
282 
286  void print (std::ostream &out,
287  const unsigned int precision = 3,
288  const bool scientific = true,
289  const bool across = true) const;
290 
295 
300 
305  int, int,
306  << "For the generation of a localized vector the map has "
307  << "to assign all elements to all vectors! "
308  << "local_size = global_size is a necessary condition, but"
309  << arg1 << " != " << arg2 << " was given!");
310 
311  };
312 
313 
314 
315  /*----------------------- Inline functions ----------------------------------*/
316 
317 
318 
319  inline
321  {}
322 
323 
324 
325  inline
326  BlockVector::BlockVector (const std::vector<Epetra_Map> &partitioning)
327  {
328  reinit (partitioning);
329  }
330 
331 
332 
333  inline
334  BlockVector::BlockVector (const std::vector<IndexSet> &partitioning,
335  const MPI_Comm &communicator)
336  {
337  reinit (partitioning, communicator);
338  }
339 
340 
341 
342  inline
343  BlockVector::BlockVector (const std::vector<size_type> &N)
344  {
345  reinit (N);
346  }
347 
348 
349 
350  template <typename InputIterator>
351  BlockVector::BlockVector (const std::vector<size_type> &n,
352  const InputIterator first,
353  const InputIterator end)
354  {
355  // first set sizes of blocks, but
356  // don't initialize them as we will
357  // copy elements soon
358  (void)end;
359  reinit (n, true);
360  InputIterator start = first;
361  for (size_type b=0; b<n.size(); ++b)
362  {
363  InputIterator end = start;
364  std::advance (end, static_cast<size_type>(n[b]));
365 
366  for (size_type i=0; i<n[b]; ++i, ++start)
367  this->block(b)(i) = *start;
368  }
370  }
371 
372 
373 
374  inline
375  BlockVector::BlockVector (const size_type num_blocks)
376  {
377  reinit (num_blocks);
378  }
379 
380 
381 
382  inline
384  {}
385 
386 
387 
388  inline
390  {
391  reinit (v);
392  }
393 
394 
395 
396  inline
398  :
400  {
401  this->components.resize (v.n_blocks());
402  this->block_indices = v.block_indices;
403 
404  for (size_type i=0; i<this->n_blocks(); ++i)
405  this->components[i] = v.components[i];
406  }
407 
408 
409  inline
410  void
412  {
413  Assert (n_blocks() == v.n_blocks(),
415 
416  for (unsigned int row=0; row<n_blocks(); ++row)
417  block(row).swap (v.block(row));
418  }
419 
420 
421  template <typename Number>
422  BlockVector &
423  BlockVector::operator = (const ::BlockVector<Number> &v)
424  {
425  if (n_blocks() != v.n_blocks())
426  {
427  std::vector<size_type> block_sizes (v.n_blocks(), 0);
428  block_indices.reinit (block_sizes);
429  if (components.size() != n_blocks())
430  components.resize(n_blocks());
431  }
432 
433  for (size_type i=0; i<this->n_blocks(); ++i)
434  this->components[i] = v.block(i);
435 
436  collect_sizes();
437 
438  return *this;
439  }
440 
441 
450  inline
451  void swap (BlockVector &u,
452  BlockVector &v)
453  {
454  u.swap (v);
455  }
456 
457 } /* namespace TrilinosWrappers */
458 
462 namespace internal
463 {
464  namespace LinearOperator
465  {
466  template <typename> class ReinitHelper;
467 
472  template<>
474  {
475  public:
476  template <typename Matrix>
477  static
478  void reinit_range_vector (const Matrix &matrix,
480  bool omit_zeroing_entries)
481  {
482  v.reinit(matrix.range_partitioner(), omit_zeroing_entries);
483  }
484 
485  template <typename Matrix>
486  static
487  void reinit_domain_vector(const Matrix &matrix,
489  bool omit_zeroing_entries)
490  {
491  v.reinit(matrix.domain_partitioner(), omit_zeroing_entries);
492  }
493  };
494 
495  } /* namespace LinearOperator */
496 } /* namespace internal */
497 
498 
504 template <>
505 struct is_serial_vector< TrilinosWrappers::BlockVector > : std_cxx11::true_type
506 {
507 };
508 
509 
515 template <>
516 struct is_serial_vector< TrilinosWrappers::MPI::BlockVector > : std_cxx11::false_type
517 {
518 };
519 
520 
521 DEAL_II_NAMESPACE_CLOSE
522 
523 #endif // DEAL_II_WITH_TRILINOS
524 
525 #endif
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:576
static ::ExceptionBase & ExcNonMatchingBlockVectors()
static void reinit_domain_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
static ::ExceptionBase & ExcNonLocalizedMap(int arg1, int arg2)
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
#define DeclException0(Exception0)
Definition: exceptions.h:541
BlockVector & operator=(const value_type s)
std::vector< Vector > components
static void reinit_range_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
virtual void swap(Vector< Number > &v)
void swap(BlockVector &u, BlockVector &v)
unsigned int n_blocks() const
void reinit(const std::vector< Epetra_Map > &partitioning, const bool omit_zeroing_entries=false)
BlockVectorBase< Vector > BaseClass
BlockType & block(const unsigned int i)
static ::ExceptionBase & ExcIteratorRangeDoesNotMatchVectorSize()
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const