28#include <boost/multi_array.hpp>
77 template<
typename GeometricSupport>
86 static std::vector<VectorZd>
complete(
size_t degree)
88 std::vector<VectorZd> powers;
90 for (
size_t l = 0; l <= degree; l++)
92 for (
size_t i = 0;
i <= l;
i++)
114 Eigen::Matrix<double, N, N> diff_x,
115 Eigen::Matrix<double, N, N> diff_y
147 Eigen::Matrix<double, N, N>
dx;
148 Eigen::Matrix<double, N, N>
dy;
195 return ( m_degree >= 0 ? (m_degree + 1) * (m_degree + 2) / 2 : 0);
226 return (x - m_xT) / m_hT;
232 std::vector<VectorZd> m_powers;
233 Eigen::Matrix2d m_rot;
285 inline double _coordinate_transform(
const VectorRd &x)
const
287 return (x - m_xE).dot(m_tE) / m_hE;
309 template <
typename BasisType>
326 static const bool hasCurl = BasisType::hasCurl;
335 const BasisType &basis,
336 const Eigen::MatrixXd &
matrix
341 assert((
size_t)
matrix.cols() == basis.dimension() ||
"Inconsistent family initialization");
347 return m_matrix.rows();
353 static_assert(
hasFunction,
"Call to function() not available");
356 for (
auto j = 1;
j < m_matrix.cols();
j++)
358 f += m_matrix(
i,
j) * m_basis.function(
j, x);
366 static_assert(
hasFunction,
"Call to function() not available");
369 for (
auto j = 1;
j < m_matrix.cols();
j++)
371 f += m_matrix(
i,
j) * ancestor_value_quad[
j][iqn];
380 static_assert(
hasGradient,
"Call to gradient() not available");
383 for (
auto j = 1;
j < m_matrix.cols();
j++)
385 G += m_matrix(
i,
j) * m_basis.gradient(
j, x);
393 static_assert(
hasGradient,
"Call to gradient() not available");
395 GradientValue G = m_matrix(
i, 0) * ancestor_gradient_quad[0][iqn];
396 for (
auto j = 1;
j < m_matrix.cols();
j++)
398 G += m_matrix(
i,
j) * ancestor_gradient_quad[
j][iqn];
406 static_assert(
hasCurl,
"Call to curl() not available");
408 CurlValue C = m_matrix(
i, 0) * m_basis.curl(0, x);
409 for (
auto j = 1;
j < m_matrix.cols();
j++)
411 C += m_matrix(
i,
j) * m_basis.curl(
j, x);
417 CurlValue curl(
size_t i,
size_t iqn,
const boost::multi_array<CurlValue, 2> &ancestor_curl_quad)
const
419 static_assert(
hasCurl,
"Call to curl() not available");
421 CurlValue C = m_matrix(
i, 0) * ancestor_curl_quad[0][iqn];
422 for (
auto j = 1;
j < m_matrix.cols();
j++)
424 C += m_matrix(
i,
j) * ancestor_curl_quad[
j][iqn];
432 static_assert(
hasDivergence,
"Call to divergence() not available");
435 for (
auto j = 1;
j < m_matrix.cols();
j++)
437 D += m_matrix(
i,
j) * m_basis.divergence(
j, x);
445 static_assert(
hasDivergence,
"Call to divergence() not available");
448 for (
auto j = 1;
j < m_matrix.cols();
j++)
450 D += m_matrix(
i,
j) * ancestor_divergence_quad[
j][iqn];
458 static_assert(
hasRotor,
"Call to rotor() not available");
460 RotorValue R = m_matrix(
i, 0) * m_basis.rotor(0, x);
461 for (
auto j = 1;
j < m_matrix.cols();
j++)
463 R += m_matrix(
i,
j) * m_basis.rotor(
j, x);
469 RotorValue rotor(
size_t i,
size_t iqn,
const boost::multi_array<RotorValue, 2> &ancestor_rotor_quad)
const
471 static_assert(
hasRotor,
"Call to rotor() not available");
473 RotorValue R = m_matrix(
i, 0) * ancestor_rotor_quad[0][iqn];
474 for (
auto j = 1;
j < m_matrix.cols();
j++)
476 R += m_matrix(
i,
j) * ancestor_rotor_quad[
j][iqn];
484 static_assert(
hasHessian,
"Call to hessian() not available");
487 for (
auto j = 1;
j < m_matrix.cols();
j++)
489 H += m_matrix(
i,
j) * m_basis.hessian(
j, x);
495 HessianValue hessian(
size_t i,
size_t iqn,
const boost::multi_array<HessianValue, 2> &ancestor_hessian_quad)
const
497 static_assert(
hasHessian,
"Call to hessian() not available");
499 HessianValue H = m_matrix(
i, 0) * ancestor_hessian_quad[0][iqn];
500 for (
auto j = 1;
j < m_matrix.cols();
j++)
502 H += m_matrix(
i,
j) * ancestor_hessian_quad[
j][iqn];
508 inline const Eigen::MatrixXd &
matrix()
const
522 return m_basis.max_degree();
527 Eigen::MatrixXd m_matrix;
557 template <
typename ScalarFamilyType,
size_t N>
585 : m_scalar_family(scalar_family)
587 static_assert(ScalarFamilyType::tensorRank ==
Scalar,
588 "Vector family can only be constructed from scalar families");
589 m_rot.row(0) << 0., 1.;
590 m_rot.row(1) << -1., 0.;
597 return m_scalar_family.dimension() * N;
603 static_assert(
hasFunction,
"Call to function() not available");
606 ek(
i / m_scalar_family.dimension()) = 1.;
607 return ek * m_scalar_family.function(
i % m_scalar_family.dimension(), x);
613 static_assert(
hasFunction,
"Call to function() not available");
616 ek(
i / m_scalar_family.dimension()) = 1.;
617 return ek * ancestor_value_quad[
i % m_scalar_family.dimension()][iqn];
623 static_assert(
hasGradient,
"Call to gradient() not available");
625 GradientValue G = Eigen::Matrix<double, N, dimspace>::Zero();
626 G.row(
i / m_scalar_family.dimension()) = m_scalar_family.gradient(
i % m_scalar_family.dimension(), x);
633 static_assert(
hasGradient,
"Call to gradient() not available");
635 GradientValue G = Eigen::Matrix<double, N, dimspace>::Zero();
636 G.row(
i / m_scalar_family.dimension()) = ancestor_gradient_quad[
i % m_scalar_family.dimension()][iqn];
643 static_assert(
hasCurl,
"Call to curl() not available");
649 CurlValue curl(
size_t i,
size_t iqn,
const boost::multi_array<VectorRd, 2> &ancestor_gradient_quad)
const
651 static_assert(
hasCurl,
"Call to curl() not available");
653 return m_rot * ancestor_gradient_quad[
i][iqn];
659 static_assert(
hasDivergence,
"Call to divergence() not available");
661 return m_scalar_family.gradient(
i % m_scalar_family.dimension(), x)(
i / m_scalar_family.dimension());
667 static_assert(
hasDivergence,
"Call to divergence() not available");
669 return ancestor_gradient_quad[
i % m_scalar_family.dimension()][iqn](
i / m_scalar_family.dimension());
677 return curl(
i, x).trace();
681 RotorValue rotor(
size_t i,
size_t iqn,
const boost::multi_array<VectorRd, 2> &ancestor_gradient_quad)
const
683 static_assert(
hasRotor,
"Call to rotor() not available");
685 return curl(
i, iqn, ancestor_gradient_quad).trace();
689 constexpr inline const ScalarFamilyType &
ancestor()
const
691 return m_scalar_family;
697 return m_scalar_family.max_degree();
701 ScalarFamilyType m_scalar_family;
702 Eigen::Matrix2d m_rot;
717 template<
typename ScalarFamilyType,
size_t N>
742 : m_scalar_family(scalar_family),
746 static_assert(ScalarFamilyType::tensorRank ==
Scalar,
747 "Vector family can only be constructed from scalar families");
750 for (
size_t j = 0;
j < N;
j++){
751 for (
size_t i = 0;
i < N;
i++){
752 m_E[
j*N +
i] = Eigen::Matrix<double, N, N>::Zero();
753 m_E[
j*N +
i](
i,
j) = 1.;
760 size_t r = m_scalar_family.dimension();
761 size_t j = ( ( (
i/r)%N )*N + (
i/r)/N ) * r + (
i%r);
762 m_transposeOperator(
i,
j) = 1.;
769 return m_scalar_family.dimension() * N * N;
775 static_assert(
hasFunction,
"Call to function() not available");
777 return m_scalar_family.function(
i % m_scalar_family.dimension(), x) * m_E[
i / m_scalar_family.dimension()];
783 static_assert(
hasFunction,
"Call to function() not available");
785 return ancestor_value_quad[
i % m_scalar_family.dimension()][iqn] * m_E[
i / m_scalar_family.dimension()];
791 static_assert(
hasDivergence,
"Call to divergence() not available");
793 Eigen::Matrix<double, N, 1> V = m_E[
i/m_scalar_family.dimension()].col( (
i/m_scalar_family.dimension()) /N);
794 return m_scalar_family.gradient(
i % m_scalar_family.dimension(), x)( (
i / m_scalar_family.dimension()) / N) * V;
800 static_assert(
hasDivergence,
"Call to divergence() not available");
802 Eigen::Matrix<double, N, 1> V = m_E[
i/m_scalar_family.dimension()].col( (
i/m_scalar_family.dimension()) /N);
803 return ancestor_gradient_quad[
i % m_scalar_family.dimension()][iqn]( (
i / m_scalar_family.dimension()) / N) * V;
809 static_assert(
hasGradient,
"Call to gradient() not available");
811 const VectorRd gradient_scalar_function = m_scalar_family.gradient(
i % m_scalar_family.dimension(), x);
812 const Eigen::Matrix<double, N, N> matrix_basis = m_E[
i / m_scalar_family.dimension()];
813 return GradientValue(gradient_scalar_function.x() * matrix_basis, gradient_scalar_function.y() * matrix_basis);
819 static_assert(
hasGradient,
"Call to gradient() not available");
821 const VectorRd gradient_scalar_function = ancestor_gradient_quad[
i % m_scalar_family.dimension()][iqn];
822 const Eigen::Matrix<double, N, N> matrix_basis = m_E[
i / m_scalar_family.dimension()];
823 return GradientValue(gradient_scalar_function.x() * matrix_basis, gradient_scalar_function.y() * matrix_basis);
829 return m_scalar_family;
841 return m_transposeOperator;
853 size_t dim_scalar = m_scalar_family.dimension();
856 Eigen::MatrixXd SB = Eigen::MatrixXd::Zero(dim_scalar * N*(N+1)/2,
dimension());
862 for (
size_t i = 0;
i<N;
i++){
864 SB.middleRows(position, (N-
i)*dim_scalar) = symOpe.middleRows(
i*(N+1)*dim_scalar, (N-
i)*dim_scalar);
866 position += (N-
i)*dim_scalar;
873 ScalarFamilyType m_scalar_family;
874 std::vector<Eigen::Matrix<double, N, N>> m_E;
875 Eigen::MatrixXd m_transposeOperator;
882 template <
typename ScalarFamilyType>
908 const ScalarFamilyType &scalar_family,
911 : m_scalar_family(scalar_family),
914 static_assert(ScalarFamilyType::hasFunction,
"Call to function() not available");
915 static_assert(std::is_same<typename ScalarFamilyType::GeometricSupport, Edge>::value,
916 "Tangent families can only be defined on edges");
922 return m_scalar_family.dimension();
928 return m_scalar_family.function(
i, x) * m_generator;
932 constexpr inline const ScalarFamilyType &
ancestor()
const
934 return m_scalar_family;
940 return m_scalar_family.max_degree();
950 ScalarFamilyType m_scalar_family;
959 template <
typename BasisType>
976 static const bool hasCurl = BasisType::hasCurl;
985 const BasisType &basis,
997 return m_basis.dimension() - m_shift;
1009 return m_basis.max_degree();
1015 static_assert(
hasFunction,
"Call to function() not available");
1017 return m_basis.function(
i + m_shift, x);
1023 static_assert(
hasGradient,
"Call to gradient() not available");
1025 return m_basis.gradient(
i + m_shift, x);
1031 static_assert(
hasCurl,
"Call to curl() not available");
1033 return m_basis.curl(
i + m_shift, x);
1039 static_assert(
hasDivergence,
"Call to divergence() not available");
1041 return m_basis.divergence(
i + m_shift, x);
1047 static_assert(
hasRotor,
"Call to rotor() not available");
1049 return m_basis.rotor(
i + m_shift, x);
1055 static_assert(
hasHessian,
"Call to hessian() not available");
1057 return m_basis.hessian(
i + m_shift, x);
1069 template <
typename BasisType>
1095 const BasisType &basis,
1132 static_assert(
hasFunction,
"Call to function() not available");
1134 return m_basis.function(
i, x);
1140 static_assert(
hasGradient,
"Call to gradient() not available");
1142 return m_basis.gradient(
i, x);
1148 static_assert(
hasCurl,
"Call to curl() not available");
1150 return m_basis.curl(
i, x);
1156 static_assert(
hasDivergence,
"Call to divergence() not available");
1158 return m_basis.divergence(
i, x);
1164 static_assert(
hasRotor,
"Call to rotor() not available");
1166 return m_basis.rotor(
i, x);
1179 template <
typename BasisType>
1207 : m_scalar_basis(basis)
1209 static_assert(BasisType::hasGradient,
1210 "Gradient basis requires gradient() for the original basis to be available");
1217 return m_scalar_basis.dimension();
1223 return m_scalar_basis.gradient(
i, x);
1229 return m_scalar_basis;
1234 BasisType m_scalar_basis;
1241 template <
typename BasisType>
1247 typedef Eigen::Matrix<double, dimspace, dimspace>
CurlValue;
1269 static_assert(BasisType::tensorRank ==
Scalar && std::is_same<typename BasisType::GeometricSupport, Cell>::value,
1270 "Curl basis can only be constructed starting from scalar bases on elements");
1271 static_assert(BasisType::hasCurl,
1272 "Curl basis requires curl() for the original basis to be available");
1278 return m_basis.dimension();
1284 return m_basis.curl(
i, x);
1304 template <
typename BasisType>
1328 : m_vector_basis(basis)
1330 static_assert(BasisType::tensorRank ==
Vector || BasisType::tensorRank ==
Matrix,
1331 "Divergence basis can only be constructed starting from vector bases");
1332 static_assert(BasisType::hasDivergence,
1333 "Divergence basis requires divergence() for the original basis to be available");
1340 return m_vector_basis.dimension();
1346 return m_vector_basis.divergence(
i, x);
1352 return m_vector_basis;
1356 BasisType m_vector_basis;
1364 template <
typename BasisType>
1390 : m_vector_basis(basis)
1392 static_assert(BasisType::tensorRank ==
Vector,
1393 "Rotor basis can only be constructed starting from vector bases");
1394 static_assert(BasisType::hasRotor,
1395 "Rotor basis requires rotor() for the original basis to be available");
1402 return m_vector_basis.dimension();
1408 return m_vector_basis.rotor(
i, x);
1414 return m_vector_basis;
1418 BasisType m_vector_basis;
1427 template <
typename BasisType>
1456 : m_scalar_basis(basis)
1458 static_assert(BasisType::tensorRank ==
Scalar,
1459 "Hessian basis can only be constructed starting from scalar bases");
1460 static_assert(BasisType::hasHessian,
1461 "Hessian basis requires hessian() for the original basis to be available");
1468 return m_scalar_basis.dimension();
1474 return m_scalar_basis.hessian(
i, x);
1480 return m_scalar_basis;
1485 BasisType m_scalar_basis;
1497 typedef Eigen::Matrix<double, dimspace, dimspace>
CurlValue;
1547 return (x - m_xT) / m_hT;
1553 std::vector<VectorZd> m_powers;
1563 typedef Eigen::Matrix<double, dimspace, dimspace>
CurlValue;
1598 inline const std::shared_ptr<RolyComplBasisCell> &
rck()
const
1605 Eigen::Matrix2d m_rot;
1606 std::shared_ptr<RolyComplBasisCell> m_Rck_basis;
1657 return (x - m_xT) / m_hT;
1663 Eigen::Matrix2d m_rot;
1664 TensorizedVectorFamily<MonomialScalarBasisCell, dimspace> m_Pkmo;
1684 template <
typename outValue,
typename inValue,
typename FunctionType>
1686 const boost::multi_array<inValue, 2> &B_quad,
1687 const FunctionType &F
1690 boost::multi_array<outValue, 2> transformed_B_quad( boost::extents[B_quad.shape()[0]][B_quad.shape()[1]] );
1692 std::transform( B_quad.origin(), B_quad.origin() + B_quad.num_elements(), transformed_B_quad.origin(), F);
1694 return transformed_B_quad;
1699 template <
typename ScalarBasisType,
size_t N>
1701 const ScalarBasisType & B,
1702 const std::vector<Eigen::VectorXd> &
v
1706 size_t k =
v.size();
1708 Eigen::MatrixXd M = Eigen::MatrixXd::Zero(r*k, r*N);
1710 for (
size_t i=0;
i < k;
i++){
1712 assert(
v[
i].size() == N);
1715 for (
size_t j=0;
j < N;
j++){
1716 M.block(
i*r,
j*r, r, r) =
v[
i](
j) * Eigen::MatrixXd::Identity(r,r);
1725 inline static std::function<Eigen::MatrixXd(
const Eigen::MatrixXd &)>
1726 symmetrise_matrix = [](
const Eigen::MatrixXd & x)->Eigen::MatrixXd {
return 0.5*(x+x.transpose());};
1729 inline static std::function<Eigen::MatrixXd(
const Eigen::MatrixXd &)>
1734 template <
typename ScalarBasisType,
size_t N>
1736 const ScalarBasisType & B
1740 Eigen::MatrixXd M = Eigen::MatrixXd::Zero(r, r*N*N);
1742 for (
size_t k=0; k < N; k++){
1743 M.block(0, k*r*(N+1), r, r) = Eigen::MatrixXd::Identity(r, r);
1756 template<
typename BasisType, BasisFunctionE BasisFunction>
1761 template<
typename BasisType>
1764 static_assert(BasisType::hasFunction,
"Call to function not available");
1768 return basis.function(
i, x);
1773 const BasisType &basis,
1776 const boost::multi_array<ReturnValue, 2> &ancestor_value_quad
1779 return basis.function(
i, iqn, ancestor_value_quad);
1785 template<
typename BasisType>
1788 static_assert(BasisType::hasGradient,
"Call to gradient not available");
1792 return basis.gradient(
i, x);
1797 const BasisType &basis,
1800 const boost::multi_array<ReturnValue, 2> &ancestor_gradient_quad
1803 return basis.gradient(
i, iqn, ancestor_gradient_quad);
1809 template<
typename BasisType>
1812 static_assert(BasisType::hasGradient,
"Call to gradient not available");
1816 return 0.5 * (basis.gradient(
i, x) + basis.gradient(
i, x).transpose());
1821 const BasisType &basis,
1824 const boost::multi_array<ReturnValue, 2> &ancestor_gradient_quad
1827 return 0.5 * (basis.gradient(
i, iqn, ancestor_gradient_quad) + basis.gradient(
i, iqn, ancestor_gradient_quad).transpose());
1833 template<
typename BasisType>
1836 static_assert(BasisType::hasGradient,
"Call to gradient not available");
1840 return 0.5 * (basis.gradient(
i, x) - basis.gradient(
i, x).transpose());
1845 const BasisType &basis,
1848 const boost::multi_array<ReturnValue, 2> &ancestor_gradient_quad
1851 return 0.5 * (basis.gradient(
i, iqn, ancestor_gradient_quad) - basis.gradient(
i, iqn, ancestor_gradient_quad).transpose());
1856 template<
typename BasisType>
1859 static_assert(BasisType::hasCurl,
"Call to curl not available");
1863 return basis.curl(
i, x);
1868 const BasisType &basis,
1871 const boost::multi_array<ReturnValue, 2> &ancestor_curl_quad
1874 return basis.curl(
i, iqn, ancestor_curl_quad);
1879 template<
typename BasisType>
1882 static_assert(BasisType::hasDivergence,
"Call to divergence not available");
1886 return basis.divergence(
i, x);
1891 const BasisType &basis,
1894 const boost::multi_array<ReturnValue, 2> &ancestor_divergence_quad
1897 return basis.divergence(
i, iqn, ancestor_divergence_quad);
1902 template<
typename BasisType>
1905 static_assert(BasisType::hasRotor,
"Call to rotor not available");
1909 return basis.rotor(
i, x);
1914 const BasisType &basis,
1917 const boost::multi_array<ReturnValue, 2> &ancestor_rotor_quad
1920 return basis.rotor(
i, iqn, ancestor_rotor_quad);
1925 template<
typename BasisType>
1928 static_assert(BasisType::hasHessian,
"Call to Hessian not available");
1932 return basis.hessian(
i, x);
1937 const BasisType &basis,
1940 const boost::multi_array<ReturnValue, 2> &ancestor_hessian_quad
1943 return basis.hessian(
i, iqn, ancestor_hessian_quad);
1949 template <
typename ScalarBasisType,
size_t N>
1973 const boost::multi_array<double, 2> &ancestor_basis_quad
1976 return basis.
function(
i, iqn, ancestor_basis_quad);
1981 template <
typename ScalarBasisType,
size_t N>
2005 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2008 return basis.
gradient(
i, iqn, ancestor_basis_quad);
2014 template <
typename ScalarBasisType,
size_t N>
2038 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2041 return 0.5 * (basis.
gradient(
i, iqn, ancestor_basis_quad) + basis.
gradient(
i, iqn, ancestor_basis_quad).transpose());
2047 template <
typename ScalarBasisType,
size_t N>
2071 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2074 return 0.5 * (basis.
gradient(
i, iqn, ancestor_basis_quad) - basis.
gradient(
i, iqn, ancestor_basis_quad).transpose());
2079 template <
typename ScalarBasisType,
size_t N>
2095 return basis.
curl(
i, x);
2103 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2106 return basis.
curl(
i, iqn, ancestor_basis_quad);
2112 template <
typename ScalarBasisType,
size_t N>
2136 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2139 return basis.
divergence(
i, iqn, ancestor_basis_quad);
2145 template <
typename ScalarBasisType,
size_t N>
2161 return basis.
rotor(
i, x);
2169 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2172 return basis.
rotor(
i, iqn, ancestor_basis_quad);
2179 template <
typename ScalarBasisType,
size_t N>
2203 const boost::multi_array<double, 2> &ancestor_basis_quad
2206 return basis.
function(
i, iqn, ancestor_basis_quad);
2211 template <
typename ScalarBasisType,
size_t N>
2235 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2238 return basis.
divergence(
i, iqn, ancestor_basis_quad);
2245 template <
typename ScalarBasisType,
size_t N>
2269 const boost::multi_array<VectorRd, 2> &ancestor_basis_quad
2272 return basis.
gradient(
i, iqn, ancestor_basis_quad);
2283 template<BasisFunctionE BasisFunction>
2286 template<
typename BasisType>
2287 static boost::multi_array<typename detail::basis_evaluation_traits<BasisType, BasisFunction>::ReturnValue, 2>
2289 const BasisType & basis,
2295 boost::multi_array<typename traits::ReturnValue, 2>
2296 basis_quad( boost::extents[basis.dimension()][quad.size()] );
2298 for (
size_t i = 0;
i < basis.dimension();
i++) {
2299 for (
size_t iqn = 0; iqn < quad.size(); iqn++) {
2300 basis_quad[
i][iqn] = traits::evaluate(basis,
i, quad[iqn].vector());
2308 template<
typename BasisType>
2309 static boost::multi_array<typename detail::basis_evaluation_traits<Family<BasisType>, BasisFunction>::ReturnValue, 2>
2317 boost::multi_array<typename traits::ReturnValue, 2>
2318 basis_quad( boost::extents[basis.
dimension()][quad.size()] );
2321 const auto & ancestor_basis = basis.
ancestor();
2322 boost::multi_array<typename traits::ReturnValue, 2>
2325 Eigen::MatrixXd M = basis.
matrix();
2326 for (
size_t i = 0;
i < size_t(M.rows());
i++) {
2327 for (
size_t iqn = 0; iqn < quad.size(); iqn++) {
2328 basis_quad[
i][iqn] = M(
i,0) * ancestor_basis_quad[0][iqn];
2330 for (
size_t j = 1;
j < size_t(M.cols());
j++) {
2331 for (
size_t iqn = 0; iqn < quad.size(); iqn++) {
2332 basis_quad[
i][iqn] += M(
i,
j) * ancestor_basis_quad[
j][iqn];
2341 template <
typename BasisType,
size_t N>
2342 static boost::multi_array<typename detail::basis_evaluation_traits<TensorizedVectorFamily<BasisType, N>, BasisFunction>::ReturnValue, 2>
2350 boost::multi_array<typename traits::ReturnValue, 2>
2351 basis_quad(boost::extents[basis.
dimension()][quad.size()]);
2353 const auto &ancestor_basis = basis.
ancestor();
2354 const boost::multi_array<typename traits::AncestorBasisFunctionValue, 2> ancestor_basis_quad
2358 for (
size_t iqn = 0; iqn < quad.size(); iqn++){
2359 basis_quad[
i][iqn] = traits::evaluate(basis,
i, iqn, ancestor_basis_quad);
2366 template <
typename BasisType,
size_t N>
2367 static boost::multi_array<typename detail::basis_evaluation_traits<MatrixFamily<BasisType, N>, BasisFunction>::ReturnValue, 2>
2375 boost::multi_array<typename traits::ReturnValue, 2>
2376 basis_quad(boost::extents[basis.
dimension()][quad.size()]);
2378 const auto &ancestor_basis = basis.
ancestor();
2379 const boost::multi_array<typename traits::AncestorBasisFunctionValue, 2> ancestor_basis_quad
2383 for (
size_t iqn = 0; iqn < quad.size(); iqn++){
2384 basis_quad[
i][iqn] = traits::evaluate(basis,
i, iqn, ancestor_basis_quad);
2403 template<
typename T>
2405 boost::multi_array<T, 2> & basis_eval,
2406 const std::function<
double(
size_t,
size_t)> & inner_product
2409 auto induced_norm = [inner_product](
size_t i)->
double
2411 return std::sqrt( inner_product(
i,
i) );
2415 size_t Nb = basis_eval.shape()[0];
2417 size_t Nqn = basis_eval.shape()[1];
2419 Eigen::MatrixXd B = Eigen::MatrixXd::Zero(Nb, Nb);
2422 double norm = induced_norm(0);
2423 for (
size_t iqn = 0; iqn < Nqn; iqn++) {
2424 basis_eval[0][iqn] /= norm;
2428 for (
size_t ib = 1; ib < Nb; ib++) {
2430 Eigen::RowVectorXd coeffs = Eigen::RowVectorXd::Zero(ib);
2431 for (
size_t pb = 0; pb < ib; pb++) {
2432 coeffs(pb) = - inner_product(ib, pb);
2436 for (
size_t pb = 0; pb < ib; pb++) {
2437 for (
size_t iqn = 0; iqn < Nqn; iqn++) {
2438 basis_eval[ib][iqn] += coeffs(pb) * basis_eval[pb][iqn];
2443 double norm = induced_norm(ib);
2444 for (
size_t iqn = 0; iqn < Nqn; iqn++) {
2445 basis_eval[ib][iqn] /= norm;
2452 B.block(ib, 0, 1, ib) = coeffs * B.topLeftCorner(ib, ib);
2453 B(ib, ib) = 1./norm;
2468 double scalar_product(
const Eigen::Matrix<double, N, N> & x,
const Eigen::Matrix<double, N, N> & y)
2470 return (x.transpose() * y).trace();
2474 template <
typename Value>
2476 const boost::multi_array<Value, 2> &basis_quad,
2480 boost::multi_array<double, 2> basis_dot_v_quad(boost::extents[basis_quad.shape()[0]][basis_quad.shape()[1]]);
2481 std::transform(basis_quad.origin(), basis_quad.origin() + basis_quad.num_elements(),
2482 basis_dot_v_quad.origin(), [&
v](
const Value &x) ->
double { return scalar_product(x,v); });
2483 return basis_dot_v_quad;
2487 template <
typename Value>
2489 const boost::multi_array<Value, 2> & basis_quad,
2490 const std::vector<Value> &
v
2493 boost::multi_array<double, 2> basis_dot_v_quad(boost::extents[basis_quad.shape()[0]][basis_quad.shape()[1]]);
2494 for (std::size_t
i = 0;
i < basis_quad.shape()[0];
i++) {
2495 for (std::size_t iqn = 0; iqn < basis_quad.shape()[1]; iqn++) {
2499 return basis_dot_v_quad;
2504 const boost::multi_array<MatrixRd, 2> & basis_quad,
2505 const std::vector<VectorRd> & v_quad
2510 const std::vector<MatrixRd> & m_quad,
2511 const boost::multi_array<VectorRd, 2> & basis_quad
2516 const std::vector<VectorRd> & v_quad,
2517 const boost::multi_array<MatrixRd, 2> & basis_quad
2521 template<
typename BasisType>
2523 const BasisType & basis,
2525 boost::multi_array<typename BasisType::FunctionValue, 2> & basis_quad
2529 assert ( basis.dimension() == basis_quad.shape()[0] && qr.size() == basis_quad.shape()[1] );
2532 std::function<double(
size_t,
size_t)> inner_product
2533 = [&basis_quad, &qr](
size_t i,
size_t j)->
double
2536 for(
size_t iqn = 0; iqn < qr.size(); iqn++) {
2542 Eigen::MatrixXd B =
gram_schmidt(basis_quad, inner_product);
2549 template <
typename BasisType>
2551 const BasisType &basis,
2552 const Eigen::MatrixXd & GM
2556 assert(basis.dimension() ==
size_t(GM.rows()) && GM.rows() == GM.cols());
2558 Eigen::MatrixXd
L = GM.llt().matrixL();
2569 template<
typename FunctionValue>
2571 const boost::multi_array<FunctionValue, 2> & B2,
2575 const std::string sym =
"nonsym"
2579 assert ( qr.size() == B1.shape()[1] && qr.size() == B2.shape()[1] );
2581 assert ( nrows <= B1.shape()[0] && ncols <= B2.shape()[0] );
2584 Eigen::ArrayXd qr_weights = Eigen::ArrayXd::Zero(qr.size());
2585 for (
size_t iqn = 0; iqn < qr.size(); iqn++){
2586 qr_weights(iqn) = qr[iqn].w;
2589 Eigen::MatrixXd M = Eigen::MatrixXd::Zero(nrows, ncols);
2590 for (
size_t i = 0;
i < nrows;
i++) {
2592 if ( sym==
"sym" ) jcut =
i;
2593 for (
size_t j = 0;
j < jcut;
j++){
2596 for(
size_t j = jcut;
j < ncols;
j++) {
2597 std::vector<double> tmp(B1.shape()[1]);
2599 auto B1i = B1[ boost::indices[
i][boost::multi_array_types::index_range(0, B1.shape()[1])] ];
2600 auto B2j = B2[ boost::indices[
j][boost::multi_array_types::index_range(0, B1.shape()[1])] ];
2602 std::transform(B1i.begin(), B1i.end(), B2j.begin(), tmp.begin(), [](FunctionValue a, FunctionValue b)->double { return scalar_product(a,b);});
2603 Eigen::ArrayXd tmp_array = Eigen::Map<Eigen::ArrayXd, Eigen::Unaligned>(tmp.data(), tmp.size());
2605 M(
i,
j) = (qr_weights * tmp_array).sum();
2621 template<
typename FunctionValue>
2623 const boost::multi_array<FunctionValue, 2> & B2,
2625 const std::string sym =
"nonsym"
2628 return compute_gram_matrix<FunctionValue>(B1, B2, qr, B1.shape()[0], B2.shape()[0], sym);
2632 template<
typename FunctionValue>
2637 return compute_gram_matrix<FunctionValue>(B, B, qr,
"sym");
2642 const boost::multi_array<double, 2> & B2,
2648 const boost::multi_array<double, 2> & B2,
2652 const std::string sym =
"nonsym"
2657 const boost::multi_array<double, 2> & B2,
2659 const std::string sym =
"nonsym"
2664 const boost::multi_array<VectorRd, 2> & B2,
2668 const std::string sym =
"nonsym"
2673 const boost::multi_array<VectorRd, 2> & B2,
2675 const std::string sym =
"nonsym"
2679 template<
typename ScalarFamilyType,
size_t N>
2681 const boost::multi_array<double, 2> & scalar_family_quad,
2685 Eigen::MatrixXd Gram = Eigen::MatrixXd::Zero(MatFam.
dimension(), MatFam.
dimension());
2688 Eigen::MatrixXd scalarGram = compute_gram_matrix<double>(scalar_family_quad, qr);
2689 for (
size_t i = 0;
i < MatFam.
dimension();
i =
i+scalarGram.rows()){
2690 Gram.block(
i,
i, scalarGram.rows(), scalarGram.cols()) = scalarGram;
2699 template <
typename T>
2710 n_rows = B.shape()[0];
2714 const size_t num_quads = qr.size();
2717 assert(num_quads == B.shape()[1]);
2720 assert(n_rows <= B.shape()[0]);
2722 Eigen::VectorXd V = Eigen::VectorXd::Zero(n_rows);
2724 for (
size_t iqn = 0; iqn < num_quads; iqn++)
2726 double qr_weight = qr[iqn].w;
2727 T f_on_qr = f(qr[iqn].vector());
2728 for (
size_t i = 0;
i < n_rows;
i++)
2738 template <
typename T,
typename U>
2746 const std::string sym =
"nonsym"
2750 if (n_rows == 0 && n_cols == 0)
2752 n_rows = B1.shape()[0];
2753 n_cols = B2.shape()[0];
2757 const size_t num_quads = qr.size();
2759 assert(num_quads == B1.shape()[1] && num_quads == B2.shape()[1]);
2761 assert(n_rows <= B1.shape()[0] && n_cols <= B2.shape()[0]);
2763 Eigen::MatrixXd M = Eigen::MatrixXd::Zero(n_rows, n_cols);
2764 for (
size_t iqn = 0; iqn < num_quads; iqn++)
2766 double qr_weight = qr[iqn].w;
2767 U f_on_qr = f(qr[iqn].vector());
2768 for (
size_t i = 0;
i < n_rows;
i++)
2770 T f_B1 = f_on_qr * B1[
i][iqn];
2774 for (
size_t j = 0;
j < jcut;
j++)
2778 for (
size_t j = jcut;
j < n_cols;
j++)
2788 template <
typename T,
typename U>
2794 const std::string sym
2802 const FType<VectorRd> &f,
2803 const BasisQuad<VectorRd> &B1,
2804 const BasisQuad<double> &B2,
2813 const FType<VectorRd> &f,
2814 const BasisQuad<double> &B1,
2815 const BasisQuad<VectorRd> &B2,
2822 template <
typename Value>
2824 const boost::multi_array<Value, 2> & f_quad,
2825 const boost::multi_array<Value, 2> & g_quad,
2830 const size_t dimf = f_quad.shape()[0];
2831 const size_t dimg = g_quad.shape()[0];
2832 Eigen::MatrixXd M = Eigen::MatrixXd::Zero(dimf, dimg);
2834 std::vector<Value> int_f(dimf);
2835 for (
size_t i = 0;
i < dimf;
i++){
2836 int_f[
i] = qr_f[0].w * f_quad[
i][0];
2837 for (
size_t iqn = 1; iqn < qr_f.size(); iqn++){
2838 int_f[
i] += qr_f[iqn].w * f_quad[
i][iqn];
2841 std::vector<Value> int_g(dimg);
2842 for (
size_t j = 0;
j < dimg;
j++){
2843 int_g[
j] = qr_g[0].w * g_quad[
j][0];
2844 for (
size_t iqn = 1; iqn < qr_g.size(); iqn++){
2845 int_g[
j] += qr_g[iqn].w * g_quad[
j][iqn];
2850 for (
size_t i = 0;
i < dimf;
i++){
2851 for (
size_t j = 0;
j < dimg;
j++){
2860 template <
typename Value>
2862 const boost::multi_array<Value, 2> & f_quad,
2863 const boost::multi_array<Value, 2> & g_quad,
2871 template <
typename Value>
2873 const boost::multi_array<Value, 2> & f_quad,
2885 template<
typename BasisType>
2887 const std::function<
typename BasisType::FunctionValue(
const VectorRd &)> & f,
2888 const BasisType & basis,
2890 const boost::multi_array<typename BasisType::FunctionValue, 2> & basis_quad,
2891 const Eigen::MatrixXd &mass_basis = Eigen::MatrixXd::Zero(1,1)
2894 Eigen::MatrixXd Mass = mass_basis;
2895 if (Mass.norm() < 1e-13){
2899 Eigen::LDLT<Eigen::MatrixXd> cholesky_mass(Mass);
2900 Eigen::VectorXd b = Eigen::VectorXd::Zero(basis.dimension());
2901 for (
size_t i = 0;
i < basis.dimension();
i++) {
2902 for (
size_t iqn = 0; iqn < quad.size(); iqn++) {
2903 b(
i) += quad[iqn].w *
scalar_product(f(quad[iqn].vector()), basis_quad[
i][iqn]);
2907 return cholesky_mass.solve(b);
2920 template <
typename BasisType>
2931 const BasisType &basis
2933 m_dim(basis.dimension()),
2946 double h = 1./std::max(1.,
double(basis.max_degree()));
2948 VectorRd xi = E.vertex(0)->coords() + double(
i) * h * (E.vertex(1)->coords() - E.vertex(0)->coords());
2949 m_nodes.emplace_back(xi.x(), xi.y(), 1.);
2961 const BasisType &basis
2963 m_dim(basis.dimension()),
2983 VectorRd xi = x0 +
T.diam()*(double(indices[
i](0)) * e1 + double(indices[
i](1)) * e2)/std::max(1.,
double(basis.max_degree()));
2984 m_nodes.emplace_back(xi.x(), xi.y(), 1.);
3001 boost::multi_array<typename BasisType::FunctionValue, 2> &values
3009 Eigen::MatrixXd
L = Gram_onbasis.ldlt().solve(Gram_P_onbasis.transpose());
Basis for the space of curls (vectorial rot) of polynomials.
Definition basis.hpp:1243
Basis (or rather family) of divergence of an existing basis.
Definition basis.hpp:1306
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:1559
Basis for the space of gradients of polynomials.
Definition basis.hpp:1181
Basis for the space of Hessians of polynomials.
Definition basis.hpp:1429
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:1611
Matrix family obtained from a scalar family.
Definition basis.hpp:719
Scalar monomial basis on a cell.
Definition basis.hpp:166
Scalar monomial basis on an edge.
Definition basis.hpp:240
Generate a basis restricted to the first "dimension" functions.
Definition basis.hpp:1071
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:1493
Basis (or rather family) of scalar rotor of an existing basis.
Definition basis.hpp:1366
Generate a basis where the function indices are shifted.
Definition basis.hpp:961
Vector family for polynomial functions that are tangent to an edge (determined by the generator)
Definition basis.hpp:884
Vector family obtained by tensorization of a scalar family.
Definition basis.hpp:559
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:1991
void RotorValue
Definition basis.hpp:1436
static constexpr const TensorRankE tensorRank
Definition basis.hpp:570
VectorRd CurlValue
Definition basis.hpp:1370
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1570
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1215
TensorizedVectorFamily< ScalarBasisType, N >::GradientValue ReturnValue
Definition basis.hpp:2019
static const bool hasRotor
Definition basis.hpp:1510
Eigen::Matrix< double, N, dimspace > GradientValue
Definition basis.hpp:562
static const bool hasGradient
Definition basis.hpp:1507
MatrixGradient(Eigen::Matrix< double, N, N > diff_x, Eigen::Matrix< double, N, N > diff_y)
Constructor.
Definition basis.hpp:113
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th function at point x.
Definition basis.hpp:378
BasisType::GradientValue GradientValue
Definition basis.hpp:314
Cell GeometricSupport
Definition basis.hpp:1621
static const bool hasFunction
Definition basis.hpp:253
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:279
BasisType AncestorType
Definition basis.hpp:1386
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:2343
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:2034
static const bool hasDivergence
Definition basis.hpp:1087
MatrixFamily< ScalarBasisType, N >::GradientValue ReturnValue
Definition basis.hpp:2250
VectorRd FunctionValue
Definition basis.hpp:1495
static const bool hasFunction
Definition basis.hpp:974
static constexpr const bool hasAncestor
Definition basis.hpp:1624
size_t dimension() const
Return the dimension of the family.
Definition basis.hpp:920
double DivergenceValue
Definition basis.hpp:1311
VectorRd CurlValue
Definition basis.hpp:244
BasisType::RotorValue ReturnValue
Definition basis.hpp:1906
HessianValue hessian(size_t i, const VectorRd &x) const
Evaluate the Hessian of the i-th basis function at point x.
Definition basis.hpp:482
static const bool hasGradient
Definition basis.hpp:573
static const bool hasFunction
Definition basis.hpp:572
const size_t matrixSize() const
Return the dimension of the matrices in the family.
Definition basis.hpp:833
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:665
BasisType m_basis
Definition basis.hpp:3016
static constexpr const bool hasAncestor
Definition basis.hpp:896
BasisType::FunctionValue FunctionValue
Definition basis.hpp:313
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1504
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1466
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th function at point x.
Definition basis.hpp:351
static std::vector< VectorZd > complete(size_t degree)
Definition basis.hpp:86
static ReturnValue evaluate(const MatrixFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2221
BasisType::FunctionValue ReturnValue
Definition basis.hpp:1765
static constexpr const bool hasAncestor
Definition basis.hpp:731
static const bool hasHessian
Definition basis.hpp:1511
CurlValue curl(size_t i, const VectorRd &x) const
Evaluate the curl of the i-th basis function at point x.
Definition basis.hpp:641
static const bool hasDivergence
Definition basis.hpp:182
size_t dimension() const
Return the dimension of the family.
Definition basis.hpp:767
void RotorValue
Definition basis.hpp:1499
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:1252
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:2823
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1130
CurlBasis(const BasisType &basis)
Constructor.
Definition basis.hpp:1266
void HessianValue
Definition basis.hpp:566
BasisType::HessianValue HessianValue
Definition basis.hpp:318
MatrixGradient()
Default constructor.
Definition basis.hpp:123
VectorRd CurlValue
Definition basis.hpp:1310
BasisFunctionE
Definition basis.hpp:1672
double FunctionValue
Definition basis.hpp:1308
void RotorValue
Definition basis.hpp:172
static const bool hasHessian
Definition basis.hpp:1630
static const bool hasCurl
Definition basis.hpp:1447
static const bool hasFunction
Definition basis.hpp:1318
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_hessian_quad)
Definition basis.hpp:1936
static const bool hasDivergence
Definition basis.hpp:900
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:320
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_gradient_quad)
Definition basis.hpp:1844
static const bool hasRotor
Definition basis.hpp:1449
static constexpr const TensorRankE tensorRank
Definition basis.hpp:895
DivergenceBasis(const BasisType &basis)
Constructor.
Definition basis.hpp:1327
static constexpr const bool hasAncestor
Definition basis.hpp:1378
static constexpr const bool hasAncestor
Definition basis.hpp:1083
VectorRd CurlValue
Definition basis.hpp:888
static const bool hasGradient
Definition basis.hpp:1257
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_divergence_quad)
Definition basis.hpp:1890
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1377
void HessianValue
Definition basis.hpp:891
Eigen::Vector2d VectorRd
Definition basis.hpp:55
Eigen::Vector2i VectorZd
Definition basis.hpp:56
BasisType::RotorValue RotorValue
Definition basis.hpp:967
void DivergenceValue
Definition basis.hpp:1564
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:1375
Eigen::Matrix< double, dimspace, dimspace > CurlValue
Definition basis.hpp:1247
static const bool hasHessian
Definition basis.hpp:902
static const bool hasGradient
Definition basis.hpp:898
Eigen::Matrix< double, dimspace, dimspace > GradientValue
Definition basis.hpp:1246
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:316
Family(const BasisType &basis, const Eigen::MatrixXd &matrix)
Constructor.
Definition basis.hpp:334
TensorizedVectorFamily< ScalarBasisType, N >::GradientValue ReturnValue
Definition basis.hpp:1986
BasisType::GradientValue GradientValue
Definition basis.hpp:1074
static ReturnValue evaluate(const MatrixFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2255
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1838
static const bool hasRotor
Definition basis.hpp:1088
static const bool hasRotor
Definition basis.hpp:735
BasisType AncestorType
Definition basis.hpp:1452
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:2310
void GradientValue
Definition basis.hpp:1615
static const bool hasFunction
Definition basis.hpp:1625
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_curl_quad)
Definition basis.hpp:1867
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:689
static const bool hasRotor
Definition basis.hpp:1260
BasisType AncestorType
Definition basis.hpp:331
static const bool hasDivergence
Definition basis.hpp:256
static const bool hasCurl
Definition basis.hpp:255
void HessianValue
Definition basis.hpp:1188
static const bool hasCurl
Definition basis.hpp:899
size_t dimension() const
Dimension of the family. This is actually the number of functions in the family, not necessarily line...
Definition basis.hpp:345
static const bool hasDivergence
Definition basis.hpp:1448
constexpr const BasisType & ancestor() const
Return the ancestor (basis that the gradient was taken of)
Definition basis.hpp:1227
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:2199
static const bool hasCurl
Definition basis.hpp:1574
void GradientValue
Definition basis.hpp:1433
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:681
void HessianValue
Definition basis.hpp:247
RestrictedBasis(const BasisType &basis, const size_t &dimension)
Constructor.
Definition basis.hpp:1094
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2155
BasisType::CurlValue CurlValue
Definition basis.hpp:965
static const bool hasFunction
Definition basis.hpp:1194
BasisType::HessianValue HessianValue
Definition basis.hpp:1078
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1282
Edge GeometricSupport
Definition basis.hpp:893
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2024
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:845
GradientBasis(const BasisType &basis)
Constructor.
Definition basis.hpp:1206
boost::multi_array< T, 2 > BasisQuad
type for a family of basis functions evaluated on quadrature nodes
Definition basis.hpp:59
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:798
VectorRd generator() const
Returns the generator of the basis functions.
Definition basis.hpp:944
static const bool hasCurl
Definition basis.hpp:577
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1254
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th basis function at point x.
Definition basis.hpp:621
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:327
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1400
TensorizedVectorFamily< ScalarBasisType, N >::FunctionValue ReturnValue
Definition basis.hpp:1954
void DivergenceValue
Definition basis.hpp:245
static const bool hasHessian
Definition basis.hpp:1322
VectorRd FunctionValue
Definition basis.hpp:1561
static const bool hasHessian
Definition basis.hpp:1201
Cell GeometricSupport
Definition basis.hpp:1568
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:839
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th basis function at point x.
Definition basis.hpp:807
MatrixGradient & operator+=(const MatrixGradient &G)
Increment.
Definition basis.hpp:135
void RotorValue
Definition basis.hpp:890
RotorBasis(const BasisType &basis)
Constructor.
Definition basis.hpp:1389
VectorRd GradientValue
Definition basis.hpp:1369
ScalarFamilyType AncestorType
Definition basis.hpp:739
size_t dimension() const
Return the dimension of the basis.
Definition basis.hpp:1107
VectorRd AncestorBasisFunctionValue
Definition basis.hpp:2252
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:631
Cell GeometricSupport
Definition basis.hpp:175
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:391
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:601
MatrixGradient< N > GradientValue
Definition basis.hpp:722
double scalar_product(const double &x, const double &y)
Scalar product between two reals.
Definition basis.cpp:163
void HessianValue
Definition basis.hpp:1437
size_t dimension() const
Dimension of the basis.
Definition basis.hpp:1586
Eigen::Matrix< double, N, 1 > FunctionValue
Definition basis.hpp:561
void DivergenceValue
Definition basis.hpp:1617
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:611
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:225
void DivergenceValue
Definition basis.hpp:171
BasisType::FunctionValue FunctionValue
Definition basis.hpp:1073
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:1439
static constexpr const bool hasAncestor
Definition basis.hpp:1193
QuadratureRule m_nodes
Nodes for the interpolation.
Definition basis.hpp:3020
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:1735
BasisType AncestorType
Definition basis.hpp:1091
BasisType::GradientValue ReturnValue
Definition basis.hpp:1789
TensorizedVectorFamily< ScalarBasisType, N >::CurlValue ReturnValue
Definition basis.hpp:2084
static constexpr const bool hasAncestor
Definition basis.hpp:571
static constexpr const TensorRankE tensorRank
Definition basis.hpp:972
MatrixRd HessianValue
Definition basis.hpp:173
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1221
Eigen::Matrix< double, dimspace, dimspace > GradientValue
Definition basis.hpp:1184
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:2132
BasisType::RotorValue RotorValue
Definition basis.hpp:317
static const bool hasGradient
Definition basis.hpp:1085
BasisType::DivergenceValue ReturnValue
Definition basis.hpp:1883
static const bool hasCurl
Definition basis.hpp:736
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:695
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_gradient_quad)
Definition basis.hpp:1820
static const bool hasHessian
Definition basis.hpp:329
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:970
static const bool hasDivergence
Definition basis.hpp:575
static const bool hasDivergence
Definition basis.hpp:1382
static const bool hasHessian
Definition basis.hpp:737
static constexpr const bool hasAncestor
Definition basis.hpp:1442
void HessianValue
Definition basis.hpp:1250
static const bool hasHessian
Definition basis.hpp:1261
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1082
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:495
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:1119
TangentFamily(const ScalarFamilyType &scalar_family, const VectorRd &generator)
Constructor.
Definition basis.hpp:907
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:926
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:673
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1639
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th basis function at point x.
Definition basis.hpp:657
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1520
MatrixFamily< ScalarBasisType, N >::DivergenceValue ReturnValue
Definition basis.hpp:2216
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:2001
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:354
BasisType::FunctionValue FunctionValue
Definition basis.hpp:963
void HessianValue
Definition basis.hpp:1500
static const bool hasHessian
Definition basis.hpp:1577
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1192
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1013
static ReturnValue evaluate(const MatrixFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2189
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:1685
static const bool hasDivergence
Definition basis.hpp:1321
static constexpr const TensorRankE tensorRank
Definition basis.hpp:251
size_t dimension() const
Return the dimension of the basis.
Definition basis.hpp:995
VectorRd FunctionValue
Definition basis.hpp:886
static const bool hasHessian
Definition basis.hpp:258
VectorRd AncestorBasisFunctionValue
Definition basis.hpp:2054
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:2368
static const bool hasFunction
Definition basis.hpp:897
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:1045
static const bool hasFunction
Definition basis.hpp:732
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:193
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:1007
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1623
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th basis function at point x.
Definition basis.hpp:1021
ScalarFamilyType::GeometricSupport GeometricSupport
Definition basis.hpp:728
CurlValue curl(size_t i, const VectorRd &x) const
Evaluate the curl of the i-th basis function at point x.
Definition basis.hpp:1029
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:2265
double AncestorBasisFunctionValue
Definition basis.hpp:2186
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1472
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:773
boost::multi_array< typename BasisType::FunctionValue, 2 > m_on_basis_nodes
Definition basis.hpp:3018
double RotorValue
Definition basis.hpp:1565
void RotorValue
Definition basis.hpp:1249
Eigen::Matrix2d MatrixRd
Definition basis.hpp:54
BasisType::DivergenceValue DivergenceValue
Definition basis.hpp:966
GradientValue CurlValue
Definition basis.hpp:563
BasisType::CurlValue ReturnValue
Definition basis.hpp:1860
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:469
static const bool hasDivergence
Definition basis.hpp:734
TensorizedVectorFamily< ScalarBasisType, N >::RotorValue ReturnValue
Definition basis.hpp:2150
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:2099
CurlValue curl(size_t i, const VectorRd &x) const
Evaluate the curl of the i-th basis function at point x.
Definition basis.hpp:1146
static const bool hasRotor
Definition basis.hpp:1200
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2057
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1959
void RotorValue
Definition basis.hpp:1618
static std::function< Eigen::MatrixXd(const Eigen::MatrixXd &)> symmetrise_matrix
Function to symmetrise a matrix (useful together with transform_values_quad)
Definition basis.hpp:1726
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th basis function at point x.
Definition basis.hpp:1037
VectorRd AncestorBasisFunctionValue
Definition basis.hpp:1988
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:217
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th function at point x.
Definition basis.hpp:430
static constexpr const bool hasAncestor
Definition basis.hpp:973
BasisType::GradientValue FunctionValue
Definition basis.hpp:1183
static const bool hasGradient
Definition basis.hpp:1573
Eigen::Matrix< double, dimspace, dimspace > CurlValue
Definition basis.hpp:1563
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:211
VectorRd FunctionValue
Definition basis.hpp:1245
void DivergenceValue
Definition basis.hpp:1186
void DivergenceValue
Definition basis.hpp:889
constexpr const BasisType & ancestor() const
Return the underlying complete basis.
Definition basis.hpp:1001
BasisType::GradientValue GradientValue
Definition basis.hpp:964
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:2067
GradientValue gradient(size_t i, const VectorRd &x) const
Evaluate the gradient of the i-th basis function at point x.
Definition basis.hpp:1138
static const bool hasHessian
Definition basis.hpp:1384
VectorRd AncestorBasisFunctionValue
Definition basis.hpp:2119
double DivergenceValue
Definition basis.hpp:1498
static constexpr const bool hasAncestor
Definition basis.hpp:178
DivergenceValue rotor(size_t i, const VectorRd &x) const
Evaluate the scalar rotor of the i-th function at point x.
Definition basis.hpp:456
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:2165
VectorRd GradientValue
Definition basis.hpp:887
constexpr const BasisType & ancestor() const
Return the underlying complete basis.
Definition basis.hpp:1113
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:2886
void CurlValue
Definition basis.hpp:1616
BasisType AncestorType
Definition basis.hpp:1263
static const bool hasCurl
Definition basis.hpp:1086
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1441
static const bool hasCurl
Definition basis.hpp:1198
static const bool hasFunction
Definition basis.hpp:1506
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:1383
void DivergenceValue
Definition basis.hpp:1435
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:781
static const bool hasGradient
Definition basis.hpp:254
Cell GeometricSupport
Definition basis.hpp:1502
VectorRd AncestorBasisFunctionValue
Definition basis.hpp:2021
MatrixFamily(const ScalarFamilyType &scalar_family)
Definition basis.hpp:741
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:65
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:154
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:168
Eigen::MatrixXd gram_schmidt(boost::multi_array< T, 2 > &basis_eval, const std::function< double(size_t, size_t)> &inner_product)
Definition basis.hpp:2404
CurlValue curl(size_t i, const VectorRd &x) const
Evaluate the curl of the i-th function at point x.
Definition basis.hpp:404
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2089
VectorRd CurlValue
Definition basis.hpp:723
VectorRd AncestorBasisFunctionValue
Definition basis.hpp:2086
static const bool hasRotor
Definition basis.hpp:183
BasisType::GradientValue ReturnValue
Definition basis.hpp:1813
static constexpr const bool hasAncestor
Definition basis.hpp:1571
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th basis function at point x.
Definition basis.hpp:1154
void RotorValue
Definition basis.hpp:1187
MatrixGradient operator+(const MatrixGradient &G)
Addition.
Definition basis.hpp:130
constexpr const BasisType & ancestor() const
Return the ancestor (basis that the gradient was taken of)
Definition basis.hpp:1288
static const bool hasCurl
Definition basis.hpp:326
BasisType::CurlValue CurlValue
Definition basis.hpp:315
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1861
static const bool hasGradient
Definition basis.hpp:1197
MatrixRd FunctionValue
Definition basis.hpp:1613
static const bool hasGradient
Definition basis.hpp:1319
Family< BasisType > m_on_basis
Definition basis.hpp:3017
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:1730
static const bool hasCurl
Definition basis.hpp:1508
void RotorValue
Definition basis.hpp:725
static constexpr const bool hasAncestor
Definition basis.hpp:252
BasisType::HessianValue HessianValue
Definition basis.hpp:968
static ReturnValue evaluate(const TensorizedVectorFamily< ScalarBasisType, N > &basis, size_t i, const VectorRd &x)
Definition basis.hpp:2122
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:443
double DivergenceValue
Definition basis.hpp:564
std::function< T(const VectorRd &)> FType
type for function of point. T is the type of value of the function
Definition basis.hpp:62
BasisType AncestorType
Definition basis.hpp:1203
static const bool hasHessian
Definition basis.hpp:184
DivergenceValue divergence(size_t i, const VectorRd &x) const
Evaluate the divergence of the i-th basis function at point x.
Definition basis.hpp:789
constexpr const BasisType & ancestor() const
Return the ancestor (basis that the gradient was taken of)
Definition basis.hpp:1350
RotorValue rotor(size_t i, const VectorRd &x) const
Evaluate the rotor of the i-th basis function at point x.
Definition basis.hpp:1162
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:417
static constexpr const TensorRankE tensorRank
Definition basis.hpp:730
static const bool hasRotor
Definition basis.hpp:576
Eigen::Matrix< double, dimspace, dimspace > GradientValue
Definition basis.hpp:1496
Eigen::Matrix< double, N, N > dx
Definition basis.hpp:147
static const bool hasGradient
Definition basis.hpp:1380
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:2522
static const bool hasDivergence
Definition basis.hpp:1628
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:1314
static constexpr const TensorRankE tensorRank
Definition basis.hpp:177
constexpr const BasisType & ancestor() const
Return the ancestor.
Definition basis.hpp:514
static const bool hasDivergence
Definition basis.hpp:1259
static const bool hasCurl
Definition basis.hpp:976
Eigen::Matrix< double, N, N > FunctionValue
Definition basis.hpp:721
VectorRd GradientValue
Definition basis.hpp:1309
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:1700
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1907
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:3000
static const bool hasRotor
Definition basis.hpp:328
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1276
static const bool hasFunction
Definition basis.hpp:1256
ScalarFamilyType AncestorType
Definition basis.hpp:904
size_t dimension() const
Return the dimension of the family.
Definition basis.hpp:595
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1930
static const bool hasCurl
Definition basis.hpp:1258
const ScalarFamilyType & ancestor() const
Return the ancestor (family that has been tensorized)
Definition basis.hpp:827
size_t dimension() const
Compute the dimension of the basis.
Definition basis.hpp:1338
static const bool hasGradient
Definition basis.hpp:1446
static const bool hasFunction
Definition basis.hpp:179
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:1913
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:1080
VectorRd AncestorBasisFunctionValue
Definition basis.hpp:2152
VectorRd CurlValue
Definition basis.hpp:1185
BasisType AncestorType
Definition basis.hpp:981
void HessianValue
Definition basis.hpp:1312
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:364
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:817
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:1532
MatrixRd FunctionValue
Definition basis.hpp:1431
Eigen::Matrix< double, dimspace, dimspace > GradientValue
Definition basis.hpp:1562
VectorRd GradientValue
Definition basis.hpp:243
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:1648
DecomposePoly(const Cell &T, const BasisType &basis)
Constructor for cell.
Definition basis.hpp:2959
static const bool hasGradient
Definition basis.hpp:733
static const bool hasRotor
Definition basis.hpp:257
static const bool hasCurl
Definition basis.hpp:1381
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1884
static const bool hasRotor
Definition basis.hpp:1629
void RotorValue
Definition basis.hpp:1372
TensorizedVectorFamily< ScalarBasisType, N >::GradientValue ReturnValue
Definition basis.hpp:2052
static const bool hasDivergence
Definition basis.hpp:1199
static const bool hasHessian
Definition basis.hpp:979
void HessianValue
Definition basis.hpp:1566
static const bool hasFunction
Definition basis.hpp:1084
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_gradient_quad)
Definition basis.hpp:1796
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:851
static const bool hasRotor
Definition basis.hpp:1576
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:938
VectorRd GradientValue
Definition basis.hpp:169
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:565
BasisType::DivergenceValue DivergenceValue
Definition basis.hpp:1076
BasisType::GradientValue ReturnValue
Definition basis.hpp:1837
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1766
ScalarFamilyType::GeometricSupport GeometricSupport
Definition basis.hpp:568
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:2700
static ReturnValue evaluate(const BasisType &basis, size_t i, size_t iqn, const boost::multi_array< ReturnValue, 2 > &ancestor_value_quad)
Definition basis.hpp:1772
VectorZd powers(size_t i) const
Returns the powers of the i-th basis function (not including the vector part)
Definition basis.hpp:1538
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1406
static const bool hasHessian
Definition basis.hpp:580
size_t max_degree() const
Returns the maximum degree of the basis functions.
Definition basis.hpp:520
HessianBasis(const BasisType &basis)
Constructor.
Definition basis.hpp:1455
double AncestorBasisFunctionValue
Definition basis.hpp:1956
static const bool hasRotor
Definition basis.hpp:901
void HessianValue
Definition basis.hpp:726
Eigen::Matrix< double, N, N > dy
Definition basis.hpp:148
BasisType::GeometricSupport GeometricSupport
Definition basis.hpp:1190
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1814
static const bool hasFunction
Definition basis.hpp:324
static const bool hasGradient
Definition basis.hpp:180
ShiftedBasis(const BasisType &basis, const int shift)
Constructor.
Definition basis.hpp:984
static constexpr const bool hasAncestor
Definition basis.hpp:1505
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th basis function at point x.
Definition basis.hpp:1344
HessianValue hessian(size_t i, const VectorRd &x) const
Evaluate the Hessian of the i-th basis function at point x.
Definition basis.hpp:1053
ScalarFamilyType AncestorType
Definition basis.hpp:582
void HessianValue
Definition basis.hpp:1619
const std::shared_ptr< RolyComplBasisCell > & rck() const
Return the Rck basis.
Definition basis.hpp:1598
void CurlValue
Definition basis.hpp:1434
static const bool hasHessian
Definition basis.hpp:1450
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:649
constexpr const BasisType & ancestor() const
Return the ancestor (basis that the gradient was taken of)
Definition basis.hpp:1412
static const bool hasGradient
Definition basis.hpp:325
static constexpr const TensorRankE tensorRank
Definition basis.hpp:322
static constexpr const TensorRankE tensorRank
Definition basis.hpp:1316
size_t dimension() const
Dimension of the basis.
Definition basis.hpp:267
static constexpr const bool hasAncestor
Definition basis.hpp:1255
void DivergenceValue
Definition basis.hpp:1248
static const bool hasFunction
Definition basis.hpp:1443
BasisType AncestorType
Definition basis.hpp:1324
static const bool hasDivergence
Definition basis.hpp:1509
QuadratureRule get_nodes() const
Return the set of nodes (useful to compute value of polynomial to decompose via evaluate_quad)
Definition basis.hpp:2994
static const bool hasFunction
Definition basis.hpp:1572
DecomposePoly(const Edge &E, const BasisType &basis)
Constructor for edge.
Definition basis.hpp:2929
void HessianValue
Definition basis.hpp:1373
VectorRd AncestorBasisFunctionValue
Definition basis.hpp:2218
size_t m_dim
Definition basis.hpp:3015
static const bool hasCurl
Definition basis.hpp:1627
static const bool hasCurl
Definition basis.hpp:181
Edge GeometricSupport
Definition basis.hpp:249
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:1969
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:932
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:2288
double FunctionValue
Definition basis.hpp:1368
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:2231
BasisType::HessianValue ReturnValue
Definition basis.hpp:1929
Eigen::Matrix< double, N, 1 > DivergenceValue
Definition basis.hpp:724
static const bool hasDivergence
Definition basis.hpp:977
const Eigen::MatrixXd & matrix() const
Return the coefficient matrix.
Definition basis.hpp:508
size_t m_nb_nodes
Definition basis.hpp:3019
constexpr const BasisType & ancestor() const
Return the ancestor (basis that the gradient was taken of)
Definition basis.hpp:1478
static const bool hasDivergence
Definition basis.hpp:1575
static const bool hasRotor
Definition basis.hpp:978
static const bool hasGradient
Definition basis.hpp:975
MatrixGradient operator-(const MatrixGradient &G)
Subtraction.
Definition basis.hpp:142
MatrixFamily< ScalarBasisType, N >::FunctionValue ReturnValue
Definition basis.hpp:2184
double FunctionValue
Definition basis.hpp:242
static const bool hasCurl
Definition basis.hpp:1320
static constexpr const bool hasAncestor
Definition basis.hpp:323
static constexpr const bool hasAncestor
Definition basis.hpp:1317
static const bool hasHessian
Definition basis.hpp:1089
static const bool hasGradient
Definition basis.hpp:1626
static ReturnValue evaluate(const BasisType &basis, size_t i, const VectorRd &x)
Definition basis.hpp:1790
TensorizedVectorFamily< ScalarBasisType, N >::DivergenceValue ReturnValue
Definition basis.hpp:2117
TensorizedVectorFamily(const ScalarFamilyType &scalar_family)
Definition basis.hpp:584
static const bool hasFunction
Definition basis.hpp:1379
BasisType::CurlValue CurlValue
Definition basis.hpp:1075
Eigen::Matrix< double, dimspace, dimspace > CurlValue
Definition basis.hpp:1497
void RotorValue
Definition basis.hpp:246
void DivergenceValue
Definition basis.hpp:1371
BasisType::RotorValue RotorValue
Definition basis.hpp:1077
VectorRd CurlValue
Definition basis.hpp:170
@ Curl
Definition basis.hpp:1677
@ Function
Definition basis.hpp:1673
@ Gradient
Definition basis.hpp:1674
@ SkewsymmetricGradient
Definition basis.hpp:1676
@ Rotor
Definition basis.hpp:1679
@ Divergence
Definition basis.hpp:1678
@ SymmetricGradient
Definition basis.hpp:1675
@ Hessian
Definition basis.hpp:1680
@ Scalar
Definition basis.hpp:66
@ Vector
Definition basis.hpp:67
@ Matrix
Definition basis.hpp:68
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:2922
Structure to store the gradient of a matrix-valued function.
Definition basis.hpp:111
Compute vectors listing the powers of monomial basis functions (for a cell, only this specialization ...
Definition basis.hpp:79
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:1757
Evaluate a basis at quadrature nodes. 'BasisFunction' (=Function, Gradient, Curl, Divergence,...
Definition basis.hpp:2284