PolyFEM
Loading...
Searching...
No Matches
ElementAssemblyValues.hpp
Go to the documentation of this file.
1#pragma once
2
5
6#include <vector>
7
8namespace polyfem
9{
10 namespace assembler
11 {
14 {
15 public:
16 // m = number of quadrature points
17
18 // vector of basis values and gradients at quadrature points for this element
19 // each element samples a single basis function at the m quadrature points
20 std::vector<AssemblyValues> basis_values;
21 // inverse transpose jacobian of geom mapping at quadrature points
22 std::vector<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3>> jac_it;
23
24 // quadrature rule to use (stores points and weights)
27
28 // img of quadrature points through the geom mapping (global pos in the mesh)
29 Eigen::MatrixXd val; // R^{m x dim}
30
31 // det(∑∇φᵢ⋅Nᵢ) det fo the jacobian of geometric mapping (constant for P1)
32 Eigen::VectorXd det; // R^{m x 1}
33
34 // only poly elements have no parameterization
36
39 void compute(const int el_index, const bool is_volume, const Eigen::MatrixXd &pts, const basis::ElementBases &basis, const basis::ElementBases &gbasis);
40
42 void compute(const int el_index, const bool is_volume, const basis::ElementBases &basis, const basis::ElementBases &gbasis);
43
45 bool is_geom_mapping_positive(const bool is_volume, const basis::ElementBases &gbasis) const;
46
47 private:
48 std::vector<AssemblyValues> g_basis_values_cache_;
49
50 void finalize_global_element(const Eigen::MatrixXd &v);
51
54
56 void finalize2d(const basis::ElementBases &gbasis, const std::vector<AssemblyValues> &gbasis_values);
57 void finalize3d(const basis::ElementBases &gbasis, const std::vector<AssemblyValues> &gbasis_values);
58
59 bool is_geom_mapping_positive(const Eigen::MatrixXd &dx, const Eigen::MatrixXd &dy, const Eigen::MatrixXd &dz) const;
60 bool is_geom_mapping_positive(const Eigen::MatrixXd &dx, const Eigen::MatrixXd &dy) const;
61 };
62 } // namespace assembler
63} // namespace polyfem
stores per element basis values at given quadrature points and geometric mapping
bool is_geom_mapping_positive(const bool is_volume, const basis::ElementBases &gbasis) const
check if the element is flipped
std::vector< AssemblyValues > g_basis_values_cache_
void compute(const int el_index, const bool is_volume, const Eigen::MatrixXd &pts, const basis::ElementBases &basis, const basis::ElementBases &gbasis)
computes the per element values at the local (ref el) points (pts) sets basis_values,...
void finalize2d(const basis::ElementBases &gbasis, const std::vector< AssemblyValues > &gbasis_values)
void finalize(const Eigen::MatrixXd &v, const Eigen::MatrixXd &dx, const Eigen::MatrixXd &dy); void f...
std::vector< Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > > jac_it
void finalize_global_element(const Eigen::MatrixXd &v)
void finalize3d(const basis::ElementBases &gbasis, const std::vector< AssemblyValues > &gbasis_values)
Stores the basis functions for a given element in a mesh (facet in 2d, cell in 3d).