Reference documentation for deal.II version GIT 5983d193e2 2023-05-27 16:50:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
trilinos_block_vector.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2021 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
17 
18 #ifdef DEAL_II_WITH_TRILINOS
19 
22 
23 
25 
26 namespace TrilinosWrappers
27 {
28  namespace MPI
29  {
30  BlockVector &
32  {
34  return *this;
35  }
36 
37 
38 
39  BlockVector &
41  {
42  // we only allow assignment to vectors with the same number of blocks
43  // or to an empty BlockVector
44  Assert(this->n_blocks() == 0 || this->n_blocks() == v.n_blocks(),
45  ExcDimensionMismatch(this->n_blocks(), v.n_blocks()));
46 
47  if (this->n_blocks() != v.n_blocks())
48  this->block_indices = v.block_indices;
49 
50  this->components.resize(this->n_blocks());
51  for (unsigned int i = 0; i < this->n_blocks(); ++i)
52  this->components[i] = v.components[i];
53 
54  this->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  // update the number of blocks
76  this->block_indices.reinit(parallel_partitioning.size(), 0);
77 
78  // initialize each block
79  this->components.resize(this->n_blocks());
80  for (unsigned int i = 0; i < this->n_blocks(); ++i)
81  this->components[i].reinit(parallel_partitioning[i],
82  communicator,
83  omit_zeroing_entries);
84 
85  // update block_indices content
86  this->collect_sizes();
87  }
88 
89 
90 
91  void
92  BlockVector::reinit(const std::vector<IndexSet> &parallel_partitioning,
93  const std::vector<IndexSet> &ghost_values,
94  const MPI_Comm communicator,
95  const bool vector_writable)
96  {
97  AssertDimension(parallel_partitioning.size(), ghost_values.size());
98 
99  // update the number of blocks
100  this->block_indices.reinit(parallel_partitioning.size(), 0);
101 
102  // initialize each block
103  this->components.resize(this->n_blocks());
104  for (unsigned int i = 0; i < this->n_blocks(); ++i)
105  this->components[i].reinit(parallel_partitioning[i],
106  ghost_values[i],
107  communicator,
108  vector_writable);
109 
110  // update block_indices content
111  this->collect_sizes();
112  }
113 
114 
115 
116  void
117  BlockVector::reinit(const BlockVector &v, const bool omit_zeroing_entries)
118  {
119  if (this->n_blocks() != v.n_blocks())
120  this->block_indices = v.get_block_indices();
121 
122  this->components.resize(this->n_blocks());
123  for (unsigned int i = 0; i < this->n_blocks(); ++i)
124  this->components[i].reinit(v.components[i],
125  omit_zeroing_entries,
126  false);
127 
128  this->collect_sizes();
129  }
130 
131 
132 
133  void
134  BlockVector::reinit(const size_type num_blocks)
135  {
136  this->block_indices.reinit(num_blocks, 0);
137 
138  this->components.resize(this->n_blocks());
139  for (auto &c : this->components)
140  c.clear();
141  }
142 
143 
144 
145  void
148  const BlockVector & v)
149  {
152 
153  if (this->n_blocks() != v.n_blocks())
154  this->block_indices = v.get_block_indices();
155 
156  this->components.resize(this->n_blocks());
157  for (unsigned int i = 0; i < this->n_blocks(); ++i)
159  v.block(i));
160 
161  this->collect_sizes();
162  }
163 
164 
165 
166  void
167  BlockVector::print(std::ostream & out,
168  const unsigned int precision,
169  const bool scientific,
170  const bool across) const
171  {
172  for (unsigned int i = 0; i < this->n_blocks(); ++i)
173  {
174  if (across)
175  out << 'C' << i << ':';
176  else
177  out << "Component " << i << std::endl;
178  this->components[i].print(out, precision, scientific, across);
179  }
180  }
181 
182  } /* end of namespace MPI */
183 } /* end of namespace TrilinosWrappers */
184 
185 
187 
188 #endif
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
BlockType & block(const unsigned int row, const unsigned int column)
unsigned int n_block_rows() const
unsigned int n_block_cols() const
const BlockIndices & get_block_indices() const
unsigned int n_blocks() const
BlockVectorBase & operator=(const value_type s)
std::vector< MPI::Vector > components
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
BlockVector & operator=(const value_type s)
void reinit(const std::vector< IndexSet > &parallel_partitioning, const MPI_Comm communicator=MPI_COMM_WORLD, const bool omit_zeroing_entries=false)
void import_nonlocal_data_for_fe(const TrilinosWrappers::BlockSparseMatrix &m, const BlockVector &v)
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:474
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:475
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define Assert(cond, exc)
Definition: exceptions.h:1614
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1787
void swap(BlockVector &u, BlockVector &v)