shards  Version of the Day
 All Classes Functions Variables Typedefs Enumerations Enumerator Groups
Shards_CellTopology.hpp
1 // @HEADER
2 // *****************************************************************************
3 // Shards : Shared Discretization Tools
4 //
5 // Copyright 2008-2011 NTESS and the Shards contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef Shards_CellTopology_hpp
11 #define Shards_CellTopology_hpp
12 
13 #ifdef HAVE_SHARDS_DEBUG
14 #define SHARDS_REQUIRE( S ) S
15 #else
16 #define SHARDS_REQUIRE( S ) /* empty */
17 #endif
18 
19 #include <string>
20 #include <vector>
21 #include <Shards_CellTopologyData.h>
22 #include <Shards_BasicTopologies.hpp>
23 
24 namespace shards {
25 
30 /*------------------------------------------------------------------------*/
31 
32 class CellTopology ;
33 
35 //std::ostream & operator << ( std::ostream & , const CellTopologyData & );
36 
37 
39 std::ostream & operator << ( std::ostream & , const CellTopology & );
40 
41 
45 enum ECellType {
46  ALL_CELLS = 0,
47  STANDARD_CELL,
48  NONSTANDARD_CELL
49 };
50 
51 inline std::string ECellTypeToString(ECellType cellType) {
52  std::string retString;
53  switch(cellType){
54  case ALL_CELLS: retString = "All"; break;
55  case STANDARD_CELL: retString = "Standard"; break;
56  case NONSTANDARD_CELL: retString = "Nonstandard"; break;
57  default: retString = "Invalid Cell";
58  }
59  return retString;
60 }
61 
62 
67  ALL_TOPOLOGIES,
68  BASE_TOPOLOGY,
69  EXTENDED_TOPOLOGY
70 };
71 
72 inline std::string ETopologyTypeToString(ETopologyType topologyType) {
73  std::string retString;
74  switch(topologyType){
75  case ALL_TOPOLOGIES: retString = "All"; break;
76  case BASE_TOPOLOGY: retString = "Base"; break;
77  case EXTENDED_TOPOLOGY: retString = "Extended"; break;
78  default: retString = "Invalid Topology";
79  }
80  return retString;
81 }
82 
83 
92 void getTopologies(std::vector<shards::CellTopology>& topologies,
93  const unsigned cellDim = 4,
94  const ECellType cellType = ALL_CELLS,
95  const ETopologyType topologyType = ALL_TOPOLOGIES);
96 
97 
98 
105 int isPredefinedCell(const CellTopology & cell);
106 
107 
108 
109 /*------------------------------------------------------------------------*/
132 private:
133 
137  void requireCell() const ;
138 
139 
144  void requireDimension( const unsigned subcellDim ) const ;
145 
146 
152  void requireSubcell( const unsigned subcellDim ,
153  const unsigned subcellOrd ) const ;
154 
155 
162  void requireNodeMap( const unsigned subcellDim ,
163  const unsigned subcellOrd ,
164  const unsigned nodeOrd ) const ;
165 
166  void requireNodePermutation( const unsigned permutationOrd ,
167  const unsigned nodeOrd ) const ;
168 
169  const CellTopologyData * m_cell ;
170 
171 public:
172 
173  /*------------------------------------------------------------------*/
180  unsigned getDimension() const
181  {
182  SHARDS_REQUIRE( requireCell() );
183  return m_cell->dimension ;
184  }
185 
186 
192  unsigned getKey() const
193  {
194  SHARDS_REQUIRE( requireCell() );
195  return m_cell->key ;
196  }
197 
198 
204  unsigned getBaseKey() const
205  {
206  SHARDS_REQUIRE( requireCell() );
207  return m_cell->base->key ;
208  }
209 
210 
211 
217  const char* getName() const
218  {
219  SHARDS_REQUIRE( requireCell() );
220  return m_cell->name ;
221  }
222 
223 
227  const char* getBaseName() const
228  {
229  SHARDS_REQUIRE( requireCell() );
230  return m_cell->base->name ;
231  }
232 
233 
235  unsigned getNodeCount() const
236  {
237  SHARDS_REQUIRE( requireCell() );
238  return m_cell->node_count ;
239  }
240 
241 
243  unsigned getVertexCount() const
244  {
245  SHARDS_REQUIRE( requireCell() );
246  return m_cell->vertex_count ;
247  }
248 
249 
251  unsigned getEdgeCount() const
252  {
253  SHARDS_REQUIRE( requireCell() );
254  return m_cell->edge_count ;
255  }
256 
258  unsigned getFaceCount() const
259  {
260  SHARDS_REQUIRE( requireCell() );
261  return m_cell->dimension == 3 ? m_cell->side_count : 0 ;
262  }
263 
264 
266  unsigned getSideCount() const
267  {
268  SHARDS_REQUIRE( requireCell() );
269  return m_cell->side_count ;
270  }
271 
272 
274  bool isValid() const
275  { return m_cell != 0 ; }
276 
277 
280  { return m_cell ; }
281 
282 
285  {
286  SHARDS_REQUIRE( requireCell() );
287  return m_cell->base ;
288  }
289 
290 
296  const CellTopologyData * getCellTopologyData( const unsigned subcell_dim ,
297  const unsigned subcell_ord ) const
298  {
299  SHARDS_REQUIRE( requireCell() );
300  SHARDS_REQUIRE( requireDimension(subcell_dim) );
301  SHARDS_REQUIRE( requireSubcell(subcell_dim,subcell_ord) );
302  return m_cell->subcell[subcell_dim][subcell_ord].topology ;
303  }
304 
305 
311  const CellTopologyData * getBaseCellTopologyData( const unsigned subcell_dim ,
312  const unsigned subcell_ord ) const
313  {
314  return getCellTopologyData(subcell_dim,subcell_ord)->base ;
315  }
316 
317 
322  unsigned getKey( const unsigned subcell_dim ,
323  const unsigned subcell_ord ) const
324  {
325  return getCellTopologyData(subcell_dim,subcell_ord)->key ;
326  }
327 
328 
329 
334  const char * getName(const unsigned subcell_dim,
335  const unsigned subcell_ord) const
336  {
337  return getCellTopologyData(subcell_dim,subcell_ord) -> name;
338  }
339 
340 
345  unsigned getNodeCount( const unsigned subcell_dim ,
346  const unsigned subcell_ord ) const
347  {
348  return getCellTopologyData(subcell_dim,subcell_ord)->node_count ;
349  }
350 
351 
356  unsigned getVertexCount( const unsigned subcell_dim ,
357  const unsigned subcell_ord ) const
358  {
359  return getCellTopologyData(subcell_dim,subcell_ord)->vertex_count ;
360  }
361 
362 
367  unsigned getEdgeCount( const unsigned subcell_dim ,
368  const unsigned subcell_ord ) const
369  {
370  return getCellTopologyData(subcell_dim,subcell_ord)->edge_count ;
371  }
372 
373 
378  unsigned getSideCount( const unsigned subcell_dim ,
379  const unsigned subcell_ord ) const
380  {
381  return getCellTopologyData(subcell_dim,subcell_ord)->side_count ;
382  }
383 
384 
388  unsigned getSubcellCount( const unsigned subcell_dim ) const
389  {
390  SHARDS_REQUIRE( requireCell() );
391  SHARDS_REQUIRE( requireDimension(subcell_dim) );
392  return m_cell->subcell_count[subcell_dim] ;
393  }
394 
395 
400  bool getSubcellHomogeneity( const unsigned subcell_dim ) const
401  {
402  SHARDS_REQUIRE( requireCell() );
403  SHARDS_REQUIRE( requireDimension(subcell_dim) );
404  return 0 != m_cell->subcell_homogeneity[subcell_dim] ;
405  }
406 
407 
414  unsigned getNodeMap( const unsigned subcell_dim ,
415  const unsigned subcell_ord ,
416  const unsigned subcell_node_ord ) const
417  {
418  SHARDS_REQUIRE( requireCell() );
419  SHARDS_REQUIRE( requireDimension(subcell_dim) );
420  SHARDS_REQUIRE( requireSubcell(subcell_dim,subcell_ord) );
421  SHARDS_REQUIRE( requireNodeMap(subcell_dim,subcell_ord,subcell_node_ord));
422  return m_cell->subcell[subcell_dim][subcell_ord].node[subcell_node_ord];
423  }
424 
425 
427  unsigned getNodePermutationCount() const
428  {
429  SHARDS_REQUIRE(requireCell());
430  return m_cell->permutation_count ;
431  }
432 
437  unsigned getNodePermutation( const unsigned permutation_ord ,
438  const unsigned node_ord ) const
439  {
440  SHARDS_REQUIRE(requireCell());
441  SHARDS_REQUIRE(requireNodePermutation(permutation_ord,node_ord));
442  return m_cell->permutation[permutation_ord].node[node_ord];
443  }
444 
449  unsigned getNodePermutationPolarity( const unsigned permutation_ord ) const
450  {
451  SHARDS_REQUIRE(requireCell());
452  SHARDS_REQUIRE(requireNodePermutation(permutation_ord,0));
453  return m_cell->permutation[permutation_ord].polarity;
454  }
455 
460  unsigned getNodePermutationInverse( const unsigned permutation_ord ,
461  const unsigned node_ord ) const
462  {
463  SHARDS_REQUIRE(requireCell());
464  SHARDS_REQUIRE(requireNodePermutation(permutation_ord,node_ord));
465  return m_cell->permutation_inverse[permutation_ord].node[node_ord];
466  }
467 
470  /*------------------------------------------------------------------*/
485  : m_cell( cell )
486  {}
487 
488 
496  CellTopology( const std::string & name,
497  const unsigned nodeCount);
498 
499 
509  CellTopology( const std::string & name,
510  const unsigned vertex_count,
511  const unsigned node_count,
512  const std::vector< const CellTopologyData * > & edges ,
513  const std::vector< unsigned > & edge_node_map ,
514  const CellTopologyData * base = NULL );
515 
516 
528  CellTopology( const std::string & name,
529  const unsigned vertex_count,
530  const unsigned node_count,
531  const std::vector< const CellTopologyData * > & edges ,
532  const std::vector< unsigned > & edge_node_map ,
533  const std::vector< const CellTopologyData * > & faces ,
534  const std::vector< unsigned > & face_node_map ,
535  const CellTopologyData * base = NULL );
536 
537 
539  CellTopology& operator = (const CellTopology& right);
540 
542  CellTopology( const CellTopology& right );
543 
545  CellTopology();
546 
548  ~CellTopology();
549 
552 }; // class CellTopology
553 
554 /*------------------------------------------------------------------------*/
555 /* \brief Find the permutation from the expected nodes to the actual nodes,
556  *
557  * Find permutation 'p' such that:
558  * actual_node[j] == expected_node[ top.permutation[p].node[j] ]
559  * for all vertices.
560  */
561 template< typename id_type >
562 int findPermutation( const CellTopologyData & top ,
563  const id_type * const expected_node ,
564  const id_type * const actual_node )
565 {
566  const int nv = top.vertex_count ;
567  const int np = top.permutation_count ;
568  int p = 0 ;
569  for ( ; p < np ; ++p ) {
570  const unsigned * const perm_node = top.permutation[p].node ;
571  int j = 0 ;
572  for ( ; j < nv && actual_node[j] == expected_node[ perm_node[j] ] ; ++j );
573  if ( nv == j ) break ;
574  }
575  if ( np == p ) p = -1 ;
576  return p ;
577 }
578 
579 template< typename id_type >
580 int findPermutation( const CellTopology & top ,
581  const id_type * const expected_node ,
582  const id_type * const actual_node )
583 {
584  return findPermutation( * top.getCellTopologyData() , expected_node , actual_node );
585 }
586 
587 /*------------------------------------------------------------------------*/
597 void badCellTopologyKey( const unsigned dimension ,
598  const unsigned face_count ,
599  const unsigned edge_count ,
600  const unsigned vertex_count ,
601  const unsigned node_count );
602 
603 
612 inline
613 unsigned cellTopologyKey( const unsigned dimension ,
614  const unsigned face_count ,
615  const unsigned edge_count ,
616  const unsigned vertex_count ,
617  const unsigned node_count )
618 {
619  const bool bad = ( dimension >> 3 ) ||
620  ( face_count >> 6 ) ||
621  ( edge_count >> 6 ) ||
622  ( vertex_count >> 6 ) ||
623  ( node_count >> 10 );
624 
625  if ( bad ) {
626  badCellTopologyKey( dimension ,
627  face_count ,
628  edge_count ,
629  vertex_count ,
630  node_count );
631  }
632 
633  const unsigned key = ( dimension << 28 ) |
634  ( face_count << 22 ) |
635  ( edge_count << 16 ) |
636  ( vertex_count << 10 ) |
637  ( node_count ) ;
638 
639  return key ;
640 }
641 
642 inline
643 bool operator==(const CellTopology &left, const CellTopology &right)
644 {
645  return left.getCellTopologyData() == right.getCellTopologyData();
646 }
647 
648 inline
649 bool operator<(const CellTopology &left, const CellTopology &right)
650 {
651  return left.getCellTopologyData() < right.getCellTopologyData();
652 // FIXME: Should is be this?
653 // return left.getKey() < right.getKey();
654 }
655 
656 inline
657 bool operator!=(const CellTopology &left, const CellTopology &right) {
658  return !(left == right);
659 }
660 
661 
664 } // namespace shards
665 
666 #undef SHARDS_REQUIRE
667 
668 #endif // Shards_CellTopology_hpp
669 
~CellTopology()
Destructor.
unsigned getVertexCount(const unsigned subcell_dim, const unsigned subcell_ord) const
Vertex count of a subcell of the given dimension and ordinal.
unsigned getNodeMap(const unsigned subcell_dim, const unsigned subcell_ord, const unsigned subcell_node_ord) const
Mapping from a subcell&#39;s node ordinal to a node ordinal of this parent cell topology.
unsigned getNodePermutation(const unsigned permutation_ord, const unsigned node_ord) const
Permutation of a cell&#39;s node ordinals.
ETopologyType
Enumeration of topology types in Shards.
unsigned getNodeCount() const
Node count of this cell topology.
struct CellTopologyData * base
Base, a.k.a. not-extended, version of this topology where vertex_count == node_count.
unsigned subcell_homogeneity[4]
Flag if the subcells of a given dimension are homogeneous.
unsigned key
Unique key for this topology.
unsigned getKey() const
Unique key for this cell topology; under certain subcell uniformity conditions.
void 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...
unsigned getNodePermutationPolarity(const unsigned permutation_ord) const
Permutation of a cell&#39;s node ordinals.
unsigned getNodeCount(const unsigned subcell_dim, const unsigned subcell_ord) const
Node count of a subcell of the given dimension and ordinal.
const CellTopologyData * getCellTopologyData() const
This cell&#39;s raw topology data.
unsigned vertex_count
Number of vertices.
const CellTopologyData * getBaseCellTopologyData(const unsigned subcell_dim, const unsigned subcell_ord) const
Raw cell topology data for the base topology of a subcell of the given dimension and ordinal...
ECellType
Enumeration of cell types in Shards.
CellTopology()
Default constructor initializes to NULL.
const char * getName() const
Unique name for this cell topology;.
unsigned subcell_count[4]
Number of subcells of each dimension.
unsigned edge_count
Number of edges (a.k.a. boundary subcells).
unsigned getEdgeCount(const unsigned subcell_dim, const unsigned subcell_ord) const
Edge count of a subcell of the given dimension and ordinal.
unsigned dimension
Topological dimension.
const CellTopologyData * getBaseCellTopologyData() const
This cell&#39;s base cell topology&#39;s raw topology data.
bool getSubcellHomogeneity(const unsigned subcell_dim) const
Query if all subcells of the given dimension have the same cell topology.
unsigned node_count
Number of nodes (a.k.a. subcells).
unsigned 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.
const unsigned * node
Subcell indexing of with respect to parent cell.
const char * name
Intuitive name for this topology.
struct CellTopologyData * topology
Subcell topology.
unsigned side_count
Number of sides (a.k.a. boundary subcells).
unsigned getSideCount(const unsigned subcell_dim, const unsigned subcell_ord) const
Side count of a subcell of the given dimension and ordinal.
unsigned getNodePermutationCount() const
Number of node permutations defined for this cell.
unsigned getKey(const unsigned subcell_dim, const unsigned subcell_ord) const
Key of a subcell of the given dimension and ordinal.
struct CellTopologyData_Subcell * subcell[4]
Array of subcells of each dimension.
unsigned getBaseKey() const
Unique key for this cell&#39;s base topology; under certain subcell uniformity conditions.
int isPredefinedCell(const CellTopology &cell)
Checks if the cell topology is predefined in shards.
unsigned getDimension() const
Dimension of this cell topology.
const CellTopologyData * getCellTopologyData(const unsigned subcell_dim, const unsigned subcell_ord) const
Raw cell topology data for a subcell of the given dimension and ordinal.
A simple &#39;C&#39; struct of cell topology attributes.
unsigned getNodePermutationInverse(const unsigned permutation_ord, const unsigned node_ord) const
Inverse permutation of a cell&#39;s node ordinals.
struct CellTopologyData_Permutation * permutation
Array of node permutations.
unsigned getFaceCount() const
Face boundary subcell count of this cell topology.
const char * getBaseName() const
Unique name for this cell&#39;s base topology.
std::ostream & operator<<(std::ostream &, const CellTopology &)
Overloaded &lt;&lt; operator for CellTopologyData objects.
CellTopology(const CellTopologyData *cell)
Wrapper for safe access to a raw cell topology data.
Provide input checked access (in debug mode) to cell topology data and a procedure to create custom c...
unsigned getEdgeCount() const
Edge boundary subcell count of this cell topology.
const char * getName(const unsigned subcell_dim, const unsigned subcell_ord) const
Name of a subcell of the given dimension and ordinal.
unsigned getSideCount() const
Side boundary subcell count of this cell topology.
CellTopology & operator=(const CellTopology &right)
Assignment operator *this = right.
bool isValid() const
This cell&#39;s raw topology data.
unsigned getSubcellCount(const unsigned subcell_dim) const
Subcell count of subcells of the given dimension.
unsigned permutation_count
Number of defined permutations.
unsigned getVertexCount() const
Vertex count of this cell topology.
void 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...