Class xlifepp::NedelecEdgeFirstHexahedronPk#

class NedelecEdgeFirstHexahedronPk : public xlifepp::NedelecEdgeHexahedron#

Inheritence diagram for xlifepp::NedelecEdgeFirstHexahedronPk:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "1" [label="xlifepp::NedelecEdgeFirstHexahedronPk" tooltip="xlifepp::NedelecEdgeFirstHexahedronPk" fillcolor="#BFBFBF"] "2" [label="xlifepp::NedelecEdgeHexahedron" tooltip="xlifepp::NedelecEdgeHexahedron"] "4" [label="xlifepp::RefElement" tooltip="xlifepp::RefElement"] "3" [label="xlifepp::RefHexahedron" tooltip="xlifepp::RefHexahedron"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "2" -> "3" [dir=forward tooltip="public-inheritance"] "3" -> "4" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for xlifepp::NedelecEdgeFirstHexahedronPk:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "17" [label="std::list< std::vector< PolynomialT< real_t > > >" tooltip="std::list< std::vector< PolynomialT< real_t > > >"] "20" [label="std::map< xlifepp::Quadrature *, std::vector< xlifepp::ShapeValues > >" tooltip="std::map< xlifepp::Quadrature *, std::vector< xlifepp::ShapeValues > >"] "18" [label="std::list< T >" tooltip="std::list< T >"] "21" [label="std::map< K, T >" tooltip="std::map< K, T >"] "7" [label="std::vector< T >" tooltip="std::vector< T >"] "8" [label="std::vector< ShapeType >" tooltip="std::vector< ShapeType >"] "6" [label="std::vector< real_t >" tooltip="std::vector< real_t >"] "10" [label="std::vector< std::vector< int_t > >" tooltip="std::vector< std::vector< int_t > >"] "9" [label="std::vector< std::vector< number_t > >" tooltip="std::vector< std::vector< number_t > >"] "11" [label="std::vector< xlifepp::GeomRefElement * >" tooltip="std::vector< xlifepp::GeomRefElement * >"] "13" [label="std::vector< xlifepp::Interpolation * >" tooltip="std::vector< xlifepp::Interpolation * >"] "19" [label="std::vector< xlifepp::PolynomialsBasisT >" tooltip="std::vector< xlifepp::PolynomialsBasisT >"] "14" [label="std::vector< xlifepp::RefDof * >" tooltip="std::vector< xlifepp::RefDof * >"] "15" [label="std::vector< xlifepp::RefElement * >" tooltip="std::vector< xlifepp::RefElement * >"] "5" [label="xlifepp::GeomRefElement" tooltip="xlifepp::GeomRefElement"] "12" [label="xlifepp::Interpolation" tooltip="xlifepp::Interpolation"] "1" [label="xlifepp::NedelecEdgeFirstHexahedronPk" tooltip="xlifepp::NedelecEdgeFirstHexahedronPk" fillcolor="#BFBFBF"] "2" [label="xlifepp::NedelecEdgeHexahedron" tooltip="xlifepp::NedelecEdgeHexahedron"] "16" [label="xlifepp::PolynomialsBasisT< K >" tooltip="xlifepp::PolynomialsBasisT< K >"] "4" [label="xlifepp::RefElement" tooltip="xlifepp::RefElement"] "3" [label="xlifepp::RefHexahedron" tooltip="xlifepp::RefHexahedron"] "17" -> "18" [dir=forward tooltip="template-instance"] "20" -> "21" [dir=forward tooltip="template-instance"] "8" -> "7" [dir=forward tooltip="template-instance"] "6" -> "7" [dir=forward tooltip="template-instance"] "10" -> "7" [dir=forward tooltip="template-instance"] "9" -> "7" [dir=forward tooltip="template-instance"] "11" -> "7" [dir=forward tooltip="template-instance"] "13" -> "7" [dir=forward tooltip="template-instance"] "19" -> "7" [dir=forward tooltip="template-instance"] "14" -> "7" [dir=forward tooltip="template-instance"] "15" -> "7" [dir=forward tooltip="template-instance"] "5" -> "6" [dir=forward tooltip="usage"] "5" -> "8" [dir=forward tooltip="usage"] "5" -> "9" [dir=forward tooltip="usage"] "5" -> "10" [dir=forward tooltip="usage"] "5" -> "11" [dir=forward tooltip="usage"] "12" -> "13" [dir=forward tooltip="usage"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "2" -> "3" [dir=forward tooltip="public-inheritance"] "16" -> "17" [dir=forward tooltip="public-inheritance"] "4" -> "5" [dir=forward tooltip="usage"] "4" -> "12" [dir=forward tooltip="usage"] "4" -> "14" [dir=forward tooltip="usage"] "4" -> "15" [dir=forward tooltip="usage"] "4" -> "9" [dir=forward tooltip="usage"] "4" -> "16" [dir=forward tooltip="usage"] "4" -> "19" [dir=forward tooltip="usage"] "4" -> "20" [dir=forward tooltip="usage"] "3" -> "4" [dir=forward tooltip="public-inheritance"] }

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)