|
HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
|
#include <boost/multi_array.hpp>#include <mesh.hpp>#include <iostream>#include <polynomialspacedimension.hpp>#include <quadraturerule.hpp>

Go to the source code of this file.
Namespaces | |
| namespace | HArDCore3D |
| namespace | HArDCore3D::detail |
Typedefs | |
| typedef Eigen::Matrix3d | HArDCore3D::MatrixRd |
| template<typename T > | |
| using | HArDCore3D::FType = std::function< T(const VectorRd &)> |
| type for function of point. T is the return type of the function | |
| template<typename T > | |
| using | HArDCore3D::CellFType = std::function< T(const VectorRd &, const Cell *)> |
| type for function of point. T is the return type of the function | |
| template<typename T > | |
| using | HArDCore3D::BasisQuad = boost::multi_array< T, 2 > |
| type for bases evaluated on quadrature nodes | |
Enumerations | |
| enum | HArDCore3D::TensorRankE { HArDCore3D::Scalar = 0 , HArDCore3D::Vector = 1 , HArDCore3D::Matrix = 2 } |
| enum | HArDCore3D::BasisFunctionE { HArDCore3D::Function , HArDCore3D::Gradient , HArDCore3D::Curl , HArDCore3D::Divergence , HArDCore3D::Hessian , HArDCore3D::CurlCurl } |
Functions | |
| template<typename ScalarFamilyType > | |
| DivergenceBasis< TangentFamily< ScalarFamilyType > > | HArDCore3D::ScalarRotFamily (const TangentFamily< ScalarFamilyType > &tangent_family, const Face &F) |
| The following function creates the "scalar rot" basis of a TangentFamily on a face. | |
| template<typename outValue , typename inValue , typename FunctionType > | |
| boost::multi_array< outValue, 2 > | HArDCore3D::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. F must take inValue and return outValue. The function must be called with outValue as template argument: transform_values_quad<outValue>(...) | |
| template<typename ScalarBasisType , size_t N> | |
| Family< TensorizedVectorFamily< ScalarBasisType, N > > | HArDCore3D::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 "TensorizedVectorFamily" (built on B, of size N) that represents the family (B_1v_1..B_rv_1 B_1v_2...B_rv_2...B_1v_k...B_rv_k). | |
| Eigen::MatrixXd | HArDCore3D::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 size b. | |
| template<typename ScalarBasisType , size_t N> | |
| Family< MatrixFamily< ScalarBasisType, N > > | HArDCore3D::IsotropicMatrixFamily (const ScalarBasisType &B) |
| From a scalar family B, constructs a "Family" of "MatrixFamily" (built on B, of size NxN) that represents the family B Id on the MatrixFamily. | |
| template<typename T > | |
| Eigen::MatrixXd | HArDCore3D::gram_schmidt (boost::multi_array< T, 2 > &basis_eval, const std::function< double(size_t, size_t)> &inner_product) |
| double | HArDCore3D::scalar_product (const double &x, const double &y) |
| Scalar product between two reals. | |
| double | HArDCore3D::scalar_product (const double &x, const Eigen::Matrix< double, 1, 1 > &y) |
| Scalar product between one real and one 1-dimension Eigen vector. | |
| double | HArDCore3D::scalar_product (const VectorRd &x, const VectorRd &y) |
| Scalar product between two vectors. | |
| template<int N> | |
| double | HArDCore3D::scalar_product (const Eigen::Matrix< double, N, N > &x, const Eigen::Matrix< double, N, N > &y) |
| Scalar product between two matrices. | |
| template<typename Value > | |
| boost::multi_array< double, 2 > | HArDCore3D::scalar_product (const boost::multi_array< Value, 2 > &basis_quad, const Value &v) |
| This overloading of the scalar_product function computes the scalar product between an evaluation of a basis and a constant value; both basis values and constant value must be of type Value. | |
| boost::multi_array< VectorRd, 2 > | HArDCore3D::vector_product (const boost::multi_array< VectorRd, 2 > &basis_quad, const VectorRd &v) |
| Compute the vector (cross) product between the evaluation of a basis and a constant vector. | |
| template<typename BasisType > | |
| Family< BasisType > | HArDCore3D::l2_orthonormalize (const BasisType &basis, const QuadratureRule &qr, boost::multi_array< typename BasisType::FunctionValue, 2 > &basis_quad) |
| \(L^2\)-orthonormalization: simply consists in using gram_schmidt() with the specific l2 inner product | |
| template<typename BasisType > | |
| Family< BasisType > | HArDCore3D::l2_orthonormalize (const BasisType &basis, const Eigen::MatrixXd &GM) |
| \(L^2\)-orthonormalization: when the Gram Matrix is passed, we use Cholesky. | |
| template<typename FunctionValue > | |
| Eigen::MatrixXd | HArDCore3D::compute_gram_matrix (const boost::multi_array< FunctionValue, 2 > &B1, const boost::multi_array< FunctionValue, 2 > &B2, const QuadratureRule &qr, const size_t nrows, const size_t ncols, const std::string sym="nonsym") |
| template<typename FunctionValue > | |
| Eigen::MatrixXd | HArDCore3D::compute_gram_matrix (const boost::multi_array< FunctionValue, 2 > &B1, const boost::multi_array< FunctionValue, 2 > &B2, const QuadratureRule &qr, const std::string sym="nonsym") |
| template<typename FunctionValue > | |
| Eigen::MatrixXd | HArDCore3D::compute_gram_matrix (const boost::multi_array< FunctionValue, 2 > &B, const QuadratureRule &qr) |
| Compute the Gram matrix given the evaluation of one family of functions at quadrature nodes. Consists in calling the generic templated version with B1=B2. | |
| Eigen::MatrixXd | HArDCore3D::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 tensorizing the latter. | |
| Eigen::MatrixXd | HArDCore3D::compute_gram_matrix (const boost::multi_array< double, 2 > &B1, const boost::multi_array< double, 2 > &B2, const QuadratureRule &qr, const size_t nrows, const size_t ncols, const std::string sym) |
| Eigen::MatrixXd | HArDCore3D::compute_gram_matrix (const boost::multi_array< double, 2 > &B1, const boost::multi_array< double, 2 > &B2, const QuadratureRule &qr, const std::string sym) |
| Eigen::MatrixXd | HArDCore3D::compute_gram_matrix (const boost::multi_array< VectorRd, 2 > &B1, const boost::multi_array< VectorRd, 2 > &B2, const QuadratureRule &qr, const size_t nrows, const size_t ncols, const std::string sym) |
| Eigen::MatrixXd | HArDCore3D::compute_gram_matrix (const boost::multi_array< VectorRd, 2 > &B1, const boost::multi_array< VectorRd, 2 > &B2, const QuadratureRule &qr, const std::string sym="nonsym") |
| Compute the Gram-like matrix given the evaluation of two families of functions at quadrature nodes. Consists in calling the Vector3d-valued version with nrows = nb of elements in B1, ncols = nb of elements in B2. | |
| template<typename ScalarFamilyType , size_t N> | |
| Eigen::MatrixXd | HArDCore3D::compute_gram_matrix (const MatrixFamily< ScalarFamilyType, N > &MatFam, const boost::multi_array< double, 2 > &scalar_family_quad, const QuadratureRule &qr) |
| Compute the Gram-like matrix for a MatrixFamily. This overload is more efficient than the generic function as it only computes the gram matrix of the underlying scalar family, and then creates the bloc-diagonal gram matrix of the MatrixFamily (which is indeed bloc diagonal given the choice of m_E elements in this class). | |
| template<typename T > | |
| Eigen::VectorXd | HArDCore3D::integrate (const FType< T > &f, const BasisQuad< T > &B, const QuadratureRule &qr, size_t n_rows=0) |
| Computes the vector of integrals (f, phi_i) | |
| template<typename T , typename U > | |
| Eigen::MatrixXd | HArDCore3D::compute_weighted_gram_matrix (const FType< U > &f, const BasisQuad< T > &B1, const BasisQuad< T > &B2, const QuadratureRule &qr, size_t n_rows=0, size_t n_cols=0, const std::string sym="nonsym") |
| Computes the Gram-like matrix of integrals (f phi_i, phi_j) | |
| template<typename T , typename U > | |
| Eigen::MatrixXd | HArDCore3D::compute_weighted_gram_matrix (const FType< U > &f, const BasisQuad< T > &B1, const BasisQuad< T > &B2, const QuadratureRule &qr, const std::string sym) |
| Computes the Gram-like matrix of integrals (f phi_i, phi_j) | |
| Eigen::MatrixXd | HArDCore3D::compute_weighted_gram_matrix (const FType< VectorRd > &f, const BasisQuad< VectorRd > &B1, const BasisQuad< double > &B2, const QuadratureRule &qr, size_t n_rows=0, size_t n_cols=0) |
| Computes the Gram-like matrix of integrals (f dot phi_i, phi_j) | |
| Eigen::MatrixXd | HArDCore3D::compute_weighted_gram_matrix (const FType< VectorRd > &f, const BasisQuad< double > &B1, const BasisQuad< VectorRd > &B2, const QuadratureRule &qr, size_t n_rows=0, size_t n_cols=0) |
| Computes the Gram-like matrix of integrals (phi_i, f dot phi_j) | |
| template<typename BasisType > | |
| Eigen::VectorXd | HArDCore3D::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. | |
Variables | |
| constexpr int | HArDCore3D::dimspace = 3 |
| Dimension, and generic types for vector in correct dimension (makes it easier to translate a code between 2D and 3D) | |
| static const std::vector< VectorRd > | HArDCore3D::basisRd = { VectorRd(1., 0., 0.), VectorRd(0., 1., 0.), VectorRd(0., 0., 1.) } |
| static std::function< Eigen::MatrixXd(const Eigen::MatrixXd &)> | HArDCore3D::symmetrise_matrix = [](const Eigen::MatrixXd & x)->Eigen::MatrixXd { return 0.5*(x+x.transpose());} |
| Function to symmetrise a matrix (useful together with transform_values_quad) | |
| static std::function< Eigen::MatrixXd(const Eigen::MatrixXd &)> | HArDCore3D::skew_symmetrise_matrix = [](const Eigen::MatrixXd & x)->Eigen::MatrixXd { return 0.5*(x-x.transpose());} |
| Function to skew-symmetrise a matrix (useful together with transform_values_quad) | |