PolyFEM
Loading...
Searching...
No Matches
polyfem::mesh::Mesh3D Class Referenceabstract

#include <Mesh3D.hpp>

Inheritance diagram for polyfem::mesh::Mesh3D:
[legend]
Collaboration diagram for polyfem::mesh::Mesh3D:
[legend]

Public Member Functions

 Mesh3D ()=default
 
virtual ~Mesh3D ()=default
 
 Mesh3D (Mesh3D &&)=default
 
Mesh3Doperator= (Mesh3D &&)=default
 
 Mesh3D (const Mesh3D &)=default
 
Mesh3Doperator= (const Mesh3D &)=default
 
bool is_volume () const override
 checks if mesh is volume
 
virtual int n_cell_edges (const int c_id) const =0
 
virtual int n_cell_faces (const int c_id) const =0
 
virtual int cell_face (const int c_id, const int lf_id) const =0
 
virtual int cell_edge (const int c_id, const int le_id) const =0
 
void elements_boxes (std::vector< std::array< Eigen::Vector3d, 2 > > &boxes) const override
 constructs a box around every element (3d cell, 2d face)
 
void barycentric_coords (const RowVectorNd &p, const int el_id, Eigen::MatrixXd &coord) const override
 constructs barycentric coodiantes for a point p.
 
void compute_cell_jacobian (const int el_id, const Eigen::MatrixXd &reference_map, Eigen::MatrixXd &jacobian) const
 
virtual RowVectorNd kernel (const int cell_id) const =0
 
double tri_area (const int gid) const override
 area of a tri face of a tet mesh
 
std::pair< RowVectorNd, int > edge_node (const Navigation3D::Index &index, const int n_new_nodes, const int i) const
 
std::pair< RowVectorNd, int > face_node (const Navigation3D::Index &index, const int n_new_nodes, const int i, const int j) const
 
std::pair< RowVectorNd, int > cell_node (const Navigation3D::Index &index, const int n_new_nodes, const int i, const int j, const int k) const
 
void get_edges (Eigen::MatrixXd &p0, Eigen::MatrixXd &p1) const override
 Get all the edges.
 
void get_edges (Eigen::MatrixXd &p0, Eigen::MatrixXd &p1, const std::vector< bool > &valid_elements) const override
 Get all the edges according to valid_elements selection.
 
virtual Navigation3D::Index get_index_from_element (int hi, int lf, int lv) const =0
 
virtual Navigation3D::Index get_index_from_element (int hi) const =0
 
virtual Navigation3D::Index get_index_from_element_edge (int hi, int v0, int v1) const =0
 
virtual Navigation3D::Index get_index_from_element_face (int hi, int v0, int v1, int v2) const =0
 
virtual std::vector< uint32_t > vertex_neighs (const int v_gid) const =0
 
virtual std::vector< uint32_t > edge_neighs (const int e_gid) const =0
 
virtual Navigation3D::Index switch_vertex (Navigation3D::Index idx) const =0
 
virtual Navigation3D::Index switch_edge (Navigation3D::Index idx) const =0
 
virtual Navigation3D::Index switch_face (Navigation3D::Index idx) const =0
 
virtual Navigation3D::Index switch_element (Navigation3D::Index idx) const =0
 
virtual Navigation3D::Index next_around_edge (Navigation3D::Index idx) const =0
 
virtual Navigation3D::Index next_around_face (Navigation3D::Index idx) const =0
 
void to_face_functions (std::array< std::function< Navigation3D::Index(Navigation3D::Index)>, 6 > &to_face) const
 
void to_vertex_functions (std::array< std::function< Navigation3D::Index(Navigation3D::Index)>, 8 > &to_vertex) const
 
void to_edge_functions (std::array< std::function< Navigation3D::Index(Navigation3D::Index)>, 12 > &to_edge) const
 
std::array< int, 8 > get_ordered_vertices_from_hex (const int element_index) const
 
