Class xlifepp::NedelecEdgeFirstHexahedronPk#
-
class NedelecEdgeFirstHexahedronPk : public xlifepp::NedelecEdgeHexahedron#
-
Inheritence diagram for xlifepp::NedelecEdgeFirstHexahedronPk:
Collaboration diagram for xlifepp::NedelecEdgeFirstHexahedronPk:
Nedelec edge first family of any order k on tetrahedron t (NCE1k) space Vk: Q_(k-1,k,k) x Q_(k,k-1,k)x Q_(k,k,k-1), Q_(k,l,m) monomials of degree at most k in x1, of degree at most l in x2, of degree at most m in x3, 3k(k+1)^2 dofs edge dofs: v-> int_e v.t q, q in P_(k-1)[e] k dofs by edge e face dofs (k>1) : v-> int_f (v x n).q, q in Q_(k-2,k-1) x Q_(k-1,k-2) 2k(k-1) dofs by face hexahedron dofs (k>1) : v-> int_t v.q, q in Q_(k-1,k-2,k-2) x Q_(k-2,k-1,k-2)x Q_(k-2,k-2,k-1) 3k(k-1)^2 dofs.
Public Functions
-
NedelecEdgeFirstHexahedronPk(const Interpolation *int_p)#
-
Nedelec edge first family of order k on hexahedron (NE1k) Degrees Of Freedom of reference element, defined as moment D.o.Fs space Vk: Q_(k-1,k,k) x Q_(k,k-1,k)x Q_(k,k,k-1), Q_(k,l,m) monomials of degree at most k in x1, at most l in x2, at most m in x3, 3k(k+1)^2 dofs edge dofs: u-> int_e u.t q, q in P_(k-1)[e] k dofs by edge e face dofs (k>1) : u-> int_f (u x n).q, q in Q_(k-2,k-1) x Q_(k-1,k-2) 2k(k-1) dofs by face tetrahedron dofs (k>1): u-> int_t u.q, q in Q_(k-1,k-2,k-2) x Q_(k-2,k-1,k-2)x Q_(k-2,k-2,k-1) 3k(k-1)^2 dofs numbering of D.oFs are given by the basis functions order of dual spaces defining D.o.Fs beginning by edges, then faces and hexahedron.
NOTE: moment D.o.Fs can be reinterpreted as punctual D.o.Fs using a well adapted quadrature formulae we do not use this representation here
-
virtual void computeShapeFunctions()#
-
compute shape functions as polynomials
compute shape functions as polynomials using general algorithm edge dofs: u-> int_e u.t q, q in P_(k-1)[e] k dofs by edge e face dofs (k>1) : u-> int_f (u x n).q, q in Q_(k-2,k-1) x Q_(k-1,k-2) 2k(k-1) dofs by face hexahedron dofs (k>1) : u-> int_t u.q, q in Q_(k-1,k-2,k-2) x Q_(k-2,k-1,k-2)x Q_(k-2,k-2,k-1) 3k(k-1)^2 dofs
-
virtual void computeShapeValues(std::vector<real_t>::const_iterator it_pt, ShapeValues &shv, const bool der1 = true, const bool der2 = false) const#
-
compute shape values at a point
-
virtual Value evalEltDof(const FeDof &dof, const Function &f, const Function &gradf, const Function &grad2f) const#
-
specific eval elt dof
-
virtual Value evalFaceDof(const FeDof &dof, const GeomElement &selt, const Function &f, const Function &gradf, const Function &grad2f) const#
-
specific eval face dof function
-
virtual void interpolationData()#
-
defines reference element interpolation data
interp defines Reference Element interpolation data
-
void pointCoordinates()#
-
builds virtual coordinates of moment dofs
-
virtual void rotateDofs(const std::vector<number_t>&, ShapeValues &shv, const bool der1 = true, const bool der2 = false) const#
-
rotation of face dofs to ensure matching of dofs on shared face
-
virtual number_t sideDofsMap(const number_t &n, const number_t &i, const number_t &j, const number_t &k = 0) const#
-
face dofs mapping
-
virtual void sideNumbering()#
-
local numbering on faces (side)
-
virtual number_t sideofsideDofsMap(const number_t &n, const number_t &i, const number_t &j) const#
-
edge dof mapping
-
void sideOfSideNumbering()#
-
local numbering on edges (side of side)
-
NedelecEdgeFirstHexahedronPk(const Interpolation *int_p)#