Reference documentation for deal.II version 9.0.0
trilinos_block_vector.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2018 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 #include <deal.II/lac/trilinos_parallel_block_vector.h>
17 
18 #ifdef DEAL_II_WITH_TRILINOS
19 
20 #include <deal.II/lac/trilinos_index_access.h>
21 
22 #include <deal.II/lac/trilinos_block_sparse_matrix.h>
23 
24 
25 DEAL_II_NAMESPACE_OPEN
26 
27 namespace TrilinosWrappers
28 {
29  namespace MPI
30  {
31  BlockVector &
33  {
35  return *this;
36  }
37 
38 
39 
40  BlockVector &
42  {
43  // we only allow assignment to vectors with the same number of blocks
44  // or to an empty BlockVector
45  Assert (n_blocks() == 0 || n_blocks() == v.n_blocks(),
47 
48  if (this->n_blocks() != v.n_blocks())
49  reinit(v.n_blocks());
50 
51  for (size_type i=0; i<this->n_blocks(); ++i)
52  this->components[i] = v.block(i);
53 
54  collect_sizes();
55 
56  return *this;
57  }
58 
59 
60 
61  BlockVector &
63  {
64  swap(v);
65  return *this;
66  }
67 
68 
69 
70  void
71  BlockVector::reinit (const std::vector<IndexSet> &parallel_partitioning,
72  const MPI_Comm &communicator,
73  const bool omit_zeroing_entries)
74  {
75  const size_type no_blocks = parallel_partitioning.size();
76  std::vector<size_type> block_sizes (no_blocks);
77 
78  for (size_type i=0; i<no_blocks; ++i)
79  {
80  block_sizes[i] = parallel_partitioning[i].size();
81  }
82 
83  this->block_indices.reinit (block_sizes);
84  if (components.size() != n_blocks())
85  components.resize(n_blocks());
86 
87  for (size_type i=0; i<n_blocks(); ++i)
88  components[i].reinit(parallel_partitioning[i], communicator, omit_zeroing_entries);
89 
90  collect_sizes();
91  }
92 
93  void
94  BlockVector::reinit (const std::vector<IndexSet> &parallel_partitioning,
95  const std::vector<IndexSet> &ghost_values,
96  const MPI_Comm &communicator,
97  const bool vector_writable)
98  {
99  const size_type no_blocks = parallel_partitioning.size();
100  std::vector<size_type> block_sizes (no_blocks);
101 
102  for (size_type i=0; i<no_blocks; ++i)
103  {
104  block_sizes[i] = parallel_partitioning[i].size();
105  }
106 
107  this->block_indices.reinit (block_sizes);
108  if (components.size() != n_blocks())
109  components.resize(n_blocks());
110 
111  for (size_type i=0; i<n_blocks(); ++i)
112  components[i].reinit(parallel_partitioning[i], ghost_values[i],
113  communicator, vector_writable);
114 
115  collect_sizes();
116  }
117 
118 
119  void
121  const bool omit_zeroing_entries)
122  {
124  if (components.size() != n_blocks())
125  components.resize(n_blocks());
126 
127  for (size_type i=0; i<n_blocks(); ++i)
128  components[i].reinit(v.block(i), omit_zeroing_entries, false);
129 
130  collect_sizes();
131  }
132 
133 
134 
135  void
136  BlockVector::reinit (const size_type num_blocks)
137  {
138  std::vector<size_type> block_sizes (num_blocks, 0);
139  this->block_indices.reinit (block_sizes);
140  if (this->components.size() != this->n_blocks())
141  this->components.resize(this->n_blocks());
142 
143  for (size_type i=0; i<this->n_blocks(); ++i)
144  components[i].clear();
145 
146  collect_sizes();
147  }
148 
149 
150 
151  void
154  const BlockVector &v)
155  {
156  Assert (m.n_block_rows() == v.n_blocks(),
158  Assert (m.n_block_cols() == v.n_blocks(),
160 
161  if (v.n_blocks() != n_blocks())
162  {
163  block_indices = v.get_block_indices();
164  components.resize(v.n_blocks());
165  }
166 
167  for (size_type i=0; i<this->n_blocks(); ++i)
168  components[i].import_nonlocal_data_for_fe(m.block(i,i), v.block(i));
169 
170  collect_sizes();
171  }
172 
173 
174 
175  void BlockVector::print (std::ostream &out,
176  const unsigned int precision,
177  const bool scientific,
178  const bool across) const
179  {
180  for (size_type i=0; i<this->n_blocks(); ++i)
181  {
182  if (across)
183  out << 'C' << i << ':';
184  else
185  out << "Component " << i << std::endl;
186  this->components[i].print(out, precision, scientific, across);
187  }
188  }
189 
190  } /* end of namespace MPI */
191 } /* end of namespace TrilinosWrappers */
192 
193 
194 DEAL_II_NAMESPACE_CLOSE
195 
196 #endif
void reinit(const std::vector< IndexSet > &parallel_partitioning, const MPI_Comm &communicator=MPI_COMM_WORLD, const bool omit_zeroing_entries=false)
unsigned int n_block_cols() const
void swap(BlockVector &u, BlockVector &v)
const BlockIndices & get_block_indices() const
#define Assert(cond, exc)
Definition: exceptions.h:1142
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)
std::vector< MPI::Vector > components
BlockVector & operator=(const value_type s)
unsigned int n_blocks() const
void import_nonlocal_data_for_fe(const TrilinosWrappers::BlockSparseMatrix &m, const BlockVector &v)
BlockType & block(const unsigned int row, const unsigned int column)
unsigned int n_block_rows() const
BlockType & block(const unsigned int i)
BlockVectorBase & operator=(const value_type s)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const