16#ifndef dealii_dof_accessor_h
17#define dealii_dof_accessor_h
27#include <boost/container/small_vector.hpp>
36template <
typename number>
38template <
typename number>
40template <
typename number>
43template <
typename Accessor>
51 namespace DoFCellAccessorImplementation
53 struct Implementation;
56 namespace DoFHandlerImplementation
58 struct Implementation;
61 struct Implementation;
67 namespace DoFHandlerImplementation
69 struct Implementation;
83 namespace DoFAccessorImplementation
101 template <
int structdim,
int dim,
int spacedim>
117 template <
int dim,
int spacedim>
127 struct Implementation;
208template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
229 using BaseClass = typename ::internal::DoFAccessorImplementation::
230 Inheritance<structdim, dimension, space_dimension>::BaseClass;
299 template <
int structdim2,
int dim2,
int spacedim2>
306 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
313 template <
bool level_dof_access2>
350 template <
bool level_dof_access2>
386 typename ::internal::DoFHandlerImplementation::
387 Iterators<dim, spacedim, level_dof_access>::line_iterator
388 line(
const unsigned int i)
const;
395 typename ::internal::DoFHandlerImplementation::
396 Iterators<dim, spacedim, level_dof_access>::quad_iterator
397 quad(
const unsigned int i)
const;
446 std::vector<types::global_dof_index> &dof_indices,
458 std::vector<types::global_dof_index> &dof_indices,
467 const std::vector<types::global_dof_index> &dof_indices,
493 const unsigned int vertex,
494 const unsigned int i,
505 const unsigned int vertex,
506 const unsigned int i,
586 std::set<types::fe_index>
619 "This accessor object has not been "
620 "associated with any DoFHandler object.");
674 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
682 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
714 const unsigned int i,
720 const unsigned int i,
726 const unsigned int vertex,
727 const unsigned int i,
735 template <
int,
int,
int,
bool>
743 friend struct ::internal::DoFHandlerImplementation::Policy::
745 friend struct ::internal::DoFHandlerImplementation::Implementation;
746 friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
747 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
748 friend struct ::internal::DoFAccessorImplementation::Implementation;
760template <
int spacedim,
bool level_dof_access>
819 const unsigned int vertex_index,
845 template <
int structdim2,
int dim2,
int spacedim2>
852 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
904 template <
bool level_dof_access2>
936 typename ::
internal::DoFHandlerImplementation::
937 Iterators<1, spacedim, level_dof_access>::line_iterator
938 line(const
unsigned int i) const;
946 typename ::
internal::DoFHandlerImplementation::
947 Iterators<1, spacedim, level_dof_access>::quad_iterator
948 quad(const
unsigned int i) const;
995 std::vector<
types::global_dof_index> &dof_indices,
996 const
types::fe_index fe_index =
numbers::invalid_fe_index) const;
1007 std::vector<
types::global_dof_index> &dof_indices,
1008 const
types::fe_index fe_index =
numbers::invalid_fe_index) const;
1028 types::global_dof_index
1030 const
unsigned int vertex,
1031 const
unsigned int i,
1032 const
types::fe_index fe_index =
numbers::invalid_fe_index) const;
1051 types::global_dof_index
1053 const
types::fe_index fe_index =
numbers::invalid_fe_index) const;
1115 "This accessor
object has not been "
1159 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
1162 const
DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1167 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
1170 const
DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1198 const
unsigned int i,
1199 const
types::global_dof_index index,
1200 const
types::fe_index fe_index =
numbers::invalid_fe_index) const;
1212 friend struct ::
internal::DoFHandlerImplementation::Policy::
1244template <
int structdim,
int dim,
int spacedim = dim>
1262 const int level = -1,
1263 const int index = -1,
1279 template <
typename OtherAccessor>
1299 const unsigned int i,
1320template <
int dimension_,
int space_dimension_,
bool level_dof_access>
1330 static const unsigned int dim = dimension_;
1335 static const unsigned int spacedim = space_dimension_;
1391 template <
int structdim2,
int dim2,
int spacedim2>
1398 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
1480 periodic_neighbor(
const unsigned int i)
const;
1488 neighbor_or_periodic_neighbor(
const unsigned int i)
const;
1501 boost::container::small_vector<
1530 neighbor_child_on_subface(
const unsigned int face_no,
1531 const unsigned int subface_no)
const;
1540 periodic_neighbor_child_on_subface(
const unsigned int face_no,
1541 const unsigned int subface_no)
const;
1571 template <
class InputVector,
typename number>
1592 template <
class InputVector,
typename ForwardIterator>
1595 ForwardIterator local_values_begin,
1596 ForwardIterator local_values_end)
const;
1618 template <
class InputVector,
typename ForwardIterator>
1622 const InputVector & values,
1623 ForwardIterator local_values_begin,
1624 ForwardIterator local_values_end)
const;
1650 template <
class OutputVector,
typename number>
1653 OutputVector & values)
const;
1686 template <
class InputVector,
typename number>
1688 get_interpolated_dof_values(
1689 const InputVector & values,
1754 template <
class OutputVector,
typename number>
1756 set_dof_values_by_interpolation(
1758 OutputVector & values,
1760 const bool perform_check =
false)
const;
1771 template <
class OutputVector,
typename number>
1773 distribute_local_to_global_by_interpolation(
1775 OutputVector & values,
1792 template <
typename number,
typename OutputVector>
1795 OutputVector & global_destination)
const;
1811 template <
typename ForwardIterator,
typename OutputVector>
1814 ForwardIterator local_source_end,
1815 OutputVector & global_destination)
const;
1830 template <
typename ForwardIterator,
typename OutputVector>
1834 ForwardIterator local_source_begin,
1835 ForwardIterator local_source_end,
1836 OutputVector & global_destination)
const;
1844 template <
typename number,
typename OutputMatrix>
1847 OutputMatrix &global_destination)
const;
1853 template <
typename number,
typename OutputMatrix,
typename OutputVector>
1857 OutputMatrix & global_matrix,
1858 OutputVector & global_vector)
const;
1886 std::vector<types::global_dof_index> &dof_indices)
const;
1920 get_dof_indices(std::vector<types::global_dof_index> &dof_indices)
const;
1927 get_mg_dof_indices(std::vector<types::global_dof_index> &dof_indices)
const;
2021 set_dof_indices(
const std::vector<types::global_dof_index> &dof_indices);
2027 set_mg_dof_indices(
const std::vector<types::global_dof_index> &dof_indices);
2113 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
2117template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
2121 return level_dof_access;
2126template <
int structdim,
int dim,
int spacedim>
2127template <
typename OtherAccessor>
2129 const OtherAccessor &)
2132 ExcMessage(
"You are attempting an illegal conversion between "
2133 "iterator/accessor types. The constructor you call "
2134 "only exists to make certain template constructs "
2135 "easier to write as dimension independent code but "
2136 "the conversion is not valid in the current context."));
2144#include <deal.II/dofs/dof_accessor.templates.h>
DoFAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
DoFAccessor< 0, 1, spacedim, level_dof_access > & operator=(const DoFAccessor< 0, 1, spacedim, level_dof_access > &da)=delete
DoFAccessor< 0, 1, spacedim, level_dof_access > & operator=(DoFAccessor< 0, 1, spacedim, level_dof_access > &&) noexcept=default
DoFAccessor(const Triangulation< 1, spacedim > *, const int=0, const int=0, const DoFHandler< 1, spacedim > *dof_handler=0)
DoFAccessor(DoFAccessor< 0, 1, spacedim, level_dof_access > &&)=default
DoFAccessor(const Triangulation< 1, spacedim > *tria, const typename TriaAccessor< 0, 1, spacedim >::VertexKind vertex_kind, const unsigned int vertex_index, const DoFHandler< 1, spacedim > *dof_handler)
DoFAccessor(const DoFAccessor< 0, 1, spacedim, level_dof_access > &)=default
DoFAccessor(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
typename::internal::DoFHandlerImplementation::Iterators< dim, spacedim, level_dof_access >::line_iterator line(const unsigned int i) const
friend struct ::internal::DoFHandlerImplementation::Policy::Implementation
static constexpr unsigned int space_dimension
const DoFHandler< dim, spacedim > & get_dof_handler() const
void set_mg_dof_index(const int level, const unsigned int i, const types::global_dof_index index) const
DoFAccessor< structdim, dim, spacedim, level_dof_access > & operator=(const DoFAccessor< structdim, dim, spacedim, level_dof_access > &da)=delete
types::fe_index nth_active_fe_index(const unsigned int n) const
DoFAccessor(DoFAccessor< structdim, dim, spacedim, level_dof_access > &&)
DoFAccessor(const DoFAccessor< structdim, dim, spacedim, level_dof_access2 > &)
TriaIterator< DoFAccessor< structdim, dim, spacedim, level_dof_access > > child(const unsigned int c) const
void set_mg_vertex_dof_index(const int level, const unsigned int vertex, const unsigned int i, const types::global_dof_index index, const types::fe_index fe_index=numbers::invalid_fe_index) const
static constexpr unsigned int dimension
void set_dof_handler(DoFHandler< dim, spacedim > *dh)
DoFAccessor(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
DoFAccessor< structdim, dim, spacedim, level_dof_access > & operator=(DoFAccessor< structdim, dim, spacedim, level_dof_access > &&)
std::set< types::fe_index > get_active_fe_indices() const
DoFAccessor(const DoFAccessor< structdim, dim, spacedim, level_dof_access > &)=default
void set_mg_dof_indices(const int level, const std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index)
void get_mg_dof_indices(const int level, std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index) const
types::global_dof_index mg_vertex_dof_index(const int level, const unsigned int vertex, const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
DoFAccessor(const Triangulation< dim, spacedim > *tria, const int level, const int index, const DoFHandler< dim, spacedim > *dof_handler)
bool operator!=(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &) const
void copy_from(const TriaAccessorBase< structdim, dim, spacedim > &da)
types::global_dof_index dof_index(const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
typename::internal::DoFHandlerImplementation::Iterators< dim, spacedim, level_dof_access >::quad_iterator quad(const unsigned int i) const
void copy_from(const DoFAccessor< structdim, dim, spacedim, level_dof_access2 > &a)
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index fe_index) const
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index) const
bool operator==(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &) const
types::global_dof_index mg_dof_index(const int level, const unsigned int i) const
void set_dof_index(const unsigned int i, const types::global_dof_index index, const types::fe_index fe_index=numbers::invalid_fe_index) const
unsigned int n_active_fe_indices() const
static bool is_level_cell()
typename ::internal::DoFAccessorImplementation::Inheritance< structdim, dimension, space_dimension >::BaseClass BaseClass
types::global_dof_index vertex_dof_index(const unsigned int vertex, const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
DoFHandler< dim, spacedim > * dof_handler
DoFAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
bool fe_index_is_active(const types::fe_index fe_index) const
boost::container::small_vector< TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > >, GeometryInfo< dimension_ >::max_children_per_cell > child_iterators() const
DoFCellAccessor(DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &&)=default
void get_active_or_mg_dof_indices(std::vector< types::global_dof_index > &dof_indices) const
void get_dof_values(const InputVector &values, Vector< number > &local_values) const
face_iterator face(const unsigned int i) const
types::fe_index future_fe_index() const
DoFCellAccessor< dimension_, space_dimension_, level_dof_access > & operator=(const DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &da)=delete
const FiniteElement< dimension_, space_dimension_ > & get_fe() const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > child(const unsigned int i) const
DoFCellAccessor(const Triangulation< dimension_, space_dimension_ > *tria, const int level, const int index, const AccessorData *local_data)
void get_dof_values(const InputVector &values, ForwardIterator local_values_begin, ForwardIterator local_values_end) const
void set_future_fe_index(const types::fe_index i) const
void distribute_local_to_global(const Vector< number > &local_source, OutputVector &global_destination) const
const FiniteElement< dimension_, space_dimension_ > & get_future_fe() const
void get_dof_values(const AffineConstraints< typename InputVector::value_type > &constraints, const InputVector &values, ForwardIterator local_values_begin, ForwardIterator local_values_end) const
void distribute_local_to_global(const FullMatrix< number > &local_source, OutputMatrix &global_destination) const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > neighbor(const unsigned int i) const
void set_active_fe_index(const types::fe_index i) const
DoFCellAccessor(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
DoFCellAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
void clear_future_fe_index() const
void distribute_local_to_global(const FullMatrix< number > &local_matrix, const Vector< number > &local_vector, OutputMatrix &global_matrix, OutputVector &global_vector) const
void distribute_local_to_global(const AffineConstraints< typename OutputVector::value_type > &constraints, ForwardIterator local_source_begin, ForwardIterator local_source_end, OutputVector &global_destination) const
void set_dof_values(const Vector< number > &local_values, OutputVector &values) const
boost::container::small_vector< face_iterator, GeometryInfo< dimension_ >::faces_per_cell > face_iterators() const
DoFCellAccessor(const DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &)=default
void distribute_local_to_global(ForwardIterator local_source_begin, ForwardIterator local_source_end, OutputVector &global_destination) const
~DoFCellAccessor()=default
bool future_fe_index_set() const
DoFCellAccessor< dimension_, space_dimension_, level_dof_access > & operator=(DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &&)=default
types::fe_index active_fe_index() const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > parent() const
DoFInvalidAccessor(const void *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
typename InvalidAccessor< structdim, dim, spacedim >::AccessorData AccessorData
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcVectorNotEmpty()
#define DeclException0(Exception0)
#define Assert(cond, exc)
static ::ExceptionBase & ExcNotActive()
static ::ExceptionBase & ExcMatrixDoesNotMatch()
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInvalidObject()
static ::ExceptionBase & ExcCantCompareIterators()
static ::ExceptionBase & ExcVectorDoesNotMatch()
static ::ExceptionBase & ExcMessage(std::string arg1)
const types::fe_index invalid_fe_index
unsigned short int fe_index
const ::Triangulation< dim, spacedim > & tria