
Classes | |
| class | InternalData |
Public Member Functions | |
| FE_Nedelec (const unsigned int p) | |
| virtual std::string | get_name () const |
| virtual bool | has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const |
| virtual bool | hp_constraints_are_implemented () const |
| virtual FiniteElementDomination::Domination | compare_for_face_domination (const FiniteElement< dim > &fe_other) const |
| virtual std::vector< std::pair < unsigned int, unsigned int > > | hp_vertex_dof_identities (const FiniteElement< dim > &fe_other) const |
| virtual std::vector< std::pair < unsigned int, unsigned int > > | hp_line_dof_identities (const FiniteElement< dim > &fe_other) const |
| virtual std::vector< std::pair < unsigned int, unsigned int > > | hp_quad_dof_identities (const FiniteElement< dim > &fe_other) const |
| virtual void | get_face_interpolation_matrix (const FiniteElement< dim > &source, FullMatrix< double > &matrix) const |
| virtual void | get_subface_interpolation_matrix (const FiniteElement< dim > &source, const unsigned int subface, FullMatrix< double > &matrix) const |
| virtual void | interpolate (std::vector< double > &local_dofs, const std::vector< double > &values) const |
| virtual void | interpolate (std::vector< double > &local_dofs, const std::vector< Vector< double > > &values, unsigned int offset=0) const |
| virtual void | interpolate (std::vector< double > &local_dofs, const VectorSlice< const std::vector< std::vector< double > > > &values) const |
| virtual std::size_t | memory_consumption () const |
| virtual FiniteElement< dim > * | clone () const |
Private Member Functions | |
| void | initialize_support_points (const unsigned int degree) |
| void | initialize_restriction () |
Static Private Member Functions | |
| static std::vector< unsigned int > | get_dpo_vector (const unsigned int degree, bool dg=false) |
Private Attributes | |
| Table< 2, double > | boundary_weights |
Friends | |
| class | FE_Nedelec |
Implementation of Nédélec elements, conforming with the space Hcurl. These elements generate vector fields with tangential components continuous between mesh cells.
We follow the convention that the degree of Nédélec elements denotes the polynomial degree of the largest complete polynomial subspace contained in the Nédélec space. This leads to the consistently numbered sequence of spaces
Consequently, approximation order of the Nedelec space equals the value degree given to the constructor. In this scheme, the lowest order element would be created by the call FE_Nedelec<dim>(0). Note that this follows the convention of Brezzi and Raviart, though not the one used in the original paper by Nedelec.
This class is not implemented for the codimension one case (spacedim != dim).
In some sense, the implementation of this element is not complete, but you will rarely notice. Here is the fact: since the element is vector-valued already on the unit cell, the Jacobian matrix (or its inverse) is needed already to generate the values of the shape functions on the cells in real space. This is in contrast to most other elements, where you only need the Jacobian for the gradients. Thus, to generate the gradients of Nedelec shape functions, one would need to have the derivatives of the inverse of the Jacobian matrix.
Basically, the Nedelec shape functions can be understood as the gradients of scalar shape functions on the real cell. They are thus the inverse Jacobian matrix times the gradients of scalar shape functions on the unit cell. The gradient of Nedelec shape functions is then, by the product rule, the sum of first the derivative (with respect to true coordinates) of the inverse Jacobian times the gradient (in unit coordinates) of the scalar shape function, plus second the inverse Jacobian times the derivative (in true coordinates) of the gradient (in unit coordinates) of the scalar shape functions. Note that each of the derivatives in true coordinates can be expressed as inverse Jacobian times gradient in unit coordinates.
The problem is the derivative of the inverse Jacobian. This rank-3 tensor can actually be computed (and we did so in very early versions of the library), but is a large task and very time consuming, so we dropped it. Since it is not available, we simply drop this first term.
What this means for the present case: first the computation of gradients of Nedelec shape functions is wrong in general. Second, in the following two cases you will not notice this:
The interpolation operators associated with the Nédélec element are constructed such that interpolation and computing the curl are commuting operations on rectangular mesh cells. We require this from interpolating arbitrary functions as well as the restriction matrices.
The node values for an element of degree k on the reference cell are:
dim-1 dimensional FE_Nedelec polynomials of degree k-1. The node values above rely on integrals, which will be computed by quadrature rules themselves. The generalized support points are a set of points such that this quadrature can be performed with sufficient accuracy. The points needed are thode of QGaussk+1 on each edge and QGaussk+2 on each face and in the interior of the cell (or none for N1).
Definition at line 154 of file fe_nedelec.h.
| FE_Nedelec< dim >::FE_Nedelec | ( | const unsigned int | p ) |
Constructor for the Nédélec element of degree p.
| virtual std::string FE_Nedelec< dim >::get_name | ( | ) | const [virtual] |
Return a string that uniquely identifies a finite element. This class returns FE_Nedelec<dim>(degree), with dim and degree replaced by appropriate values.
Implements FiniteElement< dim, dim >.
| virtual bool FE_Nedelec< dim >::has_support_on_face | ( | const unsigned int | shape_index, |
| const unsigned int | face_index | ||
| ) | const [virtual] |
Check whether a shape function may be non-zero on a face.
Reimplemented from FiniteElement< dim, dim >.
| virtual bool FE_Nedelec< dim >::hp_constraints_are_implemented | ( | ) | const [virtual] |
Return whether this element implements its hanging node constraints in the new way, which has to be used to make elements "hp compatible".
For the FE_Nedelec class the result is always true (independent of the degree of the element), as it implements the complete set of functions necessary for hp capability.
Reimplemented from FiniteElement< dim, dim >.
| virtual FiniteElementDomination::Domination FE_Nedelec< dim >::compare_for_face_domination | ( | const FiniteElement< dim > & | fe_other ) | const [virtual] |
Return whether this element dominates the one, which is given as argument.
| virtual std::vector<std::pair<unsigned int, unsigned int> > FE_Nedelec< dim >::hp_vertex_dof_identities | ( | const FiniteElement< dim > & | fe_other ) | const [virtual] |
If, on a vertex, several finite elements are active, the hp code first assigns the degrees of freedom of each of these FEs different global indices. It then calls this function to find out which of them should get identical values, and consequently can receive the same global DoF index. This function therefore returns a list of identities between DoFs of the present finite element object with the DoFs of fe_other, which is a reference to a finite element object representing one of the other finite elements active on this particular vertex. The function computes which of the degrees of freedom of the two finite element objects are equivalent, and returns a list of pairs of global dof indices in identities. The first index of each pair denotes one of the vertex dofs of the present element, whereas the second is the corresponding index of the other finite element.
| virtual std::vector<std::pair<unsigned int, unsigned int> > FE_Nedelec< dim >::hp_line_dof_identities | ( | const FiniteElement< dim > & | fe_other ) | const [virtual] |
Same as hp_vertex_dof_indices(), except that the function treats degrees of freedom on lines.
| virtual std::vector<std::pair<unsigned int, unsigned int> > FE_Nedelec< dim >::hp_quad_dof_identities | ( | const FiniteElement< dim > & | fe_other ) | const [virtual] |
Same as hp_vertex_dof_indices(), except that the function treats degrees of freedom on lines.
| virtual void FE_Nedelec< dim >::get_face_interpolation_matrix | ( | const FiniteElement< dim > & | source, |
| FullMatrix< double > & | matrix | ||
| ) | const [virtual] |
Return the matrix interpolating from a face of one element to the face of the neighboring element. The size of the matrix is then source.dofs_per_face times this->dofs_per_face.
Derived elements will have to implement this function. They may only provide interpolation matrices for certain source finite elements, for example those from the same family. If they don't implement interpolation from a given element, then they must throw an exception of type FiniteElement<dim>::ExcInterpolationNotImplemented.
| virtual void FE_Nedelec< dim >::get_subface_interpolation_matrix | ( | const FiniteElement< dim > & | source, |
| const unsigned int | subface, | ||
| FullMatrix< double > & | matrix | ||
| ) | const [virtual] |
Return the matrix interpolating from a face of one element to the subface of the neighboring element. The size of the matrix is then source.dofs_per_face times this->dofs_per_face.
Derived elements will have to implement this function. They may only provide interpolation matrices for certain source finite elements, for example those from the same family. If they don't implement interpolation from a given element, then they must throw an exception of type ExcInterpolationNotImplemented.
| virtual void FE_Nedelec< dim >::interpolate | ( | std::vector< double > & | local_dofs, |
| const std::vector< double > & | values | ||
| ) | const [virtual] |
Interpolate a set of scalar values, computed in the generalized support points.
Reimplemented from FiniteElement< dim, dim >.
| virtual void FE_Nedelec< dim >::interpolate | ( | std::vector< double > & | local_dofs, |
| const std::vector< Vector< double > > & | values, | ||
| unsigned int | offset = 0 |
||
| ) | const [virtual] |
Interpolate a set of vector values, computed in the generalized support points.
Since a finite element often only interpolates part of a vector, offset is used to determine the first component of the vector to be interpolated. Maybe consider changing your data structures to use the next function.
Reimplemented from FiniteElement< dim, dim >.
| virtual void FE_Nedelec< dim >::interpolate | ( | std::vector< double > & | local_dofs, |
| const VectorSlice< const std::vector< std::vector< double > > > & | values | ||
| ) | const [virtual] |
Interpolate a set of vector values, computed in the generalized support points.
Reimplemented from FiniteElement< dim, dim >.
| virtual std::size_t FE_Nedelec< dim >::memory_consumption | ( | ) | const [virtual] |
Determine an estimate for the memory consumption (in bytes) of this object.
This function is made virtual, since finite element objects are usually accessed through pointers to their base class, rather than the class itself.
Reimplemented from FiniteElement< dim, dim >.
| virtual FiniteElement<dim>* FE_Nedelec< dim >::clone | ( | ) | const [virtual] |
A sort of virtual copy constructor. Some places in the library, for example the constructors of FESystem as well as the hp::FECollection class, need to make copied of finite elements without knowing their exact type. They do so through this function.
Implements FiniteElement< dim, dim >.
| static std::vector<unsigned int> FE_Nedelec< dim >::get_dpo_vector | ( | const unsigned int | degree, |
| bool | dg = false |
||
| ) | [static, private] |
Only for internal use. Its full name is get_dofs_per_object_vector function and it creates the dofs_per_object vector that is needed within the constructor to be passed to the constructor of FiniteElementData.
If the optional argument dg is true, the vector returned will have all degrees of freedom assigned to the cell, none on the faces and edges.
| void FE_Nedelec< dim >::initialize_support_points | ( | const unsigned int | degree ) | [private] |
Initialize the generalized_support_points field of the FiniteElement class and fill the tables with interpolation weights (boundary_weights and interior_weights). Called from the constructor.
| void FE_Nedelec< dim >::initialize_restriction | ( | ) | [private] |
Initialize the interpolation from functions on refined mesh cells onto the father cell. According to the philosophy of the Nédélec element, this restriction operator preserves the curl of a function weakly.
friend class FE_Nedelec [friend] |
Allow access from other dimensions.
Definition at line 402 of file fe_nedelec.h.
Table<2, double> FE_Nedelec< dim >::boundary_weights [private] |
These are the factors multiplied to a function in the generalized_face_support_points when computing the integration.
Definition at line 396 of file fe_nedelec.h.
documentation generated on Fri Feb 3 2012 06:04:07 by
doxygen
1.7.2