virtual std::array< int, 4 > get_ordered_vertices_from_tet (const int element_index) const
 
virtual void get_vertex_elements_neighs (const int v_id, std::vector< int > &ids) const =0
 
virtual void get_edge_elements_neighs (const int e_id, std::vector< int > &ids) const =0
 
void compute_element_barycenters (Eigen::MatrixXd &barycenters) const override
 utility for 2d/3d.
 
- Public Member Functions inherited from polyfem::mesh::Mesh
virtual std::unique_ptr< Meshcopy () const =0
 Create a copy of the mesh.
 
virtual ~Mesh ()=default
 Destroy the Mesh object.
 
 Mesh (Mesh &&)=default
 Construct a new Mesh object.
 
Meshoperator= (Mesh &&)=default
 Copy constructor.
 
 Mesh (const Mesh &)=default
 Construct a new Mesh object.
 
Meshoperator= (const Mesh &)=default
 
virtual void refine (const int n_refinement, const double t)=0
 refine the mesh
 
int dimension () const
 utily for dimension
 
virtual bool is_conforming () const =0
 if the mesh is conforming
 
int n_elements () const
 utitlity to return the number of elements, cells or faces in 3d and 2d
 
int n_boundary_elements () const
 utitlity to return the number of boundary elements, faces or edges in 3d and 2d
 
virtual int n_cells () const =0
 number of cells
 
virtual int n_faces () const =0
 number of faces
 
virtual int n_edges () const =0
 number of edges
 
virtual int n_vertices () const =0
 number of vertices
 
virtual int n_face_vertices (const int f_id) const =0
 number of vertices of a face
 
virtual int n_cell_vertices (const int c_id) const =0
 number of vertices of a cell
 
virtual int edge_vertex (const int e_id, const int lv_id) const =0
 id of the edge vertex
 
virtual int face_vertex (const int f_id, const int lv_id) const =0
 id of the face vertex
 
virtual int cell_vertex (const int f_id, const int lv_id) const =0
 id of the vertex of a cell
 
int element_vertex (const int el_id, const int lv_id) const
 id of the vertex of a element
 
std::vector< int > element_vertices (const int el_id) const
 list of vids of an element
 
int boundary_element_vertex (const int primitive_id, const int lv_id) const
 
virtual bool is_boundary_vertex (const int vertex_global_id) const =0
 is vertex boundary
 
virtual bool is_boundary_edge (const int edge_global_id) const =0
 is edge boundary
 
virtual bool is_boundary_face (const int face_global_id) const =0
 is face boundary
 
virtual bool is_boundary_element (const int element_global_id) const =0
 is cell boundary
 
virtual bool save (const std::string &path) const =0
 
virtual void attach_higher_order_nodes (const Eigen::MatrixXd &V, const std::vector< std::vector< int > > &nodes)=0
 attach high order nodes
 
const Eigen::MatrixXi & orders () const
 order of each element
 
bool is_rational () const
 check if curved mesh has rational polynomials elements
 
void set_is_rational (const bool in_is_rational)
 Set the is rational object.
 
virtual void normalize ()=0
 normalize the mesh
 
virtual void compute_elements_tag ()=0
 compute element types, see ElementType
 
virtual void update_elements_tag ()
 Update elements types.
 
virtual double edge_length (const int gid) const
 edge length
 
virtual double quad_area (const int gid) const
 area of a quad face of an hex mesh
 
virtual RowVectorNd point (const int global_index) const =0
 point coordinates
 
virtual void set_point (const int global_index, const RowVectorNd &p)=0
 Set the point.
 
virtual RowVectorNd edge_barycenter (const int e) const =0
 edge barycenter
 
virtual RowVectorNd face_barycenter (const int f) const =0
 face barycenter
 
virtual RowVectorNd cell_barycenter (const int c) const =0
 cell barycenter
 
