DOLFIN
DOLFIN C++ interface
Public Member Functions | List of all members
dolfin::Cell Class Reference

A Cell is a MeshEntity of topological codimension 0. More...

#include <Cell.h>

Inheritance diagram for dolfin::Cell:
Inheritance graph
[legend]
Collaboration diagram for dolfin::Cell:
Collaboration graph
[legend]

Public Member Functions

 Cell ()
 Create empty cell.
 
 Cell (const Mesh &mesh, std::size_t index)
 
 ~Cell ()
 Destructor.
 
CellType::Type type () const
 Return type of cell.
 
std::size_t num_vertices () const
 Return number of vertices of cell.
 
std::size_t orientation () const
 
std::size_t orientation (const Point &up) const
 
double volume () const
 
double h () const
 
double circumradius () const
 
double inradius () const
 
double radius_ratio () const
 
double squared_distance (const Point &point) const
 
double distance (const Point &point) const
 
double normal (std::size_t facet, std::size_t i) const
 
Point normal (std::size_t facet) const
 
Point cell_normal () const
 
double facet_area (std::size_t facet) const
 
void order (const std::vector< std::int64_t > &local_to_global_vertex_indices)
 
bool ordered (const std::vector< std::int64_t > &local_to_global_vertex_indices) const
 
bool contains (const Point &point) const
 
bool collides (const Point &point) const
 
bool collides (const MeshEntity &entity) const
 
std::vector< Pointintersection (const MeshEntity &entity) const
 
void get_coordinate_dofs (std::vector< double > &coordinates) const
 Get cell coordinate dofs (not vertex coordinates)
 
void get_vertex_coordinates (std::vector< double > &coordinates) const
 Get cell vertex coordinates (not coordinate dofs)
 
void get_cell_data (ufc::cell &ufc_cell, int local_facet=-1) const
 Fill UFC cell with miscellaneous data.
 
void get_cell_topology (ufc::cell &ufc_cell) const
 Fill UFC cell with topology data.
 
- Public Member Functions inherited from dolfin::MeshEntity
 MeshEntity ()
 Default Constructor.
 
 MeshEntity (const Mesh &mesh, std::size_t dim, std::size_t index)
 
virtual ~MeshEntity ()
 Destructor.
 
void init (const Mesh &mesh, std::size_t dim, std::size_t index)
 
bool operator== (const MeshEntity &e) const
 
bool operator!= (const MeshEntity &e) const
 
const Meshmesh () const
 
std::size_t dim () const
 
std::size_t index () const
 
std::int64_t global_index () const
 
std::size_t num_entities (std::size_t dim) const
 
std::size_t num_global_entities (std::size_t dim) const
 
const unsigned int * entities (std::size_t dim) const
 
std::size_t mesh_id () const
 
bool incident (const MeshEntity &entity) const
 
std::size_t index (const MeshEntity &entity) const
 
Point midpoint () const
 
bool is_ghost () const
 
std::set< unsigned int > sharing_processes () const
 
bool is_shared () const
 
unsigned int owner () const
 
std::string str (bool verbose) const
 

Additional Inherited Members

- Protected Attributes inherited from dolfin::MeshEntity
Mesh const * _mesh
 
std::size_t _dim
 
std::size_t _local_index
 

Detailed Description

A Cell is a MeshEntity of topological codimension 0.

Constructor & Destructor Documentation

◆ Cell()

dolfin::Cell::Cell ( const Mesh mesh,
std::size_t  index 
)
inline

Create cell on given mesh with given index

Parameters
meshThe mesh.
indexThe index.

Member Function Documentation

◆ cell_normal()

Point dolfin::Cell::cell_normal ( ) const
inline

Compute normal to cell itself (viewed as embedded in 3D)

Returns
Point Normal of the cell

◆ circumradius()

double dolfin::Cell::circumradius ( ) const
inline

Compute circumradius of cell

Returns
double The circumradius of the cell.
UnitSquareMesh mesh(1, 1);
Cell cell(mesh, 0);
info("%g", cell.circumradius());

◆ collides() [1/2]

bool Cell::collides ( const Point point) const

Check whether given point collides with cell

Parameters
pointThe point to be checked.
Returns
bool True iff point collides with cell.

