Intrepid2
Public Types | Public Member Functions | Static Public Member Functions | Protected Attributes | List of all members
Intrepid2::CellTopology Class Reference

Implements arbitrary-dimensional extrusion of a base shards::CellTopology. More...

#include <Intrepid2_CellTopology.hpp>

Public Types

using CellTopoPtr = Teuchos::RCP< CellTopology >
 
using CellTopologyKey = std::pair< ordinal_type, ordinal_type >
 

Public Member Functions

 CellTopology (const shards::CellTopology &baseTopo, ordinal_type tensorialDegree)
 Constructor. More...
 
const shards::CellTopology & getBaseTopology () const
 Returns the underlying shards CellTopology.
 
ordinal_type getTensorialDegree () const
 The number of times we have taken a tensor product between a line topology and the shards topology to form this cell topology.
 
ordinal_type getDimension () const
 Dimension of this tensor topology.
 
ordinal_type getNodeCount () const
 Node count of this cell topology.
 
ordinal_type getVertexCount () const
 Vertex count of this cell topology.
 
ordinal_type getEdgeCount () const
 Edge (dimension 1) subcell count of this cell topology.
 
ordinal_type getFaceCount () const
 Face (dimension 2) subcell count of this cell topology.
 
ordinal_type getSideCount () const
 Side (dimension N-1) subcell count of this cell topology.
 
std::pair< ordinal_type,
ordinal_type > 
getKey () const
 Key that's unique for standard shards topologies and any tensorial degree.
 
ordinal_type getNodeCount (const ordinal_type subcell_dim, const ordinal_type subcell_ord) const
 Node count of a subcell of the given dimension and ordinal. More...
 
std::string getName () const
 Human-readable name of the CellTopology.
 
ordinal_type getVertexCount (const ordinal_type subcell_dim, const ordinal_type subcell_ord) const
 Vertex count of a subcell of the given dimension and ordinal. More...
 
ordinal_type getEdgeCount (const ordinal_type subcell_dim, const ordinal_type subcell_ord) const
 Edge count of a subcell of the given dimension and ordinal. More...
 
ordinal_type getExtrudedSubcellOrdinal (const ordinal_type subcell_dim_in_component_topo, const ordinal_type subcell_ord_in_component_topo) const
 Mapping from the tensorial component CellTopology's subcell ordinal to the corresponding subcell ordinal of the extruded subcell in the tensor product topology; that is if this = (shardsTopo x Line_2 x Line_2 ...) x Line_2, the mapping takes the subcell of dimension subcell_dim_in_component_topo and ordinal subcell_ord_in_component_topo in (shardsTopo x Line_2 x Line_2 ...) and returns the ordinal of that subcell extruded in the final Line_2 dimension.
 
ordinal_type getSideCount (const ordinal_type subcell_dim, const ordinal_type subcell_ord) const
 Side count of a subcell of the given dimension and ordinal. More...
 
ordinal_type getSubcellCount (const ordinal_type subcell_dim) const
 Subcell count of subcells of the given dimension. More...
 
ordinal_type getNodeFromTensorialComponentNodes (const std::vector< ordinal_type > &tensorComponentNodes) const
 Mapping from the tensorial component node ordinals to the node ordinal of this tensor cell topology. More...
 
ordinal_type getNodeMap (const ordinal_type subcell_dim, const ordinal_type subcell_ord, const ordinal_type subcell_node_ord) const
 Mapping from a subcell's node ordinal to a node ordinal of this parent cell topology. More...
 
ordinal_type getNodePermutationCount () const
 Number of node permutations defined for this cell.
 
ordinal_type getNodePermutation (const ordinal_type permutation_ord, const ordinal_type node_ord) const
 Permutation of a cell's node ordinals. More...
 
ordinal_type getNodePermutationInverse (const ordinal_type permutation_ord, const ordinal_type node_ord) const
 Inverse permutation of a cell's node ordinals. More...
 
CellTopoPtr getSide (ordinal_type sideOrdinal) const
 Returns a CellTopoPtr for the specified side.
 
CellTopoPtr getSubcell (ordinal_type scdim, ordinal_type scord) const
 Get the subcell of dimension scdim with ordinal scord. More...
 
CellTopoPtr getTensorialComponent () const
 For cell topologies of positive tensorial degree, returns the cell topology of tensorial degree one less. For cell topologies of tensorial degree zero, returns Teuchos::null.
 
