deal.II version 9.7.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_sparse_matrix.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2008 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
16
17#ifdef DEAL_II_WITH_TRILINOS
18
21
23
24namespace TrilinosWrappers
25{
27 {
28 // delete previous content of
29 // the subobjects array
30 try
31 {
32 clear();
33 }
34 catch (...)
35 {}
36 }
37
38
39
40# ifndef DOXYGEN
41 void
42 BlockSparseMatrix::reinit(const size_type n_block_rows,
43 const size_type n_block_columns)
44 {
45 // first delete previous content of
46 // the subobjects array
47 clear();
48
49 // then resize. set sizes of blocks to
50 // zero. user will later have to call
51 // collect_sizes for this
52 this->sub_objects.reinit(n_block_rows, n_block_columns);
53 this->row_block_indices.reinit(n_block_rows, 0);
54 this->column_block_indices.reinit(n_block_columns, 0);
55
56 // and reinitialize the blocks
57 for (size_type r = 0; r < this->n_block_rows(); ++r)
58 for (size_type c = 0; c < this->n_block_cols(); ++c)
59 {
60 BlockType *p = new BlockType();
61
62 Assert(this->sub_objects[r][c] == nullptr, ExcInternalError());
63 this->sub_objects[r][c] = p;
64 }
65 }
66# endif
67
68
69
70 template <typename BlockSparsityPatternType>
71 void
73 const std::vector<IndexSet> &parallel_partitioning,
74 const BlockSparsityPatternType &block_sparsity_pattern,
75 const MPI_Comm communicator,
76 const bool exchange_data)
77 {
78 std::vector<Epetra_Map> epetra_maps;
79 for (size_type i = 0; i < block_sparsity_pattern.n_block_rows(); ++i)
80 epetra_maps.push_back(
81 parallel_partitioning[i].make_trilinos_map(communicator, false));
82
83 Assert(epetra_maps.size() == block_sparsity_pattern.n_block_rows(),
84 ExcDimensionMismatch(epetra_maps.size(),
85 block_sparsity_pattern.n_block_rows()));
86 Assert(epetra_maps.size() == block_sparsity_pattern.n_block_cols(),
87 ExcDimensionMismatch(epetra_maps.size(),
88 block_sparsity_pattern.n_block_cols()));
89
90 const size_type n_block_rows = epetra_maps.size();
91 Assert(n_block_rows == block_sparsity_pattern.n_block_rows(),
93 block_sparsity_pattern.n_block_rows()));
94 Assert(n_block_rows == block_sparsity_pattern.n_block_cols(),
96 block_sparsity_pattern.n_block_cols()));
97
98 // Call the other basic reinit function, ...
99 reinit(block_sparsity_pattern.n_block_rows(),
100 block_sparsity_pattern.n_block_cols());
101
102 // ... set the correct sizes, ...
103 this->row_block_indices = block_sparsity_pattern.get_row_indices();
104 this->column_block_indices = block_sparsity_pattern.get_column_indices();
105
106 // ... and then assign the correct
107 // data to the blocks.
108 for (size_type r = 0; r < this->n_block_rows(); ++r)
109 for (size_type c = 0; c < this->n_block_cols(); ++c)
110 {
111 this->sub_objects[r][c]->reinit(parallel_partitioning[r],
112 parallel_partitioning[c],
113 block_sparsity_pattern.block(r, c),
114 communicator,
115 exchange_data);
116 }
117 }
118
119
120
121 template <typename BlockSparsityPatternType>
122 void
124 const BlockSparsityPatternType &block_sparsity_pattern)
125 {
126 std::vector<IndexSet> parallel_partitioning;
127 for (size_type i = 0; i < block_sparsity_pattern.n_block_rows(); ++i)
128 parallel_partitioning.emplace_back(
129 complete_index_set(block_sparsity_pattern.block(i, 0).n_rows()));
130
131 reinit(parallel_partitioning, block_sparsity_pattern);
132 }
133
134
135
136 template <>
137 void
138 BlockSparseMatrix::reinit(const BlockSparsityPattern &block_sparsity_pattern)
139 {
140 // Call the other basic reinit function, ...
141 reinit(block_sparsity_pattern.n_block_rows(),
142 block_sparsity_pattern.n_block_cols());
143
144 // ... set the correct sizes, ...
145 this->row_block_indices = block_sparsity_pattern.get_row_indices();
146 this->column_block_indices = block_sparsity_pattern.get_column_indices();
147
148 // ... and then assign the correct
149 // data to the blocks.
150 for (size_type r = 0; r < this->n_block_rows(); ++r)
151 for (size_type c = 0; c < this->n_block_cols(); ++c)
152 {
153 this->sub_objects[r][c]->reinit(block_sparsity_pattern.block(r, c));
154 }
155 }
156
157
158
159 void
161 const std::vector<IndexSet> &parallel_partitioning,
162 const ::BlockSparseMatrix<double> &dealii_block_sparse_matrix,
163 const MPI_Comm communicator,
164 const double drop_tolerance)
165 {
166 const size_type n_block_rows = parallel_partitioning.size();
167
168 Assert(n_block_rows == dealii_block_sparse_matrix.n_block_rows(),
170 dealii_block_sparse_matrix.n_block_rows()));
171 Assert(n_block_rows == dealii_block_sparse_matrix.n_block_cols(),
173 dealii_block_sparse_matrix.n_block_cols()));
174
175 // Call the other basic reinit function ...
177
178 // ... and then assign the correct
179 // data to the blocks.
180 for (size_type r = 0; r < this->n_block_rows(); ++r)
181 for (size_type c = 0; c < this->n_block_cols(); ++c)
182 {
183 this->sub_objects[r][c]->reinit(parallel_partitioning[r],
184 parallel_partitioning[c],
185 dealii_block_sparse_matrix.block(r,
186 c),
187 communicator,
188 drop_tolerance);
189 }
190
192 }
193
194
195
196 void
198 const ::BlockSparseMatrix<double> &dealii_block_sparse_matrix,
199 const double drop_tolerance)
200 {
201 Assert(dealii_block_sparse_matrix.n_block_rows() ==
202 dealii_block_sparse_matrix.n_block_cols(),
203 ExcDimensionMismatch(dealii_block_sparse_matrix.n_block_rows(),
204 dealii_block_sparse_matrix.n_block_cols()));
205 Assert(dealii_block_sparse_matrix.m() == dealii_block_sparse_matrix.n(),
206 ExcDimensionMismatch(dealii_block_sparse_matrix.m(),
207 dealii_block_sparse_matrix.n()));
208
209 std::vector<IndexSet> parallel_partitioning;
210 for (size_type i = 0; i < dealii_block_sparse_matrix.n_block_rows(); ++i)
211 parallel_partitioning.emplace_back(
212 complete_index_set(dealii_block_sparse_matrix.block(i, 0).m()));
213
214 reinit(parallel_partitioning,
215 dealii_block_sparse_matrix,
216 MPI_COMM_SELF,
217 drop_tolerance);
218 }
219
220
221
222 void
224 {
225 // simply forward to the (non-public) function of the base class
227 }
228
229
230
231 std::uint64_t
233 {
234 std::uint64_t n_nonzero = 0;
235 for (size_type rows = 0; rows < this->n_block_rows(); ++rows)
236 for (size_type cols = 0; cols < this->n_block_cols(); ++cols)
237 n_nonzero += this->block(rows, cols).n_nonzero_elements();
238
239 return n_nonzero;
240 }
241
242
243
246 const MPI::BlockVector &x,
247 const MPI::BlockVector &b) const
248 {
249 vmult(dst, x);
250 dst -= b;
251 dst *= -1.;
252
253 return dst.l2_norm();
254 }
255
256
257
258 // TODO: In the following we
259 // use the same code as just
260 // above three more times. Use
261 // templates.
264 const MPI::Vector &x,
265 const MPI::BlockVector &b) const
266 {
267 vmult(dst, x);
268 dst -= b;
269 dst *= -1.;
270
271 return dst.l2_norm();
272 }
273
274
275
278 const MPI::BlockVector &x,
279 const MPI::Vector &b) const
280 {
281 vmult(dst, x);
282 dst -= b;
283 dst *= -1.;
284
285 return dst.l2_norm();
286 }
287
288
289
292 const MPI::Vector &x,
293 const MPI::Vector &b) const
294 {
295 vmult(dst, x);
296 dst -= b;
297 dst *= -1.;
298
299 return dst.l2_norm();
300 }
301
302
303
306 {
307 Assert(this->n_block_cols() != 0, ExcNotInitialized());
308 Assert(this->n_block_rows() != 0, ExcNotInitialized());
309 return this->sub_objects[0][0]->get_mpi_communicator();
310 }
311
312
313
314# ifndef DOXYGEN
315 // -------------------- explicit instantiations -----------------------
316 //
317 template void
318 BlockSparseMatrix::reinit(const ::BlockSparsityPattern &);
319 template void
320 BlockSparseMatrix::reinit(const ::BlockDynamicSparsityPattern &);
321
322 template void
323 BlockSparseMatrix::reinit(const std::vector<IndexSet> &,
324 const ::BlockDynamicSparsityPattern &,
325 const MPI_Comm,
326 const bool);
327# endif // DOXYGEN
328
329} // namespace TrilinosWrappers
330
331
333
334#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)
Table< 2, ObserverPointer< BlockType, BlockMatrixBase< SparseMatrix > > > sub_objects
SparsityPatternType & block(const size_type row, const size_type column)
const BlockIndices & get_column_indices() const
const BlockIndices & get_row_indices() const
real_type l2_norm() const
std::size_t n_nonzero_elements() const
TrilinosScalar residual(MPI::BlockVector &dst, const MPI::BlockVector &x, const MPI::BlockVector &b) const
void reinit(const size_type n_block_rows, const size_type n_block_columns)
void vmult(VectorType1 &dst, const VectorType2 &src) const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcNotInitialized()
IndexSet complete_index_set(const IndexSet::size_type N)
Definition index_set.h:1219
double TrilinosScalar
Definition types.h:198