void edge_barycenters (Eigen::MatrixXd &barycenters) const
 all edges barycenters
 
void face_barycenters (Eigen::MatrixXd &barycenters) const
 all face barycenters
 
void cell_barycenters (Eigen::MatrixXd &barycenters) const
 all cells barycenters
 
virtual void bounding_box (RowVectorNd &min, RowVectorNd &max) const =0
 computes the bbox of the mesh
 
bool is_spline_compatible (const int el_id) const
 checks if element is spline compatible
 
bool is_cube (const int el_id) const
 checks if element is cube compatible
 
bool is_polytope (const int el_id) const
 checks if element is polygon compatible
 
bool is_simplex (const int el_id) const
 checks if element is simples compatible
 
const std::vector< ElementType > & elements_tag () const
 Returns the elements types.
 
void set_tag (const int el, const ElementType type)
 changes the element type
 
void compute_node_ids (const std::function< int(const size_t, const RowVectorNd &, bool)> &marker)
 computes boundary selections based on a function
 
virtual void load_boundary_ids (const std::string &path)
 loads the boundary selections for a file
 
virtual void compute_boundary_ids (const std::function< int(const size_t, const std::vector< int > &, const RowVectorNd &, bool)> &marker)=0
 computes boundary selections based on a function
 
virtual void compute_body_ids (const std::function< int(const size_t, const std::vector< int > &, const RowVectorNd &)> &marker)=0
 computes boundary selections based on a function
 
virtual void set_boundary_ids (const std::vector< int > &boundary_ids)
 Set the boundary selection from a vector.
 
virtual void set_body_ids (const std::vector< int > &body_ids)
 Set the volume sections.
 
virtual int get_default_boundary_id (const int primitive) const
 Get the default boundary selection of an element (face in 3d, edge in 2d)
 
virtual int get_boundary_id (const int primitive) const
 Get the boundary selection of an element (face in 3d, edge in 2d)
 
virtual int get_node_id (const int node_id) const
 Get the boundary selection of a node.
 
void update_nodes (const Eigen::VectorXi &in_node_to_node)
 Update the node ids to reorder them.
 
virtual int get_body_id (const int primitive) const
 Get the volume selection of an element (cell in 3d, face in 2d)
 
virtual const std::vector< int > & get_body_ids () const
 Get the volume selection of all elements (cells in 3d, faces in 2d)
 
bool has_node_ids () const
 checks if points selections are available
 
bool has_boundary_ids () const
 checks if surface selections are available
 
virtual bool has_body_ids () const
 checks if volumes selections are available
 
const std::vector< double > & cell_weights (const int cell_index) const
 weights for rational polynomial meshes
 
void set_cell_weights (const std::vector< std::vector< double > > &in_cell_weights)
 Set the cell weights for rational polynomial meshes.
 
virtual void prepare_mesh ()
 method used to finalize the mesh.
 
bool has_poly () const
 checks if the mesh has polytopes
 
bool is_simplicial () const
 checks if the mesh is simplicial
 
bool is_linear () const
 check if the mesh is linear
 
std::vector< std::pair< int, int > > edges () const
 list of sorted edges.
 
std::vector< std::vector< int > > faces () const
 list of sorted faces.
 
std::unordered_map< std::pair< int, int >, size_t, polyfem::utils::HashPairedges_to_ids () const
 map from edge (pair of v id) to the id of the edge
 
std::unordered_map< std::vector< int >, size_t, polyfem::utils::HashVectorfaces_to_ids () const
 map from face (tuple of v id) to the id of the face
 
const Eigen::VectorXi & in_ordered_vertices () const
 Order of the input vertices.
 
const Eigen::MatrixXi & in_ordered_edges () const
 Order of the input edges.
 
const Eigen::MatrixXi & in_ordered_faces () const
 Order of the input edges.
 
virtual void append (const Mesh &mesh)
 appends a new mesh to the end of this
 
