15#ifndef dealii_tpetra_trilinos_block_sparse_matrix_h
16#define dealii_tpetra_trilinos_block_sparse_matrix_h
21#ifdef DEAL_II_TRILINOS_WITH_TPETRA
39template <
typename number>
49 namespace TpetraWrappers
74 template <
typename Number,
typename MemorySpace = ::MemorySpace::Host>
159 template <
typename BlockSparsityPatternType>
161 reinit(
const std::vector<IndexSet> &input_maps,
162 const BlockSparsityPatternType &block_sparsity_pattern,
163 const MPI_Comm communicator = MPI_COMM_WORLD,
164 const bool exchange_data =
false);
171 template <
typename BlockSparsityPatternType>
173 reinit(
const BlockSparsityPatternType &block_sparsity_pattern);
183 const std::vector<IndexSet> ¶llel_partitioning,
184 const ::BlockSparseMatrix<double> &dealii_block_sparse_matrix,
185 const MPI_Comm communicator = MPI_COMM_WORLD,
186 const double drop_tolerance = 1e-13);
196 reinit(const ::BlockSparseMatrix<double> &deal_ii_sparse_matrix,
197 const double drop_tolerance = 1e-13);
241 std::vector<IndexSet>
249 std::vector<IndexSet>
258 template <
typename VectorType1,
typename VectorType2>
260 vmult(VectorType1 &dst,
const VectorType2 &src)
const;
267 template <
typename VectorType1,
typename VectorType2>
269 Tvmult(VectorType1 &dst,
const VectorType2 &src)
const;
283 template <
typename VectorType1,
284 typename VectorType2,
285 typename VectorType3>
288 const VectorType2 &x,
289 const VectorType3 &b)
const;
301 template <
typename VectorType1,
typename VectorType2>
303 vmult(VectorType1 &dst,
304 const VectorType2 &src,
306 const std::bool_constant<true>,
307 const std::bool_constant<true>)
const;
313 template <
typename VectorType1,
typename VectorType2>
315 vmult(VectorType1 &dst,
316 const VectorType2 &src,
318 const std::bool_constant<false>,
319 const std::bool_constant<true>)
const;
325 template <
typename VectorType1,
typename VectorType2>
327 vmult(VectorType1 &dst,
328 const VectorType2 &src,
330 const std::bool_constant<true>,
331 const std::bool_constant<false>)
const;
338 template <
typename VectorType1,
typename VectorType2>
340 vmult(VectorType1 &dst,
341 const VectorType2 &src,
343 const std::bool_constant<false>,
344 const std::bool_constant<false>)
const;
352 namespace TpetraWrappers
354 template <
typename Number,
typename MemorySpace>
355 template <
typename VectorType1,
typename VectorType2>
358 const VectorType2 &src)
const
369 template <
typename Number,
typename MemorySpace>
370 template <
typename VectorType1,
typename VectorType2>
373 const VectorType2 &src)
const
384 template <
typename Number,
typename MemorySpace>
385 template <
typename VectorType1,
typename VectorType2>
389 const VectorType2 &src,
391 std::bool_constant<true>,
392 std::bool_constant<true>)
const
395 BaseClass::Tvmult_block_block(dst, src);
397 BaseClass::vmult_block_block(dst, src);
402 template <
typename Number,
typename MemorySpace>
403 template <
typename VectorType1,
typename VectorType2>
407 const VectorType2 &src,
409 std::bool_constant<false>,
410 std::bool_constant<true>)
const
413 BaseClass::Tvmult_nonblock_block(dst, src);
415 BaseClass::vmult_nonblock_block(dst, src);
420 template <
typename Number,
typename MemorySpace>
421 template <
typename VectorType1,
typename VectorType2>
425 const VectorType2 &src,
427 std::bool_constant<true>,
428 std::bool_constant<false>)
const
431 BaseClass::Tvmult_block_nonblock(dst, src);
433 BaseClass::vmult_block_nonblock(dst, src);
438 template <
typename Number,
typename MemorySpace>
439 template <
typename VectorType1,
typename VectorType2>
443 const VectorType2 &src,
445 std::bool_constant<false>,
446 std::bool_constant<false>)
const
449 BaseClass::Tvmult_nonblock_nonblock(dst, src);
451 BaseClass::vmult_nonblock_nonblock(dst, src);
const value_type & const_reference
unsigned int n_block_rows() const
types::global_dof_index size_type
typename BlockType::value_type value_type
MatrixIterator< BlockMatrixIterators::Accessor< BlockMatrixBase, true > > const_iterator
const value_type * const_pointer
MatrixIterator< BlockMatrixIterators::Accessor< BlockMatrixBase, false > > iterator
typename BaseClass::value_type value_type
void reinit(const std::vector< IndexSet > ¶llel_partitioning, const ::BlockSparseMatrix< double > &dealii_block_sparse_matrix, const MPI_Comm communicator=MPI_COMM_WORLD, const double drop_tolerance=1e-13)
void Tvmult(VectorType1 &dst, const VectorType2 &src) const
BlockSparseMatrix< Number, MemorySpace > & operator=(const BlockSparseMatrix< Number, MemorySpace > &)=default
BlockSparseMatrix()=default
void vmult(VectorType1 &dst, const VectorType2 &src) const
std::vector< IndexSet > locally_owned_domain_indices() const
std::vector< IndexSet > locally_owned_range_indices() const
typename BaseClass::const_iterator const_iterator
void reinit(const ::BlockSparseMatrix< double > &deal_ii_sparse_matrix, const double drop_tolerance=1e-13)
typename BaseClass::const_pointer const_pointer
Number residual(VectorType1 &dst, const VectorType2 &x, const VectorType3 &b) const
bool is_compressed() const
std::uint64_t n_nonzero_elements() const
typename BaseClass::pointer pointer
~BlockSparseMatrix() override
typename BaseClass::BlockType BlockType
void reinit(const BlockSparsityPatternType &block_sparsity_pattern)
typename BaseClass::const_reference const_reference
typename BaseClass::iterator iterator
BlockSparseMatrix< Number, MemorySpace > & operator=(const Number d)
typename BaseClass::size_type size_type
void reinit(const size_type n_block_rows, const size_type n_block_columns)
typename BaseClass::reference reference
MPI_Comm get_mpi_communicator() const
void reinit(const std::vector< IndexSet > &input_maps, const BlockSparsityPatternType &block_sparsity_pattern, const MPI_Comm communicator=MPI_COMM_WORLD, const bool exchange_data=false)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE