Reference documentation for deal.II version 8.1.0
Public Types | Public Member Functions | Static Public Attributes | Private Member Functions | Friends | List of all members
DoFCellAccessor< DH, level_dof_access > Class Template Reference

#include <dof_accessor.h>

Inheritance diagram for DoFCellAccessor< DH, level_dof_access >:
[legend]

Public Types

typedef DH AccessorData
 
typedef DoFAccessor< DH::dimension, DH, level_dof_access > BaseClass
 
typedef DH Container
 
- Public Types inherited from DoFAccessor< DH::dimension, DH, level_dof_access >
typedef typename::internal::DoFAccessor::Inheritance< structdim, dimension, space_dimension >::BaseClass BaseClass
 
typedef DH 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< DH, 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< DH::dimension, DH::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 DH2 , bool level_dof_access2>
 DoFCellAccessor (const DoFAccessor< dim2, DH2, level_dof_access2 > &)
 
Accessing sub-objects and neighbors
TriaIterator< DoFCellAccessor< DH, level_dof_access > > neighbor (const unsigned int) const
 
TriaIterator< DoFCellAccessor< DH, level_dof_access > > child (const unsigned int) const
 
TriaIterator< DoFAccessor< DH::dimension-1, DH, level_dof_access > > face (const unsigned int i) const
 
TriaIterator< DoFCellAccessor< DH, level_dof_access > > 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
 
template<class OutputVector , typename number >
void set_dof_values_by_interpolation (const Vector< number > &local_values, OutputVector &values) 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< DH::dimension, DH::space_dimension > & get_fe () const
 
unsigned int active_fe_index () const
 
void set_active_fe_index (const unsigned int i)
 
- Public Member Functions inherited from DoFAccessor< DH::dimension, DH, level_dof_access >
const DH & get_dof_handler () const
 
void copy_from (const DoFAccessor< structdim, DH, level_dof_access2 > &a)
 
void copy_from (const TriaAccessorBase< structdim, DH::dimension, DH::space_dimension > &da)
 
TriaIterator< DoFAccessor< structdim, DH, level_dof_access > > parent () const
 
 DeclException0 (ExcInvalidObject)
 
 DeclException0 (ExcVectorNotEmpty)
 
 DeclException0 (ExcVectorDoesNotMatch)
 
 DeclException0 (ExcMatrixDoesNotMatch)
 
 DeclException0 (ExcNotActive)
 
 DeclException0 (ExcCantCompareIterators)
 
bool operator== (const DoFAccessor< dim2, DH2, level_dof_access2 > &) const
 
bool operator!= (const DoFAccessor< dim2, DH2, level_dof_access2 > &) const
 
 DoFAccessor ()
 
 DoFAccessor (const Triangulation< DH::dimension, DH::space_dimension > *tria, const int level, const int index, const DH *local_data)
 
 DoFAccessor (const InvalidAccessor< structdim2, dim2, spacedim2 > &)
 
 DoFAccessor (const DoFAccessor< dim2, DH2, level_dof_access2 > &)
 
 DoFAccessor (const DoFAccessor< structdim, DH, level_dof_access2 > &)
 
TriaIterator< DoFAccessor< structdim, DH, level_dof_access > > child (const unsigned int c) const
 
typename::internal::DoFHandler::Iterators< DH, level_dof_access >::line_iterator line (const unsigned int i) const
 
typename::internal::DoFHandler::Iterators< DH, 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=DH::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=DH::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=DH::default_fe_index)
 
types::global_dof_index vertex_dof_index (const unsigned int vertex, const unsigned int i, const unsigned int fe_index=DH::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=DH::default_fe_index) const
 
types::global_dof_index dof_index (const unsigned int i, const unsigned int fe_index=DH::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< DH::dimension, DH::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=0, const int level=-1, const int index=-1, const AccessorData *local_data=0)
 
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 > &)
 
bool used () const
 
int parent_index () const
 
unsigned int vertex_index (const unsigned int i) const
 
Point< spacedim > & vertex (const unsigned int i) const
 
typename::internal::Triangulation::Iterators< dim, spacedim >::line_iterator line (const unsigned int i) const
 
unsigned int line_index (const unsigned int i) const
 
typename::internal::Triangulation::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_indicator () const
 
void set_boundary_indicator (const types::boundary_id) const
 
void set_all_boundary_indicators (const types::boundary_id) const
 
bool at_boundary () const
 
const Boundary< dim, spacedim > & get_boundary () 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
 
double extent_in_direction (const unsigned int axis) const
 
double minimum_vertex_distance () const
 
Point< spacedim > center () 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 = DH::dimension
 
static const unsigned int spacedim = DH::space_dimension
 
- Static Public Attributes inherited from DoFAccessor< DH::dimension, DH, 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
 

Private Member Functions

DoFCellAccessor< DH, level_dof_access > & operator= (const DoFCellAccessor< DH, level_dof_access > &da)
 

Friends

template<int dim, int spacedim>
class DoFHandler
 
struct ::internal::DoFCellAccessor::Implementation
 

Additional Inherited Members

- Static Public Member Functions inherited from DoFAccessor< DH::dimension, DH, level_dof_access >
static bool is_level_cell ()
 
- Protected Types inherited from TriaAccessorBase< structdim, dim, spacedim >
typedef void AccessorData
 
- Protected Member Functions inherited from DoFAccessor< DH::dimension, DH, level_dof_access >
void set_dof_handler (DH *dh)
 
void set_dof_index (const unsigned int i, const types::global_dof_index index, const unsigned int fe_index=DH::default_fe_index) const
 
void set_mg_dof_index (const int level, const unsigned int i, const types::global_dof_index 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=DH::default_fe_index) const
 
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=DH::default_fe_index) const
 
- Protected Member Functions inherited from TriaAccessorBase< structdim, dim, spacedim >
 TriaAccessorBase (const Triangulation< dim, spacedim > *parent=0, const int level=-1, const int index=-1, const AccessorData *=0)
 
 TriaAccessorBase (const TriaAccessorBase &)
 
void copy_from (const TriaAccessorBase &)
 
TriaAccessorBaseoperator= (const TriaAccessorBase &)
 
bool operator< (const TriaAccessorBase &other) const
 
void operator= (const TriaAccessorBase *)
 
bool operator== (const TriaAccessorBase &) const
 
bool operator!= (const TriaAccessorBase &) const
 
::internal::Triangulation::TriaObjects<::internal::Triangulation::TriaObject< structdim > > & objects () const
 
void operator++ ()
 
void operator-- ()
 
- Protected Attributes inherited from DoFAccessor< DH::dimension, DH, level_dof_access >
DH * dof_handler
 
- Protected Attributes inherited from TriaAccessorBase< structdim, dim, spacedim >
typename::internal::TriaAccessor::PresentLevelType< structdim, dim >::type present_level
 
int present_index
 
const Triangulation< dim, spacedim > * tria
 

Detailed Description

template<class DH, bool level_dof_access>
class DoFCellAccessor< DH, level_dof_access >

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.

Author
Wolfgang Bangerth, 1998, Timo Heister, Guido Kanschat, 2012

Definition at line 1505 of file dof_accessor.h.

Member Typedef Documentation

template<class DH, bool level_dof_access>
typedef DH DoFCellAccessor< DH, level_dof_access >::AccessorData

Data type passed by the iterator class.

Definition at line 1522 of file dof_accessor.h.

template<class DH, bool level_dof_access>
typedef DoFAccessor<DH::dimension,DH, level_dof_access> DoFCellAccessor< DH, level_dof_access >::BaseClass

Declare a typedef to the base class to make accessing some of the exception classes simpler.

Definition at line 1530 of file dof_accessor.h.

template<class DH, bool level_dof_access>
typedef DH DoFCellAccessor< DH, level_dof_access >::Container

Define the type of the container this is part of.

Definition at line 1536 of file dof_accessor.h.

Constructor & Destructor Documentation

template<class DH , bool lda>
DoFCellAccessor< DH, lda >::DoFCellAccessor ( const Triangulation< DH::dimension, DH::space_dimension > *  tria,
const int  level,
const int  index,
const AccessorData local_data 
)
inline

Constructor

Definition at line 3130 of file dof_accessor.templates.h.

template<class DH , bool lda>
template<int structdim2, int dim2, int spacedim2>
DoFCellAccessor< DH, lda >::DoFCellAccessor ( const InvalidAccessor< structdim2, dim2, spacedim2 > &  )
inline

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).

Definition at line 3143 of file dof_accessor.templates.h.

template<class DH , bool lda>
template<int dim2, class DH2, bool lda2>
DoFCellAccessor< DH, lda >::DoFCellAccessor ( const DoFAccessor< dim2, DH2, lda2 > &  other)
inlineexplicit

Another conversion operator between objects that don't make sense, just like the previous one.

Definition at line 3153 of file dof_accessor.templates.h.

Member Function Documentation

template<class DH , bool lda>
TriaIterator< DoFCellAccessor< DH, lda > > DoFCellAccessor< DH, lda >::parent ( ) const
inline

Return the parent as a DoF cell iterator. This function is needed since the parent function of the base class returns a cell accessor without access to the DoF data.

Definition at line 3200 of file dof_accessor.templates.h.

template<class DH , bool lda>
TriaIterator< DoFCellAccessor< DH, lda > > DoFCellAccessor< DH, lda >::neighbor ( const unsigned int  i) const
inline

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.

Definition at line 3162 of file dof_accessor.templates.h.

template<class DH , bool lda>
TriaIterator< DoFCellAccessor< DH, lda > > DoFCellAccessor< DH, lda >::child ( const unsigned int  i) const
inline

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.

Definition at line 3181 of file dof_accessor.templates.h.

template<class DH , bool lda>
TriaIterator< DoFAccessor< DH::dimension-1, DH, lda > > DoFCellAccessor< DH, lda >::face ( const unsigned int  i) const
inline

Return an iterator to the ith face of this cell.

This function is not implemented in 1D, and maps to DoFAccessor::line in 2D.

Definition at line 3267 of file dof_accessor.templates.h.

template<class DH, bool level_dof_access>
TriaIterator<DoFCellAccessor<DH, level_dof_access> > DoFCellAccessor< DH, level_dof_access >::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).

template<class DH , bool lda>
template<class InputVector , typename number >
void DoFCellAccessor< DH, lda >::get_dof_values ( const InputVector &  values,
Vector< number > &  local_values 
) const
inline

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.

Definition at line 3332 of file dof_accessor.templates.h.

template<class DH , bool lda>
template<class InputVector , typename ForwardIterator >
void DoFCellAccessor< DH, lda >::get_dof_values ( const InputVector &  values,
ForwardIterator  local_values_begin,
ForwardIterator  local_values_end 
) const
inline

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.

Definition at line 3344 of file dof_accessor.templates.h.

template<class DH , bool lda>
template<class InputVector , typename ForwardIterator >
void DoFCellAccessor< DH, lda >::get_dof_values ( const ConstraintMatrix constraints,
const InputVector &  values,
ForwardIterator  local_values_begin,
ForwardIterator  local_values_end 
) const
inline

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.

Definition at line 3374 of file dof_accessor.templates.h.

template<class DH , bool lda>
template<class OutputVector , typename number >
void DoFCellAccessor< DH, lda >::set_dof_values ( const Vector< number > &  local_values,
OutputVector &  values 
) const
inline

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.

Definition at line 3405 of file dof_accessor.templates.h.

template<class DH, bool level_dof_access>
template<class InputVector , typename number >
void DoFCellAccessor< DH, level_dof_access >::get_interpolated_dof_values ( const InputVector &  values,
Vector< number > &  interpolated_values 
) 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. has no children; then, the returned value is the vector of nodal values on that cell. You could then 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.

It is assumed that both vectors already have the right size beforehand.

Unlike the get_dof_values() function, this function works on cells rather than to lines, quads, and hexes, since interpolation is presently only provided for cells by the finite element classes.

template<class DH, bool level_dof_access>
template<class OutputVector , typename number >
void DoFCellAccessor< DH, level_dof_access >::set_dof_values_by_interpolation ( const Vector< number > &  local_values,
OutputVector &  values 
) const

This, again, 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, 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.

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.

Unlike the set_dof_values() function, this function is associated to cells rather than to lines, quads, and hexes, since interpolation is presently only provided for cells by the finite element objects.

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.

template<class DH , bool lda>
template<typename number , typename OutputVector >
void DoFCellAccessor< DH, lda >::distribute_local_to_global ( const Vector< number > &  local_source,
OutputVector &  global_destination 
) const
inline

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.

Definition at line 3465 of file dof_accessor.templates.h.

template<class DH , bool lda>
template<typename ForwardIterator , typename OutputVector >
void DoFCellAccessor< DH, lda >::distribute_local_to_global ( ForwardIterator  local_source_begin,
ForwardIterator  local_source_end,
OutputVector &  global_destination 
) const
inline

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.

Definition at line 3480 of file dof_accessor.templates.h.

template<class DH , bool lda>
template<typename ForwardIterator , typename OutputVector >
void DoFCellAccessor< DH, lda >::distribute_local_to_global ( const ConstraintMatrix constraints,
ForwardIterator  local_source_begin,
ForwardIterator  local_source_end,
OutputVector &  global_destination 
) const
inline

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.

Definition at line 3496 of file dof_accessor.templates.h.

template<class DH , bool lda>
template<typename number , typename OutputMatrix >
void DoFCellAccessor< DH, lda >::distribute_local_to_global ( const FullMatrix< number > &  local_source,
OutputMatrix &  global_destination 
) const
inline

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.

Definition at line 3513 of file dof_accessor.templates.h.

template<class DH , bool lda>
template<typename number , typename OutputMatrix , typename OutputVector >
void DoFCellAccessor< DH, lda >::distribute_local_to_global ( const FullMatrix< number > &  local_matrix,
const Vector< number > &  local_vector,
OutputMatrix &  global_matrix,
OutputVector &  global_vector 
) const
inline

This function does what the two distribute_local_to_global functions with vector and matrix argument do, but all at once.

Definition at line 3527 of file dof_accessor.templates.h.

template<class DH , bool lda>
void DoFCellAccessor< DH, lda >::get_active_or_mg_dof_indices ( std::vector< types::global_dof_index > &  dof_indices) const
inline

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.

Definition at line 3318 of file dof_accessor.templates.h.

template<class DH , bool lda>
void DoFCellAccessor< DH, lda >::get_dof_indices ( std::vector< types::global_dof_index > &  dof_indices) const
inline

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).

Parameters
[out]dof_indicesThe 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().

Note
In many places in the tutorial and elsewhere in the library, the argument to this function is called 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.
Deprecated:
Currently, this function can also be called for non-active cells, if all degrees of freedom of the FiniteElement are located in vertices. This functionality will vanish in a future release.

Definition at line 3282 of file dof_accessor.templates.h.

template<class DH , bool lda>
void DoFCellAccessor< DH, lda >::get_mg_dof_indices ( std::vector< types::global_dof_index > &  dof_indices) const
inline
Deprecated:
Use get_active_or_mg_dof_indices() with level_cell_iterator returned from begin_mg().

Retrieve the global indices of the degrees of freedom on this cell in the level vector associated to the level of the cell.

Definition at line 3300 of file dof_accessor.templates.h.

template<class DH , bool lda>
const FiniteElement< DH::dimension, DH::space_dimension > & DoFCellAccessor< DH, lda >::get_fe ( ) const
inline

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.

Definition at line 3433 of file dof_accessor.templates.h.

template<class DH , bool lda>
unsigned int DoFCellAccessor< DH, lda >::active_fe_index ( ) const
inline

Returns the index inside the hp::FECollection of the FiniteElement used for this cell.

Definition at line 3443 of file dof_accessor.templates.h.

template<class DH , bool lda>
void DoFCellAccessor< DH, lda >::set_active_fe_index ( const unsigned int  i)
inline

Sets the index of the FiniteElement used for this cell.

Definition at line 3453 of file dof_accessor.templates.h.

template<class DH, bool level_dof_access>
void DoFCellAccessor< DH, level_dof_access >::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.

template<class DH , bool lda>
void DoFCellAccessor< DH, lda >::set_mg_dof_indices ( const std::vector< types::global_dof_index > &  dof_indices)
inline

Set the Level DoF indices of this cell to the given values.

Definition at line 3309 of file dof_accessor.templates.h.

template<class DH, bool level_dof_access>
void DoFCellAccessor< DH, level_dof_access >::update_cell_dof_indices_cache ( ) const

Update the cache in which we store the dof indices of this cell.

template<class DH, bool level_dof_access>
DoFCellAccessor<DH, level_dof_access>& DoFCellAccessor< DH, level_dof_access >::operator= ( const DoFCellAccessor< DH, level_dof_access > &  da)
private

Copy operator. This is normally 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 triangulations. We declare this function here private, thus it may not be used from outside. Furthermore it is not implemented and will give a linker error if used anyway.

Friends And Related Function Documentation

template<class DH, bool level_dof_access>
template<int dim, int spacedim>
friend class DoFHandler
friend

Make the DoFHandler class a friend so that it can call the update_cell_dof_indices_cache() function

Definition at line 2197 of file dof_accessor.h.

Member Data Documentation

template<class DH, bool level_dof_access>
const unsigned int DoFCellAccessor< DH, level_dof_access >::dim = DH::dimension
static

Extract dimension from DH.

Definition at line 1511 of file dof_accessor.h.

template<class DH, bool level_dof_access>
const unsigned int DoFCellAccessor< DH, level_dof_access >::spacedim = DH::space_dimension
static

Extract space dimension from DH.

Definition at line 1516 of file dof_accessor.h.


The documentation for this class was generated from the following files: