Reference documentation for deal.II version 9.0.0
|
#include <deal.II/dofs/dof_accessor.h>
Public Types | |
typedef DoFHandlerType | AccessorData |
typedef DoFAccessor< DoFHandlerType::dimension, DoFHandlerType, level_dof_access > | BaseClass |
typedef DoFHandlerType | Container |
typedef TriaIterator< DoFAccessor< DoFHandlerType::dimension-1, DoFHandlerType, level_dof_access > > | face_iterator |
Public Types inherited from DoFAccessor< DoFHandlerType::dimension, DoFHandlerType, level_dof_access > | |
typedef typename ::internal::DoFAccessorImplementation::Inheritance< structdim, dimension, space_dimension >::BaseClass | BaseClass |
typedef DoFHandlerType | AccessorData |
Public Types inherited from TriaAccessor< structdim, dim, spacedim > | |
typedef TriaAccessorBase< structdim, dim, spacedim >::AccessorData | AccessorData |
Public Types inherited from TriaAccessorBase< structdim, dim, spacedim > | |
typedef void * | LocalData |
Public Member Functions | |
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > | parent () const |
void | set_dof_indices (const std::vector< types::global_dof_index > &dof_indices) |
void | set_mg_dof_indices (const std::vector< types::global_dof_index > &dof_indices) |
void | update_cell_dof_indices_cache () const |
Constructors and initialization | |
DoFCellAccessor (const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > *tria, const int level, const int index, const AccessorData *local_data) | |
template<int structdim2, int dim2, int spacedim2> | |
DoFCellAccessor (const InvalidAccessor< structdim2, dim2, spacedim2 > &) | |
template<int dim2, class DoFHandlerType2 , bool level_dof_access2> | |
DoFCellAccessor (const DoFAccessor< dim2, DoFHandlerType2, level_dof_access2 > &) | |
DoFCellAccessor< DoFHandlerType, level_dof_access > & | operator= (const DoFCellAccessor< DoFHandlerType, level_dof_access > &da)=delete |
Accessing sub-objects and neighbors | |
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > | neighbor (const unsigned int i) const |
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > | periodic_neighbor (const unsigned int i) const |
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > | neighbor_or_periodic_neighbor (const unsigned int i) const |
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > | child (const unsigned int i) const |
face_iterator | face (const unsigned int i) const |
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > | neighbor_child_on_subface (const unsigned int face_no, const unsigned int subface_no) const |
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > | periodic_neighbor_child_on_subface (const unsigned int face_no, const unsigned int subface_no) const |
Extracting values from global vectors | |
template<class InputVector , typename number > | |
void | get_dof_values (const InputVector &values, Vector< number > &local_values) const |
template<class InputVector , typename ForwardIterator > | |
void | get_dof_values (const InputVector &values, ForwardIterator local_values_begin, ForwardIterator local_values_end) const |
template<class InputVector , typename ForwardIterator > | |
void | get_dof_values (const ConstraintMatrix &constraints, const InputVector &values, ForwardIterator local_values_begin, ForwardIterator local_values_end) const |
template<class OutputVector , typename number > | |
void | set_dof_values (const Vector< number > &local_values, OutputVector &values) const |
template<class InputVector , typename number > | |
void | get_interpolated_dof_values (const InputVector &values, Vector< number > &interpolated_values, const unsigned int fe_index=DoFHandlerType::default_fe_index) const |
template<class OutputVector , typename number > | |
void | set_dof_values_by_interpolation (const Vector< number > &local_values, OutputVector &values, const unsigned int fe_index=DoFHandlerType::default_fe_index) const |
template<typename number , typename OutputVector > | |
void | distribute_local_to_global (const Vector< number > &local_source, OutputVector &global_destination) const |
template<typename ForwardIterator , typename OutputVector > | |
void | distribute_local_to_global (ForwardIterator local_source_begin, ForwardIterator local_source_end, OutputVector &global_destination) const |
template<typename ForwardIterator , typename OutputVector > | |
void | distribute_local_to_global (const ConstraintMatrix &constraints, ForwardIterator local_source_begin, ForwardIterator local_source_end, OutputVector &global_destination) const |
template<typename number , typename OutputMatrix > | |
void | distribute_local_to_global (const FullMatrix< number > &local_source, OutputMatrix &global_destination) const |
template<typename number , typename OutputMatrix , typename OutputVector > | |
void | distribute_local_to_global (const FullMatrix< number > &local_matrix, const Vector< number > &local_vector, OutputMatrix &global_matrix, OutputVector &global_vector) const |
Accessing the DoF indices of this object | |
void | get_active_or_mg_dof_indices (std::vector< types::global_dof_index > &dof_indices) const |
void | get_dof_indices (std::vector< types::global_dof_index > &dof_indices) const |
void | get_mg_dof_indices (std::vector< types::global_dof_index > &dof_indices) const |
Accessing the finite element associated with this object | |
const FiniteElement< DoFHandlerType::dimension, DoFHandlerType::space_dimension > & | get_fe () const |
unsigned int | active_fe_index () const |
void | set_active_fe_index (const unsigned int i) const |
Public Member Functions inherited from DoFAccessor< DoFHandlerType::dimension, DoFHandlerType, level_dof_access > | |
const DoFHandlerType & | get_dof_handler () const |
void | copy_from (const DoFAccessor< structdim, DoFHandlerType, level_dof_access2 > &a) |
void | copy_from (const TriaAccessorBase< structdim, DoFHandlerType::dimension, DoFHandlerType::space_dimension > &da) |
bool | operator== (const DoFAccessor< dim2, DoFHandlerType2, level_dof_access2 > &) const |
bool | operator!= (const DoFAccessor< dim2, DoFHandlerType2, level_dof_access2 > &) const |
DoFAccessor () | |
DoFAccessor (const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > *tria, const int level, const int index, const DoFHandlerType *dof_handler) | |
DoFAccessor (const InvalidAccessor< structdim2, dim2, spacedim2 > &) | |
DoFAccessor (const DoFAccessor< dim2, DoFHandlerType2, level_dof_access2 > &) | |
DoFAccessor (const DoFAccessor< structdim, DoFHandlerType, level_dof_access2 > &) | |
DoFAccessor< structdim, DoFHandlerType, level_dof_access > & | operator= (const DoFAccessor< structdim, DoFHandlerType, level_dof_access > &da)=delete |
TriaIterator< DoFAccessor< structdim, DoFHandlerType, level_dof_access > > | child (const unsigned int c) const |
typename ::internal::DoFHandlerImplementation::Iterators< DoFHandlerType, level_dof_access >::line_iterator | line (const unsigned int i) const |
typename ::internal::DoFHandlerImplementation::Iterators< DoFHandlerType, level_dof_access >::quad_iterator | quad (const unsigned int i) const |
void | get_dof_indices (std::vector< types::global_dof_index > &dof_indices, const unsigned int fe_index=DoFHandlerType::default_fe_index) const |
void | get_mg_dof_indices (const int level, std::vector< types::global_dof_index > &dof_indices, const unsigned int fe_index=DoFHandlerType::default_fe_index) const |
void | set_mg_dof_indices (const int level, const std::vector< types::global_dof_index > &dof_indices, const unsigned int fe_index=DoFHandlerType::default_fe_index) |
types::global_dof_index | vertex_dof_index (const unsigned int vertex, const unsigned int i, const unsigned int fe_index=DoFHandlerType::default_fe_index) 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=DoFHandlerType::default_fe_index) const |
types::global_dof_index | dof_index (const unsigned int i, const unsigned int fe_index=DoFHandlerType::default_fe_index) const |
types::global_dof_index | mg_dof_index (const int level, const unsigned int i) const |
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< DoFHandlerType::dimension, DoFHandlerType::space_dimension > & | get_fe (const unsigned int fe_index) const |
Public Member Functions inherited from TriaAccessor< structdim, dim, spacedim > | |
TriaAccessor (const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr) | |
template<int structdim2, int dim2, int spacedim2> | |
TriaAccessor (const InvalidAccessor< structdim2, dim2, spacedim2 > &) | |
template<int structdim2, int dim2, int spacedim2> | |
TriaAccessor (const TriaAccessor< structdim2, dim2, spacedim2 > &) | |
void | operator= (const TriaAccessor &)=delete |
bool | used () const |
template<> | |
double | extent_in_direction (const unsigned int axis) const |
template<> | |
double | extent_in_direction (const unsigned int axis) const |
template<> | |
double | extent_in_direction (const unsigned int axis) const |
TriaIterator< TriaAccessor< 0, dim, spacedim > > | vertex_iterator (const unsigned int i) const |
unsigned int | vertex_index (const unsigned int i) const |
Point< spacedim > & | vertex (const unsigned int i) const |
typename ::internal::TriangulationImplementation::Iterators< dim, spacedim >::line_iterator | line (const unsigned int i) const |
unsigned int | line_index (const unsigned int i) const |
typename ::internal::TriangulationImplementation::Iterators< dim, spacedim >::quad_iterator | quad (const unsigned int i) const |
unsigned int | quad_index (const unsigned int i) const |
bool | face_orientation (const unsigned int face) const |
bool | face_flip (const unsigned int face) const |
bool | face_rotation (const unsigned int face) const |
bool | line_orientation (const unsigned int line) const |
bool | has_children () const |
unsigned int | n_children () const |
unsigned int | number_of_children () const |
unsigned int | max_refinement_depth () const |
TriaIterator< TriaAccessor< structdim, dim, spacedim > > | child (const unsigned int i) const |
TriaIterator< TriaAccessor< structdim, dim, spacedim > > | isotropic_child (const unsigned int i) const |
RefinementCase< structdim > | refinement_case () const |
int | child_index (const unsigned int i) const |
int | isotropic_child_index (const unsigned int i) const |
types::boundary_id | boundary_id () const |
void | set_boundary_id (const types::boundary_id) const |
void | set_all_boundary_ids (const types::boundary_id) const |
bool | at_boundary () const |
const Boundary< dim, spacedim > & | get_boundary () const |
const Manifold< dim, spacedim > & | get_manifold () const |
types::manifold_id | manifold_id () const |
void | set_manifold_id (const types::manifold_id) const |
void | set_all_manifold_ids (const types::manifold_id) const |
bool | user_flag_set () const |
void | set_user_flag () const |
void | clear_user_flag () const |
void | recursively_set_user_flag () const |
void | recursively_clear_user_flag () const |
void | clear_user_data () const |
void | set_user_pointer (void *p) const |
void | clear_user_pointer () const |
void * | user_pointer () const |
void | recursively_set_user_pointer (void *p) const |
void | recursively_clear_user_pointer () const |
void | set_user_index (const unsigned int p) const |
void | clear_user_index () const |
unsigned int | user_index () const |
void | recursively_set_user_index (const unsigned int p) const |
void | recursively_clear_user_index () const |
double | diameter () const |
std::pair< Point< spacedim >, double > | enclosing_ball () const |
BoundingBox< spacedim > | bounding_box () const |
double | extent_in_direction (const unsigned int axis) const |
double | minimum_vertex_distance () const |
Point< spacedim > | intermediate_point (const Point< structdim > &coordinates) const |
Point< structdim > | real_to_unit_cell_affine_approximation (const Point< spacedim > &point) const |
Point< spacedim > | center (const bool respect_manifold=false, const bool interpolate_from_surrounding=false) const |
Point< spacedim > | barycenter () const |
double | measure () const |
bool | is_translation_of (const TriaIterator< TriaAccessor< structdim, dim, spacedim > > &o) const |
Public Member Functions inherited from TriaAccessorBase< structdim, dim, spacedim > | |
int | level () const |
int | index () const |
IteratorState::IteratorStates | state () const |
const Triangulation< dim, spacedim > & | get_triangulation () const |
Static Public Attributes | |
static const unsigned int | dim = DoFHandlerType::dimension |
static const unsigned int | spacedim = DoFHandlerType::space_dimension |
Static Public Attributes inherited from DoFAccessor< DoFHandlerType::dimension, DoFHandlerType, level_dof_access > | |
static const unsigned int | dimension |
static const unsigned int | space_dimension |
Static Public Attributes inherited from TriaAccessorBase< structdim, dim, spacedim > | |
static const unsigned int | space_dimension = spacedim |
static const unsigned int | dimension = dim |
static const unsigned int | structure_dimension = structdim |
Friends | |
template<int dim, int spacedim> | |
class | DoFHandler |
Grant access to the degrees of freedom on a cell.
Note that since for the class we derive from, i.e. DoFAccessor<dim>
, the two template parameters are equal, the base class is actually derived from CellAccessor, which makes the functions of this class available to the DoFCellAccessor class as well.
Definition at line 1236 of file dof_accessor.h.
typedef DoFHandlerType DoFCellAccessor< DoFHandlerType, level_dof_access >::AccessorData |
Data type passed by the iterator class.
Definition at line 1253 of file dof_accessor.h.
typedef DoFAccessor<DoFHandlerType::dimension,DoFHandlerType, level_dof_access> DoFCellAccessor< DoFHandlerType, level_dof_access >::BaseClass |
Declare a typedef to the base class to make accessing some of the exception classes simpler.
Definition at line 1259 of file dof_accessor.h.
typedef DoFHandlerType DoFCellAccessor< DoFHandlerType, level_dof_access >::Container |
Define the type of the container this is part of.
Definition at line 1264 of file dof_accessor.h.
typedef TriaIterator<DoFAccessor<DoFHandlerType::dimension-1, DoFHandlerType, level_dof_access> > DoFCellAccessor< DoFHandlerType, level_dof_access >::face_iterator |
A type for an iterator over the faces of a cell. This is what the face() function returns.
Definition at line 1272 of file dof_accessor.h.
DoFCellAccessor< DoFHandlerType, level_dof_access >::DoFCellAccessor | ( | const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > * | tria, |
const int | level, | ||
const int | index, | ||
const AccessorData * | local_data | ||
) |
Constructor
DoFCellAccessor< DoFHandlerType, level_dof_access >::DoFCellAccessor | ( | 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).
|
explicit |
Another conversion operator between objects that don't make sense, just like the previous one.
|
delete |
Copy operator. These operators are usually used in a context like iterator a,b; *a=*b;
. Presumably, the intent here is to copy the object pointed to by b
to the object pointed to by a
. However, the result of dereferencing an iterator is not an object but an accessor; consequently, this operation is not useful for iterators on DoF handler objects. Consequently, this operator is declared as deleted and can not be used.
TriaIterator<DoFCellAccessor<DoFHandlerType, level_dof_access> > DoFCellAccessor< DoFHandlerType, level_dof_access >::parent | ( | ) | const |
Return the parent of this cell as a DoF cell iterator. If the parent does not exist (i.e., if the object is at the coarsest level of the mesh hierarchy), an exception is generated.
This function is needed since the parent function of the base class CellAccessor returns a triangulation cell accessor without access to the DoF data.
TriaIterator<DoFCellAccessor<DoFHandlerType, level_dof_access> > DoFCellAccessor< DoFHandlerType, level_dof_access >::neighbor | ( | const unsigned int | i | ) | const |
Return the ith
neighbor as a DoF cell iterator. This function is needed since the neighbor function of the base class returns a cell accessor without access to the DoF data.
|
inline |
Return the ith
periodic neighbor as a DoF cell iterator. This function is needed since the neighbor function of the base class returns a cell accessor without access to the DoF data.
Definition at line 149 of file dof_accessor.cc.
|
inline |
Return the ith
neighbor or periodic neighbor as a DoF cell iterator. This function is needed since the neighbor function of the base class returns a cell accessor without access to the DoF data.
Definition at line 161 of file dof_accessor.cc.
TriaIterator<DoFCellAccessor<DoFHandlerType, level_dof_access> > DoFCellAccessor< DoFHandlerType, level_dof_access >::child | ( | const unsigned int | i | ) | const |
Return the ith
child as a DoF cell iterator. This function is needed since the child function of the base class returns a cell accessor without access to the DoF data.
face_iterator DoFCellAccessor< DoFHandlerType, level_dof_access >::face | ( | const unsigned int | i | ) | const |
Return an iterator to the ith
face of this cell.
This function returns a DoFAccessor with structdim == 0
in 1D, a DoFAccessor::line in 2D, and a DoFAccessor::quad in 3d.
TriaIterator< DoFCellAccessor< DoFHandlerType, lda > > DoFCellAccessor< DoFHandlerType, lda >::neighbor_child_on_subface | ( | const unsigned int | face_no, |
const unsigned int | subface_no | ||
) | const |
Return the result of the neighbor_child_on_subface
function of the base class, but convert it so that one can also access the DoF data (the function in the base class only returns an iterator with access to the triangulation data).
Definition at line 124 of file dof_accessor.cc.
TriaIterator< DoFCellAccessor< DoFHandlerType, lda > > DoFCellAccessor< DoFHandlerType, lda >::periodic_neighbor_child_on_subface | ( | const unsigned int | face_no, |
const unsigned int | subface_no | ||
) | const |
Return the result of the periodic_neighbor_child_on_subface
function of the base class, but convert it so that one can also access the DoF data (the function in the base class only returns an iterator with access to the triangulation data).
Definition at line 136 of file dof_accessor.cc.
void DoFCellAccessor< DoFHandlerType, level_dof_access >::get_dof_values | ( | const InputVector & | values, |
Vector< number > & | local_values | ||
) | const |
Return the values of the given vector restricted to the dofs of this cell in the standard ordering: dofs on vertex 0, dofs on vertex 1, etc, dofs on line 0, dofs on line 1, etc, dofs on quad 0, etc.
The vector has to have the right size before being passed to this function. This function is only callable for active cells.
The input vector may be either a Vector<float>
, Vector<double>, or a BlockVector<double>, or a PETSc or Trilinos vector if deal.II is compiled to support these libraries. It is in the responsibility of the caller to assure that the types of the numbers stored in input and output vectors are compatible and with similar accuracy.
void DoFCellAccessor< DoFHandlerType, level_dof_access >::get_dof_values | ( | const InputVector & | values, |
ForwardIterator | local_values_begin, | ||
ForwardIterator | local_values_end | ||
) | const |
Return the values of the given vector restricted to the dofs of this cell in the standard ordering: dofs on vertex 0, dofs on vertex 1, etc, dofs on line 0, dofs on line 1, etc, dofs on quad 0, etc.
The vector has to have the right size before being passed to this function. This function is only callable for active cells.
The input vector may be either a Vector<float>
, Vector<double>, or a BlockVector<double>, or a PETSc or Trilinos vector if deal.II is compiled to support these libraries. It is in the responsibility of the caller to assure that the types of the numbers stored in input and output vectors are compatible and with similar accuracy.
void DoFCellAccessor< DoFHandlerType, level_dof_access >::get_dof_values | ( | const ConstraintMatrix & | constraints, |
const InputVector & | values, | ||
ForwardIterator | local_values_begin, | ||
ForwardIterator | local_values_end | ||
) | const |
Return the values of the given vector restricted to the dofs of this cell in the standard ordering: dofs on vertex 0, dofs on vertex 1, etc, dofs on line 0, dofs on line 1, etc, dofs on quad 0, etc.
The vector has to have the right size before being passed to this function. This function is only callable for active cells.
The input vector may be either a Vector<float>
, Vector<double>, or a BlockVector<double>, or a PETSc or Trilinos vector if deal.II is compiled to support these libraries. It is in the responsibility of the caller to assure that the types of the numbers stored in input and output vectors are compatible and with similar accuracy. The ConstraintMatrix passed as an argument to this function makes sure that constraints are correctly distributed when the dof values are calculated.
void DoFCellAccessor< DoFHandlerType, level_dof_access >::set_dof_values | ( | const Vector< number > & | local_values, |
OutputVector & | values | ||
) | const |
This function is the counterpart to get_dof_values(): it takes a vector of values for the degrees of freedom of the cell pointed to by this iterator and writes these values into the global data vector values
. This function is only callable for active cells.
Note that for continuous finite elements, calling this function affects the dof values on neighboring cells as well. It may also violate continuity requirements for hanging nodes, if neighboring cells are less refined than the present one. These requirements are not taken care of and must be enforced by the user afterwards.
The vector has to have the right size before being passed to this function.
The output vector may be either a Vector<float>, Vector<double>, or a BlockVector<double>, or a PETSc vector if deal.II is compiled to support these libraries. It is in the responsibility of the caller to assure that the types of the numbers stored in input and output vectors are compatible and with similar accuracy.
void DoFCellAccessor< DoFHandlerType, lda >::get_interpolated_dof_values | ( | const InputVector & | values, |
Vector< number > & | interpolated_values, | ||
const unsigned int | fe_index = DoFHandlerType::default_fe_index |
||
) | const |
Return the interpolation of the given finite element function to the present cell. In the simplest case, the cell is a terminal one, i.e., it has no children; then, the returned value is the vector of nodal values on that cell. You could as well get the desired values through the get_dof_values
function. In the other case, when the cell has children, we use the restriction matrices provided by the finite element class to compute the interpolation from the children to the present cell.
If the cell is part of a hp::DoFHandler object, cells only have an associated finite element space if they are active. However, this function is supposed to also provide information on inactive cells with children. Consequently, it carries a third argument that can be used in the hp context that denotes the finite element space we are supposed to interpolate onto. If the cell is active, this function then obtains the finite element function from the values
vector on this cell and interpolates it onto the space described by the fe_index
th element of the hp::FECollection associated with the hp::DoFHandler of which this cell is a part of. If the cell is not active, then we first perform this interpolation on all of its terminal children and then interpolate this function down to the cell requested keeping the function space the same.
It is assumed that both input vectors already have the right size beforehand.
Definition at line 45 of file dof_accessor_get.cc.
void DoFCellAccessor< DoFHandlerType, lda >::set_dof_values_by_interpolation | ( | const Vector< number > & | local_values, |
OutputVector & | values, | ||
const unsigned int | fe_index = DoFHandlerType::default_fe_index |
||
) | const |
This function is the counterpart to get_interpolated_dof_values(): you specify the dof values on a cell and these are interpolated to the children of the present cell and set on the terminal cells.
In principle, it works as follows: if the cell pointed to by this object is terminal (i.e., has no children), then the dof values are set in the global data vector by calling the set_dof_values() function; otherwise, the values are prolonged to each of the children and this function is called for each of them.
Using the get_interpolated_dof_values() and this function, you can compute the interpolation of a finite element function to a coarser grid by first getting the interpolated solution on a cell of the coarse grid and afterwards redistributing it using this function.
Note that for continuous finite elements, calling this function affects the dof values on neighboring cells as well. It may also violate continuity requirements for hanging nodes, if neighboring cells are less refined than the present one, or if their children are less refined than the children of this cell. These requirements are not taken care of and must be enforced by the user afterward.
If the cell is part of a hp::DoFHandler object, cells only have an associated finite element space if they are active. However, this function is supposed to also work on inactive cells with children. Consequently, it carries a third argument that can be used in the hp context that denotes the finite element space we are supposed to interpret the input vector of this function in. If the cell is active, this function then interpolates the input vector interpreted as an element of the space described by the fe_index
th element of the hp::FECollection associated with the hp::DoFHandler of which this cell is a part of, and interpolates it into the space that is associated with this cell. On the other hand, if the cell is not active, then we first perform this interpolation from this cell to its children using the given fe_index
until we end up on an active cell, at which point we follow the procedure outlined at the beginning of the paragraph.
It is assumed that both vectors already have the right size beforehand. This function relies on the existence of a natural interpolation property of finite element spaces of a cell to its children, denoted by the prolongation matrices of finite element classes. For some elements, the spaces on coarse and fine grids are not nested, in which case the interpolation to a child is not the identity; refer to the documentation of the respective finite element class for a description of what the prolongation matrices represent in this case.
Definition at line 45 of file dof_accessor_set.cc.
void DoFCellAccessor< DoFHandlerType, level_dof_access >::distribute_local_to_global | ( | const Vector< number > & | local_source, |
OutputVector & | global_destination | ||
) | const |
Distribute a local (cell based) vector to a global one by mapping the local numbering of the degrees of freedom to the global one and entering the local values into the global vector.
The elements are added up to the elements in the global vector, rather than just set, since this is usually what one wants.
void DoFCellAccessor< DoFHandlerType, level_dof_access >::distribute_local_to_global | ( | ForwardIterator | local_source_begin, |
ForwardIterator | local_source_end, | ||
OutputVector & | global_destination | ||
) | const |
Distribute a local (cell based) vector in iterator format to a global one by mapping the local numbering of the degrees of freedom to the global one and entering the local values into the global vector.
The elements are added up to the elements in the global vector, rather than just set, since this is usually what one wants.
void DoFCellAccessor< DoFHandlerType, level_dof_access >::distribute_local_to_global | ( | const ConstraintMatrix & | constraints, |
ForwardIterator | local_source_begin, | ||
ForwardIterator | local_source_end, | ||
OutputVector & | global_destination | ||
) | const |
Distribute a local (cell based) vector in iterator format to a global one by mapping the local numbering of the degrees of freedom to the global one and entering the local values into the global vector.
The elements are added up to the elements in the global vector, rather than just set, since this is usually what one wants. Moreover, the ConstraintMatrix passed to this function makes sure that also constraints are eliminated in this process.
void DoFCellAccessor< DoFHandlerType, level_dof_access >::distribute_local_to_global | ( | const FullMatrix< number > & | local_source, |
OutputMatrix & | global_destination | ||
) | const |
This function does much the same as the distribute_local_to_global(Vector,Vector)
function, but operates on matrices instead of vectors. If the matrix type is a sparse matrix then it is supposed to have non-zero entry slots where required.
void DoFCellAccessor< DoFHandlerType, level_dof_access >::distribute_local_to_global | ( | const FullMatrix< number > & | local_matrix, |
const Vector< number > & | local_vector, | ||
OutputMatrix & | global_matrix, | ||
OutputVector & | global_vector | ||
) | const |
This function does what the two distribute_local_to_global
functions with vector and matrix argument do, but all at once.
void DoFCellAccessor< DoFHandlerType, level_dof_access >::get_active_or_mg_dof_indices | ( | std::vector< types::global_dof_index > & | dof_indices | ) | const |
Obtain the global indices of the local degrees of freedom on this cell.
If this object accesses a level cell (indicated by the third template argument or is_level_cell), then return the result of get_mg_dof_indices(), else return get_dof_indices().
You will get a level_cell_iterator when calling begin_mg() and a normal one otherwise.
Examples for this use are in the implementation of DoFRenumbering.
void DoFCellAccessor< DoFHandlerType, level_dof_access >::get_dof_indices | ( | std::vector< types::global_dof_index > & | dof_indices | ) | 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).
[out] | dof_indices | The vector into which the indices will be written. It has to have the right size (namely, fe.dofs_per_cell , fe.dofs_per_face , or fe.dofs_per_line , depending on which kind of object this function is called) before being passed to this function. |
This function reimplements the same function in the base class. In contrast to the function in the base class, we do not need the fe_index
here because there is always a unique finite element index on cells.
This is a function which requires that the cell is active.
Also see get_active_or_mg_dof_indices().
local_dof_indices
by convention. The name is not meant to indicate the local numbers of degrees of freedom (which are always between zero and fe.dofs_per_cell
) but instead that the returned values are the global indices of those degrees of freedom that are located locally on the current cell.void DoFCellAccessor< DoFHandlerType, level_dof_access >::get_mg_dof_indices | ( | std::vector< types::global_dof_index > & | dof_indices | ) | const |
Retrieve the global indices of the degrees of freedom on this cell in the level vector associated to the level of the cell.
const FiniteElement<DoFHandlerType::dimension,DoFHandlerType::space_dimension>& DoFCellAccessor< DoFHandlerType, level_dof_access >::get_fe | ( | ) | const |
Return the finite element that is used on the cell pointed to by this iterator. For non-hp DoF handlers, this is of course always the same element, independent of the cell we are presently on, but for hp DoF handlers, this may change from cell to cell.
unsigned int DoFCellAccessor< DoFHandlerType, level_dof_access >::active_fe_index | ( | ) | const |
Return the index inside the hp::FECollection of the FiniteElement used for this cell. This function is only useful if the DoF handler object associated with the current cell is an hp::DoFHandler.
active_fe_index
information is exchanged from locally owned cells on one processor to other processors where they may be ghost cells, during the call to hp::DoFHandler::distribute_dofs(). See the documentation of hp::DoFHandler for more information. void DoFCellAccessor< DoFHandlerType, level_dof_access >::set_active_fe_index | ( | const unsigned int | i | ) | const |
Set the index of the FiniteElement used for this cell. This determines which element in an hp::FECollection to use. This function is only useful if the DoF handler object associated with the current cell is an hp::DoFHandler.
active_fe_index
on a ghost cell than the processor that actually owns the cell does. To avoid this mistake, one can only set active_fe_index
information on locally owned cells, and this information is then mirrored to all processors that have this cell as a ghost cell – see the documentation of the hp::DoFHandler class. void DoFCellAccessor< DoFHandlerType, lda >::set_dof_indices | ( | const std::vector< types::global_dof_index > & | dof_indices | ) |
Set the DoF indices of this cell to the given values. This function bypasses the DoF cache, if one exists for the given DoF handler class.
Definition at line 108 of file dof_accessor.cc.
void DoFCellAccessor< DoFHandlerType, level_dof_access >::set_mg_dof_indices | ( | const std::vector< types::global_dof_index > & | dof_indices | ) |
Set the Level DoF indices of this cell to the given values.
void DoFCellAccessor< DoFHandlerType, lda >::update_cell_dof_indices_cache | ( | ) | const |
Update the cache in which we store the dof indices of this cell.
Definition at line 93 of file dof_accessor.cc.
|
friend |
Make the DoFHandler class a friend so that it can call the update_cell_dof_indices_cache() function
Definition at line 1839 of file dof_accessor.h.
|
static |
Extract dimension from DoFHandlerType.
Definition at line 1242 of file dof_accessor.h.
|
static |
Extract space dimension from DoFHandlerType.
Definition at line 1247 of file dof_accessor.h.