shards  Version of the Day
 All Classes Functions Variables Typedefs Enumerations Enumerator Groups
Polytope Cell Topology Descriptions

Topological description of polytope cells (polygons and polyhedrons) in terms of their subcells. More...

Classes

struct  shards::Node
 Topological traits: Dimension = 0, Vertices = 0, Nodes = 0. More...
 
struct  shards::Particle
 Topological traits: Dimension = 0, Vertices = 1, Nodes = 1. More...
 
struct  shards::Line< NodeCount >
 Topological traits: Dimension = 1, Vertices = 2, Nodes = 2 or 3. More...
 
struct  shards::Beam< NodeCount >
 Topological traits: Dimension = 2, Edges = 1, Vertices = 2, and Nodes = 2 or 3. More...
 
struct  shards::ShellLine< NodeCount >
 Topological traits: Dimension = 2, Edges = 2, Vertices = 2, and Nodes = 2 or 3. More...
 
struct  shards::Triangle< NodeCount >
 Topological traits: Dimension = 2, Edges = 3, Vertices = 3, and Nodes = 3 or 6. More...
 
struct  shards::ShellTriangle< NodeCount >
 Topological traits: Dimension = 3, Sides = 2, Edges = 3, Vertices = 3, and Nodes = 3 or 6. More...
 
struct  shards::Quadrilateral< NodeCount >
 Topological traits: Dimension = 2, Edges = 4, Vertices = 4, and Nodes = 4, 8, or 9. More...
 
struct  shards::ShellQuadrilateral< NodeCount >
 Topological traits: Dimension = 2, Sides = 2, Edges = 4, Vertices = 4, and Nodes = 4, 8, or 9. More...
 
struct  shards::Tetrahedron< NodeCount >
 Topological traits: Dimension = 3, Sides = 4, Edges = 6, Vertices = 4, and Nodes = 4 or 10. More...
 
struct  shards::Pyramid< NodeCount >
 Topological traits: Dimension = 3, Sides = 5, Edges = 8, Vertices = 5, and Nodes = 5, 13, or 14. More...
 
struct  shards::Wedge< NodeCount >
 Topological traits: Dimension = 3, Sides = 5, Edges = 9, Vertices = 6, and Nodes = 6, 15, or 18. More...
 
struct  shards::Hexahedron< NodeCount >
 Topological traits: Dimension = 3, Sides = 6, Edges = 12, Vertices = 8, and Nodes = 8, 20, or 27. More...
 
class  shards::CellTopology
 Provide input checked access (in debug mode) to cell topology data and a procedure to create custom cell topologies. More...
 
struct  CellTopologyData
 A simple 'C' struct of cell topology attributes. More...
 
struct  CellTopologyData_Subcell
 Subcell information. More...
 
struct  CellTopologyData_Permutation
 Array of node permutations. More...
 
struct  shards::CellTopologyTraits< Dimension, Number_Vertex, Number_Node, EdgeList, EdgeMaps, FaceList, FaceMaps, PermutationMaps, PermutationPolarity >
 Compile-time traits for a cell topology. More...
 

Typedefs

typedef struct CellTopologyData CellTopologyData
 Self-typedef.
 

Enumerations

enum  shards::ECellType {
  ALL_CELLS = 0,
  STANDARD_CELL,
  NONSTANDARD_CELL
}
 Enumeration of cell types in Shards.
 
enum  shards::ETopologyType {
  ALL_TOPOLOGIES,
  BASE_TOPOLOGY,
  EXTENDED_TOPOLOGY
}
 Enumeration of topology types in Shards.
 
enum  {
  CELL_PERMUTATION_POLARITY_IRRELEVANT = 0,
  CELL_PERMUTATION_POLARITY_POSITIVE = 1,
  CELL_PERMUTATION_POLARITY_NEGATIVE = 2
}
 Values for the CellTopologyData_Permutation polarity.
 

Functions

template<>
const CellTopologyDatashards::getCellTopologyData< Node > ()
 Singleton for Node topology.
 
template<>
const CellTopologyDatashards::getCellTopologyData< Particle > ()
 Singleton for Particle topology.
 
template<>
const CellTopologyDatashards::getCellTopologyData< Line< 2 > > ()
 Singleton for line topology with two nodes.
 
template<>
const CellTopologyDatashards::getCellTopologyData< Line< 3 > > ()
 Singleton for line topology with three nodes.
 
template<>
const CellTopologyDatashards::getCellTopologyData< Beam< 2 > > ()
 Singleton for beam topology with two nodes.
 
