16#ifndef dealii_trilinos_block_sparse_matrix_h
17#define dealii_trilinos_block_sparse_matrix_h
22#ifdef DEAL_II_WITH_TRILINOS
39template <
typename number>
155 template <
typename BlockSparsityPatternType>
157 reinit(
const std::vector<IndexSet> & input_maps,
158 const BlockSparsityPatternType &block_sparsity_pattern,
159 const MPI_Comm communicator = MPI_COMM_WORLD,
160 const bool exchange_data =
false);
167 template <
typename BlockSparsityPatternType>
169 reinit(
const BlockSparsityPatternType &block_sparsity_pattern);
179 const std::vector<IndexSet> & parallel_partitioning,
180 const ::BlockSparseMatrix<double> &dealii_block_sparse_matrix,
181 const MPI_Comm communicator = MPI_COMM_WORLD,
182 const double drop_tolerance = 1e-13);
192 reinit(const ::BlockSparseMatrix<double> &deal_ii_sparse_matrix,
193 const double drop_tolerance = 1e-13);
236 std::vector<IndexSet>
244 std::vector<IndexSet>
253 template <
typename VectorType1,
typename VectorType2>
255 vmult(VectorType1 &dst,
const VectorType2 &src)
const;
262 template <
typename VectorType1,
typename VectorType2>
264 Tvmult(VectorType1 &dst,
const VectorType2 &src)
const;
338 <<
"The blocks [" << arg1 <<
',' << arg2 <<
"] and [" << arg3
339 <<
',' << arg4 <<
"] have differing row numbers.");
349 <<
"The blocks [" << arg1 <<
',' << arg2 <<
"] and [" << arg3
350 <<
',' << arg4 <<
"] have differing column numbers.");
357 template <
typename VectorType1,
typename VectorType2>
359 vmult(VectorType1 & dst,
360 const VectorType2 &src,
362 const std::integral_constant<bool, true>,
363 const std::integral_constant<bool, true>)
const;
369 template <
typename VectorType1,
typename VectorType2>
371 vmult(VectorType1 & dst,
372 const VectorType2 &src,
374 const std::integral_constant<bool, false>,
375 const std::integral_constant<bool, true>)
const;
381 template <
typename VectorType1,
typename VectorType2>
383 vmult(VectorType1 & dst,
384 const VectorType2 &src,
386 const std::integral_constant<bool, true>,
387 const std::integral_constant<bool, false>)
const;
394 template <
typename VectorType1,
typename VectorType2>
396 vmult(VectorType1 & dst,
397 const VectorType2 &src,
399 const std::integral_constant<bool, false>,
400 const std::integral_constant<bool, false>)
const;
418 this->
block(r, c) = d;
428 bool compressed =
true;
431 if (
block(row, col).is_compressed() ==
false)
442 template <
typename VectorType1,
typename VectorType2>
455 template <
typename VectorType1,
typename VectorType2>
468 template <
typename VectorType1,
typename VectorType2>
471 const VectorType2 &src,
473 std::integral_constant<bool, true>,
474 std::integral_constant<bool, true>)
const
484 template <
typename VectorType1,
typename VectorType2>
487 const VectorType2 &src,
489 std::integral_constant<bool, false>,
490 std::integral_constant<bool, true>)
const
500 template <
typename VectorType1,
typename VectorType2>
503 const VectorType2 &src,
505 std::integral_constant<bool, true>,
506 std::integral_constant<bool, false>)
const
516 template <
typename VectorType1,
typename VectorType2>
519 const VectorType2 &src,
521 std::integral_constant<bool, false>,
522 std::integral_constant<bool, false>)
const
532 inline std::vector<IndexSet>
538 std::vector<IndexSet> domain_indices;
540 domain_indices.push_back(
541 this->sub_objects[0][c]->locally_owned_domain_indices());
543 return domain_indices;
548 inline std::vector<IndexSet>
554 std::vector<IndexSet> range_indices;
556 range_indices.push_back(
557 this->sub_objects[r][0]->locally_owned_range_indices());
559 return range_indices;
566 namespace BlockLinearOperatorImplementation
582 template <
typename PayloadBlockType>
600 template <
typename... Args>
607 "TrilinosBlockPayload can only accept a payload of type TrilinosPayload.");
void Tvmult_nonblock_nonblock(VectorType &dst, const VectorType &src) const
void vmult_nonblock_nonblock(VectorType &dst, const VectorType &src) const
void Tvmult_block_block(BlockVectorType &dst, const BlockVectorType &src) const
const value_type & const_reference
void Tvmult_nonblock_block(VectorType &dst, const BlockVectorType &src) const
unsigned int n_block_rows() const
MatrixIterator< BlockMatrixIterators::Accessor< BlockMatrixBase, true > > const_iterator
void vmult_block_block(BlockVectorType &dst, const BlockVectorType &src) const
types::global_dof_index size_type
typename BlockType::value_type value_type
void vmult_nonblock_block(VectorType &dst, const BlockVectorType &src) const
unsigned int n_block_cols() const
const value_type * const_pointer
MatrixIterator< BlockMatrixIterators::Accessor< BlockMatrixBase, false > > iterator
BlockType & block(const unsigned int row, const unsigned int column)
void vmult_block_nonblock(BlockVectorType &dst, const VectorType &src) const
void Tvmult_block_nonblock(BlockVectorType &dst, const VectorType &src) const
BaseClass::value_type value_type
BlockSparseMatrix & operator=(const BlockSparseMatrix &)=default
BaseClass::size_type size_type
BlockSparseMatrix()=default
bool is_compressed() const
~BlockSparseMatrix() override
BaseClass::const_reference const_reference
TrilinosScalar residual(MPI::BlockVector &dst, const MPI::BlockVector &x, const MPI::BlockVector &b) const
MPI_Comm get_mpi_communicator() const
void Tvmult(VectorType1 &dst, const VectorType2 &src) const
BaseClass::pointer pointer
std::uint64_t n_nonzero_elements() const
std::vector< IndexSet > locally_owned_domain_indices() const
BaseClass::const_pointer const_pointer
void reinit(const size_type n_block_rows, const size_type n_block_columns)
std::vector< IndexSet > locally_owned_range_indices() const
void vmult(VectorType1 &dst, const VectorType2 &src) const
BaseClass::reference reference
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcIncompatibleRowNumbers(int arg1, int arg2, int arg3, int arg4)
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
TrilinosBlockPayload(const Args &...)
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
#define Assert(cond, exc)
static ::ExceptionBase & ExcIncompatibleColNumbers(int arg1, int arg2, int arg3, int arg4)
PayloadBlockType BlockType
static ::ExceptionBase & ExcNotInitialized()