16 #include <deal.II/lac/trilinos_block_sparse_matrix.h> 18 #ifdef DEAL_II_WITH_TRILINOS 20 # include <deal.II/lac/block_sparse_matrix.h> 21 # include <deal.II/lac/block_sparsity_pattern.h> 23 DEAL_II_NAMESPACE_OPEN
43 const size_type n_block_columns)
54 this->column_block_indices.
reinit(n_block_columns, 0);
69 template <
typename BlockSparsityPatternType>
72 const std::vector<Epetra_Map> & parallel_partitioning,
73 const BlockSparsityPatternType &block_sparsity_pattern,
74 const bool exchange_data)
76 Assert(parallel_partitioning.size() ==
77 block_sparsity_pattern.n_block_rows(),
79 block_sparsity_pattern.n_block_rows()));
80 Assert(parallel_partitioning.size() ==
81 block_sparsity_pattern.n_block_cols(),
83 block_sparsity_pattern.n_block_cols()));
85 const size_type
n_block_rows = parallel_partitioning.size();
90 block_sparsity_pattern.n_block_rows()));
93 block_sparsity_pattern.n_block_cols()));
97 reinit(block_sparsity_pattern.n_block_rows(),
98 block_sparsity_pattern.n_block_cols());
102 this->column_block_indices = block_sparsity_pattern.get_column_indices();
110 parallel_partitioning[c],
111 block_sparsity_pattern.block(r, c),
118 template <
typename BlockSparsityPatternType>
121 const std::vector<IndexSet> & parallel_partitioning,
122 const BlockSparsityPatternType &block_sparsity_pattern,
123 const MPI_Comm & communicator,
124 const bool exchange_data)
126 std::vector<Epetra_Map> epetra_maps;
127 for (size_type i = 0; i < block_sparsity_pattern.n_block_rows(); ++i)
128 epetra_maps.push_back(
129 parallel_partitioning[i].make_trilinos_map(communicator,
false));
131 reinit(epetra_maps, block_sparsity_pattern, exchange_data);
136 template <
typename BlockSparsityPatternType>
139 const BlockSparsityPatternType &block_sparsity_pattern)
141 std::vector<Epetra_Map> parallel_partitioning;
142 for (size_type i = 0; i < block_sparsity_pattern.n_block_rows(); ++i)
143 parallel_partitioning.emplace_back(
144 static_cast<TrilinosWrappers::types::int_type>(
145 block_sparsity_pattern.block(i, 0).n_rows()),
149 reinit(parallel_partitioning, block_sparsity_pattern);
179 const std::vector<Epetra_Map> & parallel_partitioning,
180 const ::BlockSparseMatrix<double> &dealii_block_sparse_matrix,
181 const double drop_tolerance)
183 const size_type
n_block_rows = parallel_partitioning.size();
187 dealii_block_sparse_matrix.n_block_rows()));
190 dealii_block_sparse_matrix.n_block_cols()));
201 parallel_partitioning[c],
202 dealii_block_sparse_matrix.block(r,
214 const ::BlockSparseMatrix<double> &dealii_block_sparse_matrix,
215 const double drop_tolerance)
217 Assert(dealii_block_sparse_matrix.n_block_rows() ==
218 dealii_block_sparse_matrix.n_block_cols(),
220 dealii_block_sparse_matrix.n_block_cols()));
221 Assert(dealii_block_sparse_matrix.m() == dealii_block_sparse_matrix.n(),
223 dealii_block_sparse_matrix.n()));
227 # ifdef DEAL_II_WITH_MPI 228 Epetra_MpiComm trilinos_communicator(MPI_COMM_SELF);
230 Epetra_SerialComm trilinos_communicator;
233 std::vector<Epetra_Map> parallel_partitioning;
234 for (size_type i = 0; i < dealii_block_sparse_matrix.n_block_rows(); ++i)
235 parallel_partitioning.emplace_back(
236 static_cast<TrilinosWrappers::types::int_type>(
237 dealii_block_sparse_matrix.block(i, 0).m()),
239 trilinos_communicator);
241 reinit(parallel_partitioning, dealii_block_sparse_matrix, drop_tolerance);
255 BlockSparseMatrix::size_type
258 size_type n_nonzero = 0;
259 for (size_type rows = 0; rows < this->
n_block_rows(); ++rows)
260 for (size_type cols = 0; cols < this->
n_block_cols(); ++cols)
328 std::vector<Epetra_Map>
337 this->sub_objects[0][c]->domain_partitioner());
344 std::vector<Epetra_Map>
364 return this->
sub_objects[0][0]->get_mpi_communicator();
378 const ::BlockSparsityPattern &,
382 const ::BlockDynamicSparsityPattern &,
387 const ::BlockDynamicSparsityPattern &,
394 DEAL_II_NAMESPACE_CLOSE
std::vector< Epetra_Map > range_partitioner() const
real_type l2_norm() const
real_type l2_norm() const
~BlockSparseMatrix() override
size_type n_block_cols() const
size_type n_block_rows() const
static ::ExceptionBase & ExcNotInitialized()
const Epetra_Comm & comm_self()
SparsityPatternType & block(const size_type row, const size_type column)
std::size_t n_nonzero_elements() const
Table< 2, SmartPointer< BlockType, BlockMatrixBase< SparseMatrix > > > sub_objects
unsigned int n_block_cols() const
std::vector< Epetra_Map > domain_partitioner() const
void vmult(VectorType1 &dst, const VectorType2 &src) const
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
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)
void reinit(const size_type n_block_rows, const size_type n_block_columns)
const BlockIndices & get_row_indices() const
MPI_Comm get_mpi_communicator() const
BlockType & block(const unsigned int row, const unsigned int column)
size_type n_nonzero_elements() const
BlockIndices row_block_indices
unsigned int n_block_rows() const
BaseClass::BlockType BlockType
const BlockIndices & get_column_indices() const
TrilinosScalar residual(MPI::BlockVector &dst, const MPI::BlockVector &x, const MPI::BlockVector &b) const
static ::ExceptionBase & ExcInternalError()