template<>
const CellTopologyDatashards::getCellTopologyData< Beam< 3 > > ()
 Singleton for beam topology with three nodes.
 
template<>
const CellTopologyDatashards::getCellTopologyData< ShellLine< 2 > > ()
 Singleton for shell-line topology with two nodes.
 
template<>
const CellTopologyDatashards::getCellTopologyData< ShellLine< 3 > > ()
 Singleton for shell-line topology with three nodes.
 
template<>
const CellTopologyDatashards::getCellTopologyData< Triangle< 3 > > ()
 Return CellTopologyData singleton for the Triangle<3>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Triangle< 6 > > ()
 Return CellTopologyData singleton for the Triangle<6>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Triangle< 4 > > ()
 Return CellTopologyData singleton for the Triangle<4> (Face of Tet8)
 
template<>
const CellTopologyDatashards::getCellTopologyData< ShellTriangle< 3 > > ()
 Return CellTopologyData singleton for the ShellTriangle<3>
 
template<>
const CellTopologyDatashards::getCellTopologyData< ShellTriangle< 6 > > ()
 Return CellTopologyData singleton for the ShellTriangle<6>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Quadrilateral< 4 > > ()
 Return CellTopologyData singleton for the Quadrilateral<4>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Quadrilateral< 8 > > ()
 Return CellTopologyData singleton for the Quadrilateral<8>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Quadrilateral< 9 > > ()
 Return CellTopologyData singleton for the Quadrilateral<9>
 
template<>
const CellTopologyDatashards::getCellTopologyData< ShellQuadrilateral< 4 > > ()
 Return CellTopologyData singleton for the ShellQuadrilateral<4>
 
template<>
const CellTopologyDatashards::getCellTopologyData< ShellQuadrilateral< 8 > > ()
 Return CellTopologyData singleton for the ShellQuadrilateral<8>
 
template<>
const CellTopologyDatashards::getCellTopologyData< ShellQuadrilateral< 9 > > ()
 Return CellTopologyData singleton for the ShellQuadrilateral<9>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Tetrahedron< 4 > > ()
 Return CellTopologyData singleton for the Tetrahedron<4>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Tetrahedron< 10 > > ()
 Return CellTopologyData singleton for the Tetrahedron<10>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Tetrahedron< 8 > > ()
 Return CellTopologyData singleton for the Tetrahedron<8>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Tetrahedron< 11 > > ()
 Return CellTopologyData singleton for the Tetrahedron<11>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Pyramid< 5 > > ()
 Return CellTopologyData singleton for the Pyramid<5>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Pyramid< 13 > > ()
 Return CellTopologyData singleton for the Pyramid<13>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Pyramid< 14 > > ()
 Return CellTopologyData singleton for the Pyramid<14>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Wedge< 6 > > ()
 Return CellTopologyData singleton for the Wedge<6>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Wedge< 15 > > ()
 Return CellTopologyData singleton for the Wedge<15>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Wedge< 18 > > ()
 Return CellTopologyData singleton for the Wedge<18>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Hexahedron< 8 > > ()
 Return CellTopologyData singleton for the Hexahedron<8>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Hexahedron< 20 > > ()
 Return CellTopologyData singleton for the Hexahedron<20>
 
template<>
const CellTopologyDatashards::getCellTopologyData< Hexahedron< 27 > > ()
 Return CellTopologyData singleton for the Hexahedron<27>
 
std::ostream & shards::operator<< (std::ostream &, const CellTopology &)
 Overloaded << operator for CellTopologyData objects. More...
 
std::string shards::ECellTypeToString (ECellType cellType)
 
std::string shards::ETopologyTypeToString (ETopologyType topologyType)
 
void shards::getTopologies (std::vector< shards::CellTopology > &topologies, const unsigned cellDim=4, const ECellType cellType=ALL_CELLS, const ETopologyType topologyType=ALL_TOPOLOGIES)
 Returns an std::vector with all cell topologies that meet the specified selection flags. More...
 
int shards::isPredefinedCell (const CellTopology &cell)
 Checks if the cell topology is predefined in shards. More...
 
template<typename id_type >
int shards::findPermutation (const CellTopologyData &top, const id_type *const expected_node, const id_type *const actual_node)
 
template<typename id_type >
int shards::findPermutation (const CellTopology &top, const id_type *const expected_node, const id_type *const actual_node)
 
void shards::badCellTopologyKey (const unsigned dimension, const unsigned face_count, const unsigned edge_count, const unsigned vertex_count, const unsigned node_count)
 Generates detailed message if one or more input parameters are out of their admissible bounds. More...
 