void append (const std::unique_ptr< Mesh > &mesh)
 appends a new mesh to the end of this, utility that takes pointer, calls other one
 
void apply_affine_transformation (const MatrixNd &A, const VectorNd &b)
 Apply an affine transformation \(Ax+b\) to the vertex positions \(x\).
 

Additional Inherited Members

- Static Public Member Functions inherited from polyfem::mesh::Mesh
static std::unique_ptr< Meshcreate (const std::string &path, const bool non_conforming=false)
 factory to build the proper mesh
 
static std::unique_ptr< Meshcreate (GEO::Mesh &M, const bool non_conforming=false)
 factory to build the proper mesh
 
static std::unique_ptr< Meshcreate (const Eigen::MatrixXd &vertices, const Eigen::MatrixXi &cells, const bool non_conforming=false)
 factory to build the proper mesh
 
static std::unique_ptr< Meshcreate (const int dim, const bool non_conforming=false)
 factory to build the proper empty mesh
 
- Protected Member Functions inherited from polyfem::mesh::Mesh
 Mesh ()=default
 Construct a new Mesh object.
 
virtual bool load (const std::string &path)=0
 loads a mesh from the path
 
virtual bool load (const GEO::Mesh &M)=0
 loads a mesh from a geo mesh
 
- Protected Attributes inherited from polyfem::mesh::Mesh
std::vector< ElementTypeelements_tag_
 list of element types
 
std::vector< int > node_ids_
 list of node labels
 
std::vector< int > boundary_ids_
 list of surface labels
 
std::vector< int > body_ids_
 list of volume labels
 
Eigen::MatrixXi orders_
 list of geometry orders, one per cell
 
bool is_rational_ = false
 stores if the mesh is rational
 
std::vector< EdgeNodesedge_nodes_
 high-order nodes associates to edges
 
std::vector< FaceNodesface_nodes_
 high-order nodes associates to faces
 
std::vector< CellNodescell_nodes_
 high-order nodes associates to cells
 
std::vector< std::vector< double > > cell_weights_
 weights associates to cells for rational polynomail meshes
 
Eigen::VectorXi in_ordered_vertices_
 Order of the input vertices.
 
Eigen::MatrixXi in_ordered_edges_
 Order of the input edges.
 
Eigen::MatrixXi in_ordered_faces_
 Order of the input faces, TODO: change to std::vector of Eigen::Vector.
 

Detailed Description

Definition at line 18 of file Mesh3D.hpp.

Constructor & Destructor Documentation

◆ Mesh3D() [1/3]

polyfem::mesh::Mesh3D::Mesh3D ( )
default

◆ ~Mesh3D()

virtual polyfem::mesh::Mesh3D::~Mesh3D ( )
virtualdefault

◆ Mesh3D() [2/3]

polyfem::mesh::Mesh3D::Mesh3D ( Mesh3D &&  )
default

◆ Mesh3D() [3/3]

polyfem::mesh::Mesh3D::Mesh3D ( const Mesh3D )
default

Member Function Documentation

◆ barycentric_coords()

void polyfem::mesh::Mesh3D::barycentric_coords ( const RowVectorNd p,
const int  el_id,
Eigen::MatrixXd &  coord 
) const
overridevirtual

constructs barycentric coodiantes for a point p.

WARNING works only for simplices

Parameters
[in]pquery point
[in]el_idelement id
[out]coordmatrix containing the barycentric coodinates

Implements polyfem::mesh::Mesh.

Definition at line 410 of file Mesh3D.cpp.

References get_ordered_vertices_from_tet(), polyfem::mesh::Mesh::is_simplex(), and polyfem::mesh::Mesh::point().

Here is the call graph for this function:

◆ cell_edge()

virtual int polyfem::mesh::Mesh3D::cell_edge ( const int  c_id,
const int  le_id 
) const
pure virtual