ordinal_type getTensorialComponentSideOrdinal (ordinal_type extrusionNodeOrdinal)
 Returns the side corresponding to the provided node in the final extrusion dimension. More...
 
bool sideIsExtrudedInFinalDimension (ordinal_type sideOrdinal) const
 Returns true if the specified side has extension in the final tensorial dimension. For topologies with zero tensorialDegree_, returns false. More...
 

Static Public Member Functions

static ordinal_type getNodeCount (const shards::CellTopology &shardsTopo)
 
static ordinal_type getSubcellOrdinalMap (CellTopoPtr cellTopo, ordinal_type subcdim, ordinal_type subcord, ordinal_type subsubcdim, ordinal_type subsubcord)
 Maps the from a subcell within a subcell of the present CellTopology to the subcell in the present CellTopology; returns the ordinal of the subcell's subcell within the cell. More...
 
static CellTopoPtr cellTopology (const shards::CellTopology &shardsCellTopo, ordinal_type tensorialDegree=0)
 static accessor that returns a CellTopoPtr; these are lazily constructed and cached. More...
 
static CellTopoPtr point ()
 
static CellTopoPtr line ()
 
static CellTopoPtr quad ()
 
static CellTopoPtr hexahedron ()
 
static CellTopoPtr triangle ()
 
static CellTopoPtr tetrahedron ()
 
static CellTopoPtr wedge ()
 

Protected Attributes

shards::CellTopology shardsBaseTopology_
 
ordinal_type tensorialDegree_
 
std::string name_
 
std::vector< std::vector
< CellTopoPtr > > 
subcells_
 

Detailed Description

Implements arbitrary-dimensional extrusion of a base shards::CellTopology.

Definition at line 62 of file Intrepid2_CellTopology.hpp.

Constructor & Destructor Documentation

Intrepid2::CellTopology::CellTopology ( const shards::CellTopology &  baseTopo,
ordinal_type  tensorialDegree 
)
inline

Constructor.

Parameters
[in]baseTopo- the base shards CellTopology.
[in]tensorialDegree- the number of dimensions to extrude into.

The dimension of the CellTopology will be baseTopo.getDimension() + tensorialDegree.

Definition at line 81 of file Intrepid2_CellTopology.hpp.

References getSubcell(), and getTensorialComponent().

Referenced by cellTopology().

Member Function Documentation

static CellTopoPtr Intrepid2::CellTopology::cellTopology ( const shards::CellTopology &  shardsCellTopo,
ordinal_type  tensorialDegree = 0 
)
inlinestatic

static accessor that returns a CellTopoPtr; these are lazily constructed and cached.

Parameters
[in]baseTopo- the base shards CellTopology.
[in]tensorialDegree- the number of dimensions to extrude into.

Definition at line 588 of file Intrepid2_CellTopology.hpp.

References CellTopology().

Referenced by Intrepid2::Basis_TensorBasis< BasisBaseClass >::Basis_TensorBasis(), getNodeFromTensorialComponentNodes(), getNodeMap(), getSubcell(), getTensorialComponent(), and Intrepid2::SerendipityBasis< FullBasis::BasisBase >::SerendipityBasis().

ordinal_type Intrepid2::CellTopology::getEdgeCount ( const ordinal_type  subcell_dim,
const ordinal_type  subcell_ord 
) const
inline

Edge count of a subcell of the given dimension and ordinal.

Parameters
subcell_dim[in] - spatial dimension of the subcell
subcell_ord[in] - subcell ordinal

Definition at line 238 of file Intrepid2_CellTopology.hpp.

ordinal_type Intrepid2::CellTopology::getNodeCount ( const ordinal_type  subcell_dim,
const ordinal_type  subcell_ord 
) const
inline

Node count of a subcell of the given dimension and ordinal.

Parameters
subcell_dim[in] - spatial dimension of the subcell
subcell_ord[in] - subcell ordinal

Definition at line 211 of file Intrepid2_CellTopology.hpp.

ordinal_type Intrepid2::CellTopology::getNodeFromTensorialComponentNodes ( const std::vector< ordinal_type > &  tensorComponentNodes) const
inline

Mapping from the tensorial component node ordinals to the node ordinal of this tensor cell topology.

Parameters
tensorComponentNodes[in] - node ordinals in the tensorial components.

Definition at line 292 of file Intrepid2_CellTopology.hpp.

References cellTopology().

