21 #ifndef ELEMENTQUAD_HPP
22 #define ELEMENTQUAD_HPP
25 #include <Eigen/Dense>
85 inline const boost::multi_array<double, 2> &
get_phiT_quadF(
size_t ilF)
const
87 return m_phiT_quadF[ilF];
93 return m_phiF_quadF[ilF];
99 return m_dphiT_quadF[ilF];
132 boost::multi_array<double, 2> m_phiT_quadT;
133 boost::multi_array<VectorRd, 2> m_dphiT_quadT;
134 boost::multi_array<VectorRd, 2> m_vec_phiT_quadT;
137 std::vector<QuadratureRule> m_quadF;
139 std::vector<boost::multi_array<double, 2>> m_phiT_quadF;
140 std::vector<boost::multi_array<double, 2>> m_phiF_quadF;
141 std::vector<boost::multi_array<VectorRd, 2>> m_dphiT_quadF;
142 std::vector<boost::multi_array<VectorRd, 2>> m_vec_phiT_quadF;
143 std::vector<boost::multi_array<VectorRd, 2>> m_vec_phiF_quadF;
Definition: elementquad.hpp:49
Definition: hybridcore.hpp:179
const boost::multi_array< double, 2 > get_phiF_quadF(size_t ilF) const
Returns values of edge basis functions at edge quadrature nodes, for edge with local number ilF.
Definition: elementquad.hpp:91
const QuadratureRule & get_quadF(size_t ilF) const
Returns quadrature rules on edge with local number ilF.
Definition: elementquad.hpp:67
ElementQuad(const HybridCore &hho, const size_t iT, const size_t doeT, const size_t doeF)
Class constructor: loads the quadrature rules and values of basis functions/gradients at these points...
Definition: elementquad.cpp:17
const QuadratureRule & get_quadT() const
Returns quadrature rules in cell.
Definition: elementquad.hpp:61
boost::multi_array< VectorRd, 2 > get_vec_phiT_quadT(size_t degree) const
Builds on the fly the values of vector cell basis functions at cell quadrature nodes....
Definition: elementquad.cpp:62
const boost::multi_array< VectorRd, 2 > & get_dphiT_quadT() const
Returns values of gradients of cell basis functions at cell quadrature nodes.
Definition: elementquad.hpp:79
boost::multi_array< VectorRd, 2 > get_vec_phiF_quadF(size_t ilF, size_t degree) const
Builds on the fly the values of vector edge basis functions at edge quadrature nodes....
Definition: elementquad.cpp:110
const boost::multi_array< double, 2 > & get_phiT_quadT() const
Returns values of cell basis functions at cell quadrature nodes.
Definition: elementquad.hpp:73
const boost::multi_array< VectorRd, 2 > & get_dphiT_quadF(size_t ilF) const
Returns values of gradients of cell basis functions at edge quadrature nodes, for edge with local num...
Definition: elementquad.hpp:97
const boost::multi_array< double, 2 > & get_phiT_quadF(size_t ilF) const
Returns values of cell basis functions at edge quadrature nodes, for edge with local number ilF.
Definition: elementquad.hpp:85
boost::multi_array< VectorRd, 2 > get_vec_phiT_quadF(size_t ilF, size_t degree) const
Builds on the fly the values of vector cell basis functions at edge quadrature nodes....
Definition: elementquad.cpp:86
Polytope< DIMENSION > Cell
Definition: Polytope2D.hpp:151
std::vector< QuadratureNode > QuadratureRule
Definition: quadraturerule.hpp:55
Definition: ddr-klplate.hpp:27