Implemented in polyfem::mesh::CMesh3D, and polyfem::mesh::NCMesh3D.

Referenced by get_edges(), and polyfem::refinement::APriori::p_refine().

Here is the caller graph for this function:

◆ cell_face()

virtual int polyfem::mesh::Mesh3D::cell_face ( const int  c_id,
const int  lf_id 
) const
pure virtual

Implemented in polyfem::mesh::CMesh3D, and polyfem::mesh::NCMesh3D.

Referenced by polyfem::io::Evaluator::interpolate_boundary_function().

Here is the caller graph for this function:

◆ cell_node()

std::pair< RowVectorNd, int > polyfem::mesh::Mesh3D::cell_node ( const Navigation3D::Index index,
const int  n_new_nodes,
const int  i,
const int  j,
const int  k 
) const

Definition at line 199 of file Mesh3D.cpp.

References polyfem::mesh::Mesh::cell_barycenter(), polyfem::mesh::Mesh::cell_nodes_, polyfem::mesh::Navigation3D::Index::element, polyfem::mesh::Mesh::is_cube(), polyfem::mesh::Mesh::is_simplex(), polyfem::mesh::Mesh::orders_, polyfem::mesh::Mesh::point(), switch_edge(), switch_face(), switch_vertex(), and polyfem::mesh::Navigation3D::Index::vertex.

Referenced by polyfem::mesh::MeshNodes::node_ids_from_cell().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_cell_jacobian()

void polyfem::mesh::Mesh3D::compute_cell_jacobian ( const int  el_id,
const Eigen::MatrixXd &  reference_map,
Eigen::MatrixXd &  jacobian 
) const

Definition at line 424 of file Mesh3D.cpp.

References get_ordered_vertices_from_tet(), polyfem::mesh::Mesh::is_simplex(), and polyfem::mesh::Mesh::point().

Referenced by polyfem::solver::AMIPSForm::AMIPSForm(), and polyfem::State::set_materials().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_element_barycenters()

void polyfem::mesh::Mesh3D::compute_element_barycenters ( Eigen::MatrixXd &  barycenters) const
inlineoverridevirtual

utility for 2d/3d.

In 2d it returns face_barycenters, in 3d it returns cell_barycenters

Parameters
[out]barycentersthe barycenters

Implements polyfem::mesh::Mesh.

Definition at line 88 of file Mesh3D.hpp.

References polyfem::mesh::Mesh::cell_barycenters().

Here is the call graph for this function:

◆ edge_neighs()

virtual std::vector< uint32_t > polyfem::mesh::Mesh3D::edge_neighs ( const int  e_gid) const
pure virtual

Implemented in polyfem::mesh::CMesh3D, and polyfem::mesh::NCMesh3D.

Referenced by polyfem::refinement::APriori::p_refine().

Here is the caller graph for this function:

◆ edge_node()

std::pair< RowVectorNd, int > polyfem::mesh::Mesh3D::edge_node ( const Navigation3D::Index index,
const int  n_new_nodes,
const int  i 
) const

Definition at line 80 of file Mesh3D.cpp.

References polyfem::mesh::Navigation3D::Index::edge, polyfem::mesh::Mesh::edge_nodes_, polyfem::mesh::Navigation3D::Index::element, polyfem::mesh::Mesh::orders_, polyfem::mesh::Mesh::point(), switch_vertex(), and polyfem::mesh::Navigation3D::Index::vertex.

Referenced by polyfem::mesh::MeshNodes::node_ids_from_edge().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ elements_boxes()

void polyfem::mesh::Mesh3D::elements_boxes ( std::vector< std::array< Eigen::Vector3d, 2 > > &  boxes) const
overridevirtual

constructs a box around every element (3d cell, 2d face)

Parameters
[out]boxesaxis aligned bounding boxes

Implements polyfem::mesh::Mesh.

Definition at line 388 of file Mesh3D.cpp.