ordinal_type Intrepid2::CellTopology::getNodeMap ( const ordinal_type  subcell_dim,
const ordinal_type  subcell_ord,
const ordinal_type  subcell_node_ord 
) const
inline

Mapping from a subcell's node ordinal to a node ordinal of this parent cell topology.

Parameters
subcell_dim[in] - spatial dimension of the subcell
subcell_ord[in] - subcell ordinal
subcell_node_ord[in] - node ordinal relative to subcell

Definition at line 327 of file Intrepid2_CellTopology.hpp.

References cellTopology(), and getDimension().

ordinal_type Intrepid2::CellTopology::getNodePermutation ( const ordinal_type  permutation_ord,
const ordinal_type  node_ord 
) const

Permutation of a cell's node ordinals.

Parameters
permutation_ordinal[in]
node_ordinal[in]
ordinal_type Intrepid2::CellTopology::getNodePermutationInverse ( const ordinal_type  permutation_ord,
const ordinal_type  node_ord 
) const

Inverse permutation of a cell's node ordinals.

Parameters
permutation_ordinal[in]
node_ordinal[in]
ordinal_type Intrepid2::CellTopology::getSideCount ( const ordinal_type  subcell_dim,
const ordinal_type  subcell_ord 
) const
inline

Side count of a subcell of the given dimension and ordinal.

Parameters
subcell_dim[in] - spatial dimension of the subcell
subcell_ord[in] - subcell ordinal

Definition at line 272 of file Intrepid2_CellTopology.hpp.

CellTopoPtr Intrepid2::CellTopology::getSubcell ( ordinal_type  scdim,
ordinal_type  scord 
) const
inline

Get the subcell of dimension scdim with ordinal scord.

Parameters
scdim[in]
scord[in] For tensor-product topologies T x L (L being the line topology), there are two "copies" of T, T0 and T1, and the enumeration of subcells of dimension d goes as follows:
  • d-dimensional subcells from T0
  • d-dimensional subcells from T1
  • ((d-1)-dimensional subcells of T) x L.

Definition at line 410 of file Intrepid2_CellTopology.hpp.

References cellTopology(), and getTensorialComponent().

Referenced by CellTopology().

ordinal_type Intrepid2::CellTopology::getSubcellCount ( const ordinal_type  subcell_dim) const
inline

Subcell count of subcells of the given dimension.

Parameters
subcell_dim[in] - spatial dimension of the subcell

Definition at line 282 of file Intrepid2_CellTopology.hpp.

Referenced by getEdgeCount(), getFaceCount(), and getSideCount().

static ordinal_type Intrepid2::CellTopology::getSubcellOrdinalMap ( CellTopoPtr  cellTopo,
ordinal_type  subcdim,
ordinal_type  subcord,
ordinal_type  subsubcdim,
ordinal_type  subsubcord 
)
inlinestatic

Maps the from a subcell within a subcell of the present CellTopology to the subcell in the present CellTopology; returns the ordinal of the subcell's subcell within the cell.

Adapted from CamelliaCellTools::getSubcellOrdinalMap().

Definition at line 437 of file Intrepid2_CellTopology.hpp.

Referenced by Intrepid2::Basis_TensorBasis< BasisBaseClass >::Basis_TensorBasis().

ordinal_type Intrepid2::CellTopology::getTensorialComponentSideOrdinal ( ordinal_type  extrusionNodeOrdinal)
inline

Returns the side corresponding to the provided node in the final extrusion dimension.

Parameters
extrusionNodeOrdinal[in] - 0 or 1, the node number for the vertex in the extrusion line topology.

Definition at line 567 of file Intrepid2_CellTopology.hpp.

ordinal_type Intrepid2::CellTopology::getVertexCount ( const ordinal_type  subcell_dim,
const ordinal_type  subcell_ord 
) const
inline

Vertex count of a subcell of the given dimension and ordinal.

Parameters
subcell_dim[in] - spatial dimension of the subcell
subcell_ord[in] - subcell ordinal

Definition at line 228 of file Intrepid2_CellTopology.hpp.

bool Intrepid2::CellTopology::sideIsExtrudedInFinalDimension ( ordinal_type  sideOrdinal) const
inline

Returns true if the specified side has extension in the final tensorial dimension. For topologies with zero tensorialDegree_, returns false.

Parameters
sideOrdinal[in] Ordinal of the side.

Definition at line 576 of file Intrepid2_CellTopology.hpp.

References getSideCount().


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