unsigned shards::cellTopologyKey (const unsigned dimension, const unsigned face_count, const unsigned edge_count, const unsigned vertex_count, const unsigned node_count)
 Generate integer key from topological dimensions. More...
 
bool shards::operator== (const CellTopology &left, const CellTopology &right)
 
bool shards::operator< (const CellTopology &left, const CellTopology &right)
 
bool shards::operator!= (const CellTopology &left, const CellTopology &right)
 
int mapCellFaceEdge (const CellTopologyData *cell_topology, unsigned face_ordinal, unsigned face_edge_ordinal)
 Map a cell->face->edge ordinal to the cell->edge ordinal. Return -1 for erroneous input.
 
template<class Traits >
const CellTopologyDatashards::getCellTopologyData ()
 Return a CellTopology singleton for the given cell topology traits.
 

Detailed Description

Topological description of polytope cells (polygons and polyhedrons) in terms of their subcells.

Polytope (2D polygon or 3D polyhedron)

A cell topology defines a polytope in terms of its vertex, edge, and side connectivity. However, a cell topology does not define a geometric realization of that polytope - it does not associate spatial coordinates with any of the vertices.

Dimension (of the polytope)

The dimension of cell topology is the smallest spatial dimension in which that cell can be realized. For example, a tetrahedron or quadrilateral shell can only be realized in a 3D space.

Base and Extended Cell Topologies

The base cell topology of a polytope associates exactly one $ {Cell}^{0} $ subcell with each vertex of the polytope. These vertices are given a canonical ordering from which vertex-ordinals are defined. For example, the vertex-ordinals for quadrilateral are defined counter-clockwise around the quadrilateral as (0,1,2,3). An extended cell topology is defined associating additional $ {Cell}^{0} $ subcells with the polytope's edges, sides, or interior. These additional $ {Cell}^{0} $ subcells are always ordered after the vertex subcells.

Node, Edge, Side, and Element Cells

Subcell Orientation

A subcell $ {Cell}^{d} $ of a parent $ {Cell}^{D} $, where $ 0 < d < D $, has vertices that are associated with the parent cell's vertices. This association defines a map from the subcell's vertex-ordinals to the parent's vertex-ordinals. For example, the four edge subcells of a quadrilateral have a vertex-ordinal map of ( (0,1), (1,2), (2,3), (3,0) ).

Function Documentation

std::ostream& shards::operator<< ( std::ostream &  ,
const CellTopology &   
)

Overloaded << operator for CellTopologyData objects.

Overloaded << operator for CellTopology objects.

void shards::getTopologies ( std::vector< shards::CellTopology > &  topologies,
const unsigned  cellDim = 4,
const ECellType  cellType = ALL_CELLS,
const ETopologyType  topologyType = ALL_TOPOLOGIES 
)

Returns an std::vector with all cell topologies that meet the specified selection flags.

Parameters
topologies[out] - vector with all topologies that
cellDim[in] - cell dimension; 0, 1, 2, 3, or 4 (default = all dimensions)
cellType[in] - cell type: default = ALL_CELLS
topologyType[in] - topology type: default = ALL_TOPOLOGIES
int shards::isPredefinedCell ( const CellTopology &  cell)

Checks if the cell topology is predefined in shards.

Parameters
cell[in] - cell topology
Returns
1 if the cell topology is defined in shards, 0 if it is a custom, user-defined cell-topology
void shards::badCellTopologyKey ( const unsigned  dimension,
const unsigned  face_count,
const unsigned  edge_count,
const unsigned  vertex_count,
const unsigned  node_count 
)

Generates detailed message if one or more input parameters are out of their admissible bounds.

Parameters
dimension[in] - maximum value = 7
face_count[in] - maximum value = 63
edge_count[in] - maximum value = 63
vertex_count[in] - maximum value = 63
node_count[in] - maximum value = 1023
unsigned shards::cellTopologyKey ( const unsigned  dimension,
const unsigned  face_count,
const unsigned  edge_count,
const unsigned  vertex_count,
const unsigned  node_count 
)
inline

Generate integer key from topological dimensions.

Parameters
dimension[in] - maximum value = 7 (3 bits)
face_count[in] - maximum value = 63 (6 bits)
edge_count[in] - maximum value = 63 (6 bits)
vertex_count[in] - maximum value = 63 (6 bits)
node_count[in] - maximum value = 1023 (10 bits) The key uses all but the first bit in a 32 bit unsigned.

Definition at line 649 of file Shards_CellTopology.hpp.