References polyfem::mesh::Mesh::cell_vertex(), polyfem::mesh::Mesh::n_cell_vertices(), polyfem::mesh::Mesh::n_elements(), and polyfem::mesh::Mesh::point().

Here is the call graph for this function:

◆ face_node()

std::pair< RowVectorNd, int > polyfem::mesh::Mesh3D::face_node ( const Navigation3D::Index index,
const int  n_new_nodes,
const int  i,
const int  j 
) const

Definition at line 99 of file Mesh3D.cpp.

References polyfem::mesh::Navigation3D::Index::element, polyfem::mesh::Navigation3D::Index::face, polyfem::mesh::Mesh::face_nodes_, polyfem::mesh::Mesh::is_cube(), polyfem::mesh::Mesh::is_simplex(), next_around_face(), polyfem::mesh::Mesh::orders_, polyfem::mesh::Mesh::point(), switch_edge(), switch_vertex(), and polyfem::mesh::Navigation3D::Index::vertex.

Referenced by polyfem::mesh::MeshNodes::node_ids_from_face().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_edge_elements_neighs()

virtual void polyfem::mesh::Mesh3D::get_edge_elements_neighs ( const int  e_id,
std::vector< int > &  ids 
) const
pure virtual

◆ get_edges() [1/2]

void polyfem::mesh::Mesh3D::get_edges ( Eigen::MatrixXd &  p0,
Eigen::MatrixXd &  p1 
) const
overridevirtual

Get all the edges.

Parameters
[out]p0edge first vertex
[out]p1edge second vertex

Implements polyfem::mesh::Mesh.

Definition at line 33 of file Mesh3D.cpp.

References polyfem::mesh::Mesh::edge_vertex(), polyfem::mesh::Mesh::n_edges(), and polyfem::mesh::Mesh::point().

Referenced by polyfem::mesh::CMesh3D::normalize(), and polyfem::refinement::APriori::p_refine().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_edges() [2/2]

void polyfem::mesh::Mesh3D::get_edges ( Eigen::MatrixXd &  p0,
Eigen::MatrixXd &  p1,
const std::vector< bool > &  valid_elements 
) const
overridevirtual

Get all the edges according to valid_elements selection.

Parameters
[out]p0edge first vertex
[out]p1edge second vertex
[in]valid_elementsflag to compute the edge

Implements polyfem::mesh::Mesh.

Definition at line 48 of file Mesh3D.cpp.

References cell_edge(), polyfem::mesh::Mesh::edge_vertex(), n_cell_edges(), and polyfem::mesh::Mesh::point().

Here is the call graph for this function:

◆ get_index_from_element() [1/2]

virtual Navigation3D::Index polyfem::mesh::Mesh3D::get_index_from_element ( int  hi) const
pure virtual

◆ get_index_from_element() [2/2]

virtual Navigation3D::Index polyfem::mesh::Mesh3D::get_index_from_element ( int  hi,
int  lf,
int  lv 
) const
pure virtual

Implemented in polyfem::mesh::CMesh3D, and polyfem::mesh::NCMesh3D.

Referenced by polyfem::mesh::extract_polyhedra(), get_ordered_vertices_from_hex(), and get_ordered_vertices_from_tet().

Here is the caller graph for this function:

◆ get_index_from_element_edge()

virtual Navigation3D::Index polyfem::mesh::Mesh3D::get_index_from_element_edge ( int  hi,
int  v0,
int  v1 
) const
pure virtual

Implemented in polyfem::mesh::CMesh3D, and polyfem::mesh::NCMesh3D.

Referenced by polyfem::basis::LagrangeBasis3d::hex_face_local_nodes(), and polyfem::basis::LagrangeBasis3d::tet_face_local_nodes().

Here is the caller graph for this function:

◆ get_index_from_element_face()

virtual Navigation3D::Index polyfem::mesh::Mesh3D::get_index_from_element_face ( int  hi,
int  v0,
int  v1,
int  v2 
) const
pure virtual