◆ collides() [2/2]

bool Cell::collides ( const MeshEntity entity) const

Check whether given entity collides with cell

Parameters
entityThe cell to be checked.
Returns
bool True iff entity collides with cell.

◆ contains()

bool Cell::contains ( const Point point) const

Check whether given point is contained in cell. This function is identical to the function collides(point).

Parameters
pointThe point to be checked.
Returns
bool True iff point is contained in cell.

◆ distance()

double dolfin::Cell::distance ( const Point point) const
inline

Compute distance to given point.

Parameters
pointThe point.
Returns
double The distance to the point.

◆ facet_area()

double dolfin::Cell::facet_area ( std::size_t  facet) const
inline

Compute the area/length of given facet with respect to the cell

Parameters
facetIndex of the facet.
Returns
double Area/length of the facet.

◆ h()

double dolfin::Cell::h ( ) const
inline

Compute greatest distance between any two vertices

Returns
double The greatest distance between any two vertices of the cell.
UnitSquareMesh mesh(1, 1);
Cell cell(mesh, 0);
info("%g", cell.h());

◆ inradius()

double dolfin::Cell::inradius ( ) const
inline

Compute inradius of cell

Returns
double Radius of the sphere inscribed in the cell.
UnitSquareMesh mesh(1, 1);
Cell cell(mesh, 0);
info("%g", cell.inradius());

◆ intersection()

std::vector< Point > Cell::intersection ( const MeshEntity entity) const

Compute triangulation of intersection with given entity

Parameters
entityThe entity with which to intersect.
Returns
std::vector<Point> A flattened array of simplices of dimension num_simplices x num_vertices x gdim = num_simplices x (tdim + 1) x gdim

◆ normal() [1/2]

double dolfin::Cell::normal ( std::size_t  facet,
std::size_t  i 
) const
inline

Compute component i of normal of given facet with respect to the cell

Parameters
facetIndex of facet.
iComponent.
Returns
double Component i of the normal of the facet.

◆ normal() [2/2]

Point dolfin::Cell::normal ( std::size_t  facet) const
inline

Compute normal of given facet with respect to the cell

Parameters
facetIndex of facet.
Returns
Point Normal of the facet.

◆ order()

void dolfin::Cell::order ( const std::vector< std::int64_t > &  local_to_global_vertex_indices)
inline

Order entities locally

Parameters
local_to_global_vertex_indicesThe global vertex indices.

◆ ordered()

bool dolfin::Cell::ordered ( const std::vector< std::int64_t > &  local_to_global_vertex_indices) const
inline

Check if entities are ordered

Parameters
local_to_global_vertex_indicesThe global vertex indices.
Returns
bool True iff ordered.

◆ orientation() [1/2]

std::size_t dolfin::Cell::orientation ( ) const
inline

Compute orientation of cell

Returns
std::size_t Orientation of the cell (0 is 'up'/'right', 1 is 'down'/'left')

◆ orientation() [2/2]

std::size_t dolfin::Cell::orientation ( const Point up) const
inline

Compute orientation of cell relative to given 'up' direction

Parameters
upThe direction defined as 'up'
Returns
std::size_t Orientation of the cell (0 is 'same', 1 is 'opposite')

◆ radius_ratio()

double dolfin::Cell::radius_ratio ( ) const
inline

Compute ratio of inradius to circumradius times dim for cell. Useful as cell quality measure. Returns 1. for equilateral and 0. for degenerate cell. See Jonathan Richard Shewchuk: What Is a Good Linear Finite Element?, online: http://www.cs.berkeley.edu/~jrs/papers/elemj.pdf

Returns
double topological_dimension * inradius / circumradius
UnitSquareMesh mesh(1, 1);
Cell cell(mesh, 0);
info("%g", cell.radius_ratio());

◆ squared_distance()

double dolfin::Cell::squared_distance ( const Point point) const
inline

Compute squared distance to given point.

Parameters
pointThe point.
Returns
double The squared distance to the point.

◆ volume()

double dolfin::Cell::volume ( ) const
inline

Compute (generalized) volume of cell

Returns
double The volume of the cell.
UnitSquare mesh(1, 1);
Cell cell(mesh, 0);
info("%g", cell.volume());

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