Intrepid2
|
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_ |
Implements arbitrary-dimensional extrusion of a base shards::CellTopology.
Definition at line 28 of file Intrepid2_CellTopology.hpp.
|
inline |
Constructor.
[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 47 of file Intrepid2_CellTopology.hpp.
References getSubcell(), and getTensorialComponent().
Referenced by cellTopology().
|
inlinestatic |
static accessor that returns a CellTopoPtr; these are lazily constructed and cached.
[in] | baseTopo | - the base shards CellTopology. |
[in] | tensorialDegree | - the number of dimensions to extrude into. |
Definition at line 554 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().
|
inline |
Edge count of a subcell of the given dimension and ordinal.
subcell_dim | [in] - spatial dimension of the subcell |
subcell_ord | [in] - subcell ordinal |
Definition at line 204 of file Intrepid2_CellTopology.hpp.
|
inline |
Node count of a subcell of the given dimension and ordinal.
subcell_dim | [in] - spatial dimension of the subcell |
subcell_ord | [in] - subcell ordinal |
Definition at line 177 of file Intrepid2_CellTopology.hpp.
|
inline |
Mapping from the tensorial component node ordinals to the node ordinal of this tensor cell topology.
tensorComponentNodes | [in] - node ordinals in the tensorial components. |
Definition at line 258 of file Intrepid2_CellTopology.hpp.
References cellTopology().
|
inline |
Mapping from a subcell's node ordinal to a node ordinal of this parent cell topology.
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 293 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.
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.
permutation_ordinal | [in] |
node_ordinal | [in] |
|
inline |
Side count of a subcell of the given dimension and ordinal.
subcell_dim | [in] - spatial dimension of the subcell |
subcell_ord | [in] - subcell ordinal |
Definition at line 238 of file Intrepid2_CellTopology.hpp.
|
inline |
Get the subcell of dimension scdim with ordinal scord.
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:
|
Definition at line 376 of file Intrepid2_CellTopology.hpp.
References cellTopology(), and getTensorialComponent().
Referenced by CellTopology().
|
inline |
Subcell count of subcells of the given dimension.
subcell_dim | [in] - spatial dimension of the subcell |
Definition at line 248 of file Intrepid2_CellTopology.hpp.
Referenced by getEdgeCount(), getFaceCount(), and getSideCount().
|
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 403 of file Intrepid2_CellTopology.hpp.
Referenced by Intrepid2::Basis_TensorBasis< BasisBaseClass >::Basis_TensorBasis().
|
inline |
Returns the side corresponding to the provided node in the final extrusion dimension.
extrusionNodeOrdinal | [in] - 0 or 1, the node number for the vertex in the extrusion line topology. |
Definition at line 533 of file Intrepid2_CellTopology.hpp.
|
inline |
Vertex count of a subcell of the given dimension and ordinal.
subcell_dim | [in] - spatial dimension of the subcell |
subcell_ord | [in] - subcell ordinal |
Definition at line 194 of file Intrepid2_CellTopology.hpp.
|
inline |
Returns true if the specified side has extension in the final tensorial dimension. For topologies with zero tensorialDegree_, returns false.
sideOrdinal | [in] Ordinal of the side. |
Definition at line 542 of file Intrepid2_CellTopology.hpp.
References getSideCount().