Reference documentation for deal.II version 9.5.0
\(\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\}}\)
Loading...
Searching...
No Matches
trilinos_block_vector.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2008 - 2023 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
26namespace TrilinosWrappers
27{
28 namespace MPI
29 {
32 {
34 return *this;
35 }
36
37
38
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(),
46
47 if (this->n_blocks() != v.n_blocks())
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
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
118 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
119 & partitioners,
120 const bool make_ghosted,
121 const bool vector_writable)
122 {
123 // update the number of blocks
124 this->block_indices.reinit(partitioners.size(), 0);
125
126 // initialize each block
127 this->components.resize(this->n_blocks());
128 for (unsigned int i = 0; i < this->n_blocks(); ++i)
129 this->components[i].reinit(partitioners[i],
130 make_ghosted,
131 vector_writable);
132
133 // update block_indices content
134 this->collect_sizes();
135 }
136
137
138
139 void
140 BlockVector::reinit(const BlockVector &v, const bool omit_zeroing_entries)
141 {
142 if (this->n_blocks() != v.n_blocks())
144
145 this->components.resize(this->n_blocks());
146 for (unsigned int i = 0; i < this->n_blocks(); ++i)
147 this->components[i].reinit(v.components[i],
148 omit_zeroing_entries,
149 false);
150
151 this->collect_sizes();
152 }
153
154
155
156 void
158 {
159 this->block_indices.reinit(num_blocks, 0);
160
161 this->components.resize(this->n_blocks());
162 for (auto &c : this->components)
163 c.clear();
164 }
165
166
167
168 void
171 const BlockVector & v)
172 {
175
176 if (this->n_blocks() != v.n_blocks())
178
179 this->components.resize(this->n_blocks());
180 for (unsigned int i = 0; i < this->n_blocks(); ++i)
182 v.block(i));
183
184 this->collect_sizes();
185 }
186
187
188
189 void
190 BlockVector::print(std::ostream & out,
191 const unsigned int precision,
192 const bool scientific,
193 const bool across) const
194 {
195 for (unsigned int i = 0; i < this->n_blocks(); ++i)
196 {
197 if (across)
198 out << 'C' << i << ':';
199 else
200 out << "Component " << i << std::endl;
201 this->components[i].print(out, precision, scientific, across);
202 }
203 }
204
205 } /* end of namespace MPI */
206} /* end of namespace TrilinosWrappers */
207
208
210
211#endif
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
unsigned int n_block_rows() const
unsigned int n_block_cols() const
BlockType & block(const unsigned int row, const unsigned int column)
unsigned int n_blocks() const
BlockVectorBase & operator=(const value_type s)
std::vector< MPI::Vector > components
BlockType & block(const unsigned int i)
const BlockIndices & get_block_indices() const
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:472
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void swap(BlockVector &u, BlockVector &v)