Implemented in polyfem::mesh::CMesh3D, and polyfem::mesh::NCMesh3D.

Referenced by polyfem::basis::LagrangeBasis3d::tet_face_local_nodes().

Here is the caller graph for this function:

◆ get_ordered_vertices_from_hex()

std::array< int, 8 > polyfem::mesh::Mesh3D::get_ordered_vertices_from_hex ( const int  element_index) const

Definition at line 322 of file Mesh3D.cpp.

References get_index_from_element(), polyfem::mesh::Mesh::is_cube(), and to_vertex_functions().

Referenced by polyfem::mesh::to_geogram_mesh().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_ordered_vertices_from_tet()

std::array< int, 4 > polyfem::mesh::Mesh3D::get_ordered_vertices_from_tet ( const int  element_index) const
virtual

Reimplemented in polyfem::mesh::NCMesh3D.

Definition at line 346 of file Mesh3D.cpp.

References get_index_from_element(), next_around_face(), switch_edge(), switch_face(), and switch_vertex().

Referenced by polyfem::mesh::CMesh3D::barycentric_coords(), barycentric_coords(), and compute_cell_jacobian().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_vertex_elements_neighs()

virtual void polyfem::mesh::Mesh3D::get_vertex_elements_neighs ( const int  v_id,
std::vector< int > &  ids 
) const
pure virtual

◆ is_volume()

bool polyfem::mesh::Mesh3D::is_volume ( ) const
inlineoverridevirtual

checks if mesh is volume

Returns
if mesh is volumetric

Implements polyfem::mesh::Mesh.

Definition at line 28 of file Mesh3D.hpp.

Referenced by polyfem::basis::PolygonalBasis3d::build_bases(), polyfem::basis::SplineBasis3d::build_bases(), polyfem::basis::PolygonalBasis3d::compute_integral_constraints(), and polyfem::compute_integral_constraints().

Here is the caller graph for this function:

◆ kernel()

virtual RowVectorNd polyfem::mesh::Mesh3D::kernel ( const int  cell_id) const
pure virtual

Implemented in polyfem::mesh::CMesh3D, and polyfem::mesh::NCMesh3D.

Referenced by polyfem::mesh::extract_polyhedra().

Here is the caller graph for this function:

◆ n_cell_edges()

virtual int polyfem::mesh::Mesh3D::n_cell_edges ( const int  c_id) const
pure virtual

Implemented in polyfem::mesh::CMesh3D, and polyfem::mesh::NCMesh3D.

Referenced by get_edges().

Here is the caller graph for this function:

◆ n_cell_faces()

virtual int polyfem::mesh::Mesh3D::n_cell_faces ( const int  c_id) const
pure virtual

Implemented in polyfem::mesh::CMesh3D, and polyfem::mesh::NCMesh3D.

Referenced by polyfem::mesh::extract_polyhedra(), and polyfem::io::Evaluator::interpolate_boundary_function().

Here is the caller graph for this function:

◆ next_around_edge()

virtual Navigation3D::Index polyfem::mesh::Mesh3D::next_around_edge ( Navigation3D::Index  idx) const
pure virtual

◆ next_around_face()

virtual Navigation3D::Index polyfem::mesh::Mesh3D::next_around_face ( Navigation3D::Index  idx) const
pure virtual

Implemented in polyfem::mesh::CMesh3D, and polyfem::mesh::NCMesh3D.

Referenced by polyfem::mesh::extract_polyhedra(), face_node(), get_ordered_vertices_from_tet(), polyfem::basis::LagrangeBasis3d::hex_face_local_nodes(), and polyfem::basis::LagrangeBasis3d::tet_face_local_nodes().

Here is the caller graph for this function:

◆ operator=() [1/2]

Mesh3D & polyfem::mesh::Mesh3D::operator= ( const Mesh3D )
default

