28#include <boost/multi_array.hpp>
79 template<
typename GeometricSupport>
88 static std::vector<VectorZd>
complete(
size_t degree)
90 std::vector<VectorZd> powers;
92 for (
size_t l = 0; l <= degree; l++)
94 for (
size_t i = 0;
i <= l;
i++)
116 Eigen::Matrix<double, N, N> diff_x,
117 Eigen::Matrix<double, N, N> diff_y
149 Eigen::Matrix<double, N, N>
dx;
150 Eigen::Matrix<double, N, N>
dy;
197 return ( m_degree >= 0 ? (m_degree + 1) * (m_degree + 2) / 2 : 0);
228 return (x - m_xT) / m_hT;
234 std::vector<VectorZd> m_powers;
235 Eigen::Matrix2d m_rot;
287 inline double _coordinate_transform(
const VectorRd &x)
const
289 return (x - m_xE).dot(m_tE) / m_hE;
311 template <
typename BasisType>
328 static const bool hasCurl = BasisType::hasCurl;
337 const BasisType &basis,
338 const Eigen::MatrixXd &
matrix
343 assert((
size_t)
matrix.cols() == basis.dimension() ||
"Inconsistent family initialization");
349 return m_matrix.rows();
355 static_assert(
hasFunction,
"Call to function() not available");
358 for (
auto j = 1;
j < m_matrix.cols();
j++)
360 f += m_matrix(
i,
j) * m_basis.function(
j, x);
368 static_assert(
hasFunction,
"Call to function() not available");
371 for (
auto j = 1;
j < m_matrix.cols();
j++)
373 f += m_matrix(
i,
j) * ancestor_value_quad[
j][iqn];
382 static_assert(
hasGradient,
"Call to gradient() not available");
385 for (
auto j = 1;
j < m_matrix.cols();
j++)
387 G += m_matrix(
i,
j) * m_basis.gradient(
j, x);
395 static_assert(
hasGradient,
"Call to gradient() not available");
397 GradientValue G = m_matrix(
i, 0) * ancestor_gradient_quad[0][iqn];
398 for (
auto j = 1;
j < m_matrix.cols();
j++)
400 G += m_matrix(
i,
j) * ancestor_gradient_quad[
j][iqn];
408 static_assert(
hasCurl,
"Call to curl() not available");
410 CurlValue C = m_matrix(
i, 0) * m_basis.curl(0, x);
411 for (
auto j = 1;
j < m_matrix.cols();
j++)
413 C += m_matrix(
i,
j) * m_basis.curl(
j, x);
419 CurlValue curl(
size_t i,
size_t iqn,
const boost::multi_array<CurlValue, 2> &ancestor_curl_quad)
const
421 static_assert(
hasCurl,
"Call to curl() not available");
423 CurlValue C = m_matrix(
i, 0) * ancestor_curl_quad[0][iqn];
424 for (
auto j = 1;
j < m_matrix.cols();
j++)
426 C += m_matrix(
i,
j) * ancestor_curl_quad[
j][iqn];
434 static_assert(
hasDivergence,
"Call to divergence() not available");
437 for (
auto j = 1;
j < m_matrix.cols();
j++)
439 D += m_matrix(
i,
j) * m_basis.divergence(
j, x);
447 static_assert(
hasDivergence,
"Call to divergence() not available");
450 for (
auto j = 1;
j < m_matrix.cols();
j++)
452 D += m_matrix(
i,
j) * ancestor_divergence_quad[
j][iqn];
460 static_assert(
hasRotor,
"Call to rotor() not available");
462 RotorValue R = m_matrix(
i, 0) * m_basis.rotor(0, x);
463 for (
auto j = 1;
j < m_matrix.cols();
j++)
465 R += m_matrix(
i,
j) * m_basis.rotor(
j, x);
471 RotorValue rotor(
size_t i,
size_t iqn,
const boost::multi_array<RotorValue, 2> &ancestor_rotor_quad)
const
473 static_assert(
hasRotor,
"Call to rotor() not available");
475 RotorValue R = m_matrix(
i, 0) * ancestor_rotor_quad[0][iqn];
476 for (
auto j = 1;
j < m_matrix.cols();
j++)
478 R += m_matrix(
i,
j) * ancestor_rotor_quad[
j][iqn];
486 static_assert(
hasHessian,
"Call to hessian() not available");
489 for (
auto j = 1;
j < m_matrix.cols();
j++)
491 H += m_matrix(
i,
j) * m_basis.hessian(
j, x);
497 HessianValue hessian(
size_t i,
size_t iqn,
const boost::multi_array<HessianValue, 2> &ancestor_hessian_quad)
const
499 static_assert(
hasHessian,
"Call to hessian() not available");
501 HessianValue H = m_matrix(
i, 0) * ancestor_hessian_quad[0][iqn];
502 for (
auto j = 1;
j < m_matrix.cols();
j++)
504 H += m_matrix(
i,
j) * ancestor_hessian_quad[
j][iqn];
510 inline const Eigen::MatrixXd &
matrix()
const
524 return m_basis.max_degree();
529 Eigen::MatrixXd m_matrix;
559 template <
typename ScalarFamilyType,
size_t N>
587 : m_scalar_family(scalar_family)
589 static_assert(ScalarFamilyType::tensorRank ==
Scalar,
590 "Vector family can only be constructed from scalar families");
591 m_rot.row(0) << 0., 1.;
592 m_rot.row(1) << -1., 0.;
599 return m_scalar_family.dimension() * N;
605 static_assert(
hasFunction,
"Call to function() not available");
608 ek(
i / m_scalar_family.dimension()) = 1.;
609 return ek * m_scalar_family.function(
i % m_scalar_family.dimension(), x);
615 static_assert(
hasFunction,
"Call to function() not available");
618 ek(
i / m_scalar_family.dimension()) = 1.;
619 return ek * ancestor_value_quad[
i % m_scalar_family.dimension()][iqn];
625 static_assert(
hasGradient,
"Call to gradient() not available");
627 GradientValue G = Eigen::Matrix<double, N, dimspace>::Zero();
628 G.row(
i / m_scalar_family.dimension()) = m_scalar_family.gradient(
i % m_scalar_family.dimension(), x);
635 static_assert(
hasGradient,
"Call to gradient() not available");
637 GradientValue G = Eigen::Matrix<double, N, dimspace>::Zero();
638 G.row(
i / m_scalar_family.dimension()) = ancestor_gradient_quad[
i % m_scalar_family.dimension()][iqn];
645 static_assert(
hasCurl,
"Call to curl() not available");
651 CurlValue curl(
size_t i,
size_t iqn,
const boost::multi_array<VectorRd, 2> &ancestor_gradient_quad)
const
653 static_assert(
hasCurl,
"Call to curl() not available");
655 return m_rot * ancestor_gradient_quad[
i][iqn];
661 static_assert(
hasDivergence,
"Call to divergence() not available");
663 return m_scalar_family.gradient(
i % m_scalar_family.dimension(), x)(
i / m_scalar_family.dimension());
669 static_assert(
hasDivergence,
"Call to divergence() not available");
671 return ancestor_gradient_quad[
i % m_scalar_family.dimension()][iqn](
i / m_scalar_family.dimension());
679 return curl(
i, x).trace();
683 RotorValue rotor(
size_t i,
size_t iqn,
const boost::multi_array<VectorRd, 2> &ancestor_gradient_quad)
const
685 static_assert(
hasRotor,
"Call to rotor() not available");
687 return curl(
i, iqn, ancestor_gradient_quad).trace();
691 constexpr inline const ScalarFamilyType &
ancestor()
const
693 return m_scalar_family;
699 return m_scalar_family.max_degree();
703 ScalarFamilyType m_scalar_family;
704 Eigen::Matrix2d m_rot;
719 template<
typename ScalarFamilyType,
size_t N>
744 : m_scalar_family(scalar_family),
748 static_assert(ScalarFamilyType::tensorRank ==
Scalar,
749 "Vector family can only be constructed from scalar families");
752 for (
size_t j = 0;
j < N;
j++){
753 for (
size_t i = 0;
i < N;
i++){
754 m_E[
j*N +
i] = Eigen::Matrix<double, N, N>::Zero();
755 m_E[
j*N +
i](
i,
j) = 1.;
762 size_t r = m_scalar_family.dimension();
763 size_t j = ( ( (
i/r)%N )*N + (
i/r)/N ) * r + (
i%r);
764 m_transposeOperator(
i,
j) = 1.;
771 return m_scalar_family.dimension() * N * N;
777 static_assert(
hasFunction,
"Call to function() not available");
779 return m_scalar_family.function(
i % m_scalar_family.dimension(), x) * m_E[
i / m_scalar_family.dimension()];
785 static_assert(
hasFunction,
"Call to function() not available");
787 return ancestor_value_quad[
i % m_scalar_family.dimension()][iqn] * m_E[
i / m_scalar_family.dimension()];
793 static_assert(
hasDivergence,
"Call to divergence() not available");
795 Eigen::Matrix<double, N, 1> V = m_E[
i/m_scalar_family.dimension()].col( (
i/m_scalar_family.dimension()) /N);
796 return m_scalar_family.gradient(
i % m_scalar_family.dimension(), x)( (
i / m_scalar_family.dimension()) / N) * V;
802 static_assert(
hasDivergence,
"Call to divergence() not available");
804 Eigen::Matrix<double, N, 1> V = m_E[
i/m_scalar_family.dimension()].col( (
i/m_scalar_family.dimension()) /N);
805 return ancestor_gradient_quad[
i % m_scalar_family.dimension()][iqn]( (
i / m_scalar_family.dimension()) / N) * V;
811 static_assert(
hasGradient,
"Call to gradient() not available");
813 const VectorRd gradient_scalar_function = m_scalar_family.gradient(
i % m_scalar_family.dimension(), x);
814 const Eigen::Matrix<double, N, N> matrix_basis = m_E[
i / m_scalar_family.dimension()];
815 return GradientValue(gradient_scalar_function.x() * matrix_basis, gradient_scalar_function.y() * matrix_basis);
821 static_assert(
hasGradient,
"Call to gradient() not available");
823 const VectorRd gradient_scalar_function = ancestor_gradient_quad[
i % m_scalar_family.dimension()][iqn];
824 const Eigen::Matrix<double, N, N> matrix_basis = m_E[
i / m_scalar_family.dimension()];
825 return GradientValue(gradient_scalar_function.x() * matrix_basis, gradient_scalar_function.y() * matrix_basis);
831 return m_scalar_family;
843 return m_transposeOperator;
855 size_t dim_scalar = m_scalar_family.dimension();
858 Eigen::MatrixXd SB = Eigen::MatrixXd::Zero(dim_scalar * N*(N+1)/2,
dimension());
864 for (
size_t i = 0;
i<N;
i++){
866 SB.middleRows(position, (N-
i)*dim_scalar) = symOpe.middleRows(
i*(N+1)*dim_scalar, (N-
i)*dim_scalar);
868 position += (N-
i)*dim_scalar;
880 size_t r = m_scalar_family.dimension();
882 Eigen::MatrixXd Tr = Eigen::MatrixXd::Zero(r,
dimension());
891 if (
size_t(l/r)%N == size_t(
size_t(l/r)/N) ){
901 ScalarFamilyType m_scalar_family;
902 std::vector<Eigen::Matrix<double, N, N>> m_E;
903 Eigen::MatrixXd m_transposeOperator;
910 template <
typename ScalarFamilyType>
936 const ScalarFamilyType &scalar_family,
939 : m_scalar_family(scalar_family),
942 static_assert(ScalarFamilyType::hasFunction,
"Call to function() not available");
943 static_assert(std::is_same<typename ScalarFamilyType::GeometricSupport, Edge>::value,
944 "Tangent families can only be defined on edges");
950 return m_scalar_family.dimension();
956 return m_scalar_family.function(
i, x) * m_generator;
960 constexpr inline const ScalarFamilyType &
ancestor()
const
962 return m_scalar_family;
968 return m_scalar_family.max_degree();
978 ScalarFamilyType m_scalar_family;
987 template <
typename BasisType>
1013 const BasisType &basis,
1025 return m_basis.dimension() - m_shift;
1037 return m_basis.max_degree();
1043 static_assert(
hasFunction,
"Call to function() not available");
1045 return m_basis.function(
i + m_shift, x);
1051 static_assert(
hasGradient,
"Call to gradient() not available");
1053 return m_basis.gradient(
i + m_shift, x);
1059 static_assert(
hasCurl,
"Call to curl() not available");
1061 return m_basis.curl(
i + m_shift, x);
1067 static_assert(
hasDivergence,
"Call to divergence() not available");
1069 return m_basis.divergence(
i + m_shift, x);
1075 static_assert(
hasRotor,
"Call to rotor() not available");
1077 return m_basis.rotor(
i + m_shift, x);
1083 static_assert(
hasHessian,
"Call to hessian() not available");
1085 return m_basis.hessian(
i + m_shift, x);
1097 template <
typename BasisType>
1123 const BasisType &basis,
1160 static_assert(
hasFunction,
"Call to function() not available");
1162 return m_basis.function(
i, x);
1168 static_assert(
hasGradient,
"Call to gradient() not available");
1170 return m_basis.gradient(
i, x);
1176 static_assert(
hasCurl,
"Call to curl() not available");
1178 return m_basis.curl(
i, x);
1184 static_assert(
hasDivergence,
"Call to divergence() not available");
1186 return m_basis.divergence(
i, x);
1192 static_assert(
hasRotor,
"Call to rotor() not available");
1194 return m_basis.rotor(
i, x);
1207 template <
typename BasisType>
1235 : m_scalar_basis(basis)
1237 static_assert(BasisType::hasGradient,
1238 "Gradient basis requires gradient() for the original basis to be available");
1245 return m_scalar_basis.dimension();
1251 return m_scalar_basis.gradient(
i, x);
1257 return m_scalar_basis;
1262 BasisType m_scalar_basis;
1269 template <
typename BasisType>
1275 typedef Eigen::Matrix<double, dimspace, dimspace>
CurlValue;
1297 static_assert(BasisType::tensorRank ==
Scalar && std::is_same<typename BasisType::GeometricSupport, Cell>::value,
1298 "Curl basis can only be constructed starting from scalar bases on elements");
1299 static_assert(BasisType::hasCurl,
1300 "Curl basis requires curl() for the original basis to be available");
1306 return m_basis.dimension();
1312 return m_basis.curl(
i, x);
1332 template <
typename BasisType>
1356 : m_vector_basis(basis)
1358 static_assert(BasisType::tensorRank ==
Vector || BasisType::tensorRank ==
Matrix,
1359 "Divergence basis can only be constructed starting from vector bases");
1360 static_assert(BasisType::hasDivergence,
1361 "Divergence basis requires divergence() for the original basis to be available");
1368 return m_vector_basis.dimension();
1374 return m_vector_basis.divergence(
i, x);
1380 return m_vector_basis;
1384 BasisType m_vector_basis;
1392 template <
typename BasisType>
1418 : m_vector_basis(basis)
1420 static_assert(BasisType::tensorRank ==
Vector,
1421 "Rotor basis can only be constructed starting from vector bases");
1422 static_assert(BasisType::hasRotor,
1423 "Rotor basis requires rotor() for the original basis to be available");
1430 return m_vector_basis.dimension();
1436 return m_vector_basis.rotor(
i, x);
1442 return m_vector_basis;
1446 BasisType m_vector_basis;
1455 template <
typename BasisType>
1484 : m_scalar_basis(basis)
1486 static_assert(BasisType::tensorRank ==
Scalar,
1487 "Hessian basis can only be constructed starting from scalar bases");
1488 static_assert(BasisType::hasHessian,
1489 "Hessian basis requires hessian() for the original basis to be available");
1496 return m_scalar_basis.dimension();
1502 return m_scalar_basis.hessian(
i, x);
1508 return m_scalar_basis;
1513 BasisType m_scalar_basis;
1525 typedef Eigen::Matrix<double, dimspace, dimspace>
CurlValue;
1575 return (x - m_xT) / m_hT;
1581 std::vector<VectorZd> m_powers;
1591 typedef Eigen::Matrix<double, dimspace, dimspace>
CurlValue;
1626 inline const std::shared_ptr<RolyComplBasisCell> &
rck()
const
1633 Eigen::Matrix2d m_rot;
1634 std::shared_ptr<RolyComplBasisCell> m_Rck_basis;
1685 return (x - m_xT) / m_hT;
1691 Eigen::Matrix2d m_rot;
1692 TensorizedVectorFamily<MonomialScalarBasisCell, dimspace> m_Pkmo;
1712 template <
typename outValue,
typename inValue,
typename FunctionType>
1714 const boost::multi_array<inValue, 2> &B_quad,
1715 const FunctionType &F
1718 boost::multi_array<outValue, 2> transformed_B_quad( boost::extents[B_quad.shape()[0]][B_quad.shape()[1]] );
1720 std::transform( B_quad.origin(), B_quad.origin() + B_quad.num_elements(), transformed_B_quad.origin(), F);
1722 return transformed_B_quad;
1727 template <
typename ScalarBasisType,
size_t N>
1729 const ScalarBasisType & B,
1730 const std::vector<Eigen::VectorXd> &
v
1734 size_t k =
v.size();
1736 Eigen::MatrixXd M = Eigen::MatrixXd::Zero(r*k, r*N);
1738 for (
size_t i=0;
i < k;
i++){
1740 assert(
v[
i].size() == N);
1743 for (
size_t j=0;
j < N;
j++){
1744 M.block(
i*r,
j*r, r, r) =
v[
i](
j) * Eigen::MatrixXd::Identity(r,r);
1761 inline static std::function<Eigen::MatrixXd(
const Eigen::MatrixXd &)>
1762 symmetrise_matrix = [](
const Eigen::MatrixXd & x)->Eigen::MatrixXd {
return 0.5*(x+x.transpose());};
1765 inline static std::function<Eigen::MatrixXd(
const Eigen::MatrixXd &)>
1770 template <
typename ScalarBasisType,
size_t N>
1772 const ScalarBasisType & B
1776 Eigen::MatrixXd M = Eigen::MatrixXd::Zero(r, r*N*N);
1778 for (
size_t k=0; k < N; k++){
1779 M.block(0, k*r*(N+1), r, r) = Eigen::MatrixXd::Identity(r, r);
1792 template<
typename BasisType, BasisFunctionE BasisFunction>
1797 template<
typename BasisType>
1800 static_assert(BasisType::hasFunction,
"Call to function not available");
1804 return basis.function(
i, x);
1809 const BasisType &basis,
1812 const boost::multi_array<ReturnValue, 2> &ancestor_value_quad
1815 return basis.function(
i, iqn, ancestor_value_quad);
1821 template<
typename BasisType>
1824 static_assert(BasisType::hasGradient,
"Call to gradient not available");
1828 return basis.gradient(
i, x);
1833 const BasisType &basis,
1836 const boost::multi_array<ReturnValue, 2> &ancestor_gradient_quad
1839 return basis.gradient(
i, iqn, ancestor_gradient_quad);
1845 template<
typename BasisType>
1848 static_assert(BasisType::hasGradient,
"Call to gradient not available");
1852 return 0.5 * (basis.gradient(
i, x) + basis.gradient(
i, x).transpose());
1857 const BasisType &basis,
1860 const boost::multi_array<ReturnValue, 2> &ancestor_gradient_quad
1863 return 0.5 * (basis.gradient(
i, iqn, ancestor_gradient_quad) + basis.gradient(
i, iqn, ancestor_gradient_quad).transpose());
1869 template<
typename BasisType>
1872 static_assert(BasisType::hasGradient,
"Call to gradient not available");
1876 return 0.5 * (basis.gradient(
i, x) - basis.gradient(
i, x).transpose());
1881 const BasisType &basis,
1884 const boost::multi_array<ReturnValue, 2> &ancestor_gradient_quad
1887 return 0.5 * (basis.gradient(
i, iqn, ancestor_gradient_quad) - basis.gradient(
i, iqn, ancestor_gradient_quad).transpose());
1892 template<
typename BasisType>
1895 static_assert(BasisType::hasCurl,
"Call to curl not available");
1899 return basis.curl(
i, x);
1904 const BasisType &basis,
1907 const boost::multi_array<ReturnValue, 2> &ancestor_curl_quad
1910 return basis.curl(
i, iqn, ancestor_curl_quad);
1915 template<
typename BasisType>
1918 static_assert(BasisType::hasDivergence,
"Call to divergence not available");
1922 return basis.divergence(
i, x);
1927 const BasisType &basis,
1930 const boost::multi_array<ReturnValue, 2> &ancestor_divergence_quad
1933 return basis.divergence(
i, iqn, ancestor_divergence_quad);
1938 template<
typename BasisType>
1941 static_assert(BasisType::hasRotor,
"Call to rotor not available");
1945 return basis.rotor(
i, x);
1950 const BasisType &basis,
1953 const boost::multi_array<ReturnValue, 2> &ancestor_rotor_quad
1956 return basis.rotor(
i, iqn, ancestor_rotor_quad);
1961 template<
typename BasisType>
1964 static_assert(BasisType::hasHessian,
"Call to Hessian not available");
1968 return basis.hessian(
i, x);
1973 const BasisType &basis,
1976 const boost::multi_array<ReturnValue, 2> &ancestor_hessian_quad
1979 return basis.hessian(
i, iqn, ancestor_hessian_quad);
1985 template <
typename ScalarBasisType,
size_t N>
2009 const boost::multi_array<double, 2> &ancestor_basis_quad
2012 return basis.
function(
i, iqn, ancestor_basis_quad);
2017 template <
typename ScalarBasisType,
size_t N>
2041 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2044 return basis.
gradient(
i, iqn, ancestor_basis_quad);
2050 template <
typename ScalarBasisType,
size_t N>
2074 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2077 return 0.5 * (basis.
gradient(
i, iqn, ancestor_basis_quad) + basis.
gradient(
i, iqn, ancestor_basis_quad).transpose());
2083 template <
typename ScalarBasisType,
size_t N>
2107 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2110 return 0.5 * (basis.
gradient(
i, iqn, ancestor_basis_quad) - basis.
gradient(
i, iqn, ancestor_basis_quad).transpose());
2115 template <
typename ScalarBasisType,
size_t N>
2131 return basis.
curl(
i, x);
2139 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2142 return basis.
curl(
i, iqn, ancestor_basis_quad);
2148 template <
typename ScalarBasisType,
size_t N>
2172 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2175 return basis.
divergence(
i, iqn, ancestor_basis_quad);
2181 template <
typename ScalarBasisType,
size_t N>
2197 return basis.
rotor(
i, x);
2205 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2208 return basis.
rotor(
i, iqn, ancestor_basis_quad);
2215 template <
typename ScalarBasisType,
size_t N>
2239 const boost::multi_array<double, 2> &ancestor_basis_quad
2242 return basis.
function(
i, iqn, ancestor_basis_quad);
2247 template <
typename ScalarBasisType,
size_t N>
2271 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2274 return basis.
divergence(
i, iqn, ancestor_basis_quad);
2281 template <
typename ScalarBasisType,
size_t N>
2305 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2308 return basis.
gradient(
i, iqn, ancestor_basis_quad);
2319 template<BasisFunctionE BasisFunction>
2322 template<
typename BasisType>
2323 static boost::multi_array<typename detail::basis_evaluation_traits<BasisType, BasisFunction>::ReturnValue, 2>
2325 const BasisType & basis,
2331 boost::multi_array<typename traits::ReturnValue, 2>
2332 basis_quad( boost::extents[basis.dimension()][quad.size()] );
2334 for (
size_t i = 0;
i < basis.dimension();
i++) {
2335 for (
size_t iqn = 0; iqn < quad.size(); iqn++) {
2336 basis_quad[
i][iqn] = traits::evaluate(basis,
i, quad[iqn].vector());
2344 template<
typename BasisType>
2345 static boost::multi_array<typename detail::basis_evaluation_traits<Family<BasisType>, BasisFunction>::ReturnValue, 2>
2353 boost::multi_array<typename traits::ReturnValue, 2>
2354 basis_quad( boost::extents[basis.
dimension()][quad.size()] );
2357 const auto & ancestor_basis = basis.
ancestor();
2358 boost::multi_array<typename traits::ReturnValue, 2>
2361 Eigen::MatrixXd M = basis.
matrix();
2362 for (
size_t i = 0;
i < size_t(M.rows());
i++) {
2363 for (
size_t iqn = 0; iqn < quad.size(); iqn++) {
2364 basis_quad[
i][iqn] = M(
i,0) * ancestor_basis_quad[0][iqn];
2366 for (
size_t j = 1;
j < size_t(M.cols());
j++) {
2367 for (
size_t iqn = 0; iqn < quad.size(); iqn++) {
2368 basis_quad[
i][iqn] += M(
i,
j) * ancestor_basis_quad[
j][iqn];
2377 template <
typename BasisType,
size_t N>
2378 static boost::multi_array<typename detail::basis_evaluation_traits<TensorizedVectorFamily<BasisType, N>, BasisFunction>::ReturnValue, 2>
2386 boost::multi_array<typename traits::ReturnValue, 2>
2387 basis_quad(boost::extents[basis.
dimension()][quad.size()]);
2389 const auto &ancestor_basis = basis.
ancestor();
2390 const boost::multi_array<typename traits::AncestorBasisFunctionValue, 2> ancestor_basis_quad
2394 for (
size_t iqn = 0; iqn < quad.size(); iqn++){
2395 basis_quad[
i][iqn] = traits::evaluate(basis,
i, iqn, ancestor_basis_quad);
2402 template <
typename BasisType,
size_t N>
2403 static boost::multi_array<typename detail::basis_evaluation_traits<MatrixFamily<BasisType, N>, BasisFunction>::ReturnValue, 2>
2411 boost::multi_array<typename traits::ReturnValue, 2>
2412 basis_quad(boost::extents[basis.
dimension()][quad.size()]);
2414 const auto &ancestor_basis = basis.
ancestor();
2415 const boost::multi_array<typename traits::AncestorBasisFunctionValue, 2> ancestor_basis_quad
2419 for (
size_t iqn = 0; iqn < quad.size(); iqn++){
2420 basis_quad[
i][iqn] = traits::evaluate(basis,
i, iqn, ancestor_basis_quad);
2439 template<
typename T>
2441 boost::multi_array<T, 2> & basis_eval,
2442 const std::function<
double(
size_t,
size_t)> & inner_product
2445 auto induced_norm = [inner_product](
size_t i)->
double
2447 return std::sqrt( inner_product(
i,
i) );
2451 size_t Nb = basis_eval.shape()[0];
2453 size_t Nqn = basis_eval.shape()[1];
2455 Eigen::MatrixXd B = Eigen::MatrixXd::Zero(Nb, Nb);
2458 double norm = induced_norm(0);
2459 for (
size_t iqn = 0; iqn < Nqn; iqn++) {
2460 basis_eval[0][iqn] /= norm;
2464 for (
size_t ib = 1; ib < Nb; ib++) {
2466 Eigen::RowVectorXd coeffs = Eigen::RowVectorXd::Zero(ib);
2467 for (
size_t pb = 0; pb < ib; pb++) {
2468 coeffs(pb) = - inner_product(ib, pb);
2472 for (
size_t pb = 0; pb < ib; pb++) {
2473 for (
size_t iqn = 0; iqn < Nqn; iqn++) {
2474 basis_eval[ib][iqn] += coeffs(pb) * basis_eval[pb][iqn];
2479 double norm = induced_norm(ib);
2480 for (
size_t iqn = 0; iqn < Nqn; iqn++) {
2481 basis_eval[ib][iqn] /= norm;
2488 B.block(ib, 0, 1, ib) = coeffs * B.topLeftCorner(ib, ib);
2489 B(ib, ib) = 1./norm;
2504 double scalar_product(
const Eigen::Matrix<double, N, N> & x,
const Eigen::Matrix<double, N, N> & y)
2506 return (x.transpose() * y).trace();
2510 template <
typename Value>
2512 const boost::multi_array<Value, 2> &basis_quad,
2516 boost::multi_array<double, 2> basis_dot_v_quad(boost::extents[basis_quad.shape()[0]][basis_quad.shape()[1]]);
2517 std::transform(basis_quad.origin(), basis_quad.origin() + basis_quad.num_elements(),
2518 basis_dot_v_quad.origin(), [&
v](
const Value &x) ->
double { return scalar_product(x,v); });
2519 return basis_dot_v_quad;
2523 template <
typename Value>
2525 const boost::multi_array<Value, 2> & basis_quad,
2526 const std::vector<Value> &
v
2529 boost::multi_array<double, 2> basis_dot_v_quad(boost::extents[basis_quad.shape()[0]][basis_quad.shape()[1]]);
2530 for (std::size_t
i = 0;
i < basis_quad.shape()[0];
i++) {
2531 for (std::size_t iqn = 0; iqn < basis_quad.shape()[1]; iqn++) {
2535 return basis_dot_v_quad;
2540 const boost::multi_array<MatrixRd, 2> & basis_quad,
2541 const std::vector<VectorRd> & v_quad
2546 const std::vector<MatrixRd> & m_quad,
2547 const boost::multi_array<VectorRd, 2> & basis_quad
2552 const std::vector<VectorRd> & v_quad,
2553 const boost::multi_array<MatrixRd, 2> & basis_quad
2557 template<
typename BasisType>
2559 const BasisType & basis,
2561 boost::multi_array<typename BasisType::FunctionValue, 2> & basis_quad
2565 assert ( basis.dimension() == basis_quad.shape()[0] && qr.size() == basis_quad.shape()[1] );
2568 std::function<double(
size_t,
size_t)> inner_product
2569 = [&basis_quad, &qr](
size_t i,
size_t j)->
double
2572 for(
size_t iqn = 0; iqn < qr.size(); iqn++) {
2578 Eigen::MatrixXd B =
gram_schmidt(basis_quad, inner_product);
2585 template <
typename BasisType>
2587 const BasisType &basis,
2588 const Eigen::MatrixXd & GM
2592 assert(basis.dimension() ==
size_t(GM.rows()) && GM.rows() == GM.cols());
2594 Eigen::MatrixXd
L = GM.llt().matrixL();
2605 template<
typename FunctionValue>
2607 const boost::multi_array<FunctionValue, 2> & B2,
2611 const std::string sym =
"nonsym"
2615 assert ( qr.size() == B1.shape()[1] && qr.size() == B2.shape()[1] );
2617 assert ( nrows <= B1.shape()[0] && ncols <= B2.shape()[0] );
2620 Eigen::ArrayXd qr_weights = Eigen::ArrayXd::Zero(qr.size());
2621 for (
size_t iqn = 0; iqn < qr.size(); iqn++){
2622 qr_weights(iqn) = qr[iqn].w;
2625 Eigen::MatrixXd M = Eigen::MatrixXd::Zero(nrows, ncols);
2626 for (
size_t i = 0;
i < nrows;
i++) {
2628 if ( sym==
"sym" ) jcut =
i;
2629 for (
size_t j = 0;
j < jcut;
j++){
2632 for(
size_t j = jcut;
j < ncols;
j++) {
2633 std::vector<double> tmp(B1.shape()[1]);
2635 auto B1i = B1[ boost::indices[
i][boost::multi_array_types::index_range(0, B1.shape()[1])] ];
2636 auto B2j = B2[ boost::indices[
j][boost::multi_array_types::index_range(0, B1.shape()[1])] ];
2638 std::transform(B1i.begin(), B1i.end(), B2j.begin(), tmp.begin(), [](FunctionValue a, FunctionValue b)->double { return scalar_product(a,b);});
2639 Eigen::ArrayXd tmp_array = Eigen::Map<Eigen::ArrayXd, Eigen::Unaligned>(tmp.data(), tmp.size());
2641 M(
i,
j) = (qr_weights * tmp_array).sum();
2657 template<
typename FunctionValue>
2659 const boost::multi_array<FunctionValue, 2> & B2,
2661 const std::string sym =
"nonsym"
2664 return compute_gram_matrix<FunctionValue>(B1, B2, qr, B1.shape()[0], B2.shape()[0], sym);
2668 template<
typename FunctionValue>
2673 return compute_gram_matrix<FunctionValue>(B, B, qr,
"sym");
2678 const boost::multi_array<double, 2> & B2,
2684 const boost::multi_array<double, 2> & B2,
2688 const std::string sym =
"nonsym"
2693 const boost::multi_array<double, 2> & B2,
2695 const std::string sym =
"nonsym"
2700 const boost::multi_array<VectorRd, 2> & B2,
2704 const std::string sym =
"nonsym"
2709 const boost::multi_array<VectorRd, 2> & B2,
2711 const std::string sym =
"nonsym"
2715 template<
typename ScalarFamilyType,
size_t N>
2717 const boost::multi_array<double, 2> & scalar_family_quad,
2721 Eigen::MatrixXd Gram = Eigen::MatrixXd::Zero(MatFam.
dimension(), MatFam.
dimension());
2724 Eigen::MatrixXd scalarGram = compute_gram_matrix<double>(scalar_family_quad, qr);
2725 for (
size_t i = 0;
i < MatFam.
dimension();
i =
i+scalarGram.rows()){
2726 Gram.block(
i,
i, scalarGram.rows(), scalarGram.cols()) = scalarGram;
2735 template <
typename T>
2746 n_rows = B.shape()[0];
2750 const size_t num_quads = qr.size();
2753 assert(num_quads == B.shape()[1]);
2756 assert(n_rows <= B.shape()[0]);
2758 Eigen::VectorXd V = Eigen::VectorXd::Zero(n_rows);
2760 for (
size_t iqn = 0; iqn < num_quads; iqn++)
2762 double qr_weight = qr[iqn].w;
2763 T f_on_qr = f(qr[iqn].vector());
2764 for (
size_t i = 0;
i < n_rows;
i++)
2774 template <
typename T,
typename U>
2782 const std::string sym =
"nonsym"
2786 if (n_rows == 0 && n_cols == 0)
2788 n_rows = B1.shape()[0];
2789 n_cols = B2.shape()[0];
2793 const size_t num_quads = qr.size();
2795 assert(num_quads == B1.shape()[1] && num_quads == B2.shape()[1]);
2797 assert(n_rows <= B1.shape()[0] && n_cols <= B2.shape()[0]);
2799 Eigen::MatrixXd M = Eigen::MatrixXd::Zero(n_rows, n_cols);
2800 for (
size_t iqn = 0; iqn < num_quads; iqn++)
2802 double qr_weight = qr[iqn].w;
2803 U f_on_qr = f(qr[iqn].vector());
2804 for (
size_t i = 0;
i < n_rows;
i++)
2806 T f_B1 = f_on_qr * B1[
i][iqn];
2810 for (
size_t j = 0;
j < jcut;
j++)
2814 for (
size_t j = jcut;
j < n_cols;
j++)
2824 template <
typename T,
typename U>
2830 const std::string sym
2838 const FType<VectorRd> &f,
2839 const BasisQuad<VectorRd> &B1,
2840 const BasisQuad<double> &B2,
2849 const FType<VectorRd> &f,
2850 const BasisQuad<double> &B1,
2851 const BasisQuad<VectorRd> &B2,
2858 template <
typename Value>
2860 const boost::multi_array<Value, 2> & f_quad,
2861 const boost::multi_array<Value, 2> & g_quad,
2866 const size_t dimf = f_quad.shape()[0];
2867 const size_t dimg = g_quad.shape()[0];
2868 Eigen::MatrixXd M = Eigen::MatrixXd::Zero(dimf, dimg);
2870 std::vector<Value> int_f(dimf);
2871 for (
size_t i = 0;
i < dimf;
i++){
2872 int_f[
i] = qr_f[0].w * f_quad[
i][0];
2873 for (
size_t iqn = 1; iqn < qr_f.size(); iqn++){
2874 int_f[
i] += qr_f[iqn].w * f_quad[
i][iqn];
2877 std::vector<Value> int_g(dimg);
2878 for (
size_t j = 0;
j < dimg;
j++){
2879 int_g[
j] = qr_g[0].w * g_quad[
j][0];
2880 for (
size_t iqn = 1; iqn < qr_g.size(); iqn++){
2881 int_g[
j] += qr_g[iqn].w * g_quad[
j][iqn];
2886 for (
size_t i = 0;
i < dimf;
i++){
2887 for (
size_t j = 0;
j < dimg;
j++){
2896 template <
typename Value>
2898 const boost::multi_array<Value, 2> & f_quad,
2899 const boost::multi_array<Value, 2> & g_quad,
2907 template <
typename Value>
2909 const boost::multi_array<Value, 2> & f_quad,
2921 template<
typename BasisType>
2923 const std::function<
typename BasisType::FunctionValue(
const VectorRd &)> & f,
2924 const BasisType & basis,
2926 const boost::multi_array<typename BasisType::FunctionValue, 2> & basis_quad,
2927 const Eigen::MatrixXd &mass_basis = Eigen::MatrixXd::Zero(1,1)
2930 Eigen::MatrixXd Mass = mass_basis;
2931 if (Mass.norm() < 1e-13){
2935 Eigen::LDLT<Eigen::MatrixXd> cholesky_mass(Mass);
2936 Eigen::VectorXd b = Eigen::VectorXd::Zero(basis.dimension());
2937 for (
size_t i = 0;
i < basis.dimension();
i++) {
2938 for (
size_t iqn = 0; iqn < quad.size(); iqn++) {
2939 b(
i) += quad[iqn].w *
scalar_product(f(quad[iqn].vector()), basis_quad[
i][iqn]);
2943 return cholesky_mass.solve(b);
2956 template <
typename BasisType>
2967 const BasisType &basis
2969 m_dim(basis.dimension()),
2982 double h = 1./std::max(1.,
double(basis.max_degree()));
2984 VectorRd xi = E.vertex(0)->coords() + double(
i) * h * (E.vertex(1)->coords() - E.vertex(0)->coords());
2985 m_nodes.emplace_back(xi.x(), xi.y(), 1.);
2997 const BasisType &basis
2999 m_dim(basis.dimension()),
3019 VectorRd xi = x0 +
T.diam()*(double(indices[
i](0)) * e1 + double(indices[
i](1)) * e2)/std::max(1.,
double(basis.max_degree()));
3020 m_nodes.emplace_back(xi.x(), xi.y(), 1.);
3037 boost::multi_array<typename BasisType::FunctionValue, 2> &values
3045 Eigen::MatrixXd
L = Gram_onbasis.ldlt().solve(Gram_P_onbasis.transpose());
Basis for the space of curls (vectorial rot) of polynomials.
Definition basis.hpp:1271
Basis (or rather family) of divergence of an existing basis.
Definition basis.hpp:1334
Basis for the complement G^{c,k}(T) in P^k(F)^2 of the range of the gradient on a face.
Definition basis.hpp:1587
Basis for the space of gradients of polynomials.
Definition basis.hpp:1209
Basis for the space of Hessians of polynomials.
Definition basis.hpp:1457
Basis for the complement H^{c,k}(T) in P^k(T)^{2x2} of the range of the Hessian on a face.
Definition basis.hpp:1639
Matrix family obtained from a scalar family.
Definition basis.hpp:721
Scalar monomial basis on a cell.
Definition basis.hpp:168
Scalar monomial basis on an edge.
Definition basis.hpp:242
Generate a basis restricted to the first "dimension" functions.
Definition basis.hpp:1099
Basis for the complement R^{c,k}(T) in P^k(T)^2 of the range of the vectorial rotational on a face.
Definition basis.hpp:1521
Basis (or rather family) of scalar rotor of an existing basis.
Definition basis.hpp:1394
Generate a basis where the function indices are shifted.
Definition basis.hpp:989
Vector family for polynomial functions that are tangent to an edge (determined by the generator)
Definition basis.hpp:912
Vector family obtained by tensorization of a scalar family.
Definition basis.hpp:561
Compute max and min eigenvalues of all matrices for i
Definition compute_eigs.m:5
function[maxeig, mineig, cond, mat]
Definition compute_eigs.m:1
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2027
void RotorValue
Definition basis.hpp:1464
static constexpr const TensorRankE tensorRank
Definition basis.hpp:572
VectorRd CurlValue
Definition basis.hpp:1398
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1598
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1243
TensorizedVectorFamily< ScalarBasisType, N >::GradientValue ReturnValue
Definition basis.hpp:2055
static const bool hasRotor
Definition basis.hpp:1538
Eigen::Matrix< double, N, dimspace > GradientValue
Definition basis.hpp:564
static const bool hasGradient
Definition basis.hpp:1535
MatrixGradient(Eigen::Matrix< double, N, N > diff_x, Eigen::Matrix< double, N, N > diff_y)
Constructor.
Definition basis.hpp:115
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th function at point x.
Definition basis.hpp:380
BasisType::GradientValue GradientValue
Definition basis.hpp:316
Cell GeometricSupport
Definition basis.hpp:1649
static const bool hasFunction
Definition basis.hpp:255
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:281
BasisType AncestorType
Definition basis.hpp:1414
static boost::multi_array< typename detail::basis_evaluation_traits< TensorizedVectorFamily< BasisType, N >, BasisFunction >::ReturnValue, 2 > compute(const TensorizedVectorFamily< BasisType, N > &basis, const QuadratureRule &quad)
Evaluate a tensorized family at quadrature nodes (optimised compared the generic basis evaluation,...
Definition basis.hpp:2379
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_basis_quad)
Definition basis.hpp:2070
static const bool hasDivergence
Definition basis.hpp:1115
MatrixFamily< ScalarBasisType, N >::GradientValue ReturnValue
Definition basis.hpp:2286
VectorRd FunctionValue
Definition basis.hpp:1523
static const bool hasFunction
Definition basis.hpp:1002
static constexpr const bool hasAncestor
Definition basis.hpp:1652
size_t dimension() const
Return the dimension of the family.
Definition basis.hpp:948
double DivergenceValue
Definition basis.hpp:1339
VectorRd CurlValue
Definition basis.hpp:246
BasisType::RotorValue ReturnValue
Definition basis.hpp:1942
HessianValue hessian(size_t i, const VectorRd &x) const
Evaluate the Hessian of the i-th basis function at point x.
Definition basis.hpp:484
static const bool hasGradient
Definition basis.hpp:575
static const bool hasFunction
Definition basis.hpp:574
const size_t matrixSize() const
Return the dimension of the matrices in the family.
Definition basis.hpp:835
DivergenceValue divergence(size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_gradient_quad) const
Evaluate the divergence of the i-th basis function at a quadrature point iqn, knowing all the gradien...
Definition basis.hpp:667
BasisType m_basis
Definition basis.hpp:3052
static constexpr const bool hasAncestor
Definition basis.hpp:924
BasisType::FunctionValue FunctionValue
Definition basis.hpp:315
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1532
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1494
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th function at point x.
Definition basis.hpp:353
static std::vector< VectorZd > complete(size_t degree)
Definition basis.hpp:88
static ReturnValue evaluate(const MatrixFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2257
BasisType::FunctionValue ReturnValue
Definition basis.hpp:1801
static constexpr const bool hasAncestor
Definition basis.hpp:733
static const bool hasHessian
Definition basis.hpp:1539
CurlValue curl(size_t i, const VectorRd &x) const
Evaluate the curl of the i-th basis function at point x.
Definition basis.hpp:643
static const bool hasDivergence
Definition basis.hpp:184
size_t dimension() const
Return the dimension of the family.
Definition basis.hpp:769
void RotorValue
Definition basis.hpp:1527
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:1280
Eigen::MatrixXd compute_closure_matrix(const boost::multi_array< Value, 2 > &f_quad, const boost::multi_array< Value, 2 > &g_quad, const QuadratureRule &qr_f, const QuadratureRule &qr_g)
Computes closure equation matrices, from evaluations at quadrature nodes. For two families of functio...
Definition basis.hpp:2859
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1158
CurlBasis(const BasisType &basis)
Constructor.
Definition basis.hpp:1294
void HessianValue
Definition basis.hpp:568
BasisType::HessianValue HessianValue
Definition basis.hpp:320
MatrixGradient()
Default constructor.
Definition basis.hpp:125
VectorRd CurlValue
Definition basis.hpp:1338
BasisFunctionE
Definition basis.hpp:1700
double FunctionValue
Definition basis.hpp:1336
void RotorValue
Definition basis.hpp:174
static const bool hasHessian
Definition basis.hpp:1658
static const bool hasCurl
Definition basis.hpp:1475
static const bool hasFunction
Definition basis.hpp:1346
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_hessian_quad)
Definition basis.hpp:1972
static const bool hasDivergence
Definition basis.hpp:928
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:322
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_gradient_quad)
Definition basis.hpp:1880
static const bool hasRotor
Definition basis.hpp:1477
static constexpr const TensorRankE tensorRank
Definition basis.hpp:923
DivergenceBasis(const BasisType &basis)
Constructor.
Definition basis.hpp:1355
static constexpr const bool hasAncestor
Definition basis.hpp:1406
static constexpr const bool hasAncestor
Definition basis.hpp:1111
VectorRd CurlValue
Definition basis.hpp:916
static const bool hasGradient
Definition basis.hpp:1285
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_divergence_quad)
Definition basis.hpp:1926
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1405
void HessianValue
Definition basis.hpp:919
Eigen::Vector2d VectorRd
Definition basis.hpp:55
Eigen::Vector2i VectorZd
Definition basis.hpp:56
BasisType::RotorValue RotorValue
Definition basis.hpp:995
void DivergenceValue
Definition basis.hpp:1592
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:1403
Eigen::Matrix< double, dimspace, dimspace > CurlValue
Definition basis.hpp:1275
static const bool hasHessian
Definition basis.hpp:930
static const bool hasGradient
Definition basis.hpp:926
Eigen::Matrix< double, dimspace, dimspace > GradientValue
Definition basis.hpp:1274
boost::multi_array< VectorRd, 2 > matrix_vector_product(const boost::multi_array< MatrixRd, 2 > &basis_quad, const std::vector< VectorRd > &v_quad)
Take the product of a matrix-valued basis with a vector.
Definition basis.cpp:175
BasisType::DivergenceValue DivergenceValue
Definition basis.hpp:318
Family(const BasisType &basis, const Eigen::MatrixXd &matrix)
Constructor.
Definition basis.hpp:336
TensorizedVectorFamily< ScalarBasisType, N >::GradientValue ReturnValue
Definition basis.hpp:2022
BasisType::GradientValue GradientValue
Definition basis.hpp:1102
static ReturnValue evaluate(const MatrixFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2291
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1874
static const bool hasRotor
Definition basis.hpp:1116
static const bool hasRotor
Definition basis.hpp:737
BasisType AncestorType
Definition basis.hpp:1480
static boost::multi_array< typename detail::basis_evaluation_traits< Family< BasisType >, BasisFunction >::ReturnValue, 2 > compute(const Family< BasisType > &basis, const QuadratureRule &quad)
Evaluate a 'Family' of functions at quadrature nodes (optimised compared the generic basis evaluation...
Definition basis.hpp:2346
void GradientValue
Definition basis.hpp:1643
static const bool hasFunction
Definition basis.hpp:1653
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_curl_quad)
Definition basis.hpp:1903
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th basis function at point x.
Definition basis.cpp:31
constexpr const ScalarFamilyType & ancestor() const
Return the ancestor (family that has been tensorized)
Definition basis.hpp:691
static const bool hasRotor
Definition basis.hpp:1288
BasisType AncestorType
Definition basis.hpp:333
static const bool hasDivergence
Definition basis.hpp:258
static const bool hasCurl
Definition basis.hpp:257
void HessianValue
Definition basis.hpp:1216
static const bool hasCurl
Definition basis.hpp:927
size_t dimension() const
Dimension of the family. This is actually the number of functions in the family, not necessarily line...
Definition basis.hpp:347
static const bool hasDivergence
Definition basis.hpp:1476
constexpr const BasisType & ancestor() const
Return the ancestor (basis that the gradient was taken of)
Definition basis.hpp:1255
static ReturnValue evaluate(const MatrixFamily< ScalarBasisType, N > &basis, size_t i, size_t iqn, const boost::multi_array< double, 2 > &ancestor_basis_quad)
Definition basis.hpp:2235
static const bool hasCurl
Definition basis.hpp:1602
void GradientValue
Definition basis.hpp:1461
RotorValue rotor(size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_gradient_quad) const
Evaluate the scalar rotor of the i-th basis function at a quadrature point iqn, knowing all the gradi...
Definition basis.hpp:683
void HessianValue
Definition basis.hpp:249
RestrictedBasis(const BasisType &basis, const size_t &dimension)
Constructor.
Definition basis.hpp:1122
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2191
BasisType::CurlValue CurlValue
Definition basis.hpp:993
static const bool hasFunction
Definition basis.hpp:1222
BasisType::HessianValue HessianValue
Definition basis.hpp:1106
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1310
Edge GeometricSupport
Definition basis.hpp:921
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2060
const Eigen::MatrixXd symmetriseOperator() const
Return the symmetrisation operator, the rN^2 square matrix that to a given vector of coefficients on ...
Definition basis.hpp:847
GradientBasis(const BasisType &basis)
Constructor.
Definition basis.hpp:1234
boost::multi_array< T, 2 > BasisQuad
type for a family of basis functions evaluated on quadrature nodes
Definition basis.hpp:61
DivergenceValue divergence(size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_gradient_quad) const
Evaluate the divergence of the i-th basis function at a quadrature point iqn, knowing all the gradien...
Definition basis.hpp:800
VectorRd generator() const
Returns the generator of the basis functions.
Definition basis.hpp:972
static const bool hasCurl
Definition basis.hpp:579
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1282
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th basis function at point x.
Definition basis.hpp:623
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th basis function at point x.
Definition basis.cpp:76
static const bool hasDivergence
Definition basis.hpp:329
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1428
TensorizedVectorFamily< ScalarBasisType, N >::FunctionValue ReturnValue
Definition basis.hpp:1990
void DivergenceValue
Definition basis.hpp:247
static const bool hasHessian
Definition basis.hpp:1350
VectorRd FunctionValue
Definition basis.hpp:1589
static const bool hasHessian
Definition basis.hpp:1229
Cell GeometricSupport
Definition basis.hpp:1596
const Eigen::MatrixXd transposeOperator() const
Return the transpose operator, the rN^2 square matrix that to a given vector of coefficients on the M...
Definition basis.hpp:841
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th basis function at point x.
Definition basis.hpp:809
MatrixGradient & operator+=(const MatrixGradient &G)
Increment.
Definition basis.hpp:137
void RotorValue
Definition basis.hpp:918
RotorBasis(const BasisType &basis)
Constructor.
Definition basis.hpp:1417
VectorRd GradientValue
Definition basis.hpp:1397
ScalarFamilyType AncestorType
Definition basis.hpp:741
size_t dimension() const
Return the dimension of the basis.
Definition basis.hpp:1135
VectorRd AncestorBasisFunctionValue
Definition basis.hpp:2288
GradientValue gradient(size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_gradient_quad) const
Evaluate the gradient of the i-th basis function at a quadrature point iqn, knowing all the gradients...
Definition basis.hpp:633
Cell GeometricSupport
Definition basis.hpp:177
GradientValue gradient(size_t i, size_t iqn, const boost::multi_array< GradientValue, 2 > &ancestor_gradient_quad) const
Evaluate the gradient of the i-th function at a quadrature point iqn, knowing all the gradients of an...
Definition basis.hpp:393
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:603
MatrixGradient< N > GradientValue
Definition basis.hpp:724
double scalar_product(const double &x, const double &y)
Scalar product between two reals.
Definition basis.cpp:163
void HessianValue
Definition basis.hpp:1465
size_t dimension() const
Dimension of the basis.
Definition basis.hpp:1614
Eigen::Matrix< double, N, 1 > FunctionValue
Definition basis.hpp:563
void DivergenceValue
Definition basis.hpp:1645
FunctionValue function(size_t i, size_t iqn, const boost::multi_array< double, 2 > &ancestor_value_quad) const
Evaluate the i-th basis function at a quadrature point iqn, knowing all the values of ancestor basis ...
Definition basis.hpp:613
Eigen::MatrixXd compute_gram_matrix(const boost::multi_array< VectorRd, 2 > &B1, const boost::multi_array< double, 2 > &B2, const QuadratureRule &qr)
Compute the Gram-like matrix given a family of vector-valued and one of scalar-valued functions by te...
Definition basis.cpp:239
void DivergenceValue
Definition basis.hpp:173
BasisType::FunctionValue FunctionValue
Definition basis.hpp:1101
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:1467
static constexpr const bool hasAncestor
Definition basis.hpp:1221
QuadratureRule m_nodes
Nodes for the interpolation.
Definition basis.hpp:3056
Family< MatrixFamily< ScalarBasisType, N > > IsotropicMatrixFamily(const ScalarBasisType &B)
From a scalar family B, constructs a "Family" of "MatrixFamily" (built on B, of size NxN) that repres...
Definition basis.hpp:1771
BasisType AncestorType
Definition basis.hpp:1119
BasisType::GradientValue ReturnValue
Definition basis.hpp:1825
TensorizedVectorFamily< ScalarBasisType, N >::CurlValue ReturnValue
Definition basis.hpp:2120
static constexpr const bool hasAncestor
Definition basis.hpp:573
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1000
MatrixRd HessianValue
Definition basis.hpp:175
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1249
Eigen::Matrix< double, dimspace, dimspace > GradientValue
Definition basis.hpp:1212
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_basis_quad)
Definition basis.hpp:2168
BasisType::RotorValue RotorValue
Definition basis.hpp:319
static const bool hasGradient
Definition basis.hpp:1113
BasisType::DivergenceValue ReturnValue
Definition basis.hpp:1919
static const bool hasCurl
Definition basis.hpp:738
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:697
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_gradient_quad)
Definition basis.hpp:1856
static const bool hasHessian
Definition basis.hpp:331
static const std::vector< VectorRd > basisRd
Definition basis.hpp:58
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:998
static const bool hasDivergence
Definition basis.hpp:577
static const bool hasDivergence
Definition basis.hpp:1410
static const bool hasHessian
Definition basis.hpp:739
static constexpr const bool hasAncestor
Definition basis.hpp:1470
void HessianValue
Definition basis.hpp:1278
static const bool hasHessian
Definition basis.hpp:1289
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1110
HessianValue hessian(size_t i, size_t iqn, const boost::multi_array< HessianValue, 2 > &ancestor_hessian_quad) const
Evaluate the hessian of the i-th function at a quadrature point iqn, knowing all the hessian of ances...
Definition basis.hpp:497
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:1147
TangentFamily(const ScalarFamilyType &scalar_family, const VectorRd &generator)
Constructor.
Definition basis.hpp:935
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:954
RotorValue rotor(size_t i, const VectorRd &x) const
Evaluate the scalar rotor of the i-th basis function at point x.
Definition basis.hpp:675
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1667
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th basis function at point x.
Definition basis.hpp:659
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1548
MatrixFamily< ScalarBasisType, N >::DivergenceValue ReturnValue
Definition basis.hpp:2252
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_basis_quad)
Definition basis.hpp:2037
Eigen::MatrixXd compute_weighted_gram_matrix(const FType< VectorRd > &f, const BasisQuad< VectorRd > &B1, const BasisQuad< double > &B2, const QuadratureRule &qr, size_t n_rows, size_t n_cols)
Computes the Gram-like matrix of integrals (f dot phi_i, phi_j)
Definition basis.cpp:368
BasisType::FunctionValue FunctionValue
Definition basis.hpp:991
void HessianValue
Definition basis.hpp:1528
static const bool hasHessian
Definition basis.hpp:1605
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1220
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1041
static ReturnValue evaluate(const MatrixFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2225
boost::multi_array< outValue, 2 > transform_values_quad(const boost::multi_array< inValue, 2 > &B_quad, const FunctionType &F)
Takes an array B_quad of values at quadrature nodes and applies the function F to all of them....
Definition basis.hpp:1713
static const bool hasDivergence
Definition basis.hpp:1349
static constexpr const TensorRankE tensorRank
Definition basis.hpp:253
size_t dimension() const
Return the dimension of the basis.
Definition basis.hpp:1023
VectorRd FunctionValue
Definition basis.hpp:914
static const bool hasHessian
Definition basis.hpp:260
VectorRd AncestorBasisFunctionValue
Definition basis.hpp:2090
static boost::multi_array< typename detail::basis_evaluation_traits< MatrixFamily< BasisType, N >, BasisFunction >::ReturnValue, 2 > compute(const MatrixFamily< BasisType, N > &basis, const QuadratureRule &quad)
Evaluate a Matrix family at quadrature nodes (optimised compared the generic basis evaluation,...
Definition basis.hpp:2404
static const bool hasFunction
Definition basis.hpp:925
RotorValue rotor(size_t i, const VectorRd &x) const
Evaluate the scalar rotor of the i-th basis function at point x.
Definition basis.hpp:1073
static const bool hasFunction
Definition basis.hpp:734
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:195
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:1035
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1651
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th basis function at point x.
Definition basis.hpp:1049
ScalarFamilyType::GeometricSupport GeometricSupport
Definition basis.hpp:730
CurlValue curl(size_t i, const VectorRd &x) const
Evaluate the curl of the i-th basis function at point x.
Definition basis.hpp:1057
static ReturnValue evaluate(const MatrixFamily< ScalarBasisType, N > &basis, size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_basis_quad)
Definition basis.hpp:2301
double AncestorBasisFunctionValue
Definition basis.hpp:2222
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1500
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:775
boost::multi_array< typename BasisType::FunctionValue, 2 > m_on_basis_nodes
Definition basis.hpp:3054
double RotorValue
Definition basis.hpp:1593
void RotorValue
Definition basis.hpp:1277
Eigen::Matrix2d MatrixRd
Definition basis.hpp:54
BasisType::DivergenceValue DivergenceValue
Definition basis.hpp:994
GradientValue CurlValue
Definition basis.hpp:565
BasisType::CurlValue ReturnValue
Definition basis.hpp:1896
RotorValue rotor(size_t i, size_t iqn, const boost::multi_array< RotorValue, 2 > &ancestor_rotor_quad) const
Evaluate the scalar rotor of the i-th function at a quadrature point iqn, knowing all the scalar roto...
Definition basis.hpp:471
static const bool hasDivergence
Definition basis.hpp:736
TensorizedVectorFamily< ScalarBasisType, N >::RotorValue ReturnValue
Definition basis.hpp:2186
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_basis_quad)
Definition basis.hpp:2135
CurlValue curl(size_t i, const VectorRd &x) const
Evaluate the curl of the i-th basis function at point x.
Definition basis.hpp:1174
static const bool hasRotor
Definition basis.hpp:1228
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2093
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1995
void RotorValue
Definition basis.hpp:1646
static std::function< Eigen::MatrixXd(const Eigen::MatrixXd &)> symmetrise_matrix
Function to symmetrise a matrix (useful together with transform_values_quad)
Definition basis.hpp:1762
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th basis function at point x.
Definition basis.hpp:1065
VectorRd AncestorBasisFunctionValue
Definition basis.hpp:2024
VectorZd powers(size_t i) const
Returns the powers of the i-th basis function (its degree can be found using powers(i)....
Definition basis.hpp:219
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th function at point x.
Definition basis.hpp:432
static constexpr const bool hasAncestor
Definition basis.hpp:1001
BasisType::GradientValue FunctionValue
Definition basis.hpp:1211
static const bool hasGradient
Definition basis.hpp:1601
Eigen::Matrix< double, dimspace, dimspace > CurlValue
Definition basis.hpp:1591
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:213
VectorRd FunctionValue
Definition basis.hpp:1273
void DivergenceValue
Definition basis.hpp:1214
void DivergenceValue
Definition basis.hpp:917
constexpr const BasisType & ancestor() const
Return the underlying complete basis.
Definition basis.hpp:1029
BasisType::GradientValue GradientValue
Definition basis.hpp:992
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_basis_quad)
Definition basis.hpp:2103
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th basis function at point x.
Definition basis.hpp:1166
static const bool hasHessian
Definition basis.hpp:1412
VectorRd AncestorBasisFunctionValue
Definition basis.hpp:2155
double DivergenceValue
Definition basis.hpp:1526
static constexpr const bool hasAncestor
Definition basis.hpp:180
DivergenceValue rotor(size_t i, const VectorRd &x) const
Evaluate the scalar rotor of the i-th function at point x.
Definition basis.hpp:458
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_basis_quad)
Definition basis.hpp:2201
VectorRd GradientValue
Definition basis.hpp:915
constexpr const BasisType & ancestor() const
Return the underlying complete basis.
Definition basis.hpp:1141
Eigen::VectorXd l2_projection(const std::function< typename BasisType::FunctionValue(const VectorRd &)> &f, const BasisType &basis, QuadratureRule &quad, const boost::multi_array< typename BasisType::FunctionValue, 2 > &basis_quad, const Eigen::MatrixXd &mass_basis=Eigen::MatrixXd::Zero(1, 1))
Compute the L2-projection of a function.
Definition basis.hpp:2922
void CurlValue
Definition basis.hpp:1644
BasisType AncestorType
Definition basis.hpp:1291
static const bool hasCurl
Definition basis.hpp:1114
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1469
static const bool hasCurl
Definition basis.hpp:1226
static const bool hasFunction
Definition basis.hpp:1534
CurlValue curl(size_t i, const VectorRd &x) const
Evaluate the curl (vector rot) of the i-th basis function at point x.
Definition basis.cpp:41
static const bool hasRotor
Definition basis.hpp:1411
void DivergenceValue
Definition basis.hpp:1463
FunctionValue function(size_t i, size_t iqn, const boost::multi_array< double, 2 > &ancestor_value_quad) const
Evaluate the i-th basis function at a quadrature point iqn, knowing all the values of ancestor basis ...
Definition basis.hpp:783
static const bool hasGradient
Definition basis.hpp:256
Eigen::MatrixXd PermuteTensorization(const size_t a, const size_t b)
Returns the matrix giving the permutation of the tensorization of a family of size a with a family of...
Definition basis.cpp:224
Cell GeometricSupport
Definition basis.hpp:1530
VectorRd AncestorBasisFunctionValue
Definition basis.hpp:2057
MatrixFamily(const ScalarFamilyType &scalar_family)
Definition basis.hpp:743
constexpr int dimspace
Dimension, and generic types for vector in correct dimension (makes it easier to translate a code bet...
Definition basis.hpp:53
TensorRankE
Definition basis.hpp:67
MatrixGradient< N > operator*(double scalar, MatrixGradient< N > const &G)
Multiplication of a MatrixGradient from the left (non-member function to be able to multiply from the...
Definition basis.hpp:156
RotorValue rotor(size_t i, const VectorRd &x) const
Evaluate the rotor of the i-th basis function at point x.
Definition basis.cpp:131
double FunctionValue
Definition basis.hpp:170
const Eigen::MatrixXd traceOperator() const
Returns the matrix corresponding to the trace, expressed as a linear operator between that MatrixFami...
Definition basis.hpp:878
Eigen::MatrixXd gram_schmidt(boost::multi_array< T, 2 > &basis_eval, const std::function< double(size_t, size_t)> &inner_product)
Definition basis.hpp:2440
CurlValue curl(size_t i, const VectorRd &x) const
Evaluate the curl of the i-th function at point x.
Definition basis.hpp:406
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2125
VectorRd CurlValue
Definition basis.hpp:725
VectorRd AncestorBasisFunctionValue
Definition basis.hpp:2122
static const bool hasRotor
Definition basis.hpp:185
BasisType::GradientValue ReturnValue
Definition basis.hpp:1849
static constexpr const bool hasAncestor
Definition basis.hpp:1599
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th basis function at point x.
Definition basis.hpp:1182
void RotorValue
Definition basis.hpp:1215
MatrixGradient operator+(const MatrixGradient &G)
Addition.
Definition basis.hpp:132
constexpr const BasisType & ancestor() const
Return the ancestor (basis that the gradient was taken of)
Definition basis.hpp:1316
static const bool hasCurl
Definition basis.hpp:328
BasisType::CurlValue CurlValue
Definition basis.hpp:317
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1897
static const bool hasGradient
Definition basis.hpp:1225
MatrixRd FunctionValue
Definition basis.hpp:1641
static const bool hasGradient
Definition basis.hpp:1347
Family< BasisType > m_on_basis
Definition basis.hpp:3053
static std::function< Eigen::MatrixXd(const Eigen::MatrixXd &)> skew_symmetrise_matrix
Function to skew-symmetrise a matrix (useful together with transform_values_quad)
Definition basis.hpp:1766
static const bool hasCurl
Definition basis.hpp:1536
void RotorValue
Definition basis.hpp:727
static constexpr const bool hasAncestor
Definition basis.hpp:254
BasisType::HessianValue HessianValue
Definition basis.hpp:996
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2158
DivergenceValue divergence(size_t i, size_t iqn, const boost::multi_array< DivergenceValue, 2 > &ancestor_divergence_quad) const
Evaluate the divergence of the i-th function at a quadrature point iqn, knowing all the divergences o...
Definition basis.hpp:445
double DivergenceValue
Definition basis.hpp:566
std::function< T(const VectorRd &)> FType
type for function of point. T is the type of value of the function
Definition basis.hpp:64
BasisType AncestorType
Definition basis.hpp:1231
static const bool hasHessian
Definition basis.hpp:186
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th basis function at point x.
Definition basis.hpp:791
constexpr const BasisType & ancestor() const
Return the ancestor (basis that the gradient was taken of)
Definition basis.hpp:1378
RotorValue rotor(size_t i, const VectorRd &x) const
Evaluate the rotor of the i-th basis function at point x.
Definition basis.hpp:1190
CurlValue curl(size_t i, size_t iqn, const boost::multi_array< CurlValue, 2 > &ancestor_curl_quad) const
Evaluate the curl of the i-th function at a quadrature point iqn, knowing all the curls of ancestor b...
Definition basis.hpp:419
static constexpr const TensorRankE tensorRank
Definition basis.hpp:732
static const bool hasRotor
Definition basis.hpp:578
Eigen::Matrix< double, dimspace, dimspace > GradientValue
Definition basis.hpp:1524
Eigen::Matrix< double, N, N > dx
Definition basis.hpp:149
static const bool hasGradient
Definition basis.hpp:1408
Family< BasisType > l2_orthonormalize(const BasisType &basis, const QuadratureRule &qr, boost::multi_array< typename BasisType::FunctionValue, 2 > &basis_quad)
-orthonormalization: simply consists in using gram_schmidt() with the specific l2 inner product
Definition basis.hpp:2558
static const bool hasDivergence
Definition basis.hpp:1656
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:1342
static constexpr const TensorRankE tensorRank
Definition basis.hpp:179
constexpr const BasisType & ancestor() const
Return the ancestor.
Definition basis.hpp:516
static const bool hasDivergence
Definition basis.hpp:1287
static const bool hasCurl
Definition basis.hpp:1004
Eigen::Matrix< double, N, N > FunctionValue
Definition basis.hpp:723
VectorRd GradientValue
Definition basis.hpp:1337
Family< TensorizedVectorFamily< ScalarBasisType, N > > GenericTensorization(const ScalarBasisType &B, const std::vector< Eigen::VectorXd > &v)
From a scalar family B=(B_1..B_r) and vectors (v_1..v_k) in R^N, constructs a "Family" of "Tensorized...
Definition basis.hpp:1728
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1943
Family< BasisType > family(boost::multi_array< typename BasisType::FunctionValue, 2 > &values)
Returns the decomposed polynomials as a Family of the provided basis.
Definition basis.hpp:3036
static const bool hasRotor
Definition basis.hpp:330
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1304
static const bool hasFunction
Definition basis.hpp:1284
ScalarFamilyType AncestorType
Definition basis.hpp:932
size_t dimension() const
Return the dimension of the family.
Definition basis.hpp:597
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1966
static const bool hasCurl
Definition basis.hpp:1286
const ScalarFamilyType & ancestor() const
Return the ancestor (family that has been tensorized)
Definition basis.hpp:829
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1366
static const bool hasGradient
Definition basis.hpp:1474
static const bool hasFunction
Definition basis.hpp:181
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th basis function at point x.
Definition basis.cpp:106
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_rotor_quad)
Definition basis.hpp:1949
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:1108
VectorRd AncestorBasisFunctionValue
Definition basis.hpp:2188
VectorRd CurlValue
Definition basis.hpp:1213
BasisType AncestorType
Definition basis.hpp:1009
void HessianValue
Definition basis.hpp:1340
FunctionValue function(size_t i, size_t iqn, const boost::multi_array< FunctionValue, 2 > &ancestor_value_quad) const
Evaluate the i-th function at a quadrature point iqn, knowing all the values of ancestor basis functi...
Definition basis.hpp:366
GradientValue gradient(size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_gradient_quad) const
Evaluate the gradient of the i-th basis function at a quadrature point iqn, knowing all the gradients...
Definition basis.hpp:819
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:1560
MatrixRd FunctionValue
Definition basis.hpp:1459
Eigen::Matrix< double, dimspace, dimspace > GradientValue
Definition basis.hpp:1590
VectorRd GradientValue
Definition basis.hpp:245
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:1676
DecomposePoly(const Cell &T, const BasisType &basis)
Constructor for cell.
Definition basis.hpp:2995
static const bool hasGradient
Definition basis.hpp:735
static const bool hasRotor
Definition basis.hpp:259
static const bool hasCurl
Definition basis.hpp:1409
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1920
static const bool hasRotor
Definition basis.hpp:1657
void RotorValue
Definition basis.hpp:1400
TensorizedVectorFamily< ScalarBasisType, N >::GradientValue ReturnValue
Definition basis.hpp:2088
static const bool hasDivergence
Definition basis.hpp:1227
static const bool hasHessian
Definition basis.hpp:1007
void HessianValue
Definition basis.hpp:1594
static const bool hasFunction
Definition basis.hpp:1112
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_gradient_quad)
Definition basis.hpp:1832
const Eigen::MatrixXd symmetricBasis() const
Returns the matrix that can be used, in a Family of this MatrixBasis, to create a basis of the subspa...
Definition basis.hpp:853
static const bool hasRotor
Definition basis.hpp:1604
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:966
VectorRd GradientValue
Definition basis.hpp:171
HessianValue hessian(size_t i, const VectorRd &x) const
Evaluate the Hessian of the i-th basis function at point x.
Definition basis.cpp:46
double RotorValue
Definition basis.hpp:567
BasisType::DivergenceValue DivergenceValue
Definition basis.hpp:1104
BasisType::GradientValue ReturnValue
Definition basis.hpp:1873
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1802
ScalarFamilyType::GeometricSupport GeometricSupport
Definition basis.hpp:570
Eigen::VectorXd integrate(const FType< T > &f, const BasisQuad< T > &B, const QuadratureRule &qr, size_t n_rows=0)
Compute the integral of a given function against all functions from a family.
Definition basis.hpp:2736
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_value_quad)
Definition basis.hpp:1808
VectorZd powers(size_t i) const
Returns the powers of the i-th basis function (not including the vector part)
Definition basis.hpp:1566
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1434
static const bool hasHessian
Definition basis.hpp:582
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:522
HessianBasis(const BasisType &basis)
Constructor.
Definition basis.hpp:1483
double AncestorBasisFunctionValue
Definition basis.hpp:1992
static const bool hasRotor
Definition basis.hpp:929
void HessianValue
Definition basis.hpp:728
Eigen::Matrix< double, N, N > dy
Definition basis.hpp:150
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:1218
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1850
static const bool hasFunction
Definition basis.hpp:326
static const bool hasGradient
Definition basis.hpp:182
ShiftedBasis(const BasisType &basis, const int shift)
Constructor.
Definition basis.hpp:1012
static constexpr const bool hasAncestor
Definition basis.hpp:1533
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1372
HessianValue hessian(size_t i, const VectorRd &x) const
Evaluate the Hessian of the i-th basis function at point x.
Definition basis.hpp:1081
ScalarFamilyType AncestorType
Definition basis.hpp:584
void HessianValue
Definition basis.hpp:1647
const std::shared_ptr< RolyComplBasisCell > & rck() const
Return the Rck basis.
Definition basis.hpp:1626
void CurlValue
Definition basis.hpp:1462
static const bool hasHessian
Definition basis.hpp:1478
CurlValue curl(size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_gradient_quad) const
Evaluate the curl of the i-th basis function at a quadrature point iqn, knowing all the gradients of ...
Definition basis.hpp:651
constexpr const BasisType & ancestor() const
Return the ancestor (basis that the gradient was taken of)
Definition basis.hpp:1440
static const bool hasGradient
Definition basis.hpp:327
static constexpr const TensorRankE tensorRank
Definition basis.hpp:324
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1344
size_t dimension() const
Dimension of the basis.
Definition basis.hpp:269
static constexpr const bool hasAncestor
Definition basis.hpp:1283
void DivergenceValue
Definition basis.hpp:1276
static const bool hasFunction
Definition basis.hpp:1471
BasisType AncestorType
Definition basis.hpp:1352
static const bool hasDivergence
Definition basis.hpp:1537
QuadratureRule get_nodes() const
Return the set of nodes (useful to compute value of polynomial to decompose via evaluate_quad)
Definition basis.hpp:3030
static const bool hasFunction
Definition basis.hpp:1600
DecomposePoly(const Edge &E, const BasisType &basis)
Constructor for edge.
Definition basis.hpp:2965
void HessianValue
Definition basis.hpp:1401
VectorRd AncestorBasisFunctionValue
Definition basis.hpp:2254
size_t m_dim
Definition basis.hpp:3051
static const bool hasCurl
Definition basis.hpp:1655
static const bool hasCurl
Definition basis.hpp:183
Edge GeometricSupport
Definition basis.hpp:251
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, size_t iqn, const boost::multi_array< double, 2 > &ancestor_basis_quad)
Definition basis.hpp:2005
boost::multi_array< VectorRd, 2 > vector_matrix_product(const std::vector< VectorRd > &v_quad, const boost::multi_array< MatrixRd, 2 > &basis_quad)
Take the product of (the transposed of) a vector with a matrix-valued basis.
Definition basis.cpp:205
constexpr const ScalarFamilyType & ancestor() const
Return the ancestor (family used for the tangent)
Definition basis.hpp:960
static boost::multi_array< typename detail::basis_evaluation_traits< BasisType, BasisFunction >::ReturnValue, 2 > compute(const BasisType &basis, const QuadratureRule &quad)
Generic basis evaluation.
Definition basis.hpp:2324
double FunctionValue
Definition basis.hpp:1396
static ReturnValue evaluate(const MatrixFamily< ScalarBasisType, N > &basis, size_t i, size_t iqn, const boost::multi_array< VectorRd, 2 > &ancestor_basis_quad)
Definition basis.hpp:2267
BasisType::HessianValue ReturnValue
Definition basis.hpp:1965
Eigen::Matrix< double, N, 1 > DivergenceValue
Definition basis.hpp:726
static const bool hasDivergence
Definition basis.hpp:1005
const Eigen::MatrixXd & matrix() const
Return the coefficient matrix.
Definition basis.hpp:510
size_t m_nb_nodes
Definition basis.hpp:3055
constexpr const BasisType & ancestor() const
Return the ancestor (basis that the gradient was taken of)
Definition basis.hpp:1506
static const bool hasDivergence
Definition basis.hpp:1603
static const bool hasRotor
Definition basis.hpp:1006
static const bool hasGradient
Definition basis.hpp:1003
MatrixGradient operator-(const MatrixGradient &G)
Subtraction.
Definition basis.hpp:144
MatrixFamily< ScalarBasisType, N >::FunctionValue ReturnValue
Definition basis.hpp:2220
double FunctionValue
Definition basis.hpp:244
static const bool hasCurl
Definition basis.hpp:1348
static constexpr const bool hasAncestor
Definition basis.hpp:325
static constexpr const bool hasAncestor
Definition basis.hpp:1345
static const bool hasHessian
Definition basis.hpp:1117
static const bool hasGradient
Definition basis.hpp:1654
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1826
TensorizedVectorFamily< ScalarBasisType, N >::DivergenceValue ReturnValue
Definition basis.hpp:2153
TensorizedVectorFamily(const ScalarFamilyType &scalar_family)
Definition basis.hpp:586
static const bool hasFunction
Definition basis.hpp:1407
BasisType::CurlValue CurlValue
Definition basis.hpp:1103
Eigen::Matrix< double, dimspace, dimspace > CurlValue
Definition basis.hpp:1525
void RotorValue
Definition basis.hpp:248
void DivergenceValue
Definition basis.hpp:1399
BasisType::RotorValue RotorValue
Definition basis.hpp:1105
VectorRd CurlValue
Definition basis.hpp:172
@ Curl
Definition basis.hpp:1705
@ Function
Definition basis.hpp:1701
@ Gradient
Definition basis.hpp:1702
@ SkewsymmetricGradient
Definition basis.hpp:1704
@ Rotor
Definition basis.hpp:1707
@ Divergence
Definition basis.hpp:1706
@ SymmetricGradient
Definition basis.hpp:1703
@ Hessian
Definition basis.hpp:1708
@ Scalar
Definition basis.hpp:68
@ Vector
Definition basis.hpp:69
@ Matrix
Definition basis.hpp:70
size_t L
Definition HHO_DiffAdvecReac.hpp:46
std::vector< QuadratureNode > QuadratureRule
Definition quadraturerule.hpp:55
if(strcmp(field, 'real')) % real valued entries T
Definition mmread.m:93
for j
Definition mmread.m:174
Definition PastixInterface.hpp:16
Definition ddr-klplate.hpp:27
static auto v
Definition ddrcore-test.hpp:32
Structure to decompose a set of polynomials on a basis on an edge or a cell.
Definition basis.hpp:2958
Structure to store the gradient of a matrix-valued function.
Definition basis.hpp:113
Compute vectors listing the powers of monomial basis functions (for a cell, only this specialization ...
Definition basis.hpp:81
Basis dimensions for various polynomial spaces on edges/faces/elements (when relevant): Pk,...
Definition polynomialspacedimension.hpp:55
Basis evaluation traits. Only specialization of 'BasisFunction' (=Function, Gradient,...
Definition basis.hpp:1793
Evaluate a basis at quadrature nodes. 'BasisFunction' (=Function, Gradient, Curl, Divergence,...
Definition basis.hpp:2320