16 #ifndef dealii_dof_accessor_h
17 #define dealii_dof_accessor_h
28 #include <boost/container/small_vector.hpp>
37 template <
typename number>
39 template <
typename number>
41 template <
typename number>
44 template <
typename Accessor>
52 namespace DoFCellAccessorImplementation
54 struct Implementation;
57 namespace DoFHandlerImplementation
59 struct Implementation;
62 struct Implementation;
68 namespace DoFHandlerImplementation
70 struct Implementation;
84 namespace DoFAccessorImplementation
102 template <
int structdim,
int dim,
int spacedim>
118 template <
int dim,
int spacedim>
128 struct Implementation;
209 template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
210 class DoFAccessor :
public ::internal::DoFAccessorImplementation::
230 using BaseClass = typename ::internal::DoFAccessorImplementation::
231 Inheritance<structdim, dimension, space_dimension>::BaseClass;
300 template <
int structdim2,
int dim2,
int spacedim2>
307 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
314 template <
bool level_dof_access2>
351 template <
bool level_dof_access2>
387 typename ::internal::DoFHandlerImplementation::
388 Iterators<dim, spacedim, level_dof_access>::line_iterator
389 line(
const unsigned int i)
const;
396 typename ::internal::DoFHandlerImplementation::
397 Iterators<dim, spacedim, level_dof_access>::quad_iterator
398 quad(
const unsigned int i)
const;
447 const unsigned int fe_index =
458 std::vector<types::global_dof_index> &dof_indices,
459 const unsigned int fe_index =
468 const std::vector<types::global_dof_index> &dof_indices,
494 const unsigned int i,
495 const unsigned int fe_index =
505 const unsigned int vertex,
506 const unsigned int i,
507 const unsigned int fe_index =
539 const unsigned int fe_index =
588 std::set<unsigned int>
609 get_fe(
const unsigned int fe_index)
const;
621 "This accessor object has not been "
622 "associated with any DoFHandler object.");
676 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
684 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
717 const unsigned int fe_index =
722 const unsigned int i,
727 const unsigned int vertex,
728 const unsigned int i,
730 const unsigned int fe_index =
737 template <
int,
int,
int,
bool>
746 friend struct ::internal::DoFHandlerImplementation::Policy::
748 friend struct ::internal::DoFHandlerImplementation::Implementation;
749 friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
750 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
751 friend struct ::internal::DoFAccessorImplementation::Implementation;
763 template <
int spacedim,
bool level_dof_access>
822 const unsigned int vertex_index,
848 template <
int structdim2,
int dim2,
int spacedim2>
855 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
907 template <
bool level_dof_access2>
939 typename ::
internal::DoFHandlerImplementation::
940 Iterators<1, spacedim, level_dof_access>::line_iterator
941 line(const
unsigned int i) const;
949 typename ::
internal::DoFHandlerImplementation::
950 Iterators<1, spacedim, level_dof_access>::quad_iterator
951 quad(const
unsigned int i) const;
999 const
unsigned int fe_index =
AccessorData::invalid_fe_index) const;
1011 const
unsigned int fe_index =
AccessorData::invalid_fe_index) const;
1033 const
unsigned int vertex,
1034 const
unsigned int i,
1035 const
unsigned int fe_index =
AccessorData::invalid_fe_index) const;
1056 const
unsigned int fe_index =
AccessorData::invalid_fe_index) const;
1118 "This accessor
object has not been "
1162 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
1165 const
DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1170 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
1173 const
DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1201 const
unsigned int i,
1203 const
unsigned int fe_index =
AccessorData::invalid_fe_index) const;
1216 friend struct ::
internal::DoFHandlerImplementation::Policy::
1218 friend struct ::
internal::DoFHandlerImplementation::Implementation;
1219 friend struct ::
internal::
hp::DoFHandlerImplementation::Implementation;
1220 friend struct ::
internal::DoFCellAccessorImplementation::Implementation;
1248 template <
int structdim,
int dim,
int spacedim = dim>
1266 const int level = -1,
1267 const int index = -1,
1283 template <
typename OtherAccessor>
1293 const unsigned int fe_index =
1304 const unsigned int fe_index =
1324 template <
int dimension_,
int space_dimension_,
bool level_dof_access>
1334 static const unsigned int dim = dimension_;
1339 static const unsigned int spacedim = space_dimension_;
1395 template <
int structdim2,
int dim2,
int spacedim2>
1402 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
1484 periodic_neighbor(
const unsigned int i)
const;
1492 neighbor_or_periodic_neighbor(
const unsigned int i)
const;
1505 boost::container::small_vector<
1534 neighbor_child_on_subface(
const unsigned int face_no,
1535 const unsigned int subface_no)
const;
1544 periodic_neighbor_child_on_subface(
const unsigned int face_no,
1545 const unsigned int subface_no)
const;
1575 template <
class InputVector,
typename number>
1596 template <
class InputVector,
typename ForwardIterator>
1599 ForwardIterator local_values_begin,
1600 ForwardIterator local_values_end)
const;
1622 template <
class InputVector,
typename ForwardIterator>
1626 const InputVector &
values,
1627 ForwardIterator local_values_begin,
1628 ForwardIterator local_values_end)
const;
1654 template <
class OutputVector,
typename number>
1657 OutputVector &
values)
const;
1690 template <
class InputVector,
typename number>
1692 get_interpolated_dof_values(
1693 const InputVector &
values,
1695 const unsigned int fe_index =
1759 template <
class OutputVector,
typename number>
1761 set_dof_values_by_interpolation(
1764 const unsigned int fe_index =
1766 const bool perform_check =
false)
const;
1777 template <
class OutputVector,
typename number>
1779 distribute_local_to_global_by_interpolation(
1782 const unsigned int fe_index =
1799 template <
typename number,
typename OutputVector>
1802 OutputVector & global_destination)
const;
1818 template <
typename ForwardIterator,
typename OutputVector>
1821 ForwardIterator local_source_end,
1822 OutputVector & global_destination)
const;
1837 template <
typename ForwardIterator,
typename OutputVector>
1841 ForwardIterator local_source_begin,
1842 ForwardIterator local_source_end,
1843 OutputVector & global_destination)
const;
1851 template <
typename number,
typename OutputMatrix>
1854 OutputMatrix &global_destination)
const;
1860 template <
typename number,
typename OutputMatrix,
typename OutputVector>
1864 OutputMatrix & global_matrix,
1865 OutputVector & global_vector)
const;
1893 std::vector<types::global_dof_index> &dof_indices)
const;
2028 set_dof_indices(
const std::vector<types::global_dof_index> &dof_indices);
2040 update_cell_dof_indices_cache()
const;
2130 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
2134 template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
2138 return level_dof_access;
2143 template <
int structdim,
int dim,
int spacedim>
2144 template <
typename OtherAccessor>
2146 const OtherAccessor &)
2149 ExcMessage(
"You are attempting an illegal conversion between "
2150 "iterator/accessor types. The constructor you call "
2151 "only exists to make certain template constructs "
2152 "easier to write as dimension independent code but "
2153 "the conversion is not valid in the current context."));
2161 #include "dof_accessor.templates.h"
DoFAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
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< 0, 1, spacedim, level_dof_access > & operator=(DoFAccessor< 0, 1, spacedim, level_dof_access > &&) noexcept=default
DoFAccessor(const DoFAccessor< 0, 1, spacedim, level_dof_access > &)=default
DoFAccessor(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
DoFAccessor< 0, 1, spacedim, level_dof_access > & operator=(const DoFAccessor< 0, 1, spacedim, level_dof_access > &da)=delete
static constexpr unsigned int space_dimension
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=(DoFAccessor< structdim, dim, spacedim, level_dof_access > &&)
DoFAccessor(DoFAccessor< structdim, dim, spacedim, level_dof_access > &&)
DoFAccessor(const DoFAccessor< structdim, dim, spacedim, level_dof_access2 > &)
void set_mg_dof_indices(const int level, const std::vector< types::global_dof_index > &dof_indices, const unsigned int fe_index=DoFHandler< dim, spacedim >::invalid_fe_index)
const DoFHandler< dim, spacedim > & get_dof_handler() const
types::global_dof_index vertex_dof_index(const unsigned int vertex, const unsigned int i, const unsigned int fe_index=DoFHandler< dim, spacedim >::invalid_fe_index) const
std::set< unsigned int > get_active_fe_indices() const
static constexpr unsigned int dimension
void set_dof_handler(DoFHandler< dim, spacedim > *dh)
DoFAccessor(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
types::global_dof_index dof_index(const unsigned int i, const unsigned int fe_index=DoFHandler< dim, spacedim >::invalid_fe_index) const
typename ::internal::DoFHandlerImplementation::Iterators< dim, spacedim, level_dof_access >::quad_iterator quad(const unsigned int i) const
DoFAccessor(const DoFAccessor< structdim, dim, spacedim, level_dof_access > &)=default
TriaIterator< DoFAccessor< structdim, dim, spacedim, level_dof_access > > child(const unsigned int c) 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 get_dof_indices(std::vector< types::global_dof_index > &dof_indices, const unsigned int fe_index=DoFHandler< dim, spacedim >::invalid_fe_index) const
void copy_from(const TriaAccessorBase< structdim, dim, spacedim > &da)
bool fe_index_is_active(const unsigned int fe_index) const
void copy_from(const DoFAccessor< structdim, dim, spacedim, level_dof_access2 > &a)
typename ::internal::DoFHandlerImplementation::Iterators< dim, spacedim, level_dof_access >::line_iterator line(const unsigned int i) const
void set_dof_index(const unsigned int i, const types::global_dof_index index, const unsigned int fe_index=DoFHandler< dim, spacedim >::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
types::global_dof_index mg_vertex_dof_index(const int level, const unsigned int vertex, const unsigned int i, const unsigned int fe_index=DoFHandler< dim, spacedim >::invalid_fe_index) const
const FiniteElement< dim, spacedim > & get_fe(const unsigned int fe_index) const
unsigned int n_active_fe_indices() const
DoFAccessor< structdim, dim, spacedim, level_dof_access > & operator=(const DoFAccessor< structdim, dim, spacedim, level_dof_access > &da)=delete
unsigned int nth_active_fe_index(const unsigned int n) const
static bool is_level_cell()
void get_mg_dof_indices(const int level, std::vector< types::global_dof_index > &dof_indices, const unsigned int fe_index=DoFHandler< dim, spacedim >::invalid_fe_index) const
typename ::internal::DoFAccessorImplementation::Inheritance< structdim, dimension, space_dimension >::BaseClass BaseClass
DoFHandler< dim, spacedim > * dof_handler
DoFAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
void set_mg_vertex_dof_index(const int level, const unsigned int vertex, const unsigned int i, const types::global_dof_index index, const unsigned int fe_index=DoFHandler< dim, spacedim >::invalid_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< dimension_, space_dimension_, level_dof_access > & operator=(DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &&)=default
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
DoFCellAccessor(const Triangulation< dimension_, space_dimension_ > *tria, const int level, const int index, const AccessorData *local_data)
void set_active_fe_index(const unsigned int i) const
void get_dof_values(const InputVector &values, ForwardIterator local_values_begin, ForwardIterator local_values_end) const
void distribute_local_to_global(const Vector< number > &local_source, OutputVector &global_destination) const
void get_dof_values(const AffineConstraints< typename InputVector::value_type > &constraints, const InputVector &values, ForwardIterator local_values_begin, ForwardIterator local_values_end) const
unsigned int active_fe_index() const
void distribute_local_to_global(const FullMatrix< number > &local_source, OutputMatrix &global_destination) const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > child(const unsigned int i) const
boost::container::small_vector< face_iterator, GeometryInfo< dimension_ >::faces_per_cell > face_iterators() const
const FiniteElement< dimension_, space_dimension_ > & get_future_fe() const
DoFCellAccessor< dimension_, space_dimension_, level_dof_access > & operator=(const DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &da)=delete
DoFCellAccessor(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > neighbor(const unsigned int i) const
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_future_fe_index(const unsigned int i) const
void set_dof_values(const Vector< number > &local_values, OutputVector &values) const
void set_mg_dof_indices(const std::vector< types::global_dof_index > &dof_indices)
DoFCellAccessor(const DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &)=default
const FiniteElement< dimension_, space_dimension_ > & get_fe() const
void get_mg_dof_indices(std::vector< types::global_dof_index > &dof_indices) const
void distribute_local_to_global(ForwardIterator local_source_begin, ForwardIterator local_source_end, OutputVector &global_destination) const
~DoFCellAccessor()=default
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > parent() const
unsigned int future_fe_index() const
bool future_fe_index_set() const
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices) const
DoFInvalidAccessor(const Triangulation< dim, spacedim > *parent=0, const int level=-1, const int index=-1, const AccessorData *local_data=0)
typename InvalidAccessor< structdim, dim, spacedim >::AccessorData AccessorData
typename TriaAccessorBase< structdim, dim, spacedim >::AccessorData AccessorData
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
#define DeclException0(Exception0)
static ::ExceptionBase & ExcVectorDoesNotMatch()
#define Assert(cond, exc)
static ::ExceptionBase & ExcNotActive()
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInvalidObject()
static ::ExceptionBase & ExcMatrixDoesNotMatch()
static ::ExceptionBase & ExcCantCompareIterators()
static ::ExceptionBase & ExcVectorNotEmpty()
static ::ExceptionBase & ExcMessage(std::string arg1)
unsigned int global_dof_index
const ::Triangulation< dim, spacedim > & tria