◆ operator=() [2/2]

Mesh3D & polyfem::mesh::Mesh3D::operator= ( Mesh3D &&  )
default

◆ switch_edge()

virtual Navigation3D::Index polyfem::mesh::Mesh3D::switch_edge ( Navigation3D::Index  idx) const
pure virtual

Implemented in polyfem::mesh::CMesh3D, and polyfem::mesh::NCMesh3D.

Referenced by cell_node(), face_node(), get_ordered_vertices_from_tet(), to_edge_functions(), to_face_functions(), and to_vertex_functions().

Here is the caller graph for this function:

◆ switch_element()

virtual Navigation3D::Index polyfem::mesh::Mesh3D::switch_element ( Navigation3D::Index  idx) const
pure virtual

Implemented in polyfem::mesh::CMesh3D, and polyfem::mesh::NCMesh3D.

Referenced by polyfem::mesh::MeshNodes::node_ids_from_face().

Here is the caller graph for this function:

◆ switch_face()

virtual Navigation3D::Index polyfem::mesh::Mesh3D::switch_face ( Navigation3D::Index  idx) const
pure virtual

Implemented in polyfem::mesh::CMesh3D, and polyfem::mesh::NCMesh3D.

Referenced by cell_node(), get_ordered_vertices_from_tet(), to_edge_functions(), to_face_functions(), and to_vertex_functions().

Here is the caller graph for this function:

◆ switch_vertex()

virtual Navigation3D::Index polyfem::mesh::Mesh3D::switch_vertex ( Navigation3D::Index  idx) const
pure virtual

Implemented in polyfem::mesh::CMesh3D, and polyfem::mesh::NCMesh3D.

Referenced by cell_node(), edge_node(), face_node(), get_ordered_vertices_from_tet(), polyfem::basis::LagrangeBasis3d::hex_face_local_nodes(), polyfem::basis::LagrangeBasis3d::tet_face_local_nodes(), to_edge_functions(), to_face_functions(), and to_vertex_functions().

Here is the caller graph for this function:

◆ to_edge_functions()

void polyfem::mesh::Mesh3D::to_edge_functions ( std::array< std::function< Navigation3D::Index(Navigation3D::Index)>, 12 > &  to_edge) const

Definition at line 298 of file Mesh3D.cpp.

References switch_edge(), switch_face(), and switch_vertex().

Here is the call graph for this function:

◆ to_face_functions()

void polyfem::mesh::Mesh3D::to_face_functions ( std::array< std::function< Navigation3D::Index(Navigation3D::Index)>, 6 > &  to_face) const

Definition at line 267 of file Mesh3D.cpp.

References switch_edge(), switch_face(), and switch_vertex().

Here is the call graph for this function:

◆ to_vertex_functions()

void polyfem::mesh::Mesh3D::to_vertex_functions ( std::array< std::function< Navigation3D::Index(Navigation3D::Index)>, 8 > &  to_vertex) const

Definition at line 285 of file Mesh3D.cpp.

References switch_edge(), switch_face(), and switch_vertex().

Referenced by get_ordered_vertices_from_hex().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ tri_area()

double polyfem::mesh::Mesh3D::tri_area ( const int  gid) const
overridevirtual

area of a tri face of a tet mesh

Parameters
[in]gidglobal face id
Returns
tet area

Reimplemented from polyfem::mesh::Mesh.

Definition at line 18 of file Mesh3D.cpp.

References polyfem::mesh::Mesh::face_vertex(), polyfem::mesh::Mesh::n_face_vertices(), polyfem::mesh::Mesh::n_vertices(), and polyfem::mesh::Mesh::point().

Here is the call graph for this function:

◆ vertex_neighs()

virtual std::vector< uint32_t > polyfem::mesh::Mesh3D::vertex_neighs ( const int  v_gid) const
pure virtual

The documentation for this class was generated from the following files: