16 #ifndef dealii_dof_accessor_h 17 #define dealii_dof_accessor_h 34 template <
typename number>
36 template <
typename number>
38 template <
typename number>
41 template <
typename Accessor>
49 namespace DoFCellAccessorImplementation
51 struct Implementation;
54 namespace DoFHandlerImplementation
56 struct Implementation;
59 struct Implementation;
65 namespace DoFHandlerImplementation
67 struct Implementation;
81 namespace DoFAccessorImplementation
100 template <
int structdim,
int dim,
int spacedim>
116 template <
int dim,
int spacedim>
210 template <
int structdim,
typename DoFHandlerType,
bool level_dof_access>
214 DoFHandlerType::dimension,
222 static const unsigned int dimension = DoFHandlerType::dimension;
228 static const unsigned int space_dimension = DoFHandlerType::space_dimension;
234 using BaseClass = typename ::internal::DoFAccessorImplementation::
235 Inheritance<structdim, dimension, space_dimension>::BaseClass;
271 DoFHandlerType::space_dimension> *tria,
274 const DoFHandlerType *dof_handler);
288 template <
int structdim2,
int dim2,
int spacedim2>
295 template <
int dim2,
class DoFHandlerType2,
bool level_dof_access2>
301 template <
bool level_dof_access2>
325 const DoFHandlerType &
326 get_dof_handler()
const;
331 template <
bool level_dof_access2>
341 DoFHandlerType::dimension,
342 DoFHandlerType::space_dimension> &da);
362 child(
const unsigned int c)
const;
369 typename ::internal::DoFHandlerImplementation::
370 Iterators<DoFHandlerType, level_dof_access>::line_iterator
371 line(
const unsigned int i)
const;
378 typename ::internal::DoFHandlerImplementation::
379 Iterators<DoFHandlerType, level_dof_access>::quad_iterator
380 quad(
const unsigned int i)
const;
429 std::vector<types::global_dof_index> &dof_indices,
430 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
441 std::vector<types::global_dof_index> &dof_indices,
442 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
450 const std::vector<types::global_dof_index> &dof_indices,
451 const unsigned int fe_index = DoFHandlerType::default_fe_index);
472 const unsigned int vertex,
473 const unsigned int i,
474 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
484 const unsigned int vertex,
485 const unsigned int i,
486 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
517 const unsigned int i,
518 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
524 mg_dof_index(
const int level,
const unsigned int i)
const;
548 n_active_fe_indices()
const;
558 nth_active_fe_index(
const unsigned int n)
const;
566 std::set<unsigned int>
578 fe_index_is_active(
const unsigned int fe_index)
const;
586 DoFHandlerType::space_dimension> &
587 get_fe(
const unsigned int fe_index)
const;
599 "This accessor object has not been " 600 "associated with any DoFHandler object.");
654 template <
int dim2,
class DoFHandlerType2,
bool level_dof_access2>
662 template <
int dim2,
class DoFHandlerType2,
bool level_dof_access2>
672 set_dof_handler(DoFHandlerType *dh);
693 const unsigned int i,
695 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
698 set_mg_dof_index(
const int level,
699 const unsigned int i,
720 set_vertex_dof_index(
721 const unsigned int vertex,
722 const unsigned int i,
724 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
727 set_mg_vertex_dof_index(
729 const unsigned int vertex,
730 const unsigned int i,
732 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
738 template <
int,
class,
bool>
744 template <
int dim,
int spacedim>
746 template <
int dim,
int spacedim>
749 friend struct ::internal::DoFHandlerImplementation::Policy::
751 friend struct ::internal::DoFHandlerImplementation::Implementation;
752 friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
753 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
754 friend struct ::internal::DoFAccessorImplementation::Implementation;
768 template <
template <
int,
int>
class DoFHandlerType,
770 bool level_dof_access>
771 class DoFAccessor<0, DoFHandlerType<1, spacedim>, level_dof_access>
779 static const unsigned int dimension = 1;
785 static const unsigned int space_dimension = spacedim;
829 const unsigned int vertex_index,
830 const DoFHandlerType<1, spacedim> * dof_handler);
840 const DoFHandlerType<1, spacedim> *dof_handler = 0);
854 template <
int structdim2,
int dim2,
int spacedim2>
861 template <
int dim2,
class DoFHandlerType2,
bool level_dof_access2>
874 operator=(
const DoFAccessor<0, DoFHandlerType<1, spacedim>, level_dof_access>
884 const DoFHandlerType<1, spacedim> &
885 get_dof_handler()
const;
890 template <
bool level_dof_access2>
893 const DoFAccessor<0, DoFHandlerType<1, spacedim>, level_dof_access2> &a);
915 child(
const unsigned int c)
const;
923 typename ::internal::DoFHandlerImplementation::
924 Iterators<DoFHandlerType<1, spacedim>, level_dof_access>::line_iterator
925 line(
const unsigned int i)
const;
933 typename ::internal::DoFHandlerImplementation::
934 Iterators<DoFHandlerType<1, spacedim>, level_dof_access>::quad_iterator
935 quad(
const unsigned int i)
const;
981 std::vector<types::global_dof_index> &dof_indices,
982 const unsigned int fe_index = AccessorData::default_fe_index)
const;
993 std::vector<types::global_dof_index> &dof_indices,
994 const unsigned int fe_index = AccessorData::default_fe_index)
const;
1015 const unsigned int vertex,
1016 const unsigned int i,
1017 const unsigned int fe_index = AccessorData::default_fe_index)
const;
1034 dof_index(
const unsigned int i,
1035 const unsigned int fe_index = AccessorData::default_fe_index)
const;
1056 n_active_fe_indices()
const;
1066 nth_active_fe_index(
const unsigned int n)
const;
1077 fe_index_is_active(
const unsigned int fe_index)
const;
1085 DoFHandlerType<1, spacedim>::space_dimension> &
1086 get_fe(
const unsigned int fe_index)
const;
1098 "This accessor object has not been " 1099 "associated with any DoFHandler object.");
1142 template <
int dim2,
class DoFHandlerType2,
bool level_dof_access2>
1150 template <
int dim2,
class DoFHandlerType2,
bool level_dof_access2>
1158 void set_dof_handler(DoFHandlerType<1, spacedim> *dh);
1179 const unsigned int i,
1181 const unsigned int fe_index = AccessorData::default_fe_index)
const;
1201 set_vertex_dof_index(
1202 const unsigned int vertex,
1203 const unsigned int i,
1205 const unsigned int fe_index = AccessorData::default_fe_index)
const;
1220 friend struct ::internal::DoFHandlerImplementation::Policy::
1222 friend struct ::internal::DoFHandlerImplementation::Implementation;
1223 friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
1224 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
1253 template <
int structdim,
int dim,
int spacedim = dim>
1271 const int level = -1,
1272 const int index = -1,
1288 template <
typename OtherAccessor>
1297 set_dof_index(
const unsigned int i,
1299 const unsigned int fe_index =
1320 template <
typename DoFHandlerType,
bool level_dof_access>
1329 static const unsigned int dim = DoFHandlerType::dimension;
1334 static const unsigned int spacedim = DoFHandlerType::space_dimension;
1373 DoFHandlerType::space_dimension> *tria,
1390 template <
int structdim2,
int dim2,
int spacedim2>
1397 template <
int dim2,
class DoFHandlerType2,
bool level_dof_access2>
1443 neighbor(
const unsigned int i)
const;
1451 periodic_neighbor(
const unsigned int i)
const;
1459 neighbor_or_periodic_neighbor(
const unsigned int i)
const;
1467 child(
const unsigned int i)
const;
1476 face(
const unsigned int i)
const;
1483 face_iterators()
const;
1492 neighbor_child_on_subface(
const unsigned int face_no,
1493 const unsigned int subface_no)
const;
1502 periodic_neighbor_child_on_subface(
const unsigned int face_no,
1503 const unsigned int subface_no)
const;
1530 template <
class InputVector,
typename number>
1532 get_dof_values(
const InputVector &values,
Vector<number> &local_values)
const;
1548 template <
class InputVector,
typename ForwardIterator>
1550 get_dof_values(
const InputVector &values,
1551 ForwardIterator local_values_begin,
1552 ForwardIterator local_values_end)
const;
1571 template <
class InputVector,
typename ForwardIterator>
1575 const InputVector & values,
1576 ForwardIterator local_values_begin,
1577 ForwardIterator local_values_end)
const;
1600 template <
class OutputVector,
typename number>
1603 OutputVector & values)
const;
1636 template <
class InputVector,
typename number>
1638 get_interpolated_dof_values(
1639 const InputVector &values,
1641 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
1695 template <
class OutputVector,
typename number>
1697 set_dof_values_by_interpolation(
1699 OutputVector & values,
1700 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
1710 template <
typename number,
typename OutputVector>
1713 OutputVector & global_destination)
const;
1723 template <
typename ForwardIterator,
typename OutputVector>
1725 distribute_local_to_global(ForwardIterator local_source_begin,
1726 ForwardIterator local_source_end,
1727 OutputVector & global_destination)
const;
1739 template <
typename ForwardIterator,
typename OutputVector>
1741 distribute_local_to_global(
1743 ForwardIterator local_source_begin,
1744 ForwardIterator local_source_end,
1745 OutputVector & global_destination)
const;
1753 template <
typename number,
typename OutputMatrix>
1756 OutputMatrix &global_destination)
const;
1762 template <
typename number,
typename OutputMatrix,
typename OutputVector>
1766 OutputMatrix & global_matrix,
1767 OutputVector & global_vector)
const;
1794 get_active_or_mg_dof_indices(
1795 std::vector<types::global_dof_index> &dof_indices)
const;
1833 get_dof_indices(std::vector<types::global_dof_index> &dof_indices)
const;
1840 get_mg_dof_indices(std::vector<types::global_dof_index> &dof_indices)
const;
1867 DoFHandlerType::space_dimension> &
1896 active_fe_index()
const;
1925 set_active_fe_index(
const unsigned int i)
const;
1935 set_dof_indices(
const std::vector<types::global_dof_index> &dof_indices);
1941 set_mg_dof_indices(
const std::vector<types::global_dof_index> &dof_indices);
1947 update_cell_dof_indices_cache()
const;
1975 DoFHandlerType::space_dimension> &
1976 get_future_fe()
const;
1997 future_fe_index()
const;
2008 set_future_fe_index(
const unsigned int i)
const;
2017 future_fe_index_set()
const;
2028 clear_future_fe_index()
const;
2036 template <
int dim,
int spacedim>
2038 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
2042 template <
int sd,
typename DoFHandlerType,
bool level_dof_access>
2046 return level_dof_access;
2051 template <
int structdim,
int dim,
int spacedim>
2052 template <
typename OtherAccessor>
2054 const OtherAccessor &)
2057 ExcMessage(
"You are attempting an illegal conversion between " 2058 "iterator/accessor types. The constructor you call " 2059 "only exists to make certain template constructs " 2060 "easier to write as dimension independent code but " 2061 "the conversion is not valid in the current context."));
2069 #include "dof_accessor.templates.h"
DoFHandlerType< 1, spacedim > AccessorData
DoFHandlerType * dof_handler
static bool is_level_cell()
bool operator!=(const AlignedVector< T > &lhs, const AlignedVector< T > &rhs)
DoFHandlerType< 1, spacedim > * dof_handler
typename ::internal::DoFAccessorImplementation::Inheritance< structdim, dimension, space_dimension >::BaseClass BaseClass
typename InvalidAccessor< structdim, dim, spacedim >::AccessorData AccessorData
bool operator==(const AlignedVector< T > &lhs, const AlignedVector< T > &rhs)
typename TriaAccessorBase< structdim, dim, spacedim >::AccessorData AccessorData
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
#define DeclExceptionMsg(Exception, defaulttext)
#define DeclException0(Exception0)
#define DEAL_II_NAMESPACE_CLOSE
DoFInvalidAccessor(const Triangulation< dim, spacedim > *parent=0, const int level=-1, const int index=-1, const AccessorData *local_data=0)
#define DEAL_II_NAMESPACE_OPEN
DoFHandlerType AccessorData
static ::ExceptionBase & ExcCantCompareIterators()