PolyFEM
|
Stores the basis functions for a given element in a mesh (facet in 2d, cell in 3d). More...
#include <ElementBases.hpp>
Public Types | |
typedef std::function< Eigen::VectorXi(const int local_index, const mesh::Mesh &mesh)> | LocalNodeFromPrimitiveFunc |
typedef std::function< void(const Eigen::MatrixXd &uv, std::vector< assembler::AssemblyValues > &basis_values)> | EvalBasesFunc |
typedef std::function< void(quadrature::Quadrature &quadrature)> | QuadratureFunction |
Public Member Functions | |
Eigen::MatrixXd | nodes () const |
Assemble the global nodal positions of the bases. | |
void | compute_quadrature (quadrature::Quadrature &quadrature) const |
quadrature points to evaluate the basis functions inside the element | |
void | compute_mass_quadrature (quadrature::Quadrature &quadrature) const |
Eigen::VectorXi | local_nodes_for_primitive (const int local_index, const mesh::Mesh &mesh) const |
void | eval_geom_mapping (const Eigen::MatrixXd &samples, Eigen::MatrixXd &mapped) const |
Map the sample positions in the parametric domain to the object domain (if the element has no parameterization, e.g. | |
void | eval_geom_mapping_grads (const Eigen::MatrixXd &samples, std::vector< Eigen::MatrixXd > &grads) const |
Evaluate the gradients of the geometric mapping defined above. | |
bool | is_complete () const |
Checks if all the bases are complete. | |
void | set_quadrature (const QuadratureFunction &fun) |
void | set_mass_quadrature (const QuadratureFunction &fun) |
void | evaluate_bases (const Eigen::MatrixXd &uv, std::vector< assembler::AssemblyValues > &basis_values) const |
evaluate stored bases at given points on the reference element saves results to basis_values | |
void | evaluate_grads (const Eigen::MatrixXd &uv, std::vector< assembler::AssemblyValues > &basis_values) const |
evaluate gradient of stored bases at given points on the reference element saves results to basis_values | |
void | set_bases_func (EvalBasesFunc fun) |
void | set_grads_func (EvalBasesFunc fun) |
void | set_local_node_from_primitive_func (LocalNodeFromPrimitiveFunc fun) |
sets mapping from local nodes to global nodes | |
Public Attributes | |
std::vector< Basis > | bases |
one basis function per node in the element | |
bool | has_parameterization = true |
Private Member Functions | |
void | evaluate_bases_default (const Eigen::MatrixXd &uv, std::vector< assembler::AssemblyValues > &basis_values) const |
default to simply calling the Basis evaluation functions | |
void | evaluate_grads_default (const Eigen::MatrixXd &uv, std::vector< assembler::AssemblyValues > &basis_values) const |
Friends | |
std::ostream & | operator<< (std::ostream &os, const ElementBases &obj) |
Stores the basis functions for a given element in a mesh (facet in 2d, cell in 3d).
Definition at line 16 of file ElementBases.hpp.
typedef std::function<void(const Eigen::MatrixXd &uv, std::vector<assembler::AssemblyValues> &basis_values)> polyfem::basis::ElementBases::EvalBasesFunc |
Definition at line 22 of file ElementBases.hpp.
typedef std::function<Eigen::VectorXi(const int local_index, const mesh::Mesh &mesh)> polyfem::basis::ElementBases::LocalNodeFromPrimitiveFunc |
Definition at line 19 of file ElementBases.hpp.
typedef std::function<void(quadrature::Quadrature &quadrature)> polyfem::basis::ElementBases::QuadratureFunction |
Definition at line 25 of file ElementBases.hpp.
|
inline |
Definition at line 34 of file ElementBases.hpp.
References mass_quadrature_builder_, and quadrature.
Referenced by polyfem::assembler::AssemblyValsCache::compute(), and polyfem::assembler::AssemblyValsCache::update().
|
inline |
quadrature points to evaluate the basis functions inside the element
Definition at line 33 of file ElementBases.hpp.
References quadrature, and quadrature_builder_.
Referenced by polyfem::assembler::ElementAssemblyValues::compute(), and polyfem::assembler::ElementAssemblyValues::is_geom_mapping_positive().
void polyfem::basis::ElementBases::eval_geom_mapping | ( | const Eigen::MatrixXd & | samples, |
Eigen::MatrixXd & | mapped | ||
) | const |
Map the sample positions in the parametric domain to the object domain (if the element has no parameterization, e.g.
harmonic bases, then the parametric domain = object domain, and the mapping is identity)
[in] | samples | #S x dim matrix of sample positions to evaluate |
[out] | mapped | #S x dim matrix of mapped positions |
Definition at line 19 of file ElementBases.cpp.
References bases, evaluate_bases(), polyfem::basis::Basis::global(), and has_parameterization.
Referenced by polyfem::assembler::MassMatrixAssembler::assemble_cross(), polyfem::solver::OperatorSplittingSolver::calculate_local_pts(), polyfem::utils::extract_nodes(), and polyfem::assembler::RhsAssembler::lsq_bc().
void polyfem::basis::ElementBases::eval_geom_mapping_grads | ( | const Eigen::MatrixXd & | samples, |
std::vector< Eigen::MatrixXd > & | grads | ||
) | const |
Evaluate the gradients of the geometric mapping defined above.
[in] | samples | #S x dim matrix of input sample positions |
[out] | grads | #S list of dim x dim matrices of gradients |
Definition at line 82 of file ElementBases.cpp.
References bases, evaluate_grads(), polyfem::basis::Basis::global(), and has_parameterization.
Referenced by polyfem::solver::OperatorSplittingSolver::calculate_local_pts().
|
inline |
evaluate stored bases at given points on the reference element saves results to basis_values
Definition at line 71 of file ElementBases.hpp.
References eval_bases_func_, and evaluate_bases_default().
Referenced by polyfem::assembler::MassMatrixAssembler::assemble_cross(), polyfem::assembler::ElementAssemblyValues::compute(), eval_geom_mapping(), polyfem::io::Evaluator::interpolate_function(), and polyfem::assembler::RhsAssembler::lsq_bc().
|
private |
default to simply calling the Basis evaluation functions
Definition at line 48 of file ElementBases.cpp.
Referenced by evaluate_bases().
|
inline |
evaluate gradient of stored bases at given points on the reference element saves results to basis_values
Definition at line 84 of file ElementBases.hpp.
References eval_grads_func_, and evaluate_grads_default().
Referenced by polyfem::assembler::ElementAssemblyValues::compute(), eval_geom_mapping_grads(), and polyfem::assembler::ElementAssemblyValues::is_geom_mapping_positive().
|
private |
Definition at line 65 of file ElementBases.cpp.
References bases.
Referenced by evaluate_grads().
bool polyfem::basis::ElementBases::is_complete | ( | ) | const |
Checks if all the bases are complete.
Definition at line 9 of file ElementBases.cpp.
References bases.
|
inline |
Definition at line 35 of file ElementBases.hpp.
References local_node_from_primitive_.
Referenced by polyfem::solver::LayerThicknessForm::build_collision_mesh(), polyfem::solver::ProxyContactForceForm::build_collision_mesh(), polyfem::assembler::RhsAssembler::compute_energy_hess(), polyfem::io::OutGeometryData::extract_boundary_mesh(), and polyfem::assembler::RhsAssembler::sample_bc().
Eigen::MatrixXd polyfem::basis::ElementBases::nodes | ( | ) | const |
Assemble the global nodal positions of the bases.
Definition at line 134 of file ElementBases.cpp.
References bases.
Referenced by polyfem::assembler::MassMatrixAssembler::assemble_cross().
|
inline |
Definition at line 96 of file ElementBases.hpp.
References eval_bases_func_.
Referenced by polyfem::basis::PolygonalBasis2d::build_bases(), polyfem::basis::PolygonalBasis3d::build_bases(), and polyfem::basis::BarycentricBasis2d::build_bases().
|
inline |
Definition at line 97 of file ElementBases.hpp.
References eval_grads_func_.
Referenced by polyfem::basis::PolygonalBasis2d::build_bases(), polyfem::basis::PolygonalBasis3d::build_bases(), and polyfem::basis::BarycentricBasis2d::build_bases().
|
inline |
sets mapping from local nodes to global nodes
Definition at line 100 of file ElementBases.hpp.
References local_node_from_primitive_.
Referenced by polyfem::basis::SplineBasis2d::build_bases(), polyfem::basis::SplineBasis3d::build_bases(), and polyfem::basis::BarycentricBasis2d::build_bases().
|
inline |
Definition at line 67 of file ElementBases.hpp.
References mass_quadrature_builder_.
Referenced by polyfem::basis::PolygonalBasis2d::build_bases(), polyfem::basis::PolygonalBasis3d::build_bases(), polyfem::basis::SplineBasis2d::build_bases(), polyfem::basis::SplineBasis3d::build_bases(), and polyfem::basis::BarycentricBasis2d::build_bases().
|
inline |
Definition at line 66 of file ElementBases.hpp.
References quadrature_builder_.
Referenced by polyfem::basis::PolygonalBasis2d::build_bases(), polyfem::basis::PolygonalBasis3d::build_bases(), polyfem::basis::SplineBasis2d::build_bases(), polyfem::basis::SplineBasis3d::build_bases(), and polyfem::basis::BarycentricBasis2d::build_bases().
|
friend |
Definition at line 57 of file ElementBases.hpp.
std::vector<Basis> polyfem::basis::ElementBases::bases |
one basis function per node in the element
Definition at line 27 of file ElementBases.hpp.
Referenced by polyfem::assembler::MassMatrixAssembler::assemble_cross(), polyfem::io::Evaluator::average_grad_based_function(), polyfem::basis::PolygonalBasis2d::build_bases(), polyfem::basis::PolygonalBasis3d::build_bases(), polyfem::basis::SplineBasis2d::build_bases(), polyfem::basis::SplineBasis3d::build_bases(), polyfem::basis::BarycentricBasis2d::build_bases(), polyfem::mesh::build_collision_proxy(), polyfem::State::build_mesh_matrices(), polyfem::assembler::ElementAssemblyValues::compute(), polyfem::compute_diplacement_grad(), polyfem::assembler::ElementAssemblyValues::eval_deformed_jacobian_determinant(), eval_geom_mapping(), eval_geom_mapping_grads(), evaluate_bases_default(), evaluate_grads_default(), polyfem::io::OutGeometryData::extract_boundary_mesh(), polyfem::assembler::ElementAssemblyValues::finalize2d(), polyfem::assembler::ElementAssemblyValues::finalize3d(), polyfem::io::Evaluator::interpolate_boundary_function(), polyfem::io::Evaluator::interpolate_function(), is_complete(), polyfem::assembler::ElementAssemblyValues::is_geom_mapping_positive(), polyfem::assembler::RhsAssembler::lsq_bc(), nodes(), polyfem::assembler::RhsAssembler::sample_bc(), and polyfem::assembler::RhsAssembler::time_bc().
|
private |
Definition at line 108 of file ElementBases.hpp.
Referenced by evaluate_bases(), and set_bases_func().
|
private |
Definition at line 109 of file ElementBases.hpp.
Referenced by evaluate_grads(), and set_grads_func().
bool polyfem::basis::ElementBases::has_parameterization = true |
Definition at line 39 of file ElementBases.hpp.
Referenced by polyfem::basis::PolygonalBasis2d::build_bases(), polyfem::basis::PolygonalBasis3d::build_bases(), polyfem::basis::BarycentricBasis2d::build_bases(), polyfem::assembler::ElementAssemblyValues::compute(), eval_geom_mapping(), eval_geom_mapping_grads(), polyfem::assembler::ElementAssemblyValues::finalize2d(), polyfem::assembler::ElementAssemblyValues::finalize3d(), and polyfem::assembler::ElementAssemblyValues::is_geom_mapping_positive().
|
private |
Definition at line 113 of file ElementBases.hpp.
Referenced by local_nodes_for_primitive(), and set_local_node_from_primitive_func().
|
private |
Definition at line 111 of file ElementBases.hpp.
Referenced by compute_mass_quadrature(), and set_mass_quadrature().
|
private |
Definition at line 110 of file ElementBases.hpp.
Referenced by compute_quadrature(), and set_quadrature().