PolyFEM
Loading...
Searching...
No Matches
ElementBases.hpp
Go to the documentation of this file.
1#pragma once
2
6
8
9#include <vector>
10
11namespace polyfem
12{
13 namespace basis
14 {
17 {
18 public:
19 typedef std::function<Eigen::VectorXi(const int local_index, const mesh::Mesh &mesh)> LocalNodeFromPrimitiveFunc;
20
21 // function type that evaluates bases at given points and saves them in basis_values
22 typedef std::function<void(const Eigen::MatrixXd &uv, std::vector<assembler::AssemblyValues> &basis_values)> EvalBasesFunc;
23
24 // function type that computes quadrature points and saves them in quadrature
25 typedef std::function<void(quadrature::Quadrature &quadrature)> QuadratureFunction;
26
27 std::vector<Basis> bases;
28
30 Eigen::MatrixXd nodes() const;
31
35 Eigen::VectorXi local_nodes_for_primitive(const int local_index, const mesh::Mesh &mesh) const { return local_node_from_primitive_(local_index, mesh); }
36
37 // whether the basis functions should be evaluated in the parametric domain (FE bases),
38 // or directly in the object domain (harmonic bases)
40
46 void eval_geom_mapping(const Eigen::MatrixXd &samples, Eigen::MatrixXd &mapped) const;
47
52 void eval_geom_mapping_grads(const Eigen::MatrixXd &samples, std::vector<Eigen::MatrixXd> &grads) const;
53
55 bool is_complete() const;
56
57 friend std::ostream &operator<<(std::ostream &os, const ElementBases &obj)
58 {
59 for (std::size_t i = 0; i < obj.bases.size(); ++i)
60 os << "local base " << i << ":\n"
61 << obj.bases[i] << "\n";
62
63 return os;
64 }
65
68
71 void evaluate_bases(const Eigen::MatrixXd &uv, std::vector<assembler::AssemblyValues> &basis_values) const
72 {
74 {
75 eval_bases_func_(uv, basis_values);
76 }
77 else
78 {
79 evaluate_bases_default(uv, basis_values);
80 }
81 }
84 void evaluate_grads(const Eigen::MatrixXd &uv, std::vector<assembler::AssemblyValues> &basis_values) const
85 {
87 {
88 eval_grads_func_(uv, basis_values);
89 }
90 else
91 {
92 evaluate_grads_default(uv, basis_values);
93 }
94 }
95
98
101
102 private:
104 void evaluate_bases_default(const Eigen::MatrixXd &uv, std::vector<assembler::AssemblyValues> &basis_values) const;
105 void evaluate_grads_default(const Eigen::MatrixXd &uv, std::vector<assembler::AssemblyValues> &basis_values) const;
106
107 private:
112
114 };
115 } // namespace basis
116} // namespace polyfem
Quadrature quadrature
Stores the basis functions for a given element in a mesh (facet in 2d, cell in 3d).
void set_bases_func(EvalBasesFunc fun)
void compute_mass_quadrature(quadrature::Quadrature &quadrature) const
void set_quadrature(const QuadratureFunction &fun)
void compute_quadrature(quadrature::Quadrature &quadrature) const
quadrature points to evaluate the basis functions inside the element
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 paramet...
void set_local_node_from_primitive_func(LocalNodeFromPrimitiveFunc fun)
sets mapping from local nodes to global nodes
QuadratureFunction quadrature_builder_
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
bool is_complete() const
Checks if all the bases are complete.
Eigen::MatrixXd nodes() const
Assemble the global nodal positions of the bases.
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_val...
std::function< Eigen::VectorXi(const int local_index, const mesh::Mesh &mesh)> LocalNodeFromPrimitiveFunc
void evaluate_bases_default(const Eigen::MatrixXd &uv, std::vector< assembler::AssemblyValues > &basis_values) const
default to simply calling the Basis evaluation functions
LocalNodeFromPrimitiveFunc local_node_from_primitive_
Eigen::VectorXi local_nodes_for_primitive(const int local_index, const mesh::Mesh &mesh) const
void set_grads_func(EvalBasesFunc fun)
void eval_geom_mapping_grads(const Eigen::MatrixXd &samples, std::vector< Eigen::MatrixXd > &grads) const
Evaluate the gradients of the geometric mapping defined above.
void evaluate_grads_default(const Eigen::MatrixXd &uv, std::vector< assembler::AssemblyValues > &basis_values) const
std::function< void(quadrature::Quadrature &quadrature)> QuadratureFunction
std::vector< Basis > bases
one basis function per node in the element
void set_mass_quadrature(const QuadratureFunction &fun)
std::function< void(const Eigen::MatrixXd &uv, std::vector< assembler::AssemblyValues > &basis_values)> EvalBasesFunc
QuadratureFunction mass_quadrature_builder_
friend std::ostream & operator<<(std::ostream &os, const ElementBases &obj)
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39