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
58 reinit (
const size_type n_block_rows,
59 const size_type n_block_columns)
71 this->column_block_indices.
reinit (n_block_columns, 0);
87 template <
typename BlockSparsityPatternType>
90 reinit (
const std::vector<Epetra_Map> ¶llel_partitioning,
91 const BlockSparsityPatternType &block_sparsity_pattern,
92 const bool exchange_data)
94 Assert (parallel_partitioning.size() == block_sparsity_pattern.n_block_rows(),
96 block_sparsity_pattern.n_block_rows()));
97 Assert (parallel_partitioning.size() == block_sparsity_pattern.n_block_cols(),
99 block_sparsity_pattern.n_block_cols()));
101 const size_type
n_block_rows = parallel_partitioning.size();
106 block_sparsity_pattern.n_block_rows()));
109 block_sparsity_pattern.n_block_cols()));
113 reinit (block_sparsity_pattern.n_block_rows(),
114 block_sparsity_pattern.n_block_cols());
118 this->column_block_indices = block_sparsity_pattern.get_column_indices();
126 parallel_partitioning[c],
127 block_sparsity_pattern.block(r,c),
134 template <
typename BlockSparsityPatternType>
137 reinit (
const std::vector<IndexSet> ¶llel_partitioning,
138 const BlockSparsityPatternType &block_sparsity_pattern,
139 const MPI_Comm &communicator,
140 const bool exchange_data)
142 std::vector<Epetra_Map> epetra_maps;
143 for (size_type i=0; i<block_sparsity_pattern.n_block_rows(); ++i)
144 epetra_maps.push_back
145 (parallel_partitioning[i].make_trilinos_map(communicator,
false));
147 reinit (epetra_maps, block_sparsity_pattern, exchange_data);
153 template <
typename BlockSparsityPatternType>
156 reinit (
const BlockSparsityPatternType &block_sparsity_pattern)
158 std::vector<Epetra_Map> parallel_partitioning;
159 for (size_type i=0; i<block_sparsity_pattern.n_block_rows(); ++i)
160 parallel_partitioning.push_back
161 (Epetra_Map(static_cast<TrilinosWrappers::types::int_type>(block_sparsity_pattern.block(i,0).n_rows()),
165 reinit (parallel_partitioning, block_sparsity_pattern);
197 reinit (
const std::vector<Epetra_Map> ¶llel_partitioning,
198 const ::BlockSparseMatrix<double> &dealii_block_sparse_matrix,
199 const double drop_tolerance)
201 const size_type
n_block_rows = parallel_partitioning.size();
205 dealii_block_sparse_matrix.n_block_rows()));
208 dealii_block_sparse_matrix.n_block_cols()));
219 parallel_partitioning[c],
220 dealii_block_sparse_matrix.block(r,c),
231 reinit (const ::BlockSparseMatrix<double> &dealii_block_sparse_matrix,
232 const double drop_tolerance)
234 Assert (dealii_block_sparse_matrix.n_block_rows() ==
235 dealii_block_sparse_matrix.n_block_cols(),
237 dealii_block_sparse_matrix.n_block_cols()));
238 Assert (dealii_block_sparse_matrix.m() ==
239 dealii_block_sparse_matrix.n(),
241 dealii_block_sparse_matrix.n()));
245 #ifdef DEAL_II_WITH_MPI 246 Epetra_MpiComm trilinos_communicator (MPI_COMM_SELF);
248 Epetra_SerialComm trilinos_communicator;
251 std::vector<Epetra_Map> parallel_partitioning;
252 for (size_type i=0; i<dealii_block_sparse_matrix.n_block_rows(); ++i)
253 parallel_partitioning.push_back (Epetra_Map(static_cast<TrilinosWrappers::types::int_type>(dealii_block_sparse_matrix.block(i,0).m()),
255 trilinos_communicator));
257 reinit (parallel_partitioning, dealii_block_sparse_matrix, drop_tolerance);
273 BlockSparseMatrix::size_type
276 size_type n_nonzero = 0;
277 for (size_type rows = 0; rows<this->
n_block_rows(); ++rows)
278 for (size_type cols = 0; cols<this->
n_block_cols(); ++cols)
388 std::vector<Epetra_Map>
403 std::vector<Epetra_Map>
431 const ::BlockSparsityPattern &,
435 const ::BlockDynamicSparsityPattern &,
440 const ::BlockDynamicSparsityPattern &,
447 DEAL_II_NAMESPACE_CLOSE
real_type l2_norm() const
Subscriptor & operator=(const Subscriptor &)
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
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)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::vector< Epetra_Map > range_partitioner() const 1
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
BaseClass::BlockType BlockType
BlockType & block(const unsigned int row, const unsigned int column)
Table< 2, SmartPointer< BlockType, BlockMatrixBase< SparseMatrix > > > sub_objects
size_type n_nonzero_elements() const
std::vector< Epetra_Map > domain_partitioner() const 1
BlockIndices row_block_indices
unsigned int n_block_rows() const
const BlockIndices & get_column_indices() const
TrilinosScalar residual(MPI::BlockVector &dst, const MPI::BlockVector &x, const MPI::BlockVector &b) const
real_type l2_norm() const
static ::ExceptionBase & ExcInternalError()