![]() |
Reference documentation for deal.II version 8.1.0
|
#include <dof_accessor.h>
Public Types | |
typedef TriaAccessor < 0, 1, spacedim > | BaseClass |
typedef DH< 1, spacedim > | AccessorData |
![]() | |
enum | VertexKind { left_vertex, interior_vertex, right_vertex } |
typedef void | AccessorData |
Public Member Functions | |
const DH< 1, spacedim > & | get_dof_handler () const |
DoFAccessor< 0, DH < 1, spacedim > , level_dof_access > & | operator= (const DoFAccessor< 0, DH< 1, spacedim >, level_dof_access > &da) |
template<bool level_dof_access2> | |
void | copy_from (const DoFAccessor< 0, DH< 1, spacedim >, level_dof_access2 > &a) |
void | copy_from (const TriaAccessorBase< 0, 1, spacedim > &da) |
TriaIterator< DoFAccessor < 0, DH< 1, spacedim > , level_dof_access > > | parent () const |
DeclException0 (ExcInvalidObject) | |
DeclException0 (ExcVectorNotEmpty) | |
DeclException0 (ExcVectorDoesNotMatch) | |
DeclException0 (ExcMatrixDoesNotMatch) | |
DeclException0 (ExcNotActive) | |
DeclException0 (ExcCantCompareIterators) | |
Constructors | |
DoFAccessor () | |
DoFAccessor (const Triangulation< 1, spacedim > *tria, const typename TriaAccessor< 0, 1, spacedim >::VertexKind vertex_kind, const unsigned int vertex_index, const DH< 1, spacedim > *dof_handler) | |
DoFAccessor (const Triangulation< 1, spacedim > *, const int=0, const int=0, const DH< 1, spacedim > *dof_handler=0) | |
template<int structdim2, int dim2, int spacedim2> | |
DoFAccessor (const InvalidAccessor< structdim2, dim2, spacedim2 > &) | |
template<int dim2, class DH2 , bool level_dof_access2> | |
DoFAccessor (const DoFAccessor< dim2, DH2, level_dof_access2 > &) | |
Accessing sub-objects | |
TriaIterator< DoFAccessor < 0, DH< 1, spacedim > , level_dof_access > > | child (const unsigned int c) const |
typename::internal::DoFHandler::Iterators < DH< 1, spacedim > , level_dof_access > ::line_iterator | line (const unsigned int i) const |
typename::internal::DoFHandler::Iterators < DH< 1, spacedim > , level_dof_access > ::quad_iterator | quad (const unsigned int i) const |
Accessing the DoF indices of this object | |
void | get_dof_indices (std::vector< types::global_dof_index > &dof_indices, const unsigned int fe_index=AccessorData::default_fe_index) const |
types::global_dof_index | vertex_dof_index (const unsigned int vertex, const unsigned int i, const unsigned int fe_index=AccessorData::default_fe_index) const |
types::global_dof_index | dof_index (const unsigned int i, const unsigned int fe_index=AccessorData::default_fe_index) const |
Accessing the finite element associated with this object | |
unsigned int | n_active_fe_indices () const |
unsigned int | nth_active_fe_index (const unsigned int n) const |
bool | fe_index_is_active (const unsigned int fe_index) const |
const FiniteElement< DH < 1, spacedim >::dimension, DH < 1, spacedim > ::space_dimension > & | get_fe (const unsigned int fe_index) const |
![]() | |
TriaAccessor (const Triangulation< 1, spacedim > *tria, const VertexKind vertex_kind, const unsigned int vertex_index) | |
TriaAccessor (const Triangulation< 1, spacedim > *tria=0, const int=0, const int=0, const AccessorData *=0) | |
template<int structdim2, int dim2, int spacedim2> | |
TriaAccessor (const TriaAccessor< structdim2, dim2, spacedim2 > &) | |
template<int structdim2, int dim2, int spacedim2> | |
TriaAccessor (const InvalidAccessor< structdim2, dim2, spacedim2 > &) | |
void | copy_from (const TriaAccessor &) |
int | index () const |
bool | at_boundary () const |
types::boundary_id | boundary_indicator () const |
bool | used () const |
void | operator++ () const |
void | operator-- () const |
bool | operator== (const TriaAccessor &) const |
bool | operator!= (const TriaAccessor &) const |
void | set_boundary_indicator (const types::boundary_id) |
void | set_all_boundary_indicators (const types::boundary_id) |
unsigned int | vertex_index (const unsigned int i=0) const |
Point< spacedim > & | vertex (const unsigned int i=0) const |
Point< spacedim > | center () const |
Static Public Attributes | |
static const unsigned int | dimension =1 |
static const unsigned int | space_dimension =spacedim |
![]() | |
static const unsigned int | space_dimension = spacedim |
static const unsigned int | dimension = 1 |
static const unsigned int | structure_dimension = 0 |
Protected Member Functions | |
template<int dim2, class DH2 , bool level_dof_access2> | |
bool | operator== (const DoFAccessor< dim2, DH2, level_dof_access2 > &) const |
template<int dim2, class DH2 , bool level_dof_access2> | |
bool | operator!= (const DoFAccessor< dim2, DH2, level_dof_access2 > &) const |
void | set_dof_handler (DH< 1, spacedim > *dh) |
void | set_dof_index (const unsigned int i, const types::global_dof_index index, const unsigned int fe_index=AccessorData::default_fe_index) const |
void | set_vertex_dof_index (const unsigned int vertex, const unsigned int i, const types::global_dof_index index, const unsigned int fe_index=AccessorData::default_fe_index) const |
Protected Attributes | |
DH< 1, spacedim > * | dof_handler |
![]() | |
const Triangulation < 1, spacedim > * | tria |
VertexKind | vertex_kind |
unsigned int | global_vertex_index |
Friends | |
template<typename > | |
class | TriaRawIterator |
template<int , int > | |
class | DoFHandler |
template<int , int > | |
class | hp::DoFHandler |
struct | ::internal::DoFHandler::Policy::Implementation |
struct | ::internal::DoFHandler::Implementation |
struct | ::internal::hp::DoFHandler::Implementation |
struct | ::internal::DoFCellAccessor::Implementation |
Additional Inherited Members | |
![]() | |
static IteratorState::IteratorStates | state () |
static int | level () |
static int | parent_index () |
static bool | face_orientation (const unsigned int face) |
Always return false. | |
static bool | face_flip (const unsigned int face) |
Always return false. | |
static bool | face_rotation (const unsigned int face) |
Always return false. | |
static bool | line_orientation (const unsigned int line) |
Always return false. | |
static bool | has_children () |
static unsigned int | n_children () |
static unsigned int | number_of_children () |
static unsigned int | max_refinement_depth () |
static TriaIterator < TriaAccessor< 0, 1, spacedim > > | child (const unsigned int) |
Return an invalid object. | |
static TriaIterator < TriaAccessor< 0, 1, spacedim > > | isotropic_child (const unsigned int) |
Return an invalid object. | |
static RefinementCase< 0 > | refinement_case () |
static int | child_index (const unsigned int i) |
Returns -1. | |
static int | isotropic_child_index (const unsigned int i) |
Returns -1. | |
static typename::internal::Triangulation::Iterators < 1, spacedim >::line_iterator | line (const unsigned int) |
static unsigned int | line_index (const unsigned int i) |
static typename::internal::Triangulation::Iterators < 1, spacedim >::quad_iterator | quad (const unsigned int i) |
static unsigned int | quad_index (const unsigned int i) |
Specialization of the general DoFAccessor class template for the case of zero-dimensional objects (a vertex) that are the face of a one-dimensional cell in spacedim space dimensions. Since vertices function differently than general faces, this class does a few things differently than the general template, but the interface should look the same.
Definition at line 858 of file dof_accessor.h.
typedef TriaAccessor<0,1,spacedim> DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::BaseClass |
Declare a typedef to the base class to make accessing some of the exception classes simpler.
Definition at line 882 of file dof_accessor.h.
typedef DH<1,spacedim> DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::AccessorData |
Data type passed by the iterator class.
Definition at line 887 of file dof_accessor.h.
DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::DoFAccessor | ( | ) |
Default constructor. Provides an accessor that can't be used.
DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::DoFAccessor | ( | const Triangulation< 1, spacedim > * | tria, |
const typename TriaAccessor< 0, 1, spacedim >::VertexKind | vertex_kind, | ||
const unsigned int | vertex_index, | ||
const DH< 1, spacedim > * | dof_handler | ||
) |
Constructor to be used if the object here refers to a vertex of a one-dimensional triangulation, i.e. a face of the triangulation.
Since there is no mapping from vertices to cells, an accessor object for a point has no way to figure out whether it is at the boundary of the domain or not. Consequently, the second argument must be passed by the object that generates this accessor – e.g. a 1d cell that can figure out whether its left or right vertex are at the boundary.
The third argument is the global index of the vertex we point to.
The fourth argument is a pointer to the DoFHandler object.
This iterator can only be called for one-dimensional triangulations.
DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::DoFAccessor | ( | const Triangulation< 1, spacedim > * | , |
const int | = 0 , |
||
const int | = 0 , |
||
const DH< 1, spacedim > * | dof_handler = 0 |
||
) |
Constructor. This constructor exists in order to maintain interface compatibility with the other accessor classes. However, it doesn't do anything useful here and so may not actually be called.
DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::DoFAccessor | ( | const InvalidAccessor< structdim2, dim2, spacedim2 > & | ) |
Conversion constructor. This constructor exists to make certain constructs simpler to write in dimension independent code. For example, it allows assigning a face iterator to a line iterator, an operation that is useful in 2d but doesn't make any sense in 3d. The constructor here exists for the purpose of making the code conform to C++ but it will unconditionally abort; in other words, assigning a face iterator to a line iterator is better put into an if-statement that checks that the dimension is two, and assign to a quad iterator in 3d (an operator that, without this constructor would be illegal if we happen to compile for 2d).
DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::DoFAccessor | ( | const DoFAccessor< dim2, DH2, level_dof_access2 > & | ) |
Another conversion operator between objects that don't make sense, just like the previous one.
const DH<1,spacedim>& DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::get_dof_handler | ( | ) | const |
Return a handle on the DoFHandler object which we are using.
DoFAccessor<0,DH<1,spacedim>, level_dof_access>& DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::operator= | ( | const DoFAccessor< 0, DH< 1, spacedim >, level_dof_access > & | da | ) |
Copy operator.
void DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::copy_from | ( | const DoFAccessor< 0, DH< 1, spacedim >, level_dof_access2 > & | a | ) |
Implement the copy operator needed for the iterator classes.
void DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::copy_from | ( | const TriaAccessorBase< 0, 1, spacedim > & | da | ) |
Copy operator used by the iterator class. Keeps the previously set dof handler, but sets the object coordinates of the TriaAccessor.
TriaIterator<DoFAccessor<0,DH<1,spacedim>, level_dof_access> > DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::parent | ( | ) | const |
Return an iterator pointing to the the parent.
TriaIterator<DoFAccessor<0,DH<1,spacedim>, level_dof_access > > DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::child | ( | const unsigned int | c | ) | const |
Return an iterator pointing to the the c-th
child.
typename ::internal::DoFHandler::Iterators<DH<1,spacedim>, level_dof_access>::line_iterator DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::line | ( | const unsigned int | i | ) | const |
Pointer to the ith
line bounding this object. If the current object is a line itself, then the only valid index is i
equals to zero, and the function returns an iterator to itself.
typename ::internal::DoFHandler::Iterators<DH<1,spacedim>, level_dof_access>::quad_iterator DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::quad | ( | const unsigned int | i | ) | const |
Pointer to the ith
quad bounding this object. If the current object is a quad itself, then the only valid index is i
equals to zero, and the function returns an iterator to itself.
void DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::get_dof_indices | ( | std::vector< types::global_dof_index > & | dof_indices, |
const unsigned int | fe_index = AccessorData::default_fe_index |
||
) | const |
Return the global indices of the degrees of freedom located on this object in the standard ordering defined by the finite element (i.e., dofs on vertex 0, dofs on vertex 1, etc, dofs on line 0, dofs on line 1, etc, dofs on quad 0, etc.) This function is only available on active objects (see this glossary entry).
The cells needs to be an active cell (and not artificial in a parallel distributed computation).
The vector has to have the right size before being passed to this function.
The last argument denotes the finite element index. For the standard DoFHandler class, this value must be equal to its default value since that class only supports the same finite element on all cells anyway.
However, for hp objects (i.e. the hp::DoFHandler class), different finite element objects may be used on different cells. On faces between two cells, as well as vertices, there may therefore be two sets of degrees of freedom, one for each of the finite elements used on the adjacent cells. In order to specify which set of degrees of freedom to work on, the last argument is used to disambiguate. Finally, if this function is called for a cell object, there can only be a single set of degrees of freedom, and fe_index has to match the result of active_fe_index().
For cells, there is only a single possible finite element index (namely the one for that cell, returned by cell->active_fe_index
. Consequently, the derived DoFCellAccessor class has an overloaded version of this function that calls the present function with cell->active_fe_index
as last argument.
types::global_dof_index DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::vertex_dof_index | ( | const unsigned int | vertex, |
const unsigned int | i, | ||
const unsigned int | fe_index = AccessorData::default_fe_index |
||
) | const |
Global DoF index of the i degree associated with the vertexth
vertex of the present cell.
The last argument denotes the finite element index. For the standard DoFHandler class, this value must be equal to its default value since that class only supports the same finite element on all cells anyway.
However, for hp objects (i.e. the hp::DoFHandler class), different finite element objects may be used on different cells. On faces between two cells, as well as vertices, there may therefore be two sets of degrees of freedom, one for each of the finite elements used on the adjacent cells. In order to specify which set of degrees of freedom to work on, the last argument is used to disambiguate. Finally, if this function is called for a cell object, there can only be a single set of degrees of freedom, and fe_index has to match the result of active_fe_index().
types::global_dof_index DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::dof_index | ( | const unsigned int | i, |
const unsigned int | fe_index = AccessorData::default_fe_index |
||
) | const |
Index of the ith degree of freedom of this object.
The last argument denotes the finite element index. For the standard DoFHandler class, this value must be equal to its default value since that class only supports the same finite element on all cells anyway.
However, for hp objects (i.e. the hp::DoFHandler class), different finite element objects may be used on different cells. On faces between two cells, as well as vertices, there may therefore be two sets of degrees of freedom, one for each of the finite elements used on the adjacent cells. In order to specify which set of degrees of freedom to work on, the last argument is used to disambiguate. Finally, if this function is called for a cell object, there can only be a single set of degrees of freedom, and fe_index has to match the result of active_fe_index().
unsigned int DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::n_active_fe_indices | ( | ) | const |
Return the number of finite elements that are active on a given object.
For non-hp DoFHandler objects, the answer is of course always one. However, for hp::DoFHandler objects, this isn't the case: If this is a cell, the answer is of course one. If it is a face, the answer may be one or two, depending on whether the two adjacent cells use the same finite element or not. If it is an edge in 3d, the possible return value may be one or any other value larger than that.
unsigned int DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::nth_active_fe_index | ( | const unsigned int | n | ) | const |
Return the n-th
active fe index on this object. For cells and all non-hp objects, there is only a single active fe index, so the argument must be equal to zero. For lower-dimensional hp objects, there are n_active_fe_indices() active finite elements, and this function can be queried for their indices.
bool DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::fe_index_is_active | ( | const unsigned int | fe_index | ) | const |
Return true if the finite element with given index is active on the present object. For non-hp DoF accessors, this is of course the case only if fe_index
equals zero. For cells, it is the case if fe_index
equals active_fe_index() of this cell. For faces and other lower-dimensional objects, there may be more than one fe_index
that are active on any given object (see n_active_fe_indices()).
const FiniteElement<DH<1,spacedim>::dimension,DH<1,spacedim>::space_dimension>& DoFAccessor< 0, DH< 1, spacedim >, level_dof_access >::get_fe | ( | const unsigned int | fe_index | ) | const |
Return a reference to the finite element used on this object with the given fe_index
. fe_index
must be used on this object, i.e. fe_index_is_active(fe_index)
must return true.
|
protected |
Compare for equality.
|
protected |
Compare for inequality.
|
protected |
Reset the DoF handler pointer.
|
protected |
Set the index of the ith degree of freedom of this object to index
.
The last argument denotes the finite element index. For the standard DoFHandler class, this value must be equal to its default value since that class only supports the same finite element on all cells anyway.
However, for hp objects (i.e. the hp::DoFHandler class), different finite element objects may be used on different cells. On faces between two cells, as well as vertices, there may therefore be two sets of degrees of freedom, one for each of the finite elements used on the adjacent cells. In order to specify which set of degrees of freedom to work on, the last argument is used to disambiguate. Finally, if this function is called for a cell object, there can only be a single set of degrees of freedom, and fe_index has to match the result of active_fe_index().
|
protected |
Set the global index of the i degree on the vertex-th
vertex of the present cell to index
.
The last argument denotes the finite element index. For the standard DoFHandler class, this value must be equal to its default value since that class only supports the same finite element on all cells anyway.
However, for hp objects (i.e. the hp::DoFHandler class), different finite element objects may be used on different cells. On faces between two cells, as well as vertices, there may therefore be two sets of degrees of freedom, one for each of the finite elements used on the adjacent cells. In order to specify which set of degrees of freedom to work on, the last argument is used to disambiguate. Finally, if this function is called for a cell object, there can only be a single set of degrees of freedom, and fe_index has to match the result of active_fe_index().
|
friend |
Iterator classes need to be friends because they need to access operator== and operator!=.
Definition at line 1471 of file dof_accessor.h.
|
friend |
Make the DoFHandler class a friend so that it can call the set_xxx() functions.
Definition at line 1479 of file dof_accessor.h.
|
static |
A static variable that allows users of this class to discover the value of the second template argument.
Definition at line 867 of file dof_accessor.h.
|
static |
A static variable that allows users of this class to discover the value of the third template argument.
Definition at line 874 of file dof_accessor.h.
|
protected |
Store the address of the DoFHandler object to be accessed.
Definition at line 1366 of file dof_accessor.h.