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
54 const size_type n_block_columns)
66 this->column_block_indices.
reinit (n_block_columns, 0);
82 template <
typename BlockSparsityPatternType>
85 reinit (
const std::vector<Epetra_Map> ¶llel_partitioning,
86 const BlockSparsityPatternType &block_sparsity_pattern,
87 const bool exchange_data)
89 Assert (parallel_partitioning.size() == block_sparsity_pattern.n_block_rows(),
90 ExcDimensionMismatch (parallel_partitioning.size(),
91 block_sparsity_pattern.n_block_rows()));
92 Assert (parallel_partitioning.size() == block_sparsity_pattern.n_block_cols(),
93 ExcDimensionMismatch (parallel_partitioning.size(),
94 block_sparsity_pattern.n_block_cols()));
96 const size_type n_block_rows = parallel_partitioning.size();
99 Assert (n_block_rows == block_sparsity_pattern.n_block_rows(),
100 ExcDimensionMismatch (n_block_rows,
101 block_sparsity_pattern.n_block_rows()));
102 Assert (n_block_rows == block_sparsity_pattern.n_block_cols(),
103 ExcDimensionMismatch (n_block_rows,
104 block_sparsity_pattern.n_block_cols()));
108 reinit (block_sparsity_pattern.n_block_rows(),
109 block_sparsity_pattern.n_block_cols());
113 this->column_block_indices = block_sparsity_pattern.get_column_indices();
121 parallel_partitioning[c],
122 block_sparsity_pattern.block(r,c),
129 template <
typename BlockSparsityPatternType>
132 reinit (
const std::vector<IndexSet> ¶llel_partitioning,
133 const BlockSparsityPatternType &block_sparsity_pattern,
134 const MPI_Comm &communicator,
135 const bool exchange_data)
137 std::vector<Epetra_Map> epetra_maps;
138 for (size_type i=0; i<block_sparsity_pattern.n_block_rows(); ++i)
139 epetra_maps.push_back
140 (parallel_partitioning[i].make_trilinos_map(communicator,
false));
142 reinit (epetra_maps, block_sparsity_pattern, exchange_data);
148 template <
typename BlockSparsityPatternType>
151 reinit (
const BlockSparsityPatternType &block_sparsity_pattern)
153 std::vector<Epetra_Map> parallel_partitioning;
154 for (size_type i=0; i<block_sparsity_pattern.n_block_rows(); ++i)
155 parallel_partitioning.push_back
156 (Epetra_Map(static_cast<TrilinosWrappers::types::int_type>(block_sparsity_pattern.block(i,0).n_rows()),
160 reinit (parallel_partitioning, block_sparsity_pattern);
168 reinit (
const BlockSparsityPattern &block_sparsity_pattern)
172 reinit (block_sparsity_pattern.n_block_rows(),
173 block_sparsity_pattern.n_block_cols());
177 this->column_block_indices = block_sparsity_pattern.get_column_indices();
192 reinit (
const std::vector<Epetra_Map> ¶llel_partitioning,
193 const ::BlockSparseMatrix<double> &dealii_block_sparse_matrix,
194 const double drop_tolerance)
196 const size_type n_block_rows = parallel_partitioning.size();
198 Assert (n_block_rows == dealii_block_sparse_matrix.n_block_rows(),
199 ExcDimensionMismatch (n_block_rows,
200 dealii_block_sparse_matrix.n_block_rows()));
201 Assert (n_block_rows == dealii_block_sparse_matrix.n_block_cols(),
202 ExcDimensionMismatch (n_block_rows,
203 dealii_block_sparse_matrix.n_block_cols()));
206 reinit (n_block_rows, n_block_rows);
214 parallel_partitioning[c],
215 dealii_block_sparse_matrix.block(r,c),
226 reinit (const ::BlockSparseMatrix<double> &dealii_block_sparse_matrix,
227 const double drop_tolerance)
229 Assert (dealii_block_sparse_matrix.n_block_rows() ==
230 dealii_block_sparse_matrix.n_block_cols(),
231 ExcDimensionMismatch (dealii_block_sparse_matrix.n_block_rows(),
232 dealii_block_sparse_matrix.n_block_cols()));
233 Assert (dealii_block_sparse_matrix.m() ==
234 dealii_block_sparse_matrix.n(),
235 ExcDimensionMismatch (dealii_block_sparse_matrix.m(),
236 dealii_block_sparse_matrix.n()));
240 #ifdef DEAL_II_WITH_MPI 241 Epetra_MpiComm trilinos_communicator (MPI_COMM_SELF);
243 Epetra_SerialComm trilinos_communicator;
246 std::vector<Epetra_Map> parallel_partitioning;
247 for (size_type i=0; i<dealii_block_sparse_matrix.n_block_rows(); ++i)
248 parallel_partitioning.push_back (Epetra_Map(static_cast<TrilinosWrappers::types::int_type>(dealii_block_sparse_matrix.block(i,0).m()),
250 trilinos_communicator));
252 reinit (parallel_partitioning, dealii_block_sparse_matrix, drop_tolerance);
268 BlockSparseMatrix::size_type
271 size_type n_nonzero = 0;
272 for (size_type rows = 0; rows<this->
n_block_rows(); ++rows)
273 for (size_type cols = 0; cols<this->
n_block_cols(); ++cols)
383 std::vector<Epetra_Map>
391 domain_partitioner.push_back(this->sub_objects[0][c]->domain_partitioner());
398 std::vector<Epetra_Map>
406 range_partitioner.push_back(this->sub_objects[r][0]->range_partitioner());
426 const ::BlockSparsityPattern &,
430 const ::BlockDynamicSparsityPattern &,
435 const ::BlockDynamicSparsityPattern &,
442 DEAL_II_NAMESPACE_CLOSE
size_type n_nonzero_elements() const
real_type l2_norm() const
Subscriptor & operator=(const Subscriptor &)
const Epetra_Comm & comm_self()
BlockSparseMatrix & operator=(const BlockSparseMatrix &)
unsigned int n_block_cols() 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)
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)
std::vector< Epetra_Map > range_partitioner() const DEAL_II_DEPRECATED
BaseClass::BlockType BlockType
std::vector< Epetra_Map > domain_partitioner() const DEAL_II_DEPRECATED
BlockType & block(const unsigned int row, const unsigned int column)
Table< 2, SmartPointer< BlockType, BlockMatrixBase< SparseMatrix > > > sub_objects
size_type n_nonzero_elements() const
BlockIndices row_block_indices
unsigned int n_block_rows() const
TrilinosScalar residual(MPI::BlockVector &dst, const MPI::BlockVector &x, const MPI::BlockVector &b) const
real_type l2_norm() const