DUNE PDELab (2.8)
#include <dune/localfunctions/utility/polynomialbasis.hh>
Public Member Functions | |
void | evaluateFunction (const typename Traits::DomainType &x, std::vector< typename Traits::RangeType > &out) const |
Evaluate all shape functions. | |
void | evaluateJacobian (const typename Traits::DomainType &x, std::vector< typename Traits::JacobianType > &out) const |
Evaluate Jacobian of all shape functions. | |
void | evaluateHessian (const typename Traits::DomainType &x, std::vector< HessianType > &out) const |
Evaluate Jacobian of all shape functions. | |
void | partial (const std::array< unsigned int, dimension > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const |
Evaluate partial derivatives of all shape functions. | |
Detailed Description
class Dune::PolynomialBasis< Eval, CM, D, R >
This is the basis class for a ''polynomial'' basis, i.e., a basis consisting of linear combiniations of a underlying second basis set. Examples are standard polynomials where the underlying basis is given by the MonomialBasis class. The basis evaluation is given by the matrix vector multiplication between the coefficient matrix and the vector filled by evaluating the underlying basis set. This class is constructed using a reference of the underlying basis and the coefficient matrix. A specialization holding an instance of the coefficient matrix is provided by the class template< class Eval, class CM = SparseCoeffMatrix<typename Eval::Field,Eval::dimRange> > class PolynomialBasisWithMatrix;
- Template Parameters
-
B Basis set with static const int dimension -> dimension of reference element typedef DomainVector -> coordinates in reference element int size(int order) const -> number of basis functions void evaluate( order, x, val ) const int order DomainVector x Container val CM stroage for coefficience with typedef Field -> field of coefficience static const int dimRange -> coeficience are of type FieldMatrix<Field,dimRange,dimRange> void mult( val, y ) Container val std::vector<RangeVector> y Container access to basis functions through forward iterator typedef value_type typedef const_iterator const_iterator begin()
The documentation for this class was generated from the following file:
- dune/localfunctions/utility/polynomialbasis.hh