16 #ifndef dealii_relaxation_block_h 17 #define dealii_relaxation_block_h 19 #include <deal.II/base/smartpointer.h> 20 #include <deal.II/base/subscriptor.h> 22 #include <deal.II/lac/precondition_block_base.h> 23 #include <deal.II/lac/sparsity_pattern.h> 24 #include <deal.II/lac/vector.h> 29 DEAL_II_NAMESPACE_OPEN
54 template <
typename MatrixType,
55 typename InverseNumberType =
typename MatrixType::value_type,
63 using number =
typename MatrixType::value_type;
173 std::vector<std::vector<unsigned int>>
order;
241 const VectorType &prev,
242 const VectorType &src,
243 const bool backward)
const;
286 template <
typename MatrixType,
287 typename InverseNumberType =
typename MatrixType::value_type,
302 using number =
typename MatrixType::value_type;
345 step(VectorType &dst,
const VectorType &rhs)
const;
351 Tstep(VectorType &dst,
const VectorType &rhs)
const;
358 vmult(VectorType &dst,
const VectorType &rhs)
const;
365 Tvmult(VectorType &dst,
const VectorType &rhs)
const;
384 template <
typename MatrixType,
385 typename InverseNumberType =
typename MatrixType::value_type,
400 using number =
typename MatrixType::value_type;
443 step(VectorType &dst,
const VectorType &rhs)
const;
449 Tstep(VectorType &dst,
const VectorType &rhs)
const;
456 vmult(VectorType &dst,
const VectorType &rhs)
const;
463 Tvmult(VectorType &dst,
const VectorType &rhs)
const;
482 template <
typename MatrixType,
483 typename InverseNumberType =
typename MatrixType::value_type,
493 using number =
typename MatrixType::value_type;
536 step(VectorType &dst,
const VectorType &rhs)
const;
542 Tstep(VectorType &dst,
const VectorType &rhs)
const;
549 vmult(VectorType &dst,
const VectorType &rhs)
const;
556 Tvmult(VectorType &dst,
const VectorType &rhs)
const;
560 DEAL_II_NAMESPACE_CLOSE
void vmult(VectorType &dst, const VectorType &rhs) const
std::size_t memory_consumption() const
void Tvmult(VectorType &dst, const VectorType &rhs) const
SmartPointer< const MatrixType, RelaxationBlock< MatrixType, InverseNumberType, VectorType > > A
SmartPointer< const AdditionalData, RelaxationBlock< MatrixType, InverseNumberType, VectorType > > additional_data
void vmult(VectorType &dst, const VectorType &rhs) const
PreconditionBlockBase< InverseNumberType >::Inversion inversion
typename MatrixType::value_type number
VectorType * temp_ghost_vector
void Tvmult(VectorType &dst, const VectorType &rhs) const
void block_kernel(const size_type block_begin, const size_type block_end)
void Tstep(VectorType &dst, const VectorType &rhs) const
AdditionalData(const double relaxation=1., const bool invert_diagonal=true, const bool same_diagonal=false, const typename PreconditionBlockBase< InverseNumberType >::Inversion inversion=PreconditionBlockBase< InverseNumberType >::gauss_jordan, const double threshold=0., VectorType *temp_ghost_vector=nullptr)
void step(VectorType &dst, const VectorType &rhs) const
InverseNumberType value_type
void vmult(VectorType &dst, const VectorType &rhs) const
void Tvmult(VectorType &dst, const VectorType &rhs) const
unsigned int global_dof_index
void Tstep(VectorType &dst, const VectorType &rhs) const
std::vector< std::vector< unsigned int > > order
void Tstep(VectorType &dst, const VectorType &rhs) const
types::global_dof_index size_type
SparsityPattern block_list
void step(VectorType &dst, const VectorType &rhs) const
Householder< InverseNumberType > & inverse_householder(size_type i)
void initialize(const MatrixType &A, const AdditionalData ¶meters)
void step(VectorType &dst, const VectorType &rhs) const
void do_step(VectorType &dst, const VectorType &prev, const VectorType &src, const bool backward) const