16 #ifndef dealii_filtered_matrix_h 17 #define dealii_filtered_matrix_h 21 #include <deal.II/base/config.h> 22 #include <deal.II/base/smartpointer.h> 23 #include <deal.II/base/thread_management.h> 24 #include <deal.II/base/memory_consumption.h> 25 #include <deal.II/lac/pointer_matrix.h> 26 #include <deal.II/lac/vector_memory.h> 30 DEAL_II_NAMESPACE_OPEN
32 template <
class VectorType>
class FilteredMatrixBlock;
194 template <
typename VectorType>
231 double value()
const;
285 const Accessor *operator-> ()
const;
345 template <
typename MatrixType>
347 const bool expect_constrained_source =
false);
363 template <
typename MatrixType>
364 void initialize (
const MatrixType &m,
365 const bool expect_constrained_source =
false);
380 void add_constraint (
const size_type i,
const double v);
397 template <
class Constra
intList>
398 void add_constraints (
const ConstraintList &new_constraints);
403 void clear_constraints ();
417 void apply_constraints (VectorType &v,
418 const bool matrix_is_symmetric)
const;
423 void apply_constraints (VectorType &v)
const;
429 void vmult (VectorType &dst,
430 const VectorType &src)
const;
437 void Tvmult (VectorType &dst,
438 const VectorType &src)
const;
447 void vmult_add (VectorType &dst,
448 const VectorType &src)
const;
457 void Tvmult_add (VectorType &dst,
458 const VectorType &src)
const;
480 std::size_t memory_consumption ()
const;
519 std::shared_ptr<PointerMatrixBase<VectorType> >
matrix;
531 void pre_filter (VectorType &v)
const;
538 void post_filter (
const VectorType &in,
539 VectorType &out)
const;
545 friend class FilteredMatrixBlock<VectorType>;
554 template <
typename VectorType>
563 Assert (index <= matrix->constraints.size(),
569 template <
typename VectorType>
579 template <
typename VectorType>
584 return matrix->constraints[index].first;
589 template <
typename VectorType>
594 return matrix->constraints[index].second;
599 template <
typename VectorType>
611 template <
typename VectorType>
622 template <
typename VectorType>
632 template <
typename number>
641 template <
typename number>
650 template <
typename number>
661 template <
typename number>
667 return ! (*
this == other);
674 template <
typename number>
683 template <
typename number>
692 template <
typename VectorType>
699 return (i1.first < i2.first);
704 template <
typename VectorType>
705 template <
typename MatrixType>
717 template <
typename VectorType>
726 template <
typename VectorType>
731 expect_constrained_source(fm.expect_constrained_source),
733 constraints (fm.constraints)
738 template <
typename VectorType>
739 template <
typename MatrixType>
745 expect_constrained_source (false)
752 template <
typename VectorType>
765 template <
typename VectorType>
776 template <
typename VectorType>
777 template <
class Constra
intList>
783 const size_type old_size = constraints.size();
784 constraints.reserve (old_size + new_constraints.size());
785 constraints.insert (constraints.end(),
786 new_constraints.begin(),
787 new_constraints.end());
790 std::inplace_merge (constraints.begin(),
791 constraints.begin()+old_size,
798 template <
typename VectorType>
804 std::vector<IndexValuePair> empty;
805 constraints.swap (empty);
810 template <
typename VectorType>
821 template <
typename VectorType>
828 apply_constraints(v);
832 template <
typename VectorType>
839 tmp_vector->reinit(v);
845 (*tmp_vector)(i->first) = -i->second;
851 matrix->vmult_add(v, *tmp_vector);
854 for (i=constraints.begin(); i!=e; ++i)
857 v(i->first) = i->second;
862 template <
typename VectorType>
877 template <
typename VectorType>
881 VectorType &out)
const 890 out(i->first) = in(i->first);
896 template <
typename VectorType>
901 if (!expect_constrained_source)
907 tmp_vector->reinit(src,
true);
909 pre_filter (*tmp_vector);
911 matrix->vmult (dst, *tmp_vector);
915 matrix->vmult (dst, src);
919 post_filter (src, dst);
924 template <
typename VectorType>
929 if (!expect_constrained_source)
935 tmp_vector->reinit(src,
true);
937 pre_filter (*tmp_vector);
939 matrix->Tvmult (dst, *tmp_vector);
943 matrix->Tvmult (dst, src);
947 post_filter (src, dst);
952 template <
typename VectorType>
957 if (!expect_constrained_source)
963 tmp_vector->reinit(src,
true);
965 pre_filter (*tmp_vector);
967 matrix->vmult_add (dst, *tmp_vector);
971 matrix->vmult_add (dst, src);
975 post_filter (src, dst);
980 template <
typename VectorType>
985 if (!expect_constrained_source)
991 tmp_vector->reinit(src,
true);
993 pre_filter (*tmp_vector);
995 matrix->Tvmult_add (dst, *tmp_vector);
999 matrix->Tvmult_add (dst, src);
1003 post_filter (src, dst);
1008 template <
typename VectorType>
1019 DEAL_II_NAMESPACE_CLOSE
void initialize(const MatrixType &m, const bool expect_constrained_source=false)
const FilteredMatrix< VectorType > * matrix
const Accessor * operator->() const
const_iterator begin() const
types::global_dof_index size_type
Subscriptor & operator=(const Subscriptor &)
void vmult_add(VectorType &dst, const VectorType &src) const
const_iterator & operator++()
void pre_filter(VectorType &v) const
std::vector< IndexValuePair >::const_iterator const_index_value_iterator
std::shared_ptr< PointerMatrixBase< VectorType > > matrix
void apply_constraints(VectorType &v, const bool matrix_is_symmetric) const
void post_filter(const VectorType &in, VectorType &out) const
std::vector< IndexValuePair > constraints
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
SymmetricTensor< rank_, dim, Number > operator*(const SymmetricTensor< rank_, dim, Number > &t, const Number &factor)
FilteredMatrix & operator=(const FilteredMatrix &fm)
Accessor(const FilteredMatrix< VectorType > *matrix, const size_type index)
void add_constraint(const size_type i, const double v)
PointerMatrixBase< VectorType > * new_pointer_matrix_base(MatrixType &matrix, const VectorType &, const char *name="PointerMatrixAux")
const Accessor & operator*() const
void vmult(VectorType &dst, const VectorType &src) const
unsigned int global_dof_index
void Tvmult(VectorType &dst, const VectorType &src) const
#define Assert(cond, exc)
bool expect_constrained_source
bool operator()(const IndexValuePair &i1, const IndexValuePair &i2) const
bool operator!=(const const_iterator &) const
bool operator==(const const_iterator &) const
void Tvmult_add(VectorType &dst, const VectorType &src) const
void add_constraints(const ConstraintList &new_constraints)
static ::ExceptionBase & ExcIteratorPastEnd()
std::pair< size_type, double > IndexValuePair
const_iterator(const FilteredMatrix< VectorType > *matrix, const size_type index)
#define AssertIsFinite(number)
std::size_t memory_consumption() const
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
const_iterator end() const