OpenMEEG
Loading...
Searching...
No Matches
OpenMEEG::Mesh Class Reference

#include <mesh.h>

Public Types

typedef std::map< const Vertex *, TrianglesRefsVertexTriangles
 

Public Member Functions

 Mesh (Geometry *geometry=nullptr)
 Default constructor or constructor using a provided geometry.
 
 Mesh (const unsigned nv, const unsigned nt, Geometry *geometry=nullptr)
 Constructor from scratch (vertices/triangles t be added)
 
 Mesh (const Mesh &)=default
 
 Mesh (Mesh &&m)=default
 
Meshoperator= (const Mesh &)=default
 
 Mesh (const std::string &filename, const bool verbose, Geometry *geometry=nullptr)
 Constructors.
 
 Mesh (const std::string &filename, Geometry *geometry=nullptr)
 
 ~Mesh ()
 Destructor.
 
std::string & name ()
 
const std::string & name () const
 
VerticesRefsvertices ()
 
const VerticesRefsvertices () const
 
Geometrygeometry () const
 
Trianglestriangles ()
 
const Trianglestriangles () const
 
TriangleIndices triangle (const Triangle &t) const
 
bool current_barrier () const
 
bool & current_barrier ()
 
bool isolated () const
 
bool & isolated ()
 
Triangleadd_triangle (const TriangleIndices inds)
 Add a triangle specified by its indices in the geometry.
 
Triangleadd_triangle (const TriangleIndices inds, const IndexMap &indmap)
 
void add (const std::vector< TriangleIndices > &trgs)
 
void add (const std::vector< TriangleIndices > &trgs, const IndexMap &indmap)
 
bool operator== (const Mesh &m) const
 
bool operator!= (const Mesh &m) const
 
void info (const bool verbose=false) const
 Print info Print to std::cout some info about the mesh.
 
bool has_self_intersection () const
 Check whether the mesh self-intersects.
 
bool intersection (const Mesh &) const
 Check whether the mesh intersects another mesh.
 
bool has_correct_orientation () const
 Check local orientation of mesh triangles.
 
void generate_indices ()
 Generate indices (if allocate).
 
void update (const bool topology_changed)
 Recompute triangles normals, area, and vertex triangles.
 
void merge (const Mesh &, const Mesh &)
 Merge two meshes.
 
Ranges vertices_ranges () const
 Get the ranges of the specific mesh in the global matrix.
 
Range triangles_range () const
 
TrianglesRefs triangles (const Vertex &V) const
 Get the triangles adjacent to vertex.
 
TrianglesRefs adjacent_triangles (const Triangle &triangle) const
 Get the triangles adjacent to.
 
void change_orientation ()
 Change mesh orientation.
 
void correct_local_orientation ()
 Correct the local orientation of the mesh triangles.
 
void correct_global_orientation ()
 Correct the global orientation (if there is one).
 
double solid_angle (const Vect3 &p) const
 Given a point p, computes the solid angle of the mesh seen from.
 
