Reference documentation for deal.II version 8.5.1
|
#include <deal.II/fe/fe_q_dg0.h>
Public Member Functions | |
FE_Q_DG0 (const unsigned int p) | |
FE_Q_DG0 (const Quadrature< 1 > &points) | |
virtual std::string | get_name () const |
virtual void | convert_generalized_support_point_values_to_nodal_values (const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const |
virtual void | interpolate (std::vector< double > &local_dofs, const std::vector< double > &values) const 1 |
virtual void | interpolate (std::vector< double > &local_dofs, const std::vector< Vector< double > > &values, const unsigned int offset=0) const 1 |
virtual void | interpolate (std::vector< double > &local_dofs, const VectorSlice< const std::vector< std::vector< double > > > &values) const 1 |
virtual void | get_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const |
virtual bool | has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const |
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > | get_constant_modes () const |
Public Member Functions inherited from FE_Q_Base< TensorProductPolynomialsConst< dim >, dim, spacedim > | |
FE_Q_Base (const TensorProductPolynomialsConst< dim > &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags) | |
virtual void | get_face_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const |
virtual void | get_subface_interpolation_matrix (const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const |
virtual const FullMatrix< double > & | get_restriction_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const |
virtual const FullMatrix< double > & | get_prolongation_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const |
virtual unsigned int | face_to_cell_index (const unsigned int face_dof_index, const unsigned int face, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false) const |
virtual bool | hp_constraints_are_implemented () const |
virtual std::vector< std::pair< unsigned int, unsigned int > > | hp_vertex_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const |
virtual std::vector< std::pair< unsigned int, unsigned int > > | hp_line_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const |
virtual std::vector< std::pair< unsigned int, unsigned int > > | hp_quad_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const |
virtual FiniteElementDomination::Domination | compare_for_face_domination (const FiniteElement< dim, spacedim > &fe_other) const |
Public Member Functions inherited from FE_Poly< TensorProductPolynomialsConst< dim >, dim, spacedim > | |
FE_Poly (const TensorProductPolynomialsConst< dim > &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components) | |
unsigned int | get_degree () const |
std::vector< unsigned int > | get_poly_space_numbering () const |
std::vector< unsigned int > | get_poly_space_numbering_inverse () const |
virtual double | shape_value (const unsigned int i, const Point< dim > &p) const |
virtual double | shape_value_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const |
virtual Tensor< 1, dim > | shape_grad (const unsigned int i, const Point< dim > &p) const |
virtual Tensor< 1, dim > | shape_grad_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const |
virtual Tensor< 2, dim > | shape_grad_grad (const unsigned int i, const Point< dim > &p) const |
virtual Tensor< 2, dim > | shape_grad_grad_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const |
virtual Tensor< 3, dim > | shape_3rd_derivative (const unsigned int i, const Point< dim > &p) const |
virtual Tensor< 3, dim > | shape_3rd_derivative_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const |
virtual Tensor< 4, dim > | shape_4th_derivative (const unsigned int i, const Point< dim > &p) const |
virtual Tensor< 4, dim > | shape_4th_derivative_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const |
Public Member Functions inherited from FiniteElement< dim, spacedim > | |
FiniteElement (const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components) | |
virtual | ~FiniteElement () |
const FiniteElement< dim, spacedim > & | operator[] (const unsigned int fe_index) const |
bool | operator== (const FiniteElement< dim, spacedim > &) const |
virtual std::size_t | memory_consumption () const |
bool | prolongation_is_implemented () const |
bool | isotropic_prolongation_is_implemented () const |
bool | restriction_is_implemented () const |
bool | isotropic_restriction_is_implemented () const |
bool | restriction_is_additive (const unsigned int index) const |
const FullMatrix< double > & | constraints (const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const |
bool | constraints_are_implemented (const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const |
std::pair< unsigned int, unsigned int > | system_to_component_index (const unsigned int index) const |
unsigned int | component_to_system_index (const unsigned int component, const unsigned int index) const |
std::pair< unsigned int, unsigned int > | face_system_to_component_index (const unsigned int index) const |
unsigned int | adjust_quad_dof_index_for_face_orientation (const unsigned int index, const bool face_orientation, const bool face_flip, const bool face_rotation) const |
unsigned int | adjust_line_dof_index_for_line_orientation (const unsigned int index, const bool line_orientation) const |
const ComponentMask & | get_nonzero_components (const unsigned int i) const |
unsigned int | n_nonzero_components (const unsigned int i) const |
bool | is_primitive () const |
bool | is_primitive (const unsigned int i) const |
unsigned int | n_base_elements () const |
virtual const FiniteElement< dim, spacedim > & | base_element (const unsigned int index) const |
unsigned int | element_multiplicity (const unsigned int index) const |
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > | system_to_base_index (const unsigned int index) const |
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > | face_system_to_base_index (const unsigned int index) const |
types::global_dof_index | first_block_of_base (const unsigned int b) const |
std::pair< unsigned int, unsigned int > | component_to_base_index (const unsigned int component) const |
std::pair< unsigned int, unsigned int > | block_to_base_index (const unsigned int block) const |
std::pair< unsigned int, types::global_dof_index > | system_to_block_index (const unsigned int component) const |
unsigned int | component_to_block_index (const unsigned int component) const |
ComponentMask | component_mask (const FEValuesExtractors::Scalar &scalar) const |
ComponentMask | component_mask (const FEValuesExtractors::Vector &vector) const |
ComponentMask | component_mask (const FEValuesExtractors::SymmetricTensor< 2 > &sym_tensor) const |
ComponentMask | component_mask (const BlockMask &block_mask) const |
BlockMask | block_mask (const FEValuesExtractors::Scalar &scalar) const |
BlockMask | block_mask (const FEValuesExtractors::Vector &vector) const |
BlockMask | block_mask (const FEValuesExtractors::SymmetricTensor< 2 > &sym_tensor) const |
BlockMask | block_mask (const ComponentMask &component_mask) const |
const std::vector< Point< dim > > & | get_unit_support_points () const |
bool | has_support_points () const |
virtual Point< dim > | unit_support_point (const unsigned int index) const |
const std::vector< Point< dim-1 > > & | get_unit_face_support_points () const |
bool | has_face_support_points () const |
virtual Point< dim-1 > | unit_face_support_point (const unsigned int index) const |
const std::vector< Point< dim > > & | get_generalized_support_points () const |
bool | has_generalized_support_points () const |
const std::vector< Point< dim-1 > > & | get_generalized_face_support_points () const |
bool | has_generalized_face_support_points () const |
GeometryPrimitive | get_associated_geometry_primitive (const unsigned int cell_dof_index) const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) |
void | subscribe (const char *identifier=0) const |
void | unsubscribe (const char *identifier=0) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Public Member Functions inherited from FiniteElementData< dim > | |
FiniteElementData (const std::vector< unsigned int > &dofs_per_object, const unsigned int n_components, const unsigned int degree, const Conformity conformity=unknown, const BlockIndices &block_indices=BlockIndices()) | |
unsigned int | n_dofs_per_vertex () const |
unsigned int | n_dofs_per_line () const |
unsigned int | n_dofs_per_quad () const |
unsigned int | n_dofs_per_hex () const |
unsigned int | n_dofs_per_face () const |
unsigned int | n_dofs_per_cell () const |
template<int structdim> | |
unsigned int | n_dofs_per_object () const |
unsigned int | n_components () const |
unsigned int | n_blocks () const |
const BlockIndices & | block_indices () const |
unsigned int | tensor_degree () const |
bool | conforms (const Conformity) const |
bool | operator== (const FiniteElementData &) const |
Protected Member Functions | |
virtual FiniteElement< dim, spacedim > * | clone () const |
Protected Member Functions inherited from FE_Q_Base< TensorProductPolynomialsConst< dim >, dim, spacedim > | |
void | initialize (const std::vector< Point< 1 > > &support_points_1d) |
void | initialize_constraints (const std::vector< Point< 1 > > &points) |
void | initialize_unit_support_points (const std::vector< Point< 1 > > &points) |
void | initialize_unit_face_support_points (const std::vector< Point< 1 > > &points) |
void | initialize_quad_dof_index_permutation () |
Protected Member Functions inherited from FE_Poly< TensorProductPolynomialsConst< dim >, dim, spacedim > | |
void | correct_third_derivatives (internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data, const internal::FEValues::MappingRelatedData< dim, spacedim > &mapping_data, const unsigned int n_q_points, const unsigned int dof) const |
Protected Member Functions inherited from FiniteElement< dim, spacedim > | |
void | reinit_restriction_and_prolongation_matrices (const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false) |
TableIndices< 2 > | interface_constraints_size () const |
virtual InternalDataBase * | get_face_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const |
virtual InternalDataBase * | get_subface_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const |
virtual void | fill_fe_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValues::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const =0 |
virtual void | fill_fe_face_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim-1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValues::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const =0 |
virtual void | fill_fe_subface_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim-1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValues::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const =0 |
Static Private Member Functions | |
static std::vector< bool > | get_riaf_vector (const unsigned int degree) |
static std::vector< unsigned int > | get_dpo_vector (const unsigned int degree) |
Additional Inherited Members | |
Public Types inherited from FiniteElementData< dim > | |
enum | Conformity { unknown = 0x00, L2 = 0x01, Hcurl = 0x02, Hdiv = 0x04, H1 = Hcurl | Hdiv, H2 = 0x0e } |
Static Public Member Functions inherited from FE_Q_Base< TensorProductPolynomialsConst< dim >, dim, spacedim > | |
static ::ExceptionBase & | ExcFEQCannotHaveDegree0 () |
Static Public Member Functions inherited from FiniteElement< dim, spacedim > | |
static ::ExceptionBase & | ExcShapeFunctionNotPrimitive (int arg1) |
static ::ExceptionBase & | ExcFENotPrimitive () |
static ::ExceptionBase & | ExcUnitShapeValuesDoNotExist () |
static ::ExceptionBase & | ExcFEHasNoSupportPoints () |
static ::ExceptionBase & | ExcEmbeddingVoid () |
static ::ExceptionBase & | ExcProjectionVoid () |
static ::ExceptionBase & | ExcWrongInterfaceMatrixSize (int arg1, int arg2) |
static ::ExceptionBase & | ExcInterpolationNotImplemented () |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, char *arg2, std::string &arg3) |
static ::ExceptionBase & | ExcNoSubscriber (char *arg1, char *arg2) |
Public Attributes inherited from FiniteElementData< dim > | |
const unsigned int | dofs_per_vertex |
const unsigned int | dofs_per_line |
const unsigned int | dofs_per_quad |
const unsigned int | dofs_per_hex |
const unsigned int | first_line_index |
const unsigned int | first_quad_index |
const unsigned int | first_hex_index |
const unsigned int | first_face_line_index |
const unsigned int | first_face_quad_index |
const unsigned int | dofs_per_face |
const unsigned int | dofs_per_cell |
const unsigned int | components |
const unsigned int | degree |
const Conformity | conforming_space |
const BlockIndices | block_indices_data |
Static Public Attributes inherited from FiniteElement< dim, spacedim > | |
static const unsigned int | space_dimension = spacedim |
Static Public Attributes inherited from FiniteElementData< dim > | |
static const unsigned int | dimension = dim |
Static Protected Member Functions inherited from FE_Q_Base< TensorProductPolynomialsConst< dim >, dim, spacedim > | |
static std::vector< unsigned int > | get_dpo_vector (const unsigned int degree) |
Static Protected Member Functions inherited from FiniteElement< dim, spacedim > | |
static std::vector< unsigned int > | compute_n_nonzero_components (const std::vector< ComponentMask > &nonzero_components) |
Protected Attributes inherited from FE_Poly< TensorProductPolynomialsConst< dim >, dim, spacedim > | |
TensorProductPolynomialsConst< dim > | poly_space |
Protected Attributes inherited from FiniteElement< dim, spacedim > | |
std::vector< std::vector< FullMatrix< double > > > | restriction |
std::vector< std::vector< FullMatrix< double > > > | prolongation |
FullMatrix< double > | interface_constraints |
std::vector< Point< dim > > | unit_support_points |
std::vector< Point< dim-1 > > | unit_face_support_points |
std::vector< Point< dim > > | generalized_support_points |
std::vector< Point< dim-1 > > | generalized_face_support_points |
Table< 2, int > | adjust_quad_dof_index_for_face_orientation_table |
std::vector< int > | adjust_line_dof_index_for_line_orientation_table |
std::vector< std::pair< unsigned int, unsigned int > > | system_to_component_table |
std::vector< std::pair< unsigned int, unsigned int > > | face_system_to_component_table |
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > | system_to_base_table |
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > | face_system_to_base_table |
BlockIndices | base_to_block_indices |
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > | component_to_base_table |
const std::vector< bool > | restriction_is_additive_flags |
const std::vector< ComponentMask > | nonzero_components |
const std::vector< unsigned int > | n_nonzero_components_table |
const bool | cached_primitivity |
Implementation of a scalar Lagrange finite element Qp+DG0
that yields the finite element space of continuous, piecewise polynomials of degree p
in each coordinate direction plus the space of locally constant functions. This class is realized using tensor product polynomials based on equidistant or given support points.
The standard constructor of this class takes the degree p
of this finite element. Alternatively, it can take a quadrature formula points
defining the support points of the Lagrange interpolation in one coordinate direction.
For more information about the spacedim
template parameter check the documentation of FiniteElement or the one of Triangulation.
For more information regarding this element see: Boffi, D., et al. "Local Mass Conservation of Stokes Finite Elements." Journal of Scientific Computing (2012): 1-18.
The constructor creates a TensorProductPolynomials object that includes the tensor product of LagrangeEquidistant
polynomials of degree p
plus the locally constant function. This TensorProductPolynomialsConst
object provides all values and derivatives of the shape functions. In case a quadrature rule is given, the constructor creates a TensorProductPolynomialsConst object that includes the tensor product of Lagrange
polynomials with the support points from points
and a locally constant function.
Furthermore the constructor fills the interface_constrains
, the prolongation
(embedding) and the restriction
matrices.
The original ordering of the shape functions represented by the TensorProductPolynomialsConst is a tensor product numbering. However, the shape functions on a cell are renumbered beginning with the shape functions whose support points are at the vertices, then on the line, on the quads, and finally (for 3d) on the hexes. Finally there is a support point for the discontinuous shape function in the middle of the cell. To be explicit, these numberings are listed in the following:
1D case:
* 0---2---1 *
2D case:
* 2-------3 * | | * | 5 | * | | * 0-------1 *
3D case:
* 6-------7 6-------7 * /| | / /| * / | | / / | * / | | / / | * 4 | 8 | 4-------5 | * | 2-------3 | | 3 * | / / | | / * | / / | | / * |/ / | |/ * 0-------1 0-------1 *
The respective coordinate values of the support points of the degrees of freedom are as follows:
[ 0, 0, 0]
; [ 1, 0, 0]
; [ 0, 1, 0]
; [ 1, 1, 0]
; [ 0, 0, 1]
; [ 1, 0, 1]
; [ 0, 1, 1]
; [ 1, 1, 1]
; [1/2, 1/2, 1/2]
; 1D case:
* 0---2---1 *
Index 3 has the same coordinates as index 2
2D case:
* 2---7---3 * | | * 4 8 5 * | | * 0---6---1 *
Index 9 has the same coordinates as index 2
3D case:
* 6--15---7 6--15---7 * /| | / /| * 12 | 19 12 1319 * / 18 | / / | * 4 | | 4---14--5 | * | 2---11--3 | | 3 * | / / | 17 / * 16 8 9 16 | 9 * |/ / | |/ * 0---10--1 0---8---1 * * *-------* *-------* * /| | / /| * / | 23 | / 25 / | * / | | / / | * * | | *-------* | * |20 *-------* | |21 * * | / / | 22 | / * | / 24 / | | / * |/ / | |/ * *-------* *-------* *
The center vertices have number 26 and 27.
The respective coordinate values of the support points of the degrees of freedom are as follows:
[0, 0, 0]
; [1, 0, 0]
; [0, 1, 0]
; [1, 1, 0]
; [0, 0, 1]
; [1, 0, 1]
; [0, 1, 1]
; [1, 1, 1]
; [0, 1/2, 0]
; [1, 1/2, 0]
; [1/2, 0, 0]
; [1/2, 1, 0]
; [0, 1/2, 1]
; [1, 1/2, 1]
; [1/2, 0, 1]
; [1/2, 1, 1]
; [0, 0, 1/2]
; [1, 0, 1/2]
; [0, 1, 1/2]
; [1, 1, 1/2]
; [0, 1/2, 1/2]
; [1, 1/2, 1/2]
; [1/2, 0, 1/2]
; [1/2, 1, 1/2]
; [1/2, 1/2, 0]
; [1/2, 1/2, 1]
; [1/2, 1/2, 1/2]
; [1/2, 1/2, 1/2]
; 1D case:
* 0--2-4-3--1 *
* 2--10-11-3 * | | * 5 14 15 7 * | 16 | * 4 12 13 6 * | | * 0--8--9--1 *
1D case:
* 0--2--3--4--1 *
Index 5 has the same coordinates as index 3
* 2--13-14-15-3 * | | * 6 22 23 24 9 * | | * 5 19 20 21 8 * | | * 4 16 17 18 7 * | | * 0--10-11-12-1 *Index 21 has the same coordinates as index 20
Definition at line 237 of file fe_q_dg0.h.
Constructor for tensor product polynomials of degree p
plus locally constant functions.
Definition at line 33 of file fe_q_dg0.cc.
FE_Q_DG0< dim, spacedim >::FE_Q_DG0 | ( | const Quadrature< 1 > & | points | ) |
Constructor for tensor product polynomials with support points points
plus locally constant functions based on a one-dimensional quadrature formula. The degree of the finite element is points.size()-1
. Note that the first point has to be 0 and the last one 1.
Definition at line 59 of file fe_q_dg0.cc.
|
virtual |
Return a string that uniquely identifies a finite element. This class returns FE_Q_DG0<dim>(degree)
, with dim
and degree
replaced by appropriate values.
Implements FiniteElement< dim, spacedim >.
Definition at line 89 of file fe_q_dg0.cc.
|
virtual |
Given the values of a function \(f(\mathbf x)\) at the (generalized) support points, this function then computes what the nodal values of the element are, i.e., \(\Psi_i[f]\), where \(\Psi_i\) are the node functionals of the element. The values \(\Psi_i[f]\) are then the expansion coefficients for the shape functions of the finite element function that interpolates the given function \(f(x)\), i.e., \( f_h(\mathbf x) = \sum_i \Psi_i[f] \varphi_i(\mathbf x) \) is the finite element interpolant of \(f\) with the current element. The operation described here is used, for example, in the FETools::compute_node_matrix() function.
In more detail, let us assume that the generalized support points (see this glossary entry ) of the current element are \(\hat{\mathbf x}_i\) and that the node functionals associated with the current element are \(\Psi_i[\cdot]\). Then, the fact that the element is based on generalized support points, implies that if we apply \(\Psi_i\) to a (possibly vector-valued) finite element function \(\varphi\), the result must have the form \(\Psi_i[\varphi] = f_i(\varphi(\hat{\mathbf x}_i))\) – in other words, the value of the node functional \(\Psi_i\) applied to \(\varphi\) only depends on the values of \(\varphi\) at \(\hat{\mathbf x}_i\) and not on values anywhere else, or integrals of \(\varphi\), or any other kind of information.
The exact form of \(f_i\) depends on the element. For example, for scalar Lagrange elements, we have that in fact \(\Psi_i[\varphi] = \varphi(\hat{\mathbf x}_i)\). If you combine multiple scalar Lagrange elements via an FESystem object, then \(\Psi_i[\varphi] = \varphi(\hat{\mathbf x}_i)_{c(i)}\) where \(c(i)\) is the result of the FiniteElement::system_to_component_index() function's return value's first component. In these two cases, \(f_i\) is therefore simply the identity (in the scalar case) or a function that selects a particular vector component of its argument. On the other hand, for Raviart-Thomas elements, one would have that \(f_i(\mathbf y) = \mathbf y \cdot \mathbf n_i\) where \(\mathbf n_i\) is the normal vector of the face at which the shape function is defined.
Given all of this, what this function does is the following: If you input a list of values of a function \(\varphi\) at all generalized support points (where each value is in fact a vector of values with as many components as the element has), then this function returns a vector of values obtained by applying the node functionals to these values. In other words, if you pass in \(\{\varphi(\hat{\mathbf x}_i)\}_{i=0}^{N-1}\) then you will get out a vector \(\{\Psi[\varphi]\}_{i=0}^{N-1}\) where \(N\) equals dofs_per_cell
.
[in] | support_point_values | An array of size dofs_per_cell (which equals the number of points the get_generalized_support_points() function will return) where each element is a vector with as many entries as the element has vector components. This array should contain the values of a function at the generalized support points of the current element. |
[out] | nodal_values | An array of size dofs_per_cell that contains the node functionals of the element applied to the given function. |
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 180 of file fe_q_dg0.cc.
|
virtual |
Interpolate a set of scalar values, computed in the generalized support points.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 206 of file fe_q_dg0.cc.
|
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, spacedim >.
Definition at line 227 of file fe_q_dg0.cc.
|
virtual |
Interpolate a set of vector values, computed in the generalized support points.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 254 of file fe_q_dg0.cc.
|
virtual |
Return the matrix interpolating from the given finite element to the present one. The size of the matrix is then dofs_per_cell
times source.dofs_per_cell
.
These matrices are only available if the source element is also a FE_Q_DG0
element. Otherwise, an exception of type FiniteElement<dim,spacedim>::ExcInterpolationNotImplemented is thrown.
Reimplemented from FE_Q_Base< TensorProductPolynomialsConst< dim >, dim, spacedim >.
Definition at line 282 of file fe_q_dg0.cc.
|
virtual |
This function returns true
, if the shape function shape_index
has non-zero function values somewhere on the face face_index
.
Reimplemented from FE_Q_Base< TensorProductPolynomialsConst< dim >, dim, spacedim >.
Definition at line 335 of file fe_q_dg0.cc.
|
virtual |
Return a list of constant modes of the element. For this element, there are two constant modes despite the element is scalar: The first constant mode is all ones for the usual FE_Q basis and the second one only using the discontinuous part.
Reimplemented from FE_Q_Base< TensorProductPolynomialsConst< dim >, dim, spacedim >.
Definition at line 349 of file fe_q_dg0.cc.
|
protectedvirtual |
clone
function instead of a copy constructor.
This function is needed by the constructors of FESystem
.
Implements FiniteElement< dim, spacedim >.
Definition at line 170 of file fe_q_dg0.cc.
|
staticprivate |
Return the restriction_is_additive flags. Only the last component is true.
Definition at line 310 of file fe_q_dg0.cc.
|
staticprivate |
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
.
Definition at line 321 of file fe_q_dg0.cc.