Reference documentation for deal.II version 8.1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
FE_PolyTensor< POLY, dim, spacedim > Class Template Reference

#include <fe_poly_tensor.h>

Inheritance diagram for FE_PolyTensor< POLY, dim, spacedim >:
[legend]

Classes

class  InternalData
 

Public Member Functions

 FE_PolyTensor (const unsigned int degree, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
 
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 UpdateFlags update_once (const UpdateFlags flags) const
 
virtual UpdateFlags update_each (const UpdateFlags flags) 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 ()
 
virtual std::string get_name () const =0
 
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
 
 DeclException1 (ExcShapeFunctionNotPrimitive, int,<< "The shape function with index "<< arg1<< " is not primitive, i.e. it is vector-valued and "<< "has more than one non-zero vector component. This "<< "function cannot be called for these shape functions. "<< "Maybe you want to use the same function with the "<< "_component suffix?")
 
 DeclException0 (ExcFENotPrimitive)
 
 DeclException0 (ExcUnitShapeValuesDoNotExist)
 
 DeclException0 (ExcFEHasNoSupportPoints)
 
 DeclException0 (ExcEmbeddingVoid)
 
 DeclException0 (ExcProjectionVoid)
 
 DeclException0 (ExcConstraintsVoid)
 
 DeclException2 (ExcWrongInterfaceMatrixSize, int, int,<< "The interface matrix has a size of "<< arg1<< "x"<< arg2<< ", which is not reasonable in the present dimension.")
 
 DeclException2 (ExcComponentIndexInvalid, int, int,<< "The component-index pair ("<< arg1<< ", "<< arg2<< ") is invalid, i.e. non-existent.")
 
 DeclException0 (ExcInterpolationNotImplemented)
 
 DeclException0 (ExcBoundaryFaceUsed)
 
 DeclException0 (ExcJacobiDeterminantHasWrongSign)
 
virtual bool has_support_on_face (const unsigned int shape_index, const unsigned int face_index) 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
 
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
 
virtual bool hp_constraints_are_implemented () const
 
virtual void get_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
 
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 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
 
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
 
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
 
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 ComponentMaskget_nonzero_components (const unsigned int i) const
 
unsigned int n_nonzero_components (const unsigned int i) 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
 
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
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
void subscribe (const char *identifier=0) const
 
void unsubscribe (const char *identifier=0) const
 
unsigned int n_subscriptions () const
 
void list_subscribers () const
 
 DeclException3 (ExcInUse, int, char *, std::string &,<< "Object of class "<< arg2<< " is still used by "<< arg1<< " other objects.\n"<< "(Additional information: "<< arg3<< ")\n"<< "Note the entry in the Frequently Asked Questions of "<< "deal.II (linked to from http://www.dealii.org/) for "<< "more information on what this error means.")
 
 DeclException2 (ExcNoSubscriber, char *, char *,<< "No subscriber with identifier \""<< arg2<< "\" did subscribe to this object of class "<< arg1)
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
- Public Member Functions inherited from FiniteElementData< dim >
 FiniteElementData ()
 
 FiniteElementData (const std::vector< unsigned int > &dofs_per_object, const unsigned int n_components, const unsigned int degree, const Conformity conformity=unknown, const unsigned int n_blocks=numbers::invalid_unsigned_int)
 
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 BlockIndicesblock_indices () const
 
bool is_primitive () const
 
unsigned int tensor_degree () const
 
bool conforms (const Conformity) const
 
bool operator== (const FiniteElementData &) const
 

Protected Member Functions

virtual Mapping< dim, spacedim >
::InternalDataBase
get_data (const UpdateFlags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature) const
 
virtual void fill_fe_values (const Mapping< dim, spacedim > &mapping, const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Quadrature< dim > &quadrature, typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, typename Mapping< dim, spacedim >::InternalDataBase &fe_internal, FEValuesData< dim, spacedim > &data, CellSimilarity::Similarity &cell_similarity) const
 
virtual void fill_fe_face_values (const Mapping< dim, spacedim > &mapping, const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim-1 > &quadrature, typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, typename Mapping< dim, spacedim >::InternalDataBase &fe_internal, FEValuesData< dim, spacedim > &data) const
 
virtual void fill_fe_subface_values (const Mapping< dim, spacedim > &mapping, const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim-1 > &quadrature, typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, typename Mapping< dim, spacedim >::InternalDataBase &fe_internal, FEValuesData< dim, spacedim > &data) 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
 
void compute_2nd (const Mapping< dim, spacedim > &mapping, const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int offset, typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, InternalDataBase &fe_internal, FEValuesData< dim, spacedim > &data) const
 
virtual FiniteElement< dim,
spacedim > * 
clone () const =0
 
- Protected Member Functions inherited from FiniteElementData< dim >
void set_primitivity (const bool value)
 

Protected Attributes

MappingType mapping_type
 
POLY poly_space
 
FullMatrix< doubleinverse_node_matrix
 
Point< dim > cached_point
 
std::vector< Tensor< 1, dim > > cached_values
 
std::vector< Tensor< 2, dim > > cached_grads
 
std::vector< Tensor< 3, dim > > cached_grad_grads
 
- Protected Attributes inherited from FiniteElement< dim, spacedim >
std::vector< std::vector
< FullMatrix< double > > > 
restriction
 
std::vector< std::vector
< FullMatrix< double > > > 
prolongation
 
FullMatrix< doubleinterface_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, intadjust_quad_dof_index_for_face_orientation_table
 
std::vector< intadjust_line_dof_index_for_line_orientation_table
 

Additional Inherited Members

- Public Types inherited from FiniteElementData< dim >
enum  Conformity {
  unknown = 0x00, L2 = 0x01, Hcurl = 0x02, Hdiv = 0x04,
  H1 = Hcurl | Hdiv, H2 = 0x0e
}
 
- 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
 
BlockIndices block_indices_data
 
- Static Public Attributes inherited from FiniteElementData< dim >
static const unsigned int dimension = dim
 
- Static Protected Member Functions inherited from FiniteElement< dim, spacedim >
static std::vector< unsigned intcompute_n_nonzero_components (const std::vector< ComponentMask > &nonzero_components)
 

Detailed Description

template<class POLY, int dim, int spacedim = dim>
class FE_PolyTensor< POLY, dim, spacedim >

This class gives a unified framework for the implementation of FiniteElement classes based on Tensor valued polynomial spaces like PolynomialsBDM and PolynomialsRaviartThomas.

Every class that implements following function can be used as template parameter POLY.

void compute (const Point<dim> &unit_point,
std::vector<Tensor<1,dim> > &values,
std::vector<Tensor<2,dim> > &grads,
std::vector<Tensor<3,dim> > &grad_grads) const;

In many cases, the node functionals depend on the shape of the mesh cell, since they evaluate normal or tangential components on the faces. In order to allow for a set of transformations, the variable mapping_type has been introduced. It should also be set in the constructor of a derived class.

This class is not a fully implemented FiniteElement class, but implements some common features of vector valued elements based on vector valued polynomial classes. What's missing here in particular is information on the topological location of the node values.

For more information on the template parameter spacedim, see the documentation for the class Triangulation.

Deriving classes

Any derived class must decide on the polynomial space to use. This polynomial space should be implemented simply as a set of vector valued polynomials like PolynomialsBDM and PolynomialsRaviartThomas. In order to facilitate this implementation, the basis of this space may be arbitrary.

Determining the correct basis

In most cases, the set of desired node values $N_i$ and the basis functions $v_j$ will not fulfill the interpolation condition $N_i(v_j) = \delta_{ij}$.

The use of the member data inverse_node_matrix allows to compute the basis $v_j$ automatically, after the node values for each original basis function of the polynomial space have been computed.

Therefore, the constructor of a derived class should have a structure like this (example for interpolation in support points):

fill_support_points();
const unsigned int n_dofs = this->dofs_per_cell;
FullMatrix<double> N(n_dofs, n_dofs);
for (unsigned int i=0;i<n_dofs;++i)
{
const Point<dim>& p = this->unit_support_point[i];
for (unsigned int j=0;j<n_dofs;++j)
for (unsigned int d=0;d<dim;++d)
N(i,j) += node_vector[i][d]
* this->shape_value_component(j, p, d);
}
this->inverse_node_matrix.reinit(n_dofs, n_dofs);
Note
The matrix inverse_node_matrix should have dimensions zero before this piece of code is executed. Only then, shape_value_component() will return the raw polynomial j as defined in the polynomial space POLY.

Setting the transformation

In most cases, vector valued basis functions must be transformed when mapped from the reference cell to the actual grid cell. These transformations can be selected from the set MappingType and stored in mapping_type. Therefore, each constructor should contain a line like:

this->mapping_type = this->mapping_none;
See also
PolynomialsBDM, PolynomialsRaviartThomas
Author
Guido Kanschat
Date
2005

Definition at line 36 of file fe.h.

Constructor & Destructor Documentation

template<class POLY, int dim, int spacedim = dim>
FE_PolyTensor< POLY, dim, spacedim >::FE_PolyTensor ( const unsigned int  degree,
const FiniteElementData< dim > &  fe_data,
const std::vector< bool > &  restriction_is_additive_flags,
const std::vector< ComponentMask > &  nonzero_components 
)

Constructor.

  • degree: constructor argument for poly. May be different from fe_data.degree.

Member Function Documentation

template<class POLY, int dim, int spacedim = dim>
virtual double FE_PolyTensor< POLY, dim, spacedim >::shape_value ( const unsigned int  i,
const Point< dim > &  p 
) const
virtual

Since these elements are vector valued, an exception is thrown.

Reimplemented from FiniteElement< dim, spacedim >.

template<class POLY, int dim, int spacedim = dim>
virtual double FE_PolyTensor< POLY, dim, spacedim >::shape_value_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const
virtual

Just like for shape_value(), but this function will be called when the shape function has more than one non-zero vector component. In that case, this function should return the value of the component-th vector component of the ith shape function at point p.

Reimplemented from FiniteElement< dim, spacedim >.

template<class POLY, int dim, int spacedim = dim>
virtual Tensor<1,dim> FE_PolyTensor< POLY, dim, spacedim >::shape_grad ( const unsigned int  i,
const Point< dim > &  p 
) const
virtual

Since these elements are vector valued, an exception is thrown.

Reimplemented from FiniteElement< dim, spacedim >.

template<class POLY, int dim, int spacedim = dim>
virtual Tensor<1,dim> FE_PolyTensor< POLY, dim, spacedim >::shape_grad_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const
virtual

Just like for shape_grad(), but this function will be called when the shape function has more than one non-zero vector component. In that case, this function should return the gradient of the component-th vector component of the ith shape function at point p.

Reimplemented from FiniteElement< dim, spacedim >.

template<class POLY, int dim, int spacedim = dim>
virtual Tensor<2,dim> FE_PolyTensor< POLY, dim, spacedim >::shape_grad_grad ( const unsigned int  i,
const Point< dim > &  p 
) const
virtual

Since these elements are vector valued, an exception is thrown.

Reimplemented from FiniteElement< dim, spacedim >.

template<class POLY, int dim, int spacedim = dim>
virtual Tensor<2,dim> FE_PolyTensor< POLY, dim, spacedim >::shape_grad_grad_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const
virtual

Just like for shape_grad_grad(), but this function will be called when the shape function has more than one non-zero vector component. In that case, this function should return the gradient of the component-th vector component of the ith shape function at point p.

Reimplemented from FiniteElement< dim, spacedim >.

template<class POLY, int dim, int spacedim = dim>
virtual UpdateFlags FE_PolyTensor< POLY, dim, spacedim >::update_once ( const UpdateFlags  flags) const
virtual

Given flags, determines the values which must be computed only for the reference cell. Make sure, that mapping_type is set by the derived class, such that this function can operate correctly.

Implements FiniteElement< dim, spacedim >.

Reimplemented in FE_ABF< dim >.

template<class POLY, int dim, int spacedim = dim>
virtual UpdateFlags FE_PolyTensor< POLY, dim, spacedim >::update_each ( const UpdateFlags  flags) const
virtual

Given flags, determines the values which must be computed in each cell cell. Make sure, that mapping_type is set by the derived class, such that this function can operate correctly.

Implements FiniteElement< dim, spacedim >.

Reimplemented in FE_ABF< dim >.

template<class POLY, int dim, int spacedim = dim>
virtual Mapping<dim,spacedim>::InternalDataBase* FE_PolyTensor< POLY, dim, spacedim >::get_data ( const UpdateFlags  flags,
const Mapping< dim, spacedim > &  mapping,
const Quadrature< dim > &  quadrature 
) const
protectedvirtual

Prepare internal data structures and fill in values independent of the cell. Returns a pointer to an object of which the caller of this function then has to assume ownership (which includes destruction when it is no more needed).

Implements FiniteElement< dim, spacedim >.

template<class POLY, int dim, int spacedim = dim>
virtual void FE_PolyTensor< POLY, dim, spacedim >::fill_fe_values ( const Mapping< dim, spacedim > &  mapping,
const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const Quadrature< dim > &  quadrature,
typename Mapping< dim, spacedim >::InternalDataBase mapping_internal,
typename Mapping< dim, spacedim >::InternalDataBase fe_internal,
FEValuesData< dim, spacedim > &  data,
CellSimilarity::Similarity &  cell_similarity 
) const
protectedvirtual

Fill the fields of FEValues. This function performs all the operations needed to compute the data of an FEValues object.

The same function in mapping must have been called for the same cell first!

Implements FiniteElement< dim, spacedim >.

template<class POLY, int dim, int spacedim = dim>
virtual void FE_PolyTensor< POLY, dim, spacedim >::fill_fe_face_values ( const Mapping< dim, spacedim > &  mapping,
const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const unsigned int  face_no,
const Quadrature< dim-1 > &  quadrature,
typename Mapping< dim, spacedim >::InternalDataBase mapping_internal,
typename Mapping< dim, spacedim >::InternalDataBase fe_internal,
FEValuesData< dim, spacedim > &  data 
) const
protectedvirtual

Fill the fields of FEFaceValues. This function performs all the operations needed to compute the data of an FEFaceValues object.

The same function in mapping must have been called for the same cell first!

Implements FiniteElement< dim, spacedim >.

template<class POLY, int dim, int spacedim = dim>
virtual void FE_PolyTensor< POLY, dim, spacedim >::fill_fe_subface_values ( const Mapping< dim, spacedim > &  mapping,
const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const unsigned int  face_no,
const unsigned int  sub_no,
const Quadrature< dim-1 > &  quadrature,
typename Mapping< dim, spacedim >::InternalDataBase mapping_internal,
typename Mapping< dim, spacedim >::InternalDataBase fe_internal,
FEValuesData< dim, spacedim > &  data 
) const
protectedvirtual

Fill the fields of FESubfaceValues. This function performs all the operations needed to compute the data of an FESubfaceValues object.

The same function in mapping must have been called for the same cell first!

Implements FiniteElement< dim, spacedim >.

Member Data Documentation

template<class POLY, int dim, int spacedim = dim>
MappingType FE_PolyTensor< POLY, dim, spacedim >::mapping_type
protected

The mapping type to be used to map shape functions from the reference cell to the mesh cell.

Definition at line 202 of file fe_poly_tensor.h.

template<class POLY, int dim, int spacedim = dim>
POLY FE_PolyTensor< POLY, dim, spacedim >::poly_space
protected

The polynomial space. Its type is given by the template parameter POLY.

Definition at line 284 of file fe_poly_tensor.h.

template<class POLY, int dim, int spacedim = dim>
FullMatrix<double> FE_PolyTensor< POLY, dim, spacedim >::inverse_node_matrix
protected

The inverse of the matrix aij of node values Ni applied to polynomial pj. This matrix is used to convert polynomials in the "raw" basis provided in poly_space to the basis dual to the node functionals on the reference cell.

This object is not filled by FE_PolyTensor, but is a chance for a derived class to allow for reorganization of the basis functions. If it is left empty, the basis in poly_space is used.

Definition at line 306 of file fe_poly_tensor.h.

template<class POLY, int dim, int spacedim = dim>
Point<dim> FE_PolyTensor< POLY, dim, spacedim >::cached_point
mutableprotected

If a shape function is computed at a single point, we must compute all of them to apply inverse_node_matrix. In order to avoid too much overhead, we cache the point and the function values for the next evaluation.

Definition at line 318 of file fe_poly_tensor.h.

template<class POLY, int dim, int spacedim = dim>
std::vector<Tensor<1,dim> > FE_PolyTensor< POLY, dim, spacedim >::cached_values
mutableprotected

Cached shape function values after call to shape_value_component().

Definition at line 325 of file fe_poly_tensor.h.

template<class POLY, int dim, int spacedim = dim>
std::vector<Tensor<2,dim> > FE_PolyTensor< POLY, dim, spacedim >::cached_grads
mutableprotected

Cached shape function gradients after call to shape_grad_component().

Definition at line 332 of file fe_poly_tensor.h.

template<class POLY, int dim, int spacedim = dim>
std::vector<Tensor<3,dim> > FE_PolyTensor< POLY, dim, spacedim >::cached_grad_grads
mutableprotected

Cached second derivatives of shape functions after call to shape_grad_grad_component().

Definition at line 339 of file fe_poly_tensor.h.


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