Normal normal (const Vertex &v) const
 Get normal at vertex.`.
 
void laplacian (SymMatrix &A) const
 Compute mesh laplacian.
 
bool & outermost ()
 
bool outermost () const
 Returns True if it is an outermost mesh.
 
void smooth (const double &smoothing_intensity, const unsigned &niter)
 Smooth Mesh.
 
void gradient_norm2 (SymMatrix &A) const
 Compute the square norm of the surfacic gradient.
 
void load (const std::string &filename, const bool verbose=true)
 Read mesh from file.
 
void save (const std::string &filename) const
 Save mesh to file.
 

Friends

class Geometry
 
class MeshIO
 

Detailed Description

Definition at line 34 of file mesh.h.

Member Typedef Documentation

◆ VertexTriangles

Definition at line 43 of file mesh.h.

Constructor & Destructor Documentation

◆ Mesh() [1/6]

OpenMEEG::Mesh::Mesh ( Geometry * geometry = nullptr)
inline

Default constructor or constructor using a provided geometry.

Parameters
geometry

Definition at line 48 of file mesh.h.

◆ Mesh() [2/6]

OpenMEEG::Mesh::Mesh ( const unsigned nv,
const unsigned nt,
Geometry * geometry = nullptr )

Constructor from scratch (vertices/triangles t be added)

Parameters
nvspace to allocate for vertices
ntspace to allocate for triangles
geometrythe geometry to use

◆ Mesh() [3/6]

OpenMEEG::Mesh::Mesh ( const Mesh & )
default

◆ Mesh() [4/6]

OpenMEEG::Mesh::Mesh ( Mesh && m)
default

◆ Mesh() [5/6]

OpenMEEG::Mesh::Mesh ( const std::string & filename,
const bool verbose,
Geometry * geometry = nullptr )
inline

Constructors.

Parameters
filenamemesh file name
verboseverbose mode
geometrygeometry

Definition at line 68 of file mesh.h.

◆ Mesh() [6/6]

OpenMEEG::Mesh::Mesh ( const std::string & filename,
Geometry * geometry = nullptr )
inline
Parameters
filenamemesh file name
geometrygeometry

Definition at line 75 of file mesh.h.

◆ ~Mesh()

OpenMEEG::Mesh::~Mesh ( )
inline

Destructor.

Definition at line 79 of file mesh.h.

Member Function Documentation

◆ operator=()

Mesh & OpenMEEG::Mesh::operator= ( const Mesh & )
default

◆ name() [1/2]

std::string & OpenMEEG::Mesh::name ( )
inline
Returns
the mesh name

Definition at line 81 of file mesh.h.

◆ name() [2/2]

const std::string & OpenMEEG::Mesh::name ( ) const
inline
Returns
the mesh name

Definition at line 82 of file mesh.h.

◆ vertices() [1/2]

VerticesRefs & OpenMEEG::Mesh::vertices ( )
inline
Returns
the vector of pointers to the mesh vertices

Definition at line 84 of file mesh.h.

◆ vertices() [2/2]

const VerticesRefs & OpenMEEG::Mesh::vertices ( ) const
inline
Returns
the vector of pointers to the mesh vertices

Definition at line 85 of file mesh.h.

◆ geometry()

Geometry & OpenMEEG::Mesh::geometry ( ) const
inline

Definition at line 87 of file mesh.h.

◆ triangles() [1/3]

Triangles & OpenMEEG::Mesh::triangles ( )
inline
Returns
the triangles of the mesh

Definition at line 89 of file mesh.h.

◆ triangles() [2/3]

const Triangles & OpenMEEG::Mesh::triangles ( ) const
inline
Returns
the triangles of the mesh

Definition at line 90 of file mesh.h.

◆ triangle()

TriangleIndices OpenMEEG::Mesh::triangle ( const Triangle & t) const

◆ current_barrier() [1/2]

bool OpenMEEG::Mesh::current_barrier ( ) const
inline

Definition at line 94 of file mesh.h.

◆ current_barrier() [2/2]

bool & OpenMEEG::Mesh::current_barrier ( )
inline

Definition at line 95 of file mesh.h.

◆ isolated() [1/2]

bool OpenMEEG::Mesh::isolated ( ) const
inline

Definition at line 96 of file mesh.h.

◆ isolated() [2/2]

bool & OpenMEEG::Mesh::isolated ( )
inline

Definition at line 97 of file mesh.h.

◆ add_triangle() [1/2]

Triangle & OpenMEEG::Mesh::add_triangle ( const TriangleIndices inds)

Add a triangle specified by its indices in the geometry.

◆ add_triangle() [2/2]

Triangle & OpenMEEG::Mesh::add_triangle ( const TriangleIndices inds,
const IndexMap & indmap )
inline

Definition at line 103 of file mesh.h.

◆ add() [1/2]

void OpenMEEG::Mesh::add ( const std::vector< TriangleIndices > & trgs)
inline

Definition at line 108 of file mesh.h.

◆ add() [2/2]

void OpenMEEG::Mesh::add ( const std::vector< TriangleIndices > & trgs,
const IndexMap & indmap )
inline

Definition at line 113 of file mesh.h.

◆ operator==()

bool OpenMEEG::Mesh::operator== ( const Mesh & m) const
inline

Definition at line 118 of file mesh.h.

◆ operator!=()

bool OpenMEEG::Mesh::operator!= ( const Mesh & m) const
inline

Definition at line 119 of file mesh.h.

◆ info()

void OpenMEEG::Mesh::info ( const bool verbose = false) const

Print info Print to std::cout some info about the mesh.

Returns
void
See also

Print mesh information.

◆ has_self_intersection()

bool OpenMEEG::Mesh::has_self_intersection ( ) const

Check whether the mesh self-intersects.

◆ intersection()

bool OpenMEEG::Mesh::intersection ( const Mesh & ) const

Check whether the mesh intersects another mesh.

◆ has_correct_orientation()

bool OpenMEEG::Mesh::has_correct_orientation ( ) const

Check local orientation of mesh triangles.

◆ generate_indices()

void OpenMEEG::Mesh::generate_indices ( )

Generate indices (if allocate).

◆ update()

void OpenMEEG::Mesh::update ( const bool topology_changed)

Recompute triangles normals, area, and vertex triangles.

◆ merge()

void OpenMEEG::Mesh::merge ( const Mesh & ,
const Mesh &  )

Merge two meshes.

◆ vertices_ranges()

Ranges OpenMEEG::Mesh::vertices_ranges ( ) const
inline

Get the ranges of the specific mesh in the global matrix.

Returns
vector of Range
See also

Definition at line 140 of file mesh.h.

◆ triangles_range()

Range OpenMEEG::Mesh::triangles_range ( ) const
inline

Definition at line 157 of file mesh.h.

◆ triangles() [3/3]

TrianglesRefs OpenMEEG::Mesh::triangles ( const Vertex & V) const
inline

Get the triangles adjacent to vertex.

Parameters
V.

Definition at line 161 of file mesh.h.

◆ adjacent_triangles()

TrianglesRefs OpenMEEG::Mesh::adjacent_triangles ( const Triangle & triangle) const
inline

Get the triangles adjacent to.

Parameters
triangle.

Definition at line 165 of file mesh.h.

◆ change_orientation()

void OpenMEEG::Mesh::change_orientation ( )
inline

Change mesh orientation.

Definition at line 177 of file mesh.h.

◆ correct_local_orientation()

void OpenMEEG::Mesh::correct_local_orientation ( )

Correct the local orientation of the mesh triangles.

◆ correct_global_orientation()

void OpenMEEG::Mesh::correct_global_orientation ( )

Correct the global orientation (if there is one).

◆ solid_angle()

double OpenMEEG::Mesh::solid_angle ( const Vect3 & p) const

Given a point p, computes the solid angle of the mesh seen from.

Parameters
p.

◆ normal()

Normal OpenMEEG::Mesh::normal ( const Vertex & v) const

Get normal at vertex.`.

◆ laplacian()

void OpenMEEG::Mesh::laplacian ( SymMatrix & A) const

Compute mesh laplacian.

◆ outermost() [1/2]

bool & OpenMEEG::Mesh::outermost ( )
inline

Definition at line 188 of file mesh.h.

◆ outermost() [2/2]

bool OpenMEEG::Mesh::outermost ( ) const
inline

Returns True if it is an outermost mesh.

Definition at line 189 of file mesh.h.

◆ smooth()

void OpenMEEG::Mesh::smooth ( const double & smoothing_intensity,
const unsigned & niter )

Smooth Mesh.

Parameters
smoothing_intensity
niter
Returns
void

◆ gradient_norm2()

void OpenMEEG::Mesh::gradient_norm2 ( SymMatrix & A) const

Compute the square norm of the surfacic gradient.

◆ load()

void OpenMEEG::Mesh::load ( const std::string & filename,
const bool verbose = true )

Read mesh from file.

Parameters
filenamecan be .vtk, .tri (ascii), .off, .bnd or .mesh. Be verbose if
verboseis true. The difference between read and load is that read just reads the file and does not update the geometry. Read has to be used when multiple meshes are used in a geometry. load reads a mesh.

◆ save()

void OpenMEEG::Mesh::save ( const std::string & filename) const

Save mesh to file.

Parameters
filenamecan be .vtk, .tri (ascii), .bnd, .off or .mesh

Friends And Related Symbol Documentation

◆ Geometry

friend class Geometry
friend

Definition at line 40 of file mesh.h.

◆ MeshIO

friend class MeshIO
friend

Definition at line 41 of file mesh.h.


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