16 #ifndef dealii__dof_accessor_h 17 #define dealii__dof_accessor_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/grid/tria_accessor.h> 22 #include <deal.II/dofs/dof_handler.h> 23 #include <deal.II/hp/dof_handler.h> 27 DEAL_II_NAMESPACE_OPEN
31 template <
typename number>
class Vector;
43 struct Implementation;
48 struct Implementation;
51 struct Implementation;
59 struct Implementation;
91 template <
int structdim,
int dim,
int spacedim>
98 typedef ::TriaAccessor<structdim,dim,spacedim>
BaseClass;
107 template <
int dim,
int spacedim>
184 template <
int structdim,
typename DoFHandlerType,
bool level_dof_access>
193 static const unsigned int dimension=DoFHandlerType::dimension;
199 static const unsigned int space_dimension=DoFHandlerType::space_dimension;
206 typename ::internal::DoFAccessor::Inheritance<structdim, dimension, space_dimension>::BaseClass
232 const DoFHandlerType *local_data);
246 template <
int structdim2,
int dim2,
int spacedim2>
253 template <
int dim2,
class DoFHandlerType2,
bool level_dof_access2>
259 template <
bool level_dof_access2>
268 const DoFHandlerType &
269 get_dof_handler ()
const;
274 template <
bool level_dof_access2>
287 static bool is_level_cell();
300 child (
const unsigned int c)
const;
307 typename ::internal::DoFHandler::Iterators<DoFHandlerType, level_dof_access>::line_iterator
308 line (
const unsigned int i)
const;
315 typename ::internal::DoFHandler::Iterators<DoFHandlerType, level_dof_access>::quad_iterator
316 quad (
const unsigned int i)
const;
363 void get_dof_indices (std::vector<types::global_dof_index> &dof_indices,
364 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
372 void get_mg_dof_indices (
const int level,
373 std::vector<types::global_dof_index> &dof_indices,
374 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
379 void set_mg_dof_indices (
const int level,
380 const std::vector<types::global_dof_index> &dof_indices,
381 const unsigned int fe_index = DoFHandlerType::default_fe_index);
401 (
const unsigned int vertex,
402 const unsigned int i,
403 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
412 const unsigned int vertex,
413 const unsigned int i,
414 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
444 (
const unsigned int i,
445 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
474 n_active_fe_indices ()
const;
484 nth_active_fe_index (
const unsigned int n)
const;
495 fe_index_is_active (
const unsigned int fe_index)
const;
503 get_fe (
const unsigned int fe_index)
const;
568 template <
int dim2,
class DoFHandlerType2,
bool level_dof_access2>
574 template <
int dim2,
class DoFHandlerType2,
bool level_dof_access2>
580 void set_dof_handler (DoFHandlerType *dh);
600 (
const unsigned int i,
602 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
623 void set_vertex_dof_index
624 (
const unsigned int vertex,
625 const unsigned int i,
627 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
629 void set_mg_vertex_dof_index
631 const unsigned int vertex,
632 const unsigned int i,
634 const unsigned int fe_index = DoFHandlerType::default_fe_index)
const;
641 template <
int,
class,
bool>
friend class DoFAccessor;
664 friend struct ::internal::DoFHandler::Policy::Implementation;
665 friend struct ::internal::DoFHandler::Implementation;
666 friend struct ::internal::hp::DoFHandler::Implementation;
667 friend struct ::internal::DoFCellAccessor::Implementation;
668 friend struct ::internal::DoFAccessor::Implementation;
682 template <
template <
int,
int>
class DoFHandlerType,
int spacedim,
bool level_dof_access>
691 static const unsigned int dimension=1;
697 static const unsigned int space_dimension=spacedim;
740 const unsigned int vertex_index,
741 const DoFHandlerType<1,spacedim> *dof_handler);
751 const DoFHandlerType<1,spacedim> *dof_handler = 0);
765 template <
int structdim2,
int dim2,
int spacedim2>
772 template <
int dim2,
class DoFHandlerType2,
bool level_dof_access2>
782 const DoFHandlerType<1,spacedim> &
783 get_dof_handler ()
const;
789 operator = (
const DoFAccessor<0,DoFHandlerType<1,spacedim>, level_dof_access> &da);
794 template <
bool level_dof_access2>
795 void copy_from (
const DoFAccessor<0, DoFHandlerType<1,spacedim>, level_dof_access2> &a);
816 child (
const unsigned int c)
const;
823 typename ::internal::DoFHandler::Iterators<DoFHandlerType<1,spacedim>, level_dof_access>::line_iterator
824 line (
const unsigned int i)
const;
831 typename ::internal::DoFHandler::Iterators<DoFHandlerType<1,spacedim>, level_dof_access>::quad_iterator
832 quad (
const unsigned int i)
const;
878 void get_dof_indices (std::vector<types::global_dof_index> &dof_indices,
879 const unsigned int fe_index = AccessorData::default_fe_index)
const;
899 const unsigned int i,
900 const unsigned int fe_index = AccessorData::default_fe_index)
const;
930 const unsigned int fe_index = AccessorData::default_fe_index)
const;
954 n_active_fe_indices ()
const;
964 nth_active_fe_index (
const unsigned int n)
const;
975 fe_index_is_active (
const unsigned int fe_index)
const;
983 get_fe (
const unsigned int fe_index)
const;
1038 template <
int dim2,
class DoFHandlerType2,
bool level_dof_access2>
1044 template <
int dim2,
class DoFHandlerType2,
bool level_dof_access2>
1050 void set_dof_handler (DoFHandlerType<1,spacedim> *dh);
1069 void set_dof_index (
const unsigned int i,
1071 const unsigned int fe_index = AccessorData::default_fe_index)
const;
1090 void set_vertex_dof_index (
const unsigned int vertex,
1091 const unsigned int i,
1093 const unsigned int fe_index = AccessorData::default_fe_index)
const;
1109 friend struct ::internal::DoFHandler::Policy::Implementation;
1110 friend struct ::internal::DoFHandler::Implementation;
1111 friend struct ::internal::hp::DoFHandler::Implementation;
1112 friend struct ::internal::DoFCellAccessor::Implementation;
1131 template <
typename DoFHandlerType,
bool level_dof_access>
1138 static const unsigned int dim = DoFHandlerType::dimension;
1143 static const unsigned int spacedim = DoFHandlerType::space_dimension;
1183 const AccessorData *local_data);
1197 template <
int structdim2,
int dim2,
int spacedim2>
1204 template <
int dim2,
class DoFHandlerType2,
bool level_dof_access2>
1237 neighbor (
const unsigned int)
const;
1245 child (
const unsigned int)
const;
1254 face (
const unsigned int i)
const;
1263 neighbor_child_on_subface (
const unsigned int face_no,
1264 const unsigned int subface_no)
const;
1291 template <
class InputVector,
typename number>
1292 void get_dof_values (
const InputVector &values,
1309 template <
class InputVector,
typename ForwardIterator>
1310 void get_dof_values (
const InputVector &values,
1311 ForwardIterator local_values_begin,
1312 ForwardIterator local_values_end)
const;
1330 template <
class InputVector,
typename ForwardIterator>
1332 const InputVector &values,
1333 ForwardIterator local_values_begin,
1334 ForwardIterator local_values_end)
const;
1357 template <
class OutputVector,
typename number>
1359 OutputVector &values)
const;
1392 template <
class InputVector,
typename number>
1393 void get_interpolated_dof_values (
const InputVector &values,
1395 const unsigned int fe_index
1396 = DoFHandlerType::default_fe_index)
const;
1450 template <
class OutputVector,
typename number>
1451 void set_dof_values_by_interpolation (
const Vector<number> &local_values,
1452 OutputVector &values,
1453 const unsigned int fe_index
1454 = DoFHandlerType::default_fe_index)
const;
1464 template <
typename number,
typename OutputVector>
1467 OutputVector &global_destination)
const;
1477 template <
typename ForwardIterator,
typename OutputVector>
1479 distribute_local_to_global (ForwardIterator local_source_begin,
1480 ForwardIterator local_source_end,
1481 OutputVector &global_destination)
const;
1493 template <
typename ForwardIterator,
typename OutputVector>
1496 ForwardIterator local_source_begin,
1497 ForwardIterator local_source_end,
1498 OutputVector &global_destination)
const;
1506 template <
typename number,
typename OutputMatrix>
1509 OutputMatrix &global_destination)
const;
1515 template <
typename number,
typename OutputMatrix,
typename OutputVector>
1519 OutputMatrix &global_matrix,
1520 OutputVector &global_vector)
const;
1546 void get_active_or_mg_dof_indices (std::vector<types::global_dof_index> &dof_indices)
const;
1583 void get_dof_indices (std::vector<types::global_dof_index> &dof_indices)
const;
1592 void get_mg_dof_indices (std::vector<types::global_dof_index> &dof_indices)
const;
1633 unsigned int active_fe_index ()
const;
1648 void set_active_fe_index (
const unsigned int i);
1657 void set_dof_indices (
const std::vector<types::global_dof_index> &dof_indices);
1662 void set_mg_dof_indices (
const std::vector<types::global_dof_index> &dof_indices);
1667 void update_cell_dof_indices_cache ()
const;
1688 friend struct ::internal::DoFCellAccessor::Implementation;
1692 template <
int sd,
typename DoFHandlerType,
bool level_dof_access>
1697 return level_dof_access;
1702 DEAL_II_NAMESPACE_CLOSE
1705 #include "dof_accessor.templates.h" DoFHandlerType * dof_handler
::TriaAccessor< structdim, dim, spacedim > BaseClass
static bool is_level_cell()
DoFHandlerType< 1, spacedim > * dof_handler
TriaIterator< DoFAccessor< DoFHandlerType::dimension-1, DoFHandlerType, level_dof_access > > face_iterator
TriaAccessor< 0, 1, spacedim > BaseClass
DoFHandlerType AccessorData
typename ::internal::DoFAccessor::Inheritance< structdim, dimension, space_dimension >::BaseClass BaseClass
unsigned int global_dof_index
#define DeclException0(Exception0)
::CellAccessor< dim, spacedim > BaseClass
DoFHandlerType< 1, spacedim > AccessorData
DoFHandlerType AccessorData
DoFAccessor< structdim, DoFHandlerType, level_dof_access > & operator=(const DoFAccessor< structdim, DoFHandlerType, level_dof_access > &da)
DoFAccessor< DoFHandlerType::dimension, DoFHandlerType, level_dof_access > BaseClass