16 #ifndef dealii__trilinos_block_sparse_matrix_h 17 #define dealii__trilinos_block_sparse_matrix_h 20 #include <deal.II/base/config.h> 22 #ifdef DEAL_II_WITH_TRILINOS 24 # include <deal.II/base/table.h> 25 # include <deal.II/base/template_constraints.h> 26 # include <deal.II/lac/block_matrix_base.h> 27 # include <deal.II/lac/trilinos_sparse_matrix.h> 28 # include <deal.II/lac/trilinos_block_vector.h> 29 # include <deal.II/lac/full_matrix.h> 30 # include <deal.II/lac/exceptions.h> 34 DEAL_II_NAMESPACE_OPEN
88 typedef BaseClass::pointer pointer;
89 typedef BaseClass::const_pointer const_pointer;
90 typedef BaseClass::reference reference;
91 typedef BaseClass::const_reference const_reference;
92 typedef BaseClass::size_type size_type;
147 const size_type n_block_columns);
154 template <
typename BlockSparsityPatternType>
155 void reinit (
const std::vector<Epetra_Map> &input_maps,
156 const BlockSparsityPatternType &block_sparsity_pattern,
157 const bool exchange_data =
false);
164 template <
typename BlockSparsityPatternType>
165 void reinit (
const std::vector<IndexSet> &input_maps,
166 const BlockSparsityPatternType &block_sparsity_pattern,
167 const MPI_Comm &communicator = MPI_COMM_WORLD,
168 const bool exchange_data =
false);
175 template <
typename BlockSparsityPatternType>
176 void reinit (
const BlockSparsityPatternType &block_sparsity_pattern);
186 void reinit (
const std::vector<Epetra_Map> &input_maps,
187 const ::BlockSparseMatrix<double> &deal_ii_sparse_matrix,
188 const double drop_tolerance=1e-13) DEAL_II_DEPRECATED;
198 const
double drop_tolerance=1e-13);
256 template <typename VectorType1, typename VectorType2>
257 void vmult (VectorType1 &dst,
258 const VectorType2 &src) const;
265 template <typename VectorType1, typename VectorType2>
266 void Tvmult (VectorType1 &dst,
267 const VectorType2 &src) const;
334 const MPI::
Vector &b) const;
374 << "The blocks [" << arg1 << ',' << arg2 << "] and ["
375 << arg3 << ',' << arg4 << "] have differing row
numbers.");
382 << "The blocks [" << arg1 << ',' << arg2 << "] and ["
383 << arg3 << ',' << arg4 << "] have differing column
numbers.");
390 template <typename VectorType1, typename VectorType2>
391 void vmult (VectorType1 &dst,
392 const VectorType2 &src,
393 const
bool transpose,
395 const ::
internal::bool2type<true>) const;
401 template <typename VectorType1, typename VectorType2>
402 void vmult (VectorType1 &dst,
403 const VectorType2 &src,
404 const
bool transpose,
406 const ::
internal::bool2type<true>) const;
412 template <typename VectorType1, typename VectorType2>
413 void vmult (VectorType1 &dst,
414 const VectorType2 &src,
415 const
bool transpose,
417 const ::
internal::bool2type<false>) const;
424 template <typename VectorType1, typename VectorType2>
425 void vmult (VectorType1 &dst,
426 const VectorType2 &src,
427 const
bool transpose,
429 const ::
internal::bool2type<false>) const;
448 this->
block(r,c) = d;
459 bool compressed =
true;
462 if (
block(row, col).is_compressed() ==
false)
473 template <
typename VectorType1,
typename VectorType2>
477 const VectorType2 &src)
const 479 vmult(dst, src,
false,
486 template <
typename VectorType1,
typename VectorType2>
490 const VectorType2 &src)
const 492 vmult(dst, src,
true,
499 template <
typename VectorType1,
typename VectorType2>
503 const VectorType2 &src,
504 const bool transpose,
508 if (transpose ==
true)
517 template <
typename VectorType1,
typename VectorType2>
521 const VectorType2 &src,
522 const bool transpose,
526 if (transpose ==
true)
534 template <
typename VectorType1,
typename VectorType2>
538 const VectorType2 &src,
539 const bool transpose,
543 if (transpose ==
true)
551 template <
typename VectorType1,
typename VectorType2>
555 const VectorType2 &src,
556 const bool transpose,
560 if (transpose ==
true)
567 #ifdef DEAL_II_WITH_CXX11 590 template<
typename PayloadBlockType>
608 template <
typename... Args>
612 "TrilinosBlockPayload can only accept a payload of type TrilinosPayload.");
619 #endif // DEAL_II_WITH_CXX11 624 DEAL_II_NAMESPACE_CLOSE
626 #endif // DEAL_II_WITH_TRILINOS 628 #endif // dealii__trilinos_block_sparse_matrix_h
void vmult_block_nonblock(BlockVectorType &dst, const VectorType &src) const
static ::ExceptionBase & ExcIncompatibleRowNumbers(int arg1, int arg2, int arg3, int arg4)
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
BlockMatrixBase< SparseMatrix > BaseClass
void vmult_nonblock_block(VectorType &dst, const BlockVectorType &src) const
void Tvmult_block_nonblock(BlockVectorType &dst, const VectorType &src) const
void Tvmult_block_block(BlockVectorType &dst, const BlockVectorType &src) const
void Tvmult_nonblock_block(VectorType &dst, const BlockVectorType &src) const
void vmult_block_block(BlockVectorType &dst, const BlockVectorType &src) const
void Tvmult_nonblock_nonblock(VectorType &dst, const VectorType &src) const
static ::ExceptionBase & ExcIncompatibleColNumbers(int arg1, int arg2, int arg3, int arg4)
BlockSparseMatrix & operator=(const BlockSparseMatrix &)
unsigned int n_block_cols() const
void vmult(VectorType1 &dst, const VectorType2 &src) const
#define Assert(cond, exc)
std::vector< Epetra_Map > range_partitioner() const 1
BaseClass::value_type value_type
void reinit(const size_type n_block_rows, const size_type n_block_columns)
BlockType::value_type value_type
bool is_compressed() const
void Tvmult(VectorType1 &dst, const VectorType2 &src) const
PayloadBlockType BlockType
TrilinosBlockPayload(const Args &...)
BaseClass::BlockType BlockType
void vmult_nonblock_nonblock(VectorType &dst, const VectorType &src) const
BlockType & block(const unsigned int row, const unsigned int column)
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
size_type n_nonzero_elements() const
std::vector< Epetra_Map > domain_partitioner() const 1
unsigned int n_block_rows() const
TrilinosScalar residual(MPI::BlockVector &dst, const MPI::BlockVector &x, const MPI::